Actual source code: test16.c
slepc-3.17.2 2022-08-09
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;
50: SlepcInitialize(&argc,&argv,(char*)0,help);
51: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
52: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n);
53: ctx.h = 1.0/(PetscReal)n;
54: ctx.kappa = 1.0;
56: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
57: Prepare nonlinear eigensolver context
58: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60: NEPCreate(PETSC_COMM_WORLD,&nep);
62: MatCreate(PETSC_COMM_WORLD,&F);
63: MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n);
64: MatSetFromOptions(F);
65: MatSeqAIJSetPreallocation(F,3,NULL);
66: MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
67: MatSetUp(F);
68: NEPSetFunction(nep,F,F,FormFunction,&ctx);
70: MatCreate(PETSC_COMM_WORLD,&J);
71: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);
72: MatSetFromOptions(J);
73: MatSeqAIJSetPreallocation(J,3,NULL);
74: MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
75: MatSetUp(J);
76: NEPSetJacobian(nep,J,FormJacobian,&ctx);
78: NEPSetType(nep,NEPNLEIGS);
79: NEPGetRG(nep,&rg);
80: RGSetType(rg,RGINTERVAL);
81: #if defined(PETSC_USE_COMPLEX)
82: RGIntervalSetEndpoints(rg,2.0,400.0,-0.001,0.001);
83: #else
84: RGIntervalSetEndpoints(rg,2.0,400.0,0,0);
85: #endif
86: NEPSetTarget(nep,25.0);
87: NEPSetEigenvalueComparison(nep,MyEigenSort,&target);
88: NEPSetTolerances(nep,PETSC_SMALL,PETSC_DEFAULT);
89: NEPSetFromOptions(nep);
90: NEPGetTarget(nep,&target);
92: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
93: Solve the eigensystem and display the solution
94: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
96: NEPSolve(nep);
98: /* show detailed info unless -terse option is given by user */
99: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
100: if (terse) NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
101: else {
102: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
103: NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
104: NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
105: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
106: }
108: NEPDestroy(&nep);
109: MatDestroy(&F);
110: MatDestroy(&J);
111: SlepcFinalize();
112: return 0;
113: }
115: /* ------------------------------------------------------------------- */
116: /*
117: FormFunction - Computes Function matrix T(lambda)
119: Input Parameters:
120: . nep - the NEP context
121: . lambda - the scalar argument
122: . ctx - optional user-defined context, as set by NEPSetFunction()
124: Output Parameters:
125: . fun - Function matrix
126: . B - optionally different preconditioning matrix
127: */
128: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
129: {
130: ApplicationCtx *user = (ApplicationCtx*)ctx;
131: PetscScalar A[3],c,d;
132: PetscReal h;
133: PetscInt i,n,j[3],Istart,Iend;
134: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
137: /*
138: Compute Function entries and insert into matrix
139: */
140: MatGetSize(fun,&n,NULL);
141: MatGetOwnershipRange(fun,&Istart,&Iend);
142: if (Istart==0) FirstBlock=PETSC_TRUE;
143: if (Iend==n) LastBlock=PETSC_TRUE;
144: h = user->h;
145: c = user->kappa/(lambda-user->kappa);
146: d = n;
148: /*
149: Interior grid points
150: */
151: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
152: j[0] = i-1; j[1] = i; j[2] = i+1;
153: A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
154: MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES);
155: }
157: /*
158: Boundary points
159: */
160: if (FirstBlock) {
161: i = 0;
162: j[0] = 0; j[1] = 1;
163: A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
164: MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
165: }
167: if (LastBlock) {
168: i = n-1;
169: j[0] = n-2; j[1] = n-1;
170: A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
171: MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
172: }
174: /*
175: Assemble matrix
176: */
177: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
178: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
179: if (fun != B) {
180: MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
181: MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
182: }
183: PetscFunctionReturn(0);
184: }
186: /* ------------------------------------------------------------------- */
187: /*
188: FormJacobian - Computes Jacobian matrix T'(lambda)
190: Input Parameters:
191: . nep - the NEP context
192: . lambda - the scalar argument
193: . ctx - optional user-defined context, as set by NEPSetJacobian()
195: Output Parameters:
196: . jac - Jacobian matrix
197: . B - optionally different preconditioning matrix
198: */
199: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
200: {
201: ApplicationCtx *user = (ApplicationCtx*)ctx;
202: PetscScalar A[3],c;
203: PetscReal h;
204: PetscInt i,n,j[3],Istart,Iend;
205: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
208: /*
209: Compute Jacobian entries and insert into matrix
210: */
211: MatGetSize(jac,&n,NULL);
212: MatGetOwnershipRange(jac,&Istart,&Iend);
213: if (Istart==0) FirstBlock=PETSC_TRUE;
214: if (Iend==n) LastBlock=PETSC_TRUE;
215: h = user->h;
216: c = user->kappa/(lambda-user->kappa);
218: /*
219: Interior grid points
220: */
221: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
222: j[0] = i-1; j[1] = i; j[2] = i+1;
223: A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
224: MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES);
225: }
227: /*
228: Boundary points
229: */
230: if (FirstBlock) {
231: i = 0;
232: j[0] = 0; j[1] = 1;
233: A[0] = -2.0*h/3.0; A[1] = -h/6.0;
234: MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
235: }
237: if (LastBlock) {
238: i = n-1;
239: j[0] = n-2; j[1] = n-1;
240: A[0] = -h/6.0; A[1] = -h/3.0-c*c;
241: MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
242: }
244: /*
245: Assemble matrix
246: */
247: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
248: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
249: PetscFunctionReturn(0);
250: }
252: /*
253: Function for user-defined eigenvalue ordering criterion.
255: Given two eigenvalues ar+i*ai and br+i*bi, the subroutine must choose
256: one of them as the preferred one according to the criterion.
257: In this example, eigenvalues are sorted with respect to the target,
258: but those on the right of the target are preferred.
259: */
260: PetscErrorCode MyEigenSort(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *r,void *ctx)
261: {
262: PetscReal a,b;
263: PetscScalar target = *(PetscScalar*)ctx;
266: if (PetscRealPart(ar-target)<0.0 && PetscRealPart(br-target)>0.0) *r = 1;
267: else {
268: a = SlepcAbsEigenvalue(ar-target,ai);
269: b = SlepcAbsEigenvalue(br-target,bi);
270: if (a>b) *r = 1;
271: else if (a<b) *r = -1;
272: else *r = 0;
273: }
274: PetscFunctionReturn(0);
275: }
277: /*TEST
279: test:
280: suffix: 1
281: args: -nep_nev 4 -nep_ncv 8 -terse
282: requires: double !complex
284: TEST*/