Actual source code: test17.c
slepc-3.21.2 2024-09-25
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[] = "Tests a user-provided preconditioner.\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"
15: " -a <a>, where <a> is the coefficient that multiplies u in the equation.\n"
16: " -split <0/1>, to select the split form in the problem definition (enabled by default).\n";
18: /* Based on ex22.c (delay) */
20: #include <slepcnep.h>
22: /*
23: User-defined application context
24: */
25: typedef struct {
26: PetscScalar tau;
27: PetscReal a;
28: } ApplicationCtx;
30: /*
31: Create problem matrices in split form
32: */
33: PetscErrorCode BuildSplitMatrices(PetscInt n,PetscReal a,Mat *Id,Mat *A,Mat *B)
34: {
35: PetscInt i,Istart,Iend;
36: PetscReal h,xi;
37: PetscScalar b;
39: PetscFunctionBeginUser;
40: h = PETSC_PI/(PetscReal)(n+1);
42: /* Identity matrix */
43: PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,Id));
44: PetscCall(MatSetOption(*Id,MAT_HERMITIAN,PETSC_TRUE));
46: /* A = 1/h^2*tridiag(1,-2,1) + a*I */
47: PetscCall(MatCreate(PETSC_COMM_WORLD,A));
48: PetscCall(MatSetSizes(*A,PETSC_DECIDE,PETSC_DECIDE,n,n));
49: PetscCall(MatSetFromOptions(*A));
50: PetscCall(MatGetOwnershipRange(*A,&Istart,&Iend));
51: for (i=Istart;i<Iend;i++) {
52: if (i>0) PetscCall(MatSetValue(*A,i,i-1,1.0/(h*h),INSERT_VALUES));
53: if (i<n-1) PetscCall(MatSetValue(*A,i,i+1,1.0/(h*h),INSERT_VALUES));
54: PetscCall(MatSetValue(*A,i,i,-2.0/(h*h)+a,INSERT_VALUES));
55: }
56: PetscCall(MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY));
57: PetscCall(MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY));
58: PetscCall(MatSetOption(*A,MAT_HERMITIAN,PETSC_TRUE));
60: /* B = diag(b(xi)) */
61: PetscCall(MatCreate(PETSC_COMM_WORLD,B));
62: PetscCall(MatSetSizes(*B,PETSC_DECIDE,PETSC_DECIDE,n,n));
63: PetscCall(MatSetFromOptions(*B));
64: PetscCall(MatGetOwnershipRange(*B,&Istart,&Iend));
65: for (i=Istart;i<Iend;i++) {
66: xi = (i+1)*h;
67: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
68: PetscCall(MatSetValue(*B,i,i,b,INSERT_VALUES));
69: }
70: PetscCall(MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY));
71: PetscCall(MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY));
72: PetscCall(MatSetOption(*B,MAT_HERMITIAN,PETSC_TRUE));
73: PetscFunctionReturn(PETSC_SUCCESS);
74: }
76: /*
77: Create preconditioner matrices (only Ap=diag(A))
78: */
79: PetscErrorCode BuildSplitPreconditioner(PetscInt n,PetscReal a,Mat *Ap)
80: {
81: PetscInt i,Istart,Iend;
82: PetscReal h;
84: PetscFunctionBeginUser;
85: h = PETSC_PI/(PetscReal)(n+1);
87: /* Ap = diag(A) */
88: PetscCall(MatCreate(PETSC_COMM_WORLD,Ap));
89: PetscCall(MatSetSizes(*Ap,PETSC_DECIDE,PETSC_DECIDE,n,n));
90: PetscCall(MatSetFromOptions(*Ap));
91: PetscCall(MatGetOwnershipRange(*Ap,&Istart,&Iend));
92: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(*Ap,i,i,-2.0/(h*h)+a,INSERT_VALUES));
93: PetscCall(MatAssemblyBegin(*Ap,MAT_FINAL_ASSEMBLY));
94: PetscCall(MatAssemblyEnd(*Ap,MAT_FINAL_ASSEMBLY));
95: PetscCall(MatSetOption(*Ap,MAT_HERMITIAN,PETSC_TRUE));
96: PetscFunctionReturn(PETSC_SUCCESS);
97: }
99: /*
100: Compute Function matrix T(lambda)
101: */
102: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
103: {
104: ApplicationCtx *user = (ApplicationCtx*)ctx;
105: PetscInt i,n,Istart,Iend;
106: PetscReal h,xi;
107: PetscScalar b;
109: PetscFunctionBeginUser;
110: PetscCall(MatGetSize(fun,&n,NULL));
111: h = PETSC_PI/(PetscReal)(n+1);
112: PetscCall(MatGetOwnershipRange(fun,&Istart,&Iend));
113: for (i=Istart;i<Iend;i++) {
114: if (i>0) PetscCall(MatSetValue(fun,i,i-1,1.0/(h*h),INSERT_VALUES));
115: if (i<n-1) PetscCall(MatSetValue(fun,i,i+1,1.0/(h*h),INSERT_VALUES));
116: xi = (i+1)*h;
117: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
118: PetscCall(MatSetValue(fun,i,i,-lambda-2.0/(h*h)+user->a+PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES));
119: if (B!=fun) PetscCall(MatSetValue(B,i,i,-lambda-2.0/(h*h)+user->a+PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES));
120: }
121: PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY));
122: PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY));
123: if (fun != B) {
124: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
125: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
126: }
127: PetscFunctionReturn(PETSC_SUCCESS);
128: }
130: /*
131: Compute Jacobian matrix T'(lambda)
132: */
133: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
134: {
135: ApplicationCtx *user = (ApplicationCtx*)ctx;
136: PetscInt i,n,Istart,Iend;
137: PetscReal h,xi;
138: PetscScalar b;
140: PetscFunctionBeginUser;
141: PetscCall(MatGetSize(jac,&n,NULL));
142: h = PETSC_PI/(PetscReal)(n+1);
143: PetscCall(MatGetOwnershipRange(jac,&Istart,&Iend));
144: for (i=Istart;i<Iend;i++) {
145: xi = (i+1)*h;
146: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
147: PetscCall(MatSetValue(jac,i,i,-1.0-user->tau*PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES));
148: }
149: PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
150: PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
151: PetscFunctionReturn(PETSC_SUCCESS);
152: }
154: int main(int argc,char **argv)
155: {
156: NEP nep; /* nonlinear eigensolver context */
157: Mat Id,A,B,Ap,J,F,P; /* problem matrices */
158: FN f1,f2,f3; /* functions to define the nonlinear operator */
159: ApplicationCtx ctx; /* user-defined context */
160: Mat mats[3],M;
161: FN funs[3];
162: PetscScalar coeffs[2];
163: PetscInt n=128,nterm;
164: PetscReal tau=0.001,a=20;
165: PetscBool split=PETSC_TRUE;
166: MatStructure mstr;
168: PetscFunctionBeginUser;
169: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
170: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
171: PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
172: PetscCall(PetscOptionsGetReal(NULL,NULL,"-a",&a,NULL));
173: PetscCall(PetscOptionsGetBool(NULL,NULL,"-split",&split,NULL));
174: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g, a=%g\n\n",n,(double)tau,(double)a));
176: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177: Create nonlinear eigensolver and solve the problem
178: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180: PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
181: if (split) {
182: PetscCall(BuildSplitMatrices(n,a,&Id,&A,&B));
183: /* f1=-lambda */
184: PetscCall(FNCreate(PETSC_COMM_WORLD,&f1));
185: PetscCall(FNSetType(f1,FNRATIONAL));
186: coeffs[0] = -1.0; coeffs[1] = 0.0;
187: PetscCall(FNRationalSetNumerator(f1,2,coeffs));
188: /* f2=1.0 */
189: PetscCall(FNCreate(PETSC_COMM_WORLD,&f2));
190: PetscCall(FNSetType(f2,FNRATIONAL));
191: coeffs[0] = 1.0;
192: PetscCall(FNRationalSetNumerator(f2,1,coeffs));
193: /* f3=exp(-tau*lambda) */
194: PetscCall(FNCreate(PETSC_COMM_WORLD,&f3));
195: PetscCall(FNSetType(f3,FNEXP));
196: PetscCall(FNSetScale(f3,-tau,1.0));
197: mats[0] = A; funs[0] = f2;
198: mats[1] = Id; funs[1] = f1;
199: mats[2] = B; funs[2] = f3;
200: PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN));
201: PetscCall(BuildSplitPreconditioner(n,a,&Ap));
202: mats[0] = Ap;
203: mats[1] = Id;
204: mats[2] = B;
205: PetscCall(NEPSetSplitPreconditioner(nep,3,mats,SAME_NONZERO_PATTERN));
206: } else {
207: /* callback form */
208: ctx.tau = tau;
209: ctx.a = a;
210: PetscCall(MatCreate(PETSC_COMM_WORLD,&F));
211: PetscCall(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n));
212: PetscCall(MatSetFromOptions(F));
213: PetscCall(MatSeqAIJSetPreallocation(F,3,NULL));
214: PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
215: PetscCall(MatDuplicate(F,MAT_DO_NOT_COPY_VALUES,&P));
216: PetscCall(NEPSetFunction(nep,F,P,FormFunction,&ctx));
217: PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
218: PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
219: PetscCall(MatSetFromOptions(J));
220: PetscCall(MatSeqAIJSetPreallocation(J,3,NULL));
221: PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
222: PetscCall(NEPSetJacobian(nep,J,FormJacobian,&ctx));
223: }
225: /* Set solver parameters at runtime */
226: PetscCall(NEPSetFromOptions(nep));
228: if (split) {
229: PetscCall(NEPGetSplitPreconditionerInfo(nep,&nterm,&mstr));
230: PetscCall(PetscPrintf(PETSC_COMM_WORLD," Nonlinear preconditioner with %" PetscInt_FMT " terms, with %s nonzero pattern\n",nterm,MatStructures[mstr]));
231: PetscCall(NEPGetSplitPreconditionerTerm(nep,0,&M));
232: }
234: /* Solve the eigensystem */
235: PetscCall(NEPSolve(nep));
236: PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
238: PetscCall(NEPDestroy(&nep));
239: if (split) {
240: PetscCall(MatDestroy(&Id));
241: PetscCall(MatDestroy(&A));
242: PetscCall(MatDestroy(&B));
243: PetscCall(MatDestroy(&Ap));
244: PetscCall(FNDestroy(&f1));
245: PetscCall(FNDestroy(&f2));
246: PetscCall(FNDestroy(&f3));
247: } else {
248: PetscCall(MatDestroy(&F));
249: PetscCall(MatDestroy(&P));
250: PetscCall(MatDestroy(&J));
251: }
252: PetscCall(SlepcFinalize());
253: return 0;
254: }
256: /*TEST
258: testset:
259: args: -a 90000 -nep_nev 2
260: requires: double !defined(PETSCTEST_VALGRIND)
261: output_file: output/test17_1.out
262: timeoutfactor: 2
263: filter: grep -v "with 3 terms, with SAME"
264: test:
265: suffix: 1
266: args: -nep_type slp -nep_two_sided {{0 1}} -split {{0 1}}
268: testset:
269: args: -nep_nev 2 -rg_type interval -rg_interval_endpoints .5,15,-.1,.1 -nep_target .7
270: requires: !single
271: output_file: output/test17_2.out
272: filter: sed -e "s/[+-]0\.0*i//g" | grep -v "with 3 terms, with SAME"
273: test:
274: suffix: 2_interpol
275: args: -nep_type interpol -nep_interpol_st_ksp_type bcgs -nep_interpol_st_pc_type sor -nep_tol 1e-6 -nep_interpol_st_ksp_rtol 1e-7
276: test:
277: suffix: 2_nleigs
278: args: -nep_type nleigs -split {{0 1}}
279: requires: complex
280: test:
281: suffix: 2_nleigs_real
282: args: -nep_type nleigs -rg_interval_endpoints .5,15 -split {{0 1}} -nep_nleigs_ksp_type tfqmr
283: requires: !complex
285: testset:
286: args: -nep_type ciss -rg_type ellipse -rg_ellipse_center 10 -rg_ellipse_radius 9.5 -rg_ellipse_vscale 0.1 -nep_ciss_ksp_type bcgs -nep_ciss_pc_type sor
287: output_file: output/test17_3.out
288: requires: complex !single !defined(PETSCTEST_VALGRIND)
289: filter: grep -v "with 3 terms, with SAME"
290: test:
291: suffix: 3
292: args: -split {{0 1}}
293: test:
294: suffix: 3_par
295: nsize: 2
296: args: -nep_ciss_partitions 2
298: TEST*/