Actual source code: test5.c
slepc-3.21.1 2024-04-26
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[] = "Test the INTERPOL solver with a user-provided PEP.\n\n"
12: "This is based on ex22.\n"
13: "The command line options are:\n"
14: " -n <n>, where <n> = number of grid subdivisions.\n"
15: " -tau <tau>, where <tau> is the delay parameter.\n\n";
17: /*
18: Solve parabolic partial differential equation with time delay tau
20: u_t = u_xx + a*u(t) + b*u(t-tau)
21: u(0,t) = u(pi,t) = 0
23: with a = 20 and b(x) = -4.1+x*(1-exp(x-pi)).
25: Discretization leads to a DDE of dimension n
27: -u' = A*u(t) + B*u(t-tau)
29: which results in the nonlinear eigenproblem
31: (-lambda*I + A + exp(-tau*lambda)*B)*u = 0
32: */
34: #include <slepcnep.h>
36: int main(int argc,char **argv)
37: {
38: NEP nep;
39: PEP pep;
40: Mat Id,A,B;
41: FN f1,f2,f3;
42: RG rg;
43: Mat mats[3];
44: FN funs[3];
45: PetscScalar coeffs[2],b;
46: PetscInt n=128,nev,Istart,Iend,i,deg;
47: PetscReal tau=0.001,h,a=20,xi,tol;
48: PetscBool terse;
50: PetscFunctionBeginUser;
51: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
52: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
53: PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
54: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau));
55: h = PETSC_PI/(PetscReal)(n+1);
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Create a standalone PEP and RG objects with appropriate settings
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: PetscCall(PEPCreate(PETSC_COMM_WORLD,&pep));
62: PetscCall(PEPSetType(pep,PEPTOAR));
63: PetscCall(PEPSetFromOptions(pep));
65: PetscCall(RGCreate(PETSC_COMM_WORLD,&rg));
66: PetscCall(RGSetType(rg,RGINTERVAL));
67: PetscCall(RGIntervalSetEndpoints(rg,5,20,-0.1,0.1));
68: PetscCall(RGSetFromOptions(rg));
70: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
71: Create nonlinear eigensolver context
72: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
74: PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
76: /* Identity matrix */
77: PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id));
78: PetscCall(MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE));
80: /* A = 1/h^2*tridiag(1,-2,1) + a*I */
81: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
82: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
83: PetscCall(MatSetFromOptions(A));
84: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
85: for (i=Istart;i<Iend;i++) {
86: if (i>0) PetscCall(MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES));
87: if (i<n-1) PetscCall(MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES));
88: PetscCall(MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES));
89: }
90: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
91: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
92: PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));
94: /* B = diag(b(xi)) */
95: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
96: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
97: PetscCall(MatSetFromOptions(B));
98: PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
99: for (i=Istart;i<Iend;i++) {
100: xi = (i+1)*h;
101: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
102: PetscCall(MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES));
103: }
104: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
105: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
106: PetscCall(MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE));
108: /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */
109: PetscCall(FNCreate(PETSC_COMM_WORLD,&f1));
110: PetscCall(FNSetType(f1,FNRATIONAL));
111: coeffs[0] = -1.0; coeffs[1] = 0.0;
112: PetscCall(FNRationalSetNumerator(f1,2,coeffs));
114: PetscCall(FNCreate(PETSC_COMM_WORLD,&f2));
115: PetscCall(FNSetType(f2,FNRATIONAL));
116: coeffs[0] = 1.0;
117: PetscCall(FNRationalSetNumerator(f2,1,coeffs));
119: PetscCall(FNCreate(PETSC_COMM_WORLD,&f3));
120: PetscCall(FNSetType(f3,FNEXP));
121: PetscCall(FNSetScale(f3,-tau,1.0));
123: /* Set the split operator */
124: mats[0] = A; funs[0] = f2;
125: mats[1] = Id; funs[1] = f1;
126: mats[2] = B; funs[2] = f3;
127: PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN));
129: /* Customize nonlinear solver; set runtime options */
130: PetscCall(NEPSetType(nep,NEPINTERPOL));
131: PetscCall(NEPSetRG(nep,rg));
132: PetscCall(NEPInterpolSetPEP(nep,pep));
133: PetscCall(NEPInterpolGetInterpolation(nep,&tol,°));
134: PetscCall(NEPInterpolSetInterpolation(nep,tol,deg+2)); /* increase degree of interpolation */
135: PetscCall(NEPSetFromOptions(nep));
137: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
138: Solve the eigensystem
139: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141: PetscCall(NEPSolve(nep));
142: PetscCall(NEPGetDimensions(nep,&nev,NULL,NULL));
143: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev));
145: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
146: Display solution and clean up
147: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149: /* show detailed info unless -terse option is given by user */
150: PetscCall(PetscOptionsHasName(NULL,NULL,"-terse",&terse));
151: if (terse) PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
152: else {
153: PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL));
154: PetscCall(NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD));
155: PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD));
156: PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
157: }
158: PetscCall(NEPDestroy(&nep));
159: PetscCall(PEPDestroy(&pep));
160: PetscCall(RGDestroy(&rg));
161: PetscCall(MatDestroy(&Id));
162: PetscCall(MatDestroy(&A));
163: PetscCall(MatDestroy(&B));
164: PetscCall(FNDestroy(&f1));
165: PetscCall(FNDestroy(&f2));
166: PetscCall(FNDestroy(&f3));
167: PetscCall(SlepcFinalize());
168: return 0;
169: }
171: /*TEST
173: testset:
174: args: -nep_nev 3 -nep_target 5 -terse
175: output_file: output/test5_1.out
176: filter: sed -e "s/[+-]0\.0*i//g"
177: requires: !single
178: test:
179: suffix: 1
180: test:
181: suffix: 2_cuda
182: args: -mat_type aijcusparse
183: requires: cuda
184: test:
185: suffix: 3
186: args: -nep_view_values draw
187: requires: x
189: TEST*/