Actual source code: ex22.c
slepc-main 2024-12-17
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: PetscFunctionBeginUser;
48: PetscCall(SlepcInitialize(&argc,&argv,NULL,help));
49: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
50: PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
51: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau));
52: h = PETSC_PI/(PetscReal)(n+1);
54: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55: Create nonlinear eigensolver context
56: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58: PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
60: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
61: Create problem matrices and coefficient functions. Pass them to NEP
62: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
64: /*
65: Identity matrix
66: */
67: PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id));
68: PetscCall(MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE));
70: /*
71: A = 1/h^2*tridiag(1,-2,1) + a*I
72: */
73: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
74: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
75: PetscCall(MatSetFromOptions(A));
76: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
77: for (i=Istart;i<Iend;i++) {
78: if (i>0) PetscCall(MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES));
79: if (i<n-1) PetscCall(MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES));
80: PetscCall(MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES));
81: }
82: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
83: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
84: PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));
86: /*
87: B = diag(b(xi))
88: */
89: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
90: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
91: PetscCall(MatSetFromOptions(B));
92: PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
93: for (i=Istart;i<Iend;i++) {
94: xi = (i+1)*h;
95: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
96: PetscCall(MatSetValue(B,i,i,b,INSERT_VALUES));
97: }
98: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
99: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
100: PetscCall(MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE));
102: /*
103: Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda)
104: */
105: PetscCall(FNCreate(PETSC_COMM_WORLD,&f1));
106: PetscCall(FNSetType(f1,FNRATIONAL));
107: coeffs[0] = -1.0; coeffs[1] = 0.0;
108: PetscCall(FNRationalSetNumerator(f1,2,coeffs));
110: PetscCall(FNCreate(PETSC_COMM_WORLD,&f2));
111: PetscCall(FNSetType(f2,FNRATIONAL));
112: coeffs[0] = 1.0;
113: PetscCall(FNRationalSetNumerator(f2,1,coeffs));
115: PetscCall(FNCreate(PETSC_COMM_WORLD,&f3));
116: PetscCall(FNSetType(f3,FNEXP));
117: PetscCall(FNSetScale(f3,-tau,1.0));
119: /*
120: Set the split operator. Note that A is passed first so that
121: SUBSET_NONZERO_PATTERN can be used
122: */
123: mats[0] = A; funs[0] = f2;
124: mats[1] = Id; funs[1] = f1;
125: mats[2] = B; funs[2] = f3;
126: PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN));
128: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129: Customize nonlinear solver; set runtime options
130: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132: PetscCall(NEPSetTolerances(nep,1e-9,PETSC_CURRENT));
133: PetscCall(NEPSetDimensions(nep,1,PETSC_DETERMINE,PETSC_DETERMINE));
134: PetscCall(NEPRIISetLagPreconditioner(nep,0));
136: /*
137: Set solver parameters at runtime
138: */
139: PetscCall(NEPSetFromOptions(nep));
141: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142: Solve the eigensystem
143: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
145: PetscCall(NEPSolve(nep));
146: PetscCall(NEPGetDimensions(nep,&nev,NULL,NULL));
147: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
149: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: Display solution and clean up
151: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153: /* show detailed info unless -terse option is given by user */
154: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
155: if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
156: else {
157: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
158: PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD));
159: PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
160: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
161: }
162: PetscCall(NEPDestroy(&nep));
163: PetscCall(MatDestroy(&Id));
164: PetscCall(MatDestroy(&A));
165: PetscCall(MatDestroy(&B));
166: PetscCall(FNDestroy(&f1));
167: PetscCall(FNDestroy(&f2));
168: PetscCall(FNDestroy(&f3));
169: PetscCall(SlepcFinalize());
170: return 0;
171: }
173: /*TEST
175: testset:
176: suffix: 1
177: args: -nep_type {{rii slp narnoldi}} -terse
178: filter: sed -e "s/[+-]0\.0*i//g"
179: requires: !single
181: test:
182: suffix: 1_ciss
183: 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
184: requires: complex !single
186: test:
187: suffix: 2
188: 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
189: filter: sed -e "s/[+-]0\.0*i//g"
190: requires: !single
192: testset:
193: args: -n 512 -nep_target 10 -nep_nev 3 -terse
194: filter: sed -e "s/[+-]0\.0*i//g"
195: requires: !single
196: output_file: output/ex22_3.out
197: test:
198: suffix: 3
199: args: -nep_type {{rii slp narnoldi}}
200: test:
201: suffix: 3_simpleu
202: args: -nep_type {{rii slp narnoldi}} -nep_deflation_simpleu
203: test:
204: suffix: 3_slp_thres
205: args: -nep_type slp -nep_slp_deflation_threshold 1e-8
206: test:
207: suffix: 3_rii_thres
208: args: -nep_type rii -nep_rii_deflation_threshold 1e-8
210: test:
211: suffix: 4
212: 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
213: filter: sed -e "s/[+-]0\.0*i//g"
214: requires: x !single
215: output_file: output/ex22_2.out
217: TEST*/