Actual source code: ex27.c
slepc-3.17.1 2022-04-11
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[] = "Simple nonlinear eigenproblem using the NLEIGS solver.\n\n"
12: "The command line options are:\n"
13: " -n <n>, where <n> = matrix dimension.\n"
14: " -split <0/1>, to select the split form in the problem definition (enabled by default)\n";
16: /*
17: Solve T(lambda)x=0 using NLEIGS solver
18: with T(lambda) = -D+sqrt(lambda)*I
19: where D is the Laplacian operator in 1 dimension
20: and with the interpolation interval [.01,16]
21: */
23: #include <slepcnep.h>
25: /*
26: User-defined routines
27: */
28: PetscErrorCode FormFunction(NEP,PetscScalar,Mat,Mat,void*);
29: PetscErrorCode FormJacobian(NEP,PetscScalar,Mat,void*);
30: PetscErrorCode ComputeSingularities(NEP,PetscInt*,PetscScalar*,void*);
32: int main(int argc,char **argv)
33: {
34: NEP nep; /* nonlinear eigensolver context */
35: Mat F,J,A[2];
36: NEPType type;
37: PetscInt n=100,nev,Istart,Iend,i;
38: PetscBool terse,split=PETSC_TRUE;
39: RG rg;
40: FN f[2];
41: PetscScalar coeffs;
43: SlepcInitialize(&argc,&argv,(char*)0,help);
44: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
45: PetscOptionsGetBool(NULL,NULL,"-split",&split,NULL);
46: PetscPrintf(PETSC_COMM_WORLD,"\nSquare root eigenproblem, n=%" PetscInt_FMT "%s\n\n",n,split?" (in split form)":"");
48: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
49: Create nonlinear eigensolver context
50: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
52: NEPCreate(PETSC_COMM_WORLD,&nep);
54: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
55: Select the NLEIGS solver and set required options for it
56: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
58: NEPSetType(nep,NEPNLEIGS);
59: NEPNLEIGSSetSingularitiesFunction(nep,ComputeSingularities,NULL);
60: NEPGetRG(nep,&rg);
61: RGSetType(rg,RGINTERVAL);
62: #if defined(PETSC_USE_COMPLEX)
63: RGIntervalSetEndpoints(rg,0.01,16.0,-0.001,0.001);
64: #else
65: RGIntervalSetEndpoints(rg,0.01,16.0,0,0);
66: #endif
67: NEPSetTarget(nep,1.1);
69: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
70: Define the nonlinear problem
71: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
73: if (split) {
74: /*
75: Create matrices for the split form
76: */
77: MatCreate(PETSC_COMM_WORLD,&A[0]);
78: MatSetSizes(A[0],PETSC_DECIDE,PETSC_DECIDE,n,n);
79: MatSetFromOptions(A[0]);
80: MatSetUp(A[0]);
81: MatGetOwnershipRange(A[0],&Istart,&Iend);
82: for (i=Istart;i<Iend;i++) {
83: if (i>0) MatSetValue(A[0],i,i-1,1.0,INSERT_VALUES);
84: if (i<n-1) MatSetValue(A[0],i,i+1,1.0,INSERT_VALUES);
85: MatSetValue(A[0],i,i,-2.0,INSERT_VALUES);
86: }
87: MatAssemblyBegin(A[0],MAT_FINAL_ASSEMBLY);
88: MatAssemblyEnd(A[0],MAT_FINAL_ASSEMBLY);
90: MatCreateConstantDiagonal(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,n,n,1.0,&A[1]);
92: /*
93: Define functions for the split form
94: */
95: FNCreate(PETSC_COMM_WORLD,&f[0]);
96: FNSetType(f[0],FNRATIONAL);
97: coeffs = 1.0;
98: FNRationalSetNumerator(f[0],1,&coeffs);
99: FNCreate(PETSC_COMM_WORLD,&f[1]);
100: FNSetType(f[1],FNSQRT);
101: NEPSetSplitOperator(nep,2,A,f,SUBSET_NONZERO_PATTERN);
103: } else {
104: /*
105: Callback form: create matrix and set Function evaluation routine
106: */
107: MatCreate(PETSC_COMM_WORLD,&F);
108: MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n);
109: MatSetFromOptions(F);
110: MatSeqAIJSetPreallocation(F,3,NULL);
111: MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
112: MatSetUp(F);
113: NEPSetFunction(nep,F,F,FormFunction,NULL);
115: MatCreate(PETSC_COMM_WORLD,&J);
116: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);
117: MatSetFromOptions(J);
118: MatSeqAIJSetPreallocation(J,1,NULL);
119: MatMPIAIJSetPreallocation(J,1,NULL,1,NULL);
120: MatSetUp(J);
121: NEPSetJacobian(nep,J,FormJacobian,NULL);
122: }
124: /*
125: Set solver parameters at runtime
126: */
127: NEPSetFromOptions(nep);
129: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
130: Solve the eigensystem
131: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
132: NEPSolve(nep);
133: NEPGetType(nep,&type);
134: PetscPrintf(PETSC_COMM_WORLD," Solution method: %s\n",type);
135: NEPGetDimensions(nep,&nev,NULL,NULL);
136: PetscPrintf(PETSC_COMM_WORLD," Number of requested eigenvalues: %" PetscInt_FMT "\n",nev);
138: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139: Display solution and clean up
140: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142: /* show detailed info unless -terse option is given by user */
143: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
144: if (terse) NEPErrorView(nep,NEP_ERROR_BACKWARD,NULL);
145: else {
146: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
147: NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
148: NEPErrorView(nep,NEP_ERROR_BACKWARD,PETSC_VIEWER_STDOUT_WORLD);
149: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
150: }
151: NEPDestroy(&nep);
152: if (split) {
153: MatDestroy(&A[0]);
154: MatDestroy(&A[1]);
155: FNDestroy(&f[0]);
156: FNDestroy(&f[1]);
157: } else {
158: MatDestroy(&F);
159: MatDestroy(&J);
160: }
161: SlepcFinalize();
162: return 0;
163: }
165: /* ------------------------------------------------------------------- */
166: /*
167: FormFunction - Computes Function matrix T(lambda)
168: */
169: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
170: {
171: PetscInt i,n,col[3],Istart,Iend;
172: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
173: PetscScalar value[3],t;
176: /*
177: Compute Function entries and insert into matrix
178: */
179: t = PetscSqrtScalar(lambda);
180: MatGetSize(fun,&n,NULL);
181: MatGetOwnershipRange(fun,&Istart,&Iend);
182: if (Istart==0) FirstBlock=PETSC_TRUE;
183: if (Iend==n) LastBlock=PETSC_TRUE;
184: value[0]=1.0; value[1]=t-2.0; value[2]=1.0;
185: for (i=(FirstBlock? Istart+1: Istart); i<(LastBlock? Iend-1: Iend); i++) {
186: col[0]=i-1; col[1]=i; col[2]=i+1;
187: MatSetValues(fun,1,&i,3,col,value,INSERT_VALUES);
188: }
189: if (LastBlock) {
190: i=n-1; col[0]=n-2; col[1]=n-1;
191: MatSetValues(fun,1,&i,2,col,value,INSERT_VALUES);
192: }
193: if (FirstBlock) {
194: i=0; col[0]=0; col[1]=1; value[0]=t-2.0; value[1]=1.0;
195: MatSetValues(fun,1,&i,2,col,value,INSERT_VALUES);
196: }
198: /*
199: Assemble matrix
200: */
201: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
202: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
203: if (fun != B) {
204: MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
205: MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
206: }
207: PetscFunctionReturn(0);
208: }
210: /* ------------------------------------------------------------------- */
211: /*
212: FormJacobian - Computes Jacobian matrix T'(lambda)
213: */
214: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
215: {
216: Vec d;
219: MatCreateVecs(jac,&d,NULL);
220: VecSet(d,0.5/PetscSqrtScalar(lambda));
221: MatDiagonalSet(jac,d,INSERT_VALUES);
222: VecDestroy(&d);
223: PetscFunctionReturn(0);
224: }
226: /* ------------------------------------------------------------------- */
227: /*
228: ComputeSingularities - Computes maxnp points (at most) in the complex plane where
229: the function T(.) is not analytic.
231: In this case, we discretize the singularity region (-inf,0)~(-10e+6,-10e-6)
232: */
233: PetscErrorCode ComputeSingularities(NEP nep,PetscInt *maxnp,PetscScalar *xi,void *pt)
234: {
235: PetscReal h;
236: PetscInt i;
239: h = 11.0/(*maxnp-1);
240: xi[0] = -1e-5; xi[*maxnp-1] = -1e+6;
241: for (i=1;i<*maxnp-1;i++) xi[i] = -PetscPowReal(10,-5+h*i);
242: PetscFunctionReturn(0);
243: }
245: /*TEST
247: testset:
248: args: -nep_nev 3 -terse
249: output_file: output/ex27_1.out
250: requires: !single
251: filter: sed -e "s/[+-]0\.0*i//g"
252: test:
253: suffix: 1
254: args: -nep_nleigs_interpolation_degree 90
255: test:
256: suffix: 3
257: args: -nep_tol 1e-8 -nep_nleigs_rk_shifts 1.06,1.1,1.12,1.15 -nep_conv_norm -nep_nleigs_interpolation_degree 20
258: test:
259: suffix: 5
260: args: -mat_type aijcusparse
261: requires: cuda
263: testset:
264: args: -split 0 -nep_nev 3 -terse
265: output_file: output/ex27_2.out
266: filter: sed -e "s/[+-]0\.0*i//g"
267: test:
268: suffix: 2
269: args: -nep_nleigs_interpolation_degree 90
270: requires: !single
271: test:
272: suffix: 4
273: args: -nep_nleigs_rk_shifts 1.06,1.1,1.12,1.15 -nep_nleigs_interpolation_degree 20
274: requires: double
275: test:
276: suffix: 6
277: args: -mat_type aijcusparse
278: requires: cuda !single
280: testset:
281: args: -split 0 -nep_type ciss -nep_ciss_extraction {{ritz hankel caa}} -rg_type ellipse -rg_ellipse_center 8 -rg_ellipse_radius .7 -nep_ciss_moments 4 -rg_ellipse_vscale 0.1 -terse
282: requires: complex !single
283: output_file: output/ex27_7.out
284: timeoutfactor: 2
285: test:
286: suffix: 7
287: test:
288: suffix: 7_par
289: nsize: 2
290: args: -nep_ciss_partitions 2
292: testset:
293: args: -nep_type ciss -rg_type ellipse -rg_ellipse_center 8 -rg_ellipse_radius .7 -rg_ellipse_vscale 0.1 -terse
294: requires: complex
295: filter: sed -e "s/ (in split form)//" | sed -e "s/56925/56924/" | sed -e "s/60753/60754/" | sed -e "s/92630/92629/"
296: output_file: output/ex27_7.out
297: timeoutfactor: 2
298: test:
299: suffix: 8
300: test:
301: suffix: 8_parallel
302: nsize: 4
303: args: -nep_ciss_partitions 4 -ds_parallel distributed
304: test:
305: suffix: 8_hpddm
306: args: -nep_ciss_ksp_type hpddm
307: requires: hpddm
309: test:
310: suffix: 9
311: args: -nep_nev 4 -n 20 -terse
312: requires: !single
313: filter: sed -e "s/[+-]0\.0*i//g"
315: TEST*/