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