Actual source code: test16.c
slepc-3.15.0 2021-03-31
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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[] = "Illustrates use of NEPSetEigenvalueComparison().\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*);
30: PetscErrorCode MyEigenSort(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*);
32: /*
33: User-defined application context
34: */
35: typedef struct {
36: PetscScalar kappa; /* ratio between stiffness of spring and attached mass */
37: PetscReal h; /* mesh spacing */
38: } ApplicationCtx;
40: int main(int argc,char **argv)
41: {
42: NEP nep; /* nonlinear eigensolver context */
43: Mat F,J; /* Function and Jacobian matrices */
44: ApplicationCtx ctx; /* user-defined context */
45: PetscScalar target;
46: RG rg;
47: PetscInt n=128;
48: PetscBool terse;
51: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
52: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
53: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%D\n\n",n);
54: ctx.h = 1.0/(PetscReal)n;
55: ctx.kappa = 1.0;
57: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58: Prepare nonlinear eigensolver context
59: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61: NEPCreate(PETSC_COMM_WORLD,&nep);
63: MatCreate(PETSC_COMM_WORLD,&F);
64: MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n);
65: MatSetFromOptions(F);
66: MatSeqAIJSetPreallocation(F,3,NULL);
67: MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
68: MatSetUp(F);
69: NEPSetFunction(nep,F,F,FormFunction,&ctx);
71: MatCreate(PETSC_COMM_WORLD,&J);
72: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);
73: MatSetFromOptions(J);
74: MatSeqAIJSetPreallocation(J,3,NULL);
75: MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
76: MatSetUp(J);
77: NEPSetJacobian(nep,J,FormJacobian,&ctx);
79: NEPSetType(nep,NEPNLEIGS);
80: NEPGetRG(nep,&rg);
81: RGSetType(rg,RGINTERVAL);
82: #if defined(PETSC_USE_COMPLEX)
83: RGIntervalSetEndpoints(rg,2.0,400.0,-0.001,0.001);
84: #else
85: RGIntervalSetEndpoints(rg,2.0,400.0,0,0);
86: #endif
87: NEPSetTarget(nep,25.0);
88: NEPSetEigenvalueComparison(nep,MyEigenSort,&target);
89: NEPSetTolerances(nep,PETSC_SMALL,PETSC_DEFAULT);
90: NEPSetFromOptions(nep);
91: NEPGetTarget(nep,&target);
93: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
94: Solve the eigensystem and display the solution
95: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
97: NEPSolve(nep);
99: /* show detailed info unless -terse option is given by user */
100: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
101: if (terse) {
102: NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
103: } else {
104: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
105: NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
106: NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
107: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
108: }
110: NEPDestroy(&nep);
111: MatDestroy(&F);
112: MatDestroy(&J);
113: SlepcFinalize();
114: return ierr;
115: }
117: /* ------------------------------------------------------------------- */
118: /*
119: FormFunction - Computes Function matrix T(lambda)
121: Input Parameters:
122: . nep - the NEP context
123: . lambda - the scalar argument
124: . ctx - optional user-defined context, as set by NEPSetFunction()
126: Output Parameters:
127: . fun - Function matrix
128: . B - optionally different preconditioning matrix
129: */
130: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
131: {
133: ApplicationCtx *user = (ApplicationCtx*)ctx;
134: PetscScalar A[3],c,d;
135: PetscReal h;
136: PetscInt i,n,j[3],Istart,Iend;
137: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
140: /*
141: Compute Function entries and insert into matrix
142: */
143: MatGetSize(fun,&n,NULL);
144: MatGetOwnershipRange(fun,&Istart,&Iend);
145: if (Istart==0) FirstBlock=PETSC_TRUE;
146: if (Iend==n) LastBlock=PETSC_TRUE;
147: h = user->h;
148: c = user->kappa/(lambda-user->kappa);
149: d = n;
151: /*
152: Interior grid points
153: */
154: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
155: j[0] = i-1; j[1] = i; j[2] = i+1;
156: A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
157: MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES);
158: }
160: /*
161: Boundary points
162: */
163: if (FirstBlock) {
164: i = 0;
165: j[0] = 0; j[1] = 1;
166: A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
167: MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
168: }
170: if (LastBlock) {
171: i = n-1;
172: j[0] = n-2; j[1] = n-1;
173: A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
174: MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
175: }
177: /*
178: Assemble matrix
179: */
180: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
181: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
182: if (fun != B) {
183: MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
184: MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
185: }
186: return(0);
187: }
189: /* ------------------------------------------------------------------- */
190: /*
191: FormJacobian - Computes Jacobian matrix T'(lambda)
193: Input Parameters:
194: . nep - the NEP context
195: . lambda - the scalar argument
196: . ctx - optional user-defined context, as set by NEPSetJacobian()
198: Output Parameters:
199: . jac - Jacobian matrix
200: . B - optionally different preconditioning matrix
201: */
202: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
203: {
205: ApplicationCtx *user = (ApplicationCtx*)ctx;
206: PetscScalar A[3],c;
207: PetscReal h;
208: PetscInt i,n,j[3],Istart,Iend;
209: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
212: /*
213: Compute Jacobian entries and insert into matrix
214: */
215: MatGetSize(jac,&n,NULL);
216: MatGetOwnershipRange(jac,&Istart,&Iend);
217: if (Istart==0) FirstBlock=PETSC_TRUE;
218: if (Iend==n) LastBlock=PETSC_TRUE;
219: h = user->h;
220: c = user->kappa/(lambda-user->kappa);
222: /*
223: Interior grid points
224: */
225: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
226: j[0] = i-1; j[1] = i; j[2] = i+1;
227: A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
228: MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES);
229: }
231: /*
232: Boundary points
233: */
234: if (FirstBlock) {
235: i = 0;
236: j[0] = 0; j[1] = 1;
237: A[0] = -2.0*h/3.0; A[1] = -h/6.0;
238: MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
239: }
241: if (LastBlock) {
242: i = n-1;
243: j[0] = n-2; j[1] = n-1;
244: A[0] = -h/6.0; A[1] = -h/3.0-c*c;
245: MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
246: }
248: /*
249: Assemble matrix
250: */
251: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
252: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
253: return(0);
254: }
256: /*
257: Function for user-defined eigenvalue ordering criterion.
259: Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
260: one of them as the preferred one according to the criterion.
261: In this example, eigenvalues are sorted with respect to the target,
262: but those on the right of the target are preferred.
263: */
264: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
265: {
266: PetscReal a,b;
267: PetscScalar target = *(PetscScalar*)ctx;
270: if (PetscRealPart(ar-target)<0.0 && PetscRealPart(br-target)>0.0) *r = 1;
271: else {
272: a = SlepcAbsEigenvalue(ar-target,ai);
273: b = SlepcAbsEigenvalue(br-target,bi);
274: if (a>b) *r = 1;
275: else if (a<b) *r = -1;
276: else *r = 0;
277: }
278: return(0);
279: }
281: /*TEST
283: test:
284: suffix: 1
285: args: -nep_nev 4 -nep_ncv 8 -terse
286: requires: double !complex
288: TEST*/