Actual source code: test16.c
slepc-3.18.3 2023-03-24
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[] = "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);
52: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
53: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%" PetscInt_FMT "\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) NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
102: else {
103: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
104: NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
105: NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
106: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
107: }
109: NEPDestroy(&nep);
110: MatDestroy(&F);
111: MatDestroy(&J);
112: SlepcFinalize();
113: return 0;
114: }
116: /* ------------------------------------------------------------------- */
117: /*
118: FormFunction - Computes Function matrix T(lambda)
120: Input Parameters:
121: . nep - the NEP context
122: . lambda - the scalar argument
123: . ctx - optional user-defined context, as set by NEPSetFunction()
125: Output Parameters:
126: . fun - Function matrix
127: . B - optionally different preconditioning matrix
128: */
129: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
130: {
131: ApplicationCtx *user = (ApplicationCtx*)ctx;
132: PetscScalar A[3],c,d;
133: PetscReal h;
134: PetscInt i,n,j[3],Istart,Iend;
135: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
138: /*
139: Compute Function entries and insert into matrix
140: */
141: MatGetSize(fun,&n,NULL);
142: MatGetOwnershipRange(fun,&Istart,&Iend);
143: if (Istart==0) FirstBlock=PETSC_TRUE;
144: if (Iend==n) LastBlock=PETSC_TRUE;
145: h = user->h;
146: c = user->kappa/(lambda-user->kappa);
147: d = n;
149: /*
150: Interior grid points
151: */
152: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
153: j[0] = i-1; j[1] = i; j[2] = i+1;
154: A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
155: MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES);
156: }
158: /*
159: Boundary points
160: */
161: if (FirstBlock) {
162: i = 0;
163: j[0] = 0; j[1] = 1;
164: A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
165: MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
166: }
168: if (LastBlock) {
169: i = n-1;
170: j[0] = n-2; j[1] = n-1;
171: A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
172: MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
173: }
175: /*
176: Assemble matrix
177: */
178: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
179: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
180: if (fun != B) {
181: MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
182: MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
183: }
184: return 0;
185: }
187: /* ------------------------------------------------------------------- */
188: /*
189: FormJacobian - Computes Jacobian matrix T'(lambda)
191: Input Parameters:
192: . nep - the NEP context
193: . lambda - the scalar argument
194: . ctx - optional user-defined context, as set by NEPSetJacobian()
196: Output Parameters:
197: . jac - Jacobian matrix
198: . B - optionally different preconditioning matrix
199: */
200: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
201: {
202: ApplicationCtx *user = (ApplicationCtx*)ctx;
203: PetscScalar A[3],c;
204: PetscReal h;
205: PetscInt i,n,j[3],Istart,Iend;
206: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
209: /*
210: Compute Jacobian entries and insert into matrix
211: */
212: MatGetSize(jac,&n,NULL);
213: MatGetOwnershipRange(jac,&Istart,&Iend);
214: if (Istart==0) FirstBlock=PETSC_TRUE;
215: if (Iend==n) LastBlock=PETSC_TRUE;
216: h = user->h;
217: c = user->kappa/(lambda-user->kappa);
219: /*
220: Interior grid points
221: */
222: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
223: j[0] = i-1; j[1] = i; j[2] = i+1;
224: A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
225: MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES);
226: }
228: /*
229: Boundary points
230: */
231: if (FirstBlock) {
232: i = 0;
233: j[0] = 0; j[1] = 1;
234: A[0] = -2.0*h/3.0; A[1] = -h/6.0;
235: MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
236: }
238: if (LastBlock) {
239: i = n-1;
240: j[0] = n-2; j[1] = n-1;
241: A[0] = -h/6.0; A[1] = -h/3.0-c*c;
242: MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
243: }
245: /*
246: Assemble matrix
247: */
248: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
249: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
250: return 0;
251: }
253: /*
254: Function for user-defined eigenvalue ordering criterion.
256: Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
257: one of them as the preferred one according to the criterion.
258: In this example, eigenvalues are sorted with respect to the target,
259: but those on the right of the target are preferred.
260: */
261: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
262: {
263: PetscReal a,b;
264: PetscScalar target = *(PetscScalar*)ctx;
267: if (PetscRealPart(ar-target)<0.0 && PetscRealPart(br-target)>0.0) *r = 1;
268: else {
269: a = SlepcAbsEigenvalue(ar-target,ai);
270: b = SlepcAbsEigenvalue(br-target,bi);
271: if (a>b) *r = 1;
272: else if (a<b) *r = -1;
273: else *r = 0;
274: }
275: return 0;
276: }
278: /*TEST
280: test:
281: suffix: 1
282: args: -nep_nev 4 -nep_ncv 8 -terse
283: requires: double !complex
285: TEST*/