Actual source code: ex22.c
slepc-3.17.0 2022-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: */
11: static char help[] = "Delay differential equation.\n\n"
12: "The command line options are:\n"
13: " -n <n>, where <n> = number of grid subdivisions.\n"
14: " -tau <tau>, where <tau> is the delay parameter.\n\n";
16: /*
17: Solve parabolic partial differential equation with time delay tau
19: u_t = u_xx + a*u(t) + b*u(t-tau)
20: u(0,t) = u(pi,t) = 0
22: with a = 20 and b(x) = -4.1+x*(1-exp(x-pi)).
24: Discretization leads to a DDE of dimension n
26: -u' = A*u(t) + B*u(t-tau)
28: which results in the nonlinear eigenproblem
30: (-lambda*I + A + exp(-tau*lambda)*B)*u = 0
31: */
33: #include <slepcnep.h>
35: int main(int argc,char **argv)
36: {
37: NEP nep; /* nonlinear eigensolver context */
38: Mat Id,A,B; /* problem matrices */
39: FN f1,f2,f3; /* functions to define the nonlinear operator */
40: Mat mats[3];
41: FN funs[3];
42: PetscScalar coeffs[2],b;
43: PetscInt n=128,nev,Istart,Iend,i;
44: PetscReal tau=0.001,h,a=20,xi;
45: PetscBool terse;
47: SlepcInitialize(&argc,&argv,(char*)0,help);
48: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
49: PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
50: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau);
51: h = PETSC_PI/(PetscReal)(n+1);
53: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
54: Create nonlinear eigensolver context
55: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
57: NEPCreate(PETSC_COMM_WORLD,&nep);
59: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
60: Create problem matrices and coefficient functions. Pass them to NEP
61: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
63: /*
64: Identity matrix
65: */
66: MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id);
67: MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE);
69: /*
70: A = 1/h^2*tridiag(1,-2,1) + a*I
71: */
72: MatCreate(PETSC_COMM_WORLD,&A);
73: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
74: MatSetFromOptions(A);
75: MatSetUp(A);
76: MatGetOwnershipRange(A,&Istart,&Iend);
77: for (i=Istart;i<Iend;i++) {
78: if (i>0) MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES);
79: if (i<n-1) MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES);
80: MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES);
81: }
82: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
83: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
84: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
86: /*
87: B = diag(b(xi))
88: */
89: MatCreate(PETSC_COMM_WORLD,&B);
90: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
91: MatSetFromOptions(B);
92: MatSetUp(B);
93: MatGetOwnershipRange(B,&Istart,&Iend);
94: for (i=Istart;i<Iend;i++) {
95: xi = (i+1)*h;
96: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
97: MatSetValue(B,i,i,b,INSERT_VALUES);
98: }
99: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
100: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
101: MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);
103: /*
104: Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda)
105: */
106: FNCreate(PETSC_COMM_WORLD,&f1);
107: FNSetType(f1,FNRATIONAL);
108: coeffs[0] = -1.0; coeffs[1] = 0.0;
109: FNRationalSetNumerator(f1,2,coeffs);
111: FNCreate(PETSC_COMM_WORLD,&f2);
112: FNSetType(f2,FNRATIONAL);
113: coeffs[0] = 1.0;
114: FNRationalSetNumerator(f2,1,coeffs);
116: FNCreate(PETSC_COMM_WORLD,&f3);
117: FNSetType(f3,FNEXP);
118: FNSetScale(f3,-tau,1.0);
120: /*
121: Set the split operator. Note that A is passed first so that
122: SUBSET_NONZERO_PATTERN can be used
123: */
124: mats[0] = A; funs[0] = f2;
125: mats[1] = Id; funs[1] = f1;
126: mats[2] = B; funs[2] = f3;
127: NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);
129: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130: Customize nonlinear solver; set runtime options
131: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133: NEPSetTolerances(nep,1e-9,PETSC_DEFAULT);
134: NEPSetDimensions(nep,1,PETSC_DEFAULT,PETSC_DEFAULT);
135: NEPRIISetLagPreconditioner(nep,0);
137: /*
138: Set solver parameters at runtime
139: */
140: NEPSetFromOptions(nep);
142: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
143: Solve the eigensystem
144: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
146: NEPSolve(nep);
147: NEPGetDimensions(nep,&nev,NULL,NULL);
148: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
150: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
151: Display solution and clean up
152: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
154: /* show detailed info unless -terse option is given by user */
155: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
156: if (terse) NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
157: else {
158: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
159: NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
160: NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
161: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
162: }
163: NEPDestroy(&nep);
164: MatDestroy(&Id);
165: MatDestroy(&A);
166: MatDestroy(&B);
167: FNDestroy(&f1);
168: FNDestroy(&f2);
169: FNDestroy(&f3);
170: SlepcFinalize();
171: return 0;
172: }
174: /*TEST
176: testset:
177: suffix: 1
178: args: -nep_type {{rii slp narnoldi}} -terse
179: filter: sed -e "s/[+-]0\.0*i//g"
180: requires: !single
182: test:
183: suffix: 1_ciss
184: args: -nep_type ciss -nep_ciss_extraction {{ritz hankel caa}} -rg_type ellipse -rg_ellipse_center 10 -rg_ellipse_radius 9.5 -nep_ncv 24 -terse
185: requires: complex !single
187: test:
188: suffix: 2
189: args: -nep_type interpol -nep_interpol_pep_extract {{none norm residual}} -rg_type interval -rg_interval_endpoints 5,20,-.1,.1 -nep_nev 3 -nep_target 5 -terse
190: filter: sed -e "s/[+-]0\.0*i//g"
191: requires: !single
193: testset:
194: args: -n 512 -nep_target 10 -nep_nev 3 -terse
195: filter: sed -e "s/[+-]0\.0*i//g"
196: requires: !single
197: output_file: output/ex22_3.out
198: test:
199: suffix: 3
200: args: -nep_type {{rii slp narnoldi}}
201: test:
202: suffix: 3_simpleu
203: args: -nep_type {{rii slp narnoldi}} -nep_deflation_simpleu
204: test:
205: suffix: 3_slp_thres
206: args: -nep_type slp -nep_slp_deflation_threshold 1e-8
207: test:
208: suffix: 3_rii_thres
209: args: -nep_type rii -nep_rii_deflation_threshold 1e-8
211: test:
212: suffix: 4
213: args: -nep_type interpol -rg_type interval -rg_interval_endpoints 5,20,-.1,.1 -nep_nev 3 -nep_target 5 -terse -nep_monitor draw::draw_lg
214: filter: sed -e "s/[+-]0\.0*i//g"
215: requires: x !single
216: output_file: output/ex22_2.out
218: TEST*/