Actual source code: test5.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[] = "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: SlepcInitialize(&argc,&argv,(char*)0,help);
51: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
52: PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL);
53: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau);
54: h = PETSC_PI/(PetscReal)(n+1);
56: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: Create a standalone PEP and RG objects with appropriate settings
58: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60: PEPCreate(PETSC_COMM_WORLD,&pep);
61: PEPSetType(pep,PEPTOAR);
62: PEPSetFromOptions(pep);
64: RGCreate(PETSC_COMM_WORLD,&rg);
65: RGSetType(rg,RGINTERVAL);
66: RGIntervalSetEndpoints(rg,5,20,-0.1,0.1);
67: RGSetFromOptions(rg);
69: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
70: Create nonlinear eigensolver context
71: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73: NEPCreate(PETSC_COMM_WORLD,&nep);
75: /* Identity matrix */
76: MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id);
77: MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE);
79: /* A = 1/h^2*tridiag(1,-2,1) + a*I */
80: MatCreate(PETSC_COMM_WORLD,&A);
81: MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n);
82: MatSetFromOptions(A);
83: MatSetUp(A);
84: MatGetOwnershipRange(A,&Istart,&Iend);
85: for (i=Istart;i<Iend;i++) {
86: if (i>0) MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES);
87: if (i<n-1) MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES);
88: MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES);
89: }
90: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
91: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
92: MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE);
94: /* B = diag(b(xi)) */
95: MatCreate(PETSC_COMM_WORLD,&B);
96: MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n);
97: MatSetFromOptions(B);
98: MatSetUp(B);
99: MatGetOwnershipRange(B,&Istart,&Iend);
100: for (i=Istart;i<Iend;i++) {
101: xi = (i+1)*h;
102: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
103: MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES);
104: }
105: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
106: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
107: MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE);
109: /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */
110: FNCreate(PETSC_COMM_WORLD,&f1);
111: FNSetType(f1,FNRATIONAL);
112: coeffs[0] = -1.0; coeffs[1] = 0.0;
113: FNRationalSetNumerator(f1,2,coeffs);
115: FNCreate(PETSC_COMM_WORLD,&f2);
116: FNSetType(f2,FNRATIONAL);
117: coeffs[0] = 1.0;
118: FNRationalSetNumerator(f2,1,coeffs);
120: FNCreate(PETSC_COMM_WORLD,&f3);
121: FNSetType(f3,FNEXP);
122: FNSetScale(f3,-tau,1.0);
124: /* Set the split operator */
125: mats[0] = A; funs[0] = f2;
126: mats[1] = Id; funs[1] = f1;
127: mats[2] = B; funs[2] = f3;
128: NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN);
130: /* Customize nonlinear solver; set runtime options */
131: NEPSetType(nep,NEPINTERPOL);
132: NEPSetRG(nep,rg);
133: NEPInterpolSetPEP(nep,pep);
134: NEPInterpolGetInterpolation(nep,&tol,°);
135: NEPInterpolSetInterpolation(nep,tol,deg+2); /* increase degree of interpolation */
136: NEPSetFromOptions(nep);
138: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139: Solve the eigensystem
140: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142: NEPSolve(nep);
143: NEPGetDimensions(nep,&nev,NULL,NULL);
144: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
146: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147: Display solution and clean up
148: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150: /* show detailed info unless -terse option is given by user */
151: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
152: if (terse) NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
153: else {
154: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
155: NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
156: NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
157: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
158: }
159: NEPDestroy(&nep);
160: PEPDestroy(&pep);
161: RGDestroy(&rg);
162: MatDestroy(&Id);
163: MatDestroy(&A);
164: MatDestroy(&B);
165: FNDestroy(&f1);
166: FNDestroy(&f2);
167: FNDestroy(&f3);
168: SlepcFinalize();
169: return 0;
170: }
172: /*TEST
174: testset:
175: args: -nep_nev 3 -nep_target 5 -terse
176: output_file: output/test5_1.out
177: filter: sed -e "s/[+-]0\.0*i//g"
178: requires: !single
179: test:
180: suffix: 1
181: test:
182: suffix: 2_cuda
183: args: -mat_type aijcusparse
184: requires: cuda
185: test:
186: suffix: 3
187: args: -nep_view_values draw
188: requires: x
190: TEST*/