Actual source code: pdde_stability.c
slepc-3.19.0 2023-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: This example implements one of the problems found at
12: NLEVP: A Collection of Nonlinear Eigenvalue Problems,
13: The University of Manchester.
14: The details of the collection can be found at:
15: [1] T. Betcke et al., "NLEVP: A Collection of Nonlinear Eigenvalue
16: Problems", ACM Trans. Math. Software 39(2), Article 7, 2013.
18: The pdde_stability problem is a complex-symmetric QEP from the stability
19: analysis of a discretized partial delay-differential equation. It requires
20: complex scalars.
21: */
23: static char help[] = "Stability analysis of a discretized partial delay-differential equation.\n\n"
24: "The command line options are:\n"
25: " -m <m>, grid size, the matrices have dimension n=m*m.\n"
26: " -c <a0,b0,a1,b1,a2,b2,phi1>, comma-separated list of 7 real parameters.\n\n";
28: #include <slepcpep.h>
30: #define NMAT 3
32: /*
33: Function for user-defined eigenvalue ordering criterion.
35: Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
36: one of them as the preferred one according to the criterion.
37: In this example, the preferred value is the one with absolute value closest to 1.
38: */
39: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
40: {
41: PetscReal aa,ab;
43: PetscFunctionBeginUser;
44: aa = PetscAbsReal(SlepcAbsEigenvalue(ar,ai)-PetscRealConstant(1.0));
45: ab = PetscAbsReal(SlepcAbsEigenvalue(br,bi)-PetscRealConstant(1.0));
46: *r = aa > ab ? 1 : (aa < ab ? -1 : 0);
47: PetscFunctionReturn(PETSC_SUCCESS);
48: }
50: int main(int argc,char **argv)
51: {
52: Mat A[NMAT]; /* problem matrices */
53: PEP pep; /* polynomial eigenproblem solver context */
54: PetscInt m=15,n,II,Istart,Iend,i,j,k;
55: PetscReal h,xi,xj,c[7] = { 2, .3, -2, .2, -2, -.3, -PETSC_PI/2 };
56: PetscScalar alpha,beta,gamma;
57: PetscBool flg,terse;
59: PetscFunctionBeginUser;
60: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
61: #if !defined(PETSC_USE_COMPLEX)
62: SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"This example requires complex scalars");
63: #endif
65: PetscCall(PetscOptionsGetInt(NULL,NULL,"-m",&m,NULL));
66: n = m*m;
67: h = PETSC_PI/(m+1);
68: gamma = PetscExpScalar(PETSC_i*c[6]);
69: gamma = gamma/PetscAbsScalar(gamma);
70: k = 7;
71: PetscCall(PetscOptionsGetRealArray(NULL,NULL,"-c",c,&k,&flg));
72: PetscCheck(!flg || k==7,PETSC_COMM_WORLD,PETSC_ERR_USER,"The number of parameters -c should be 7, you provided %" PetscInt_FMT,k);
73: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\nPDDE stability, n=%" PetscInt_FMT " (m=%" PetscInt_FMT ")\n\n",n,m));
75: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
76: Compute the polynomial matrices
77: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
79: /* initialize matrices */
80: for (i=0;i<NMAT;i++) {
81: PetscCall(MatCreate(PETSC_COMM_WORLD,&A[i]));
82: PetscCall(MatSetSizes(A[i],PETSC_DECIDE,PETSC_DECIDE,n,n));
83: PetscCall(MatSetFromOptions(A[i]));
84: PetscCall(MatSetUp(A[i]));
85: }
86: PetscCall(MatGetOwnershipRange(A[0],&Istart,&Iend));
88: /* A[1] has a pattern similar to the 2D Laplacian */
89: for (II=Istart;II<Iend;II++) {
90: i = II/m; j = II-i*m;
91: xi = (i+1)*h; xj = (j+1)*h;
92: alpha = c[0]+c[1]*PetscSinReal(xi)+gamma*(c[2]+c[3]*xi*(1.0-PetscExpReal(xi-PETSC_PI)));
93: beta = c[0]+c[1]*PetscSinReal(xj)-gamma*(c[2]+c[3]*xj*(1.0-PetscExpReal(xj-PETSC_PI)));
94: PetscCall(MatSetValue(A[1],II,II,alpha+beta-4.0/(h*h),INSERT_VALUES));
95: if (j>0) PetscCall(MatSetValue(A[1],II,II-1,1.0/(h*h),INSERT_VALUES));
96: if (j<m-1) PetscCall(MatSetValue(A[1],II,II+1,1.0/(h*h),INSERT_VALUES));
97: if (i>0) PetscCall(MatSetValue(A[1],II,II-m,1.0/(h*h),INSERT_VALUES));
98: if (i<m-1) PetscCall(MatSetValue(A[1],II,II+m,1.0/(h*h),INSERT_VALUES));
99: }
101: /* A[0] and A[2] are diagonal */
102: for (II=Istart;II<Iend;II++) {
103: i = II/m; j = II-i*m;
104: xi = (i+1)*h; xj = (j+1)*h;
105: alpha = c[4]+c[5]*xi*(PETSC_PI-xi);
106: beta = c[4]+c[5]*xj*(PETSC_PI-xj);
107: PetscCall(MatSetValue(A[0],II,II,alpha,INSERT_VALUES));
108: PetscCall(MatSetValue(A[2],II,II,beta,INSERT_VALUES));
109: }
111: /* assemble matrices */
112: for (i=0;i<NMAT;i++) PetscCall(MatAssemblyBegin(A[i],MAT_FINAL_ASSEMBLY));
113: for (i=0;i<NMAT;i++) PetscCall(MatAssemblyEnd(A[i],MAT_FINAL_ASSEMBLY));
115: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116: Create the eigensolver and solve the problem
117: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
120: PetscCall(PEPSetOperators(pep,NMAT,A));
121: PetscCall(PEPSetEigenvalueComparison(pep,MyEigenSort,NULL));
122: PetscCall(PEPSetDimensions(pep,4,PETSC_DEFAULT,PETSC_DEFAULT));
123: PetscCall(PEPSetFromOptions(pep));
124: PetscCall(PEPSolve(pep));
126: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127: Display solution and clean up
128: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
130: /* show detailed info unless -terse option is given by user */
131: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
132: if (terse) PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,NULL));
133: else {
134: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
135: PetscCall(PEPConvergedReasonView(pep,PETSC_VIEWER_STDOUT_WORLD));
136: PetscCall(PEPErrorView(pep,PEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD));
137: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
138: }
139: PetscCall(PEPDestroy(&pep));
140: for (i=0;i<NMAT;i++) PetscCall(MatDestroy(&A[i]));
141: PetscCall(SlepcFinalize());
142: return 0;
143: }
145: /*TEST
147: build:
148: requires: complex
150: test:
151: suffix: 1
152: args: -pep_type {{toar qarnoldi linear}} -pep_ncv 25 -terse
153: requires: complex double
155: TEST*/