Actual source code: test3.c
slepc-3.18.2 2023-01-26
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[] = "Test the SLP solver with a user-provided EPS.\n\n"
12: "This is a simplified version of ex20.\n"
13: "The command line options are:\n"
14: " -n <n>, where <n> = number of grid subdivisions.\n";
16: /*
17: Solve 1-D PDE
18: -u'' = lambda*u
19: on [0,1] subject to
20: u(0)=0, u'(1)=u(1)*lambda*kappa/(kappa-lambda)
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*);
31: /*
32: User-defined application context
33: */
34: typedef struct {
35: PetscScalar kappa; /* ratio between stiffness of spring and attached mass */
36: PetscReal h; /* mesh spacing */
37: } ApplicationCtx;
39: int main(int argc,char **argv)
40: {
41: NEP nep;
42: EPS eps;
43: KSP ksp;
44: PC pc;
45: Mat F,J;
46: ApplicationCtx ctx;
47: PetscInt n=128;
48: PetscReal deftol;
49: PetscBool terse,flag,ts;
52: SlepcInitialize(&argc,&argv,(char*)0,help);
53: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
54: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n);
55: ctx.h = 1.0/(PetscReal)n;
56: ctx.kappa = 1.0;
58: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59: Create a standalone EPS with appropriate settings
60: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
62: EPSCreate(PETSC_COMM_WORLD,&eps);
63: EPSSetType(eps,EPSGD);
64: EPSSetFromOptions(eps);
66: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67: Create a standalone KSP with appropriate settings
68: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
70: KSPCreate(PETSC_COMM_WORLD,&ksp);
71: KSPSetType(ksp,KSPBCGS);
72: KSPGetPC(ksp,&pc);
73: PCSetType(pc,PCBJACOBI);
74: KSPSetFromOptions(ksp);
76: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
77: Prepare nonlinear eigensolver context
78: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
80: NEPCreate(PETSC_COMM_WORLD,&nep);
82: /* Create Function and Jacobian matrices; set evaluation routines */
83: MatCreate(PETSC_COMM_WORLD,&F);
84: MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n);
85: MatSetFromOptions(F);
86: MatSeqAIJSetPreallocation(F,3,NULL);
87: MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
88: MatSetUp(F);
89: NEPSetFunction(nep,F,F,FormFunction,&ctx);
91: MatCreate(PETSC_COMM_WORLD,&J);
92: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);
93: MatSetFromOptions(J);
94: MatSeqAIJSetPreallocation(J,3,NULL);
95: MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
96: MatSetUp(J);
97: NEPSetJacobian(nep,J,FormJacobian,&ctx);
99: /* Set options */
100: NEPSetType(nep,NEPSLP);
101: NEPSLPSetEPS(nep,eps);
102: NEPSLPSetKSP(nep,ksp);
103: NEPSetFromOptions(nep);
105: /* Print some options */
106: PetscObjectTypeCompare((PetscObject)nep,NEPSLP,&flag);
107: if (flag) {
108: NEPGetTwoSided(nep,&ts);
109: if (ts) {
110: NEPSLPGetDeflationThreshold(nep,&deftol);
111: PetscPrintf(PETSC_COMM_WORLD," Two-sided solve with deflation threshold=%g\n",(double)deftol);
112: }
113: }
115: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
116: Solve the eigensystem
117: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
119: NEPSolve(nep);
121: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
122: Display solution and clean up
123: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125: /* show detailed info unless -terse option is given by user */
126: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
127: if (terse) NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
128: else {
129: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
130: NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
131: NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
132: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
133: }
135: NEPDestroy(&nep);
136: EPSDestroy(&eps);
137: KSPDestroy(&ksp);
138: MatDestroy(&F);
139: MatDestroy(&J);
140: SlepcFinalize();
141: return 0;
142: }
144: /* ------------------------------------------------------------------- */
145: /*
146: FormFunction - Computes Function matrix T(lambda)
148: Input Parameters:
149: . nep - the NEP context
150: . lambda - the scalar argument
151: . ctx - optional user-defined context, as set by NEPSetFunction()
153: Output Parameters:
154: . fun - Function matrix
155: . B - optionally different preconditioning matrix
156: */
157: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
158: {
159: ApplicationCtx *user = (ApplicationCtx*)ctx;
160: PetscScalar A[3],c,d;
161: PetscReal h;
162: PetscInt i,n,j[3],Istart,Iend;
163: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
166: /*
167: Compute Function entries and insert into matrix
168: */
169: MatGetSize(fun,&n,NULL);
170: MatGetOwnershipRange(fun,&Istart,&Iend);
171: if (Istart==0) FirstBlock=PETSC_TRUE;
172: if (Iend==n) LastBlock=PETSC_TRUE;
173: h = user->h;
174: c = user->kappa/(lambda-user->kappa);
175: d = n;
177: /*
178: Interior grid points
179: */
180: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
181: j[0] = i-1; j[1] = i; j[2] = i+1;
182: A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
183: MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES);
184: }
186: /*
187: Boundary points
188: */
189: if (FirstBlock) {
190: i = 0;
191: j[0] = 0; j[1] = 1;
192: A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
193: MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
194: }
196: if (LastBlock) {
197: i = n-1;
198: j[0] = n-2; j[1] = n-1;
199: A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
200: MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
201: }
203: /*
204: Assemble matrix
205: */
206: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
207: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
208: if (fun != B) {
209: MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
210: MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
211: }
212: return 0;
213: }
215: /* ------------------------------------------------------------------- */
216: /*
217: FormJacobian - Computes Jacobian matrix T'(lambda)
219: Input Parameters:
220: . nep - the NEP context
221: . lambda - the scalar argument
222: . ctx - optional user-defined context, as set by NEPSetJacobian()
224: Output Parameters:
225: . jac - Jacobian matrix
226: . B - optionally different preconditioning matrix
227: */
228: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
229: {
230: ApplicationCtx *user = (ApplicationCtx*)ctx;
231: PetscScalar A[3],c;
232: PetscReal h;
233: PetscInt i,n,j[3],Istart,Iend;
234: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
237: /*
238: Compute Jacobian entries and insert into matrix
239: */
240: MatGetSize(jac,&n,NULL);
241: MatGetOwnershipRange(jac,&Istart,&Iend);
242: if (Istart==0) FirstBlock=PETSC_TRUE;
243: if (Iend==n) LastBlock=PETSC_TRUE;
244: h = user->h;
245: c = user->kappa/(lambda-user->kappa);
247: /*
248: Interior grid points
249: */
250: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
251: j[0] = i-1; j[1] = i; j[2] = i+1;
252: A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
253: MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES);
254: }
256: /*
257: Boundary points
258: */
259: if (FirstBlock) {
260: i = 0;
261: j[0] = 0; j[1] = 1;
262: A[0] = -2.0*h/3.0; A[1] = -h/6.0;
263: MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
264: }
266: if (LastBlock) {
267: i = n-1;
268: j[0] = n-2; j[1] = n-1;
269: A[0] = -h/6.0; A[1] = -h/3.0-c*c;
270: MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
271: }
273: /*
274: Assemble matrix
275: */
276: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
277: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
278: return 0;
279: }
281: /*TEST
283: test:
284: args: -nep_target 21 -terse
285: requires: !single
286: test:
287: suffix: 1
288: test:
289: suffix: 1_ts
290: args: -nep_two_sided
292: TEST*/