Actual source code: test11.c
slepc-3.20.2 2024-03-15
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 CISS solver with the problem of ex22.\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;
38: Mat Id,A,B,mats[3];
39: FN f1,f2,f3,funs[3];
40: RG rg;
41: KSP *ksp;
42: PC pc;
43: NEPCISSExtraction ext;
44: PetscScalar coeffs[2],b;
45: PetscInt n=128,Istart,Iend,i,nsolve;
46: PetscReal tau=0.001,h,a=20,xi;
47: PetscBool flg,terse;
49: PetscFunctionBeginUser;
50: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
51: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
52: PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
53: PetscCall(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 nonlinear eigensolver context
58: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60: PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
62: /* Identity matrix */
63: PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&Id));
64: PetscCall(MatSetOption(Id,MAT_HERMITIAN,PETSC_TRUE));
66: /* A = 1/h^2*tridiag(1,-2,1) + a*I */
67: PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
68: PetscCall(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,n,n));
69: PetscCall(MatSetFromOptions(A));
70: PetscCall(MatSetUp(A));
71: PetscCall(MatGetOwnershipRange(A,&Istart,&Iend));
72: for (i=Istart;i<Iend;i++) {
73: if (i>0) PetscCall(MatSetValue(A,i,i-1,1.0/(h*h),INSERT_VALUES));
74: if (i<n-1) PetscCall(MatSetValue(A,i,i+1,1.0/(h*h),INSERT_VALUES));
75: PetscCall(MatSetValue(A,i,i,-2.0/(h*h)+a,INSERT_VALUES));
76: }
77: PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
78: PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
79: PetscCall(MatSetOption(A,MAT_HERMITIAN,PETSC_TRUE));
81: /* B = diag(b(xi)) */
82: PetscCall(MatCreate(PETSC_COMM_WORLD,&B));
83: PetscCall(MatSetSizes(B,PETSC_DECIDE,PETSC_DECIDE,n,n));
84: PetscCall(MatSetFromOptions(B));
85: PetscCall(MatSetUp(B));
86: PetscCall(MatGetOwnershipRange(B,&Istart,&Iend));
87: for (i=Istart;i<Iend;i++) {
88: xi = (i+1)*h;
89: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
90: PetscCall(MatSetValues(B,1,&i,1,&i,&b,INSERT_VALUES));
91: }
92: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
93: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
94: PetscCall(MatSetOption(B,MAT_HERMITIAN,PETSC_TRUE));
96: /* Functions: f1=-lambda, f2=1.0, f3=exp(-tau*lambda) */
97: PetscCall(FNCreate(PETSC_COMM_WORLD,&f1));
98: PetscCall(FNSetType(f1,FNRATIONAL));
99: coeffs[0] = -1.0; coeffs[1] = 0.0;
100: PetscCall(FNRationalSetNumerator(f1,2,coeffs));
102: PetscCall(FNCreate(PETSC_COMM_WORLD,&f2));
103: PetscCall(FNSetType(f2,FNRATIONAL));
104: coeffs[0] = 1.0;
105: PetscCall(FNRationalSetNumerator(f2,1,coeffs));
107: PetscCall(FNCreate(PETSC_COMM_WORLD,&f3));
108: PetscCall(FNSetType(f3,FNEXP));
109: PetscCall(FNSetScale(f3,-tau,1.0));
111: /* Set the split operator */
112: mats[0] = A; funs[0] = f2;
113: mats[1] = Id; funs[1] = f1;
114: mats[2] = B; funs[2] = f3;
115: PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN));
117: /* Customize nonlinear solver; set runtime options */
118: PetscCall(NEPSetType(nep,NEPCISS));
119: PetscCall(NEPSetDimensions(nep,1,24,PETSC_DEFAULT));
120: PetscCall(NEPSetTolerances(nep,1e-9,PETSC_DEFAULT));
121: PetscCall(NEPGetRG(nep,&rg));
122: PetscCall(RGSetType(rg,RGELLIPSE));
123: PetscCall(RGEllipseSetParameters(rg,10.0,9.5,0.1));
124: PetscCall(NEPCISSSetSizes(nep,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,1,PETSC_DEFAULT,PETSC_TRUE));
125: PetscCall(NEPCISSGetKSPs(nep,&nsolve,&ksp));
126: for (i=0;i<nsolve;i++) {
127: PetscCall(KSPSetTolerances(ksp[i],1e-12,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
128: PetscCall(KSPSetType(ksp[i],KSPBCGS));
129: PetscCall(KSPGetPC(ksp[i],&pc));
130: PetscCall(PCSetType(pc,PCSOR));
131: }
132: PetscCall(NEPSetFromOptions(nep));
134: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
135: Solve the eigensystem
136: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138: PetscCall(PetscObjectTypeCompare((PetscObject)nep,NEPCISS,&flg));
139: if (flg) {
140: PetscCall(NEPCISSGetExtraction(nep,&ext));
141: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Running CISS with %" PetscInt_FMT " KSP solvers (%s extraction)\n",nsolve,NEPCISSExtractions[ext]));
142: }
143: PetscCall(NEPSolve(nep));
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(MatDestroy(&Id));
160: PetscCall(MatDestroy(&A));
161: PetscCall(MatDestroy(&B));
162: PetscCall(FNDestroy(&f1));
163: PetscCall(FNDestroy(&f2));
164: PetscCall(FNDestroy(&f3));
165: PetscCall(SlepcFinalize());
166: return 0;
167: }
169: /*TEST
171: build:
172: requires: complex
174: testset:
175: args: -nep_ciss_extraction {{ritz hankel caa}} -terse
176: requires: complex !single
177: output_file: output/test11_1.out
178: filter: sed -e "s/([A-Z]* extraction)//"
179: test:
180: suffix: 1
181: args: -nep_ciss_delta 1e-10
182: test:
183: suffix: 2
184: nsize: 2
185: args: -nep_ciss_partitions 2
186: test:
187: suffix: 3
188: args: -nep_ciss_moments 6 -nep_ciss_blocksize 4 -nep_ciss_refine_inner 1 -nep_ciss_refine_blocksize 2
190: test:
191: suffix: 4
192: args: -terse -nep_view
193: requires: complex !single
194: filter: grep -v tolerance | grep -v threshold
196: TEST*/