Actual source code: test17.c
slepc-3.20.1 2023-11-27
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(MatSetUp(*A));
51: PetscCall(MatGetOwnershipRange(*A,&Istart,&Iend));
52: for (i=Istart;i<Iend;i++) {
53: if (i>0) PetscCall(MatSetValue(*A,i,i-1,1.0/(h*h),INSERT_VALUES));
54: if (i<n-1) PetscCall(MatSetValue(*A,i,i+1,1.0/(h*h),INSERT_VALUES));
55: PetscCall(MatSetValue(*A,i,i,-2.0/(h*h)+a,INSERT_VALUES));
56: }
57: PetscCall(MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY));
58: PetscCall(MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY));
59: PetscCall(MatSetOption(*A,MAT_HERMITIAN,PETSC_TRUE));
61: /* B = diag(b(xi)) */
62: PetscCall(MatCreate(PETSC_COMM_WORLD,B));
63: PetscCall(MatSetSizes(*B,PETSC_DECIDE,PETSC_DECIDE,n,n));
64: PetscCall(MatSetFromOptions(*B));
65: PetscCall(MatSetUp(*B));
66: PetscCall(MatGetOwnershipRange(*B,&Istart,&Iend));
67: for (i=Istart;i<Iend;i++) {
68: xi = (i+1)*h;
69: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
70: PetscCall(MatSetValue(*B,i,i,b,INSERT_VALUES));
71: }
72: PetscCall(MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY));
73: PetscCall(MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY));
74: PetscCall(MatSetOption(*B,MAT_HERMITIAN,PETSC_TRUE));
75: PetscFunctionReturn(PETSC_SUCCESS);
76: }
78: /*
79: Create preconditioner matrices (only Ap=diag(A))
80: */
81: PetscErrorCode BuildSplitPreconditioner(PetscInt n,PetscReal a,Mat *Ap)
82: {
83: PetscInt i,Istart,Iend;
84: PetscReal h;
86: PetscFunctionBeginUser;
87: h = PETSC_PI/(PetscReal)(n+1);
89: /* Ap = diag(A) */
90: PetscCall(MatCreate(PETSC_COMM_WORLD,Ap));
91: PetscCall(MatSetSizes(*Ap,PETSC_DECIDE,PETSC_DECIDE,n,n));
92: PetscCall(MatSetFromOptions(*Ap));
93: PetscCall(MatSetUp(*Ap));
94: PetscCall(MatGetOwnershipRange(*Ap,&Istart,&Iend));
95: for (i=Istart;i<Iend;i++) PetscCall(MatSetValue(*Ap,i,i,-2.0/(h*h)+a,INSERT_VALUES));
96: PetscCall(MatAssemblyBegin(*Ap,MAT_FINAL_ASSEMBLY));
97: PetscCall(MatAssemblyEnd(*Ap,MAT_FINAL_ASSEMBLY));
98: PetscCall(MatSetOption(*Ap,MAT_HERMITIAN,PETSC_TRUE));
99: PetscFunctionReturn(PETSC_SUCCESS);
100: }
102: /*
103: Compute Function matrix T(lambda)
104: */
105: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
106: {
107: ApplicationCtx *user = (ApplicationCtx*)ctx;
108: PetscInt i,n,Istart,Iend;
109: PetscReal h,xi;
110: PetscScalar b;
112: PetscFunctionBeginUser;
113: PetscCall(MatGetSize(fun,&n,NULL));
114: h = PETSC_PI/(PetscReal)(n+1);
115: PetscCall(MatGetOwnershipRange(fun,&Istart,&Iend));
116: for (i=Istart;i<Iend;i++) {
117: if (i>0) PetscCall(MatSetValue(fun,i,i-1,1.0/(h*h),INSERT_VALUES));
118: if (i<n-1) PetscCall(MatSetValue(fun,i,i+1,1.0/(h*h),INSERT_VALUES));
119: xi = (i+1)*h;
120: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
121: PetscCall(MatSetValue(fun,i,i,-lambda-2.0/(h*h)+user->a+PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES));
122: if (B!=fun) PetscCall(MatSetValue(B,i,i,-lambda-2.0/(h*h)+user->a+PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES));
123: }
124: PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY));
125: PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY));
126: if (fun != B) {
127: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
128: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
129: }
130: PetscFunctionReturn(PETSC_SUCCESS);
131: }
133: /*
134: Compute Jacobian matrix T'(lambda)
135: */
136: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
137: {
138: ApplicationCtx *user = (ApplicationCtx*)ctx;
139: PetscInt i,n,Istart,Iend;
140: PetscReal h,xi;
141: PetscScalar b;
143: PetscFunctionBeginUser;
144: PetscCall(MatGetSize(jac,&n,NULL));
145: h = PETSC_PI/(PetscReal)(n+1);
146: PetscCall(MatGetOwnershipRange(jac,&Istart,&Iend));
147: for (i=Istart;i<Iend;i++) {
148: xi = (i+1)*h;
149: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
150: PetscCall(MatSetValue(jac,i,i,-1.0-user->tau*PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES));
151: }
152: PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
153: PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
154: PetscFunctionReturn(PETSC_SUCCESS);
155: }
157: int main(int argc,char **argv)
158: {
159: NEP nep; /* nonlinear eigensolver context */
160: Mat Id,A,B,Ap,J,F,P; /* problem matrices */
161: FN f1,f2,f3; /* functions to define the nonlinear operator */
162: ApplicationCtx ctx; /* user-defined context */
163: Mat mats[3];
164: FN funs[3];
165: PetscScalar coeffs[2];
166: PetscInt n=128;
167: PetscReal tau=0.001,a=20;
168: PetscBool split=PETSC_TRUE;
170: PetscFunctionBeginUser;
171: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
172: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
173: PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
174: PetscCall(PetscOptionsGetReal(NULL,NULL,"-a",&a,NULL));
175: PetscCall(PetscOptionsGetBool(NULL,NULL,"-split",&split,NULL));
176: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g, a=%g\n\n",n,(double)tau,(double)a));
178: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179: Create nonlinear eigensolver and solve the problem
180: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
182: PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
183: if (split) {
184: PetscCall(BuildSplitMatrices(n,a,&Id,&A,&B));
185: /* f1=-lambda */
186: PetscCall(FNCreate(PETSC_COMM_WORLD,&f1));
187: PetscCall(FNSetType(f1,FNRATIONAL));
188: coeffs[0] = -1.0; coeffs[1] = 0.0;
189: PetscCall(FNRationalSetNumerator(f1,2,coeffs));
190: /* f2=1.0 */
191: PetscCall(FNCreate(PETSC_COMM_WORLD,&f2));
192: PetscCall(FNSetType(f2,FNRATIONAL));
193: coeffs[0] = 1.0;
194: PetscCall(FNRationalSetNumerator(f2,1,coeffs));
195: /* f3=exp(-tau*lambda) */
196: PetscCall(FNCreate(PETSC_COMM_WORLD,&f3));
197: PetscCall(FNSetType(f3,FNEXP));
198: PetscCall(FNSetScale(f3,-tau,1.0));
199: mats[0] = A; funs[0] = f2;
200: mats[1] = Id; funs[1] = f1;
201: mats[2] = B; funs[2] = f3;
202: PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN));
203: PetscCall(BuildSplitPreconditioner(n,a,&Ap));
204: mats[0] = Ap;
205: mats[1] = Id;
206: mats[2] = B;
207: PetscCall(NEPSetSplitPreconditioner(nep,3,mats,SAME_NONZERO_PATTERN));
208: } else {
209: /* callback form */
210: ctx.tau = tau;
211: ctx.a = a;
212: PetscCall(MatCreate(PETSC_COMM_WORLD,&F));
213: PetscCall(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n));
214: PetscCall(MatSetFromOptions(F));
215: PetscCall(MatSeqAIJSetPreallocation(F,3,NULL));
216: PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
217: PetscCall(MatSetUp(F));
218: PetscCall(MatDuplicate(F,MAT_DO_NOT_COPY_VALUES,&P));
219: PetscCall(NEPSetFunction(nep,F,P,FormFunction,&ctx));
220: PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
221: PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
222: PetscCall(MatSetFromOptions(J));
223: PetscCall(MatSeqAIJSetPreallocation(J,3,NULL));
224: PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
225: PetscCall(MatSetUp(J));
226: PetscCall(NEPSetJacobian(nep,J,FormJacobian,&ctx));
227: }
229: /* Set solver parameters at runtime */
230: PetscCall(NEPSetFromOptions(nep));
232: /* Solve the eigensystem */
233: PetscCall(NEPSolve(nep));
234: PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
236: PetscCall(NEPDestroy(&nep));
237: if (split) {
238: PetscCall(MatDestroy(&Id));
239: PetscCall(MatDestroy(&A));
240: PetscCall(MatDestroy(&B));
241: PetscCall(MatDestroy(&Ap));
242: PetscCall(FNDestroy(&f1));
243: PetscCall(FNDestroy(&f2));
244: PetscCall(FNDestroy(&f3));
245: } else {
246: PetscCall(MatDestroy(&F));
247: PetscCall(MatDestroy(&P));
248: PetscCall(MatDestroy(&J));
249: }
250: PetscCall(SlepcFinalize());
251: return 0;
252: }
254: /*TEST
256: testset:
257: args: -a 90000 -nep_nev 2
258: requires: double !defined(PETSCTEST_VALGRIND)
259: output_file: output/test17_1.out
260: timeoutfactor: 2
261: test:
262: suffix: 1
263: args: -nep_type slp -nep_two_sided {{0 1}} -split {{0 1}}
265: testset:
266: args: -nep_nev 2 -rg_type interval -rg_interval_endpoints .5,15,-.1,.1 -nep_target .7
267: requires: !single
268: output_file: output/test17_2.out
269: filter: sed -e "s/[+-]0\.0*i//g"
270: test:
271: suffix: 2_interpol
272: 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
273: test:
274: suffix: 2_nleigs
275: args: -nep_type nleigs -split {{0 1}}
276: requires: complex
277: test:
278: suffix: 2_nleigs_real
279: args: -nep_type nleigs -rg_interval_endpoints .5,15 -split {{0 1}} -nep_nleigs_ksp_type tfqmr
280: requires: !complex
282: testset:
283: 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
284: output_file: output/test17_3.out
285: requires: complex !single !defined(PETSCTEST_VALGRIND)
286: test:
287: suffix: 3
288: args: -split {{0 1}}
289: test:
290: suffix: 3_par
291: nsize: 2
292: args: -nep_ciss_partitions 2
294: TEST*/