Actual source code: test10.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 multiple calls to NEPSolve() with different matrix size.\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: " -split <0/1>, to select the split form in the problem definition (enabled by default).\n";
17: /* Based on ex22.c (delay) */
19: #include <slepcnep.h>
21: /*
22: User-defined application context
23: */
24: typedef struct {
25: PetscScalar tau;
26: PetscReal a;
27: } ApplicationCtx;
29: /*
30: Create problem matrices in split form
31: */
32: PetscErrorCode BuildSplitMatrices(PetscInt n,PetscReal a,Mat *Id,Mat *A,Mat *B)
33: {
34: PetscInt i,Istart,Iend;
35: PetscReal h,xi;
36: PetscScalar b;
38: PetscFunctionBeginUser;
39: h = PETSC_PI/(PetscReal)(n+1);
41: /* Identity matrix */
42: PetscCall(MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,Id));
43: PetscCall(MatSetOption(*Id,MAT_HERMITIAN,PETSC_TRUE));
45: /* A = 1/h^2*tridiag(1,-2,1) + a*I */
46: PetscCall(MatCreate(PETSC_COMM_WORLD,A));
47: PetscCall(MatSetSizes(*A,PETSC_DECIDE,PETSC_DECIDE,n,n));
48: PetscCall(MatSetFromOptions(*A));
49: PetscCall(MatGetOwnershipRange(*A,&Istart,&Iend));
50: for (i=Istart;i<Iend;i++) {
51: if (i>0) PetscCall(MatSetValue(*A,i,i-1,1.0/(h*h),INSERT_VALUES));
52: if (i<n-1) PetscCall(MatSetValue(*A,i,i+1,1.0/(h*h),INSERT_VALUES));
53: PetscCall(MatSetValue(*A,i,i,-2.0/(h*h)+a,INSERT_VALUES));
54: }
55: PetscCall(MatAssemblyBegin(*A,MAT_FINAL_ASSEMBLY));
56: PetscCall(MatAssemblyEnd(*A,MAT_FINAL_ASSEMBLY));
57: PetscCall(MatSetOption(*A,MAT_HERMITIAN,PETSC_TRUE));
59: /* B = diag(b(xi)) */
60: PetscCall(MatCreate(PETSC_COMM_WORLD,B));
61: PetscCall(MatSetSizes(*B,PETSC_DECIDE,PETSC_DECIDE,n,n));
62: PetscCall(MatSetFromOptions(*B));
63: PetscCall(MatGetOwnershipRange(*B,&Istart,&Iend));
64: for (i=Istart;i<Iend;i++) {
65: xi = (i+1)*h;
66: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
67: PetscCall(MatSetValue(*B,i,i,b,INSERT_VALUES));
68: }
69: PetscCall(MatAssemblyBegin(*B,MAT_FINAL_ASSEMBLY));
70: PetscCall(MatAssemblyEnd(*B,MAT_FINAL_ASSEMBLY));
71: PetscCall(MatSetOption(*B,MAT_HERMITIAN,PETSC_TRUE));
72: PetscFunctionReturn(PETSC_SUCCESS);
73: }
75: /*
76: Compute Function matrix T(lambda)
77: */
78: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
79: {
80: ApplicationCtx *user = (ApplicationCtx*)ctx;
81: PetscInt i,n,Istart,Iend;
82: PetscReal h,xi;
83: PetscScalar b;
85: PetscFunctionBeginUser;
86: PetscCall(MatGetSize(fun,&n,NULL));
87: h = PETSC_PI/(PetscReal)(n+1);
88: PetscCall(MatGetOwnershipRange(fun,&Istart,&Iend));
89: for (i=Istart;i<Iend;i++) {
90: if (i>0) PetscCall(MatSetValue(fun,i,i-1,1.0/(h*h),INSERT_VALUES));
91: if (i<n-1) PetscCall(MatSetValue(fun,i,i+1,1.0/(h*h),INSERT_VALUES));
92: xi = (i+1)*h;
93: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
94: PetscCall(MatSetValue(fun,i,i,-lambda-2.0/(h*h)+user->a+PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES));
95: }
96: PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
97: PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
98: if (fun != B) {
99: PetscCall(MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY));
100: PetscCall(MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY));
101: }
102: PetscFunctionReturn(PETSC_SUCCESS);
103: }
105: /*
106: Compute Jacobian matrix T'(lambda)
107: */
108: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
109: {
110: ApplicationCtx *user = (ApplicationCtx*)ctx;
111: PetscInt i,n,Istart,Iend;
112: PetscReal h,xi;
113: PetscScalar b;
115: PetscFunctionBeginUser;
116: PetscCall(MatGetSize(jac,&n,NULL));
117: h = PETSC_PI/(PetscReal)(n+1);
118: PetscCall(MatGetOwnershipRange(jac,&Istart,&Iend));
119: for (i=Istart;i<Iend;i++) {
120: xi = (i+1)*h;
121: b = -4.1+xi*(1.0-PetscExpReal(xi-PETSC_PI));
122: PetscCall(MatSetValue(jac,i,i,-1.0-user->tau*PetscExpScalar(-user->tau*lambda)*b,INSERT_VALUES));
123: }
124: PetscCall(MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY));
125: PetscCall(MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY));
126: PetscFunctionReturn(PETSC_SUCCESS);
127: }
129: int main(int argc,char **argv)
130: {
131: NEP nep; /* nonlinear eigensolver context */
132: Mat Id,A,B,J,F; /* problem matrices */
133: FN f1,f2,f3; /* functions to define the nonlinear operator */
134: ApplicationCtx ctx; /* user-defined context */
135: Mat mats[3];
136: FN funs[3];
137: PetscScalar coeffs[2];
138: PetscInt n=128;
139: PetscReal tau=0.001,a=20;
140: PetscBool split=PETSC_TRUE;
142: PetscFunctionBeginUser;
143: PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
144: PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
145: PetscCall(PetscOptionsGetReal(NULL,NULL,"-tau",&tau,NULL));
146: PetscCall(PetscOptionsGetBool(NULL,NULL,"-split",&split,NULL));
147: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau));
149: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150: Create nonlinear eigensolver and set options
151: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
153: PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));
154: PetscCall(NEPSetTolerances(nep,1e-9,PETSC_DEFAULT));
156: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157: First solve
158: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
160: if (split) {
161: PetscCall(BuildSplitMatrices(n,a,&Id,&A,&B));
162: /* f1=-lambda */
163: PetscCall(FNCreate(PETSC_COMM_WORLD,&f1));
164: PetscCall(FNSetType(f1,FNRATIONAL));
165: coeffs[0] = -1.0; coeffs[1] = 0.0;
166: PetscCall(FNRationalSetNumerator(f1,2,coeffs));
167: /* f2=1.0 */
168: PetscCall(FNCreate(PETSC_COMM_WORLD,&f2));
169: PetscCall(FNSetType(f2,FNRATIONAL));
170: coeffs[0] = 1.0;
171: PetscCall(FNRationalSetNumerator(f2,1,coeffs));
172: /* f3=exp(-tau*lambda) */
173: PetscCall(FNCreate(PETSC_COMM_WORLD,&f3));
174: PetscCall(FNSetType(f3,FNEXP));
175: PetscCall(FNSetScale(f3,-tau,1.0));
176: mats[0] = A; funs[0] = f2;
177: mats[1] = Id; funs[1] = f1;
178: mats[2] = B; funs[2] = f3;
179: PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN));
180: } else {
181: /* callback form */
182: ctx.tau = tau;
183: ctx.a = a;
184: PetscCall(MatCreate(PETSC_COMM_WORLD,&F));
185: PetscCall(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n));
186: PetscCall(MatSetFromOptions(F));
187: PetscCall(MatSeqAIJSetPreallocation(F,3,NULL));
188: PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
189: PetscCall(NEPSetFunction(nep,F,F,FormFunction,&ctx));
190: PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
191: PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
192: PetscCall(MatSetFromOptions(J));
193: PetscCall(MatSeqAIJSetPreallocation(J,3,NULL));
194: PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
195: PetscCall(NEPSetJacobian(nep,J,FormJacobian,&ctx));
196: }
198: /* Set solver parameters at runtime */
199: PetscCall(NEPSetFromOptions(nep));
201: /* Solve the eigensystem */
202: PetscCall(NEPSolve(nep));
203: PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
205: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206: Second solve, with problem matrices of size 2*n
207: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
209: n *= 2;
210: PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Delay Eigenproblem, n=%" PetscInt_FMT ", tau=%g\n\n",n,(double)tau));
211: if (split) {
212: PetscCall(MatDestroy(&Id));
213: PetscCall(MatDestroy(&A));
214: PetscCall(MatDestroy(&B));
215: PetscCall(BuildSplitMatrices(n,a,&Id,&A,&B));
216: mats[0] = A;
217: mats[1] = Id;
218: mats[2] = B;
219: PetscCall(NEPSetSplitOperator(nep,3,mats,funs,SUBSET_NONZERO_PATTERN));
220: } else {
221: /* callback form */
222: PetscCall(MatDestroy(&F));
223: PetscCall(MatDestroy(&J));
224: PetscCall(MatCreate(PETSC_COMM_WORLD,&F));
225: PetscCall(MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n));
226: PetscCall(MatSetFromOptions(F));
227: PetscCall(MatSeqAIJSetPreallocation(F,3,NULL));
228: PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
229: PetscCall(NEPSetFunction(nep,F,F,FormFunction,&ctx));
230: PetscCall(MatCreate(PETSC_COMM_WORLD,&J));
231: PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n));
232: PetscCall(MatSetFromOptions(J));
233: PetscCall(MatSeqAIJSetPreallocation(J,3,NULL));
234: PetscCall(MatMPIAIJSetPreallocation(F,3,NULL,1,NULL));
235: PetscCall(NEPSetJacobian(nep,J,FormJacobian,&ctx));
236: }
238: /* Solve the eigensystem */
239: PetscCall(NEPSolve(nep));
240: PetscCall(NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL));
242: PetscCall(NEPDestroy(&nep));
243: if (split) {
244: PetscCall(MatDestroy(&Id));
245: PetscCall(MatDestroy(&A));
246: PetscCall(MatDestroy(&B));
247: PetscCall(FNDestroy(&f1));
248: PetscCall(FNDestroy(&f2));
249: PetscCall(FNDestroy(&f3));
250: } else {
251: PetscCall(MatDestroy(&F));
252: PetscCall(MatDestroy(&J));
253: }
254: PetscCall(SlepcFinalize());
255: return 0;
256: }
258: /*TEST
260: testset:
261: nsize: 2
262: requires: !single
263: output_file: output/test10_1.out
264: test:
265: suffix: 1
266: args: -nep_type narnoldi -nep_target 0.55
267: test:
268: suffix: 1_rii
269: args: -nep_type rii -nep_target 0.55 -nep_rii_hermitian -split {{0 1}}
270: test:
271: suffix: 1_narnoldi
272: args: -nep_type narnoldi -nep_target 0.55 -nep_narnoldi_lag_preconditioner 2
273: test:
274: suffix: 1_slp
275: args: -nep_type slp -nep_slp_st_pc_type redundant -split {{0 1}}
276: test:
277: suffix: 1_interpol
278: args: -nep_type interpol -rg_type interval -rg_interval_endpoints .5,1,-.1,.1 -nep_target .7 -nep_interpol_st_pc_type redundant
279: test:
280: suffix: 1_narnoldi_sync
281: args: -nep_type narnoldi -ds_parallel synchronized
283: testset:
284: args: -nep_nev 2 -rg_type interval -rg_interval_endpoints .5,15,-.1,.1 -nep_target .7
285: requires: !single
286: output_file: output/test10_2.out
287: filter: sed -e "s/[+-]0\.0*i//g"
288: test:
289: suffix: 2_interpol
290: args: -nep_type interpol -nep_interpol_pep_type jd -nep_interpol_st_pc_type sor
291: test:
292: suffix: 2_nleigs
293: args: -nep_type nleigs -split {{0 1}}
294: requires: complex
295: test:
296: suffix: 2_nleigs_real
297: args: -nep_type nleigs -rg_interval_endpoints .5,15 -split {{0 1}}
298: requires: !complex
300: test:
301: suffix: 3
302: requires: complex !single
303: args: -nep_type ciss -rg_type ellipse -rg_ellipse_center 10 -rg_ellipse_radius 9.5 -rg_ellipse_vscale 0.1 -split {{0 1}}
304: timeoutfactor: 2
306: TEST*/