Actual source code: test3.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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: Mat F,J;
44: ApplicationCtx ctx;
45: PetscInt n=128;
46: PetscBool terse;
49: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
50: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
51: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%D\n\n",n);
52: ctx.h = 1.0/(PetscReal)n;
53: ctx.kappa = 1.0;
55: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
56: Create a standalone EPS with appropriate settings
57: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
59: EPSCreate(PETSC_COMM_WORLD,&eps);
60: EPSSetType(eps,EPSGD);
61: EPSSetFromOptions(eps);
63: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
64: Prepare nonlinear eigensolver context
65: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
67: NEPCreate(PETSC_COMM_WORLD,&nep);
69: /* Create Function and Jacobian matrices; set evaluation routines */
70: MatCreate(PETSC_COMM_WORLD,&F);
71: MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n);
72: MatSetFromOptions(F);
73: MatSeqAIJSetPreallocation(F,3,NULL);
74: MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
75: MatSetUp(F);
76: NEPSetFunction(nep,F,F,FormFunction,&ctx);
78: MatCreate(PETSC_COMM_WORLD,&J);
79: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);
80: MatSetFromOptions(J);
81: MatSeqAIJSetPreallocation(J,3,NULL);
82: MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
83: MatSetUp(J);
84: NEPSetJacobian(nep,J,FormJacobian,&ctx);
86: NEPSetType(nep,NEPSLP);
87: NEPSLPSetEPS(nep,eps);
88: NEPSetFromOptions(nep);
90: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
91: Solve the eigensystem
92: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
94: NEPSolve(nep);
96: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
97: Display solution and clean up
98: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
100: /* show detailed info unless -terse option is given by user */
101: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
102: if (terse) {
103: NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
104: } else {
105: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
106: NEPReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
107: NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
108: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
109: }
111: NEPDestroy(&nep);
112: EPSDestroy(&eps);
113: MatDestroy(&F);
114: MatDestroy(&J);
115: SlepcFinalize();
116: return ierr;
117: }
119: /* ------------------------------------------------------------------- */
120: /*
121: FormFunction - Computes Function matrix T(lambda)
123: Input Parameters:
124: . nep - the NEP context
125: . lambda - the scalar argument
126: . ctx - optional user-defined context, as set by NEPSetFunction()
128: Output Parameters:
129: . fun - Function matrix
130: . B - optionally different preconditioning matrix
131: */
132: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
133: {
135: ApplicationCtx *user = (ApplicationCtx*)ctx;
136: PetscScalar A[3],c,d;
137: PetscReal h;
138: PetscInt i,n,j[3],Istart,Iend;
139: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
142: /*
143: Compute Function entries and insert into matrix
144: */
145: MatGetSize(fun,&n,NULL);
146: MatGetOwnershipRange(fun,&Istart,&Iend);
147: if (Istart==0) FirstBlock=PETSC_TRUE;
148: if (Iend==n) LastBlock=PETSC_TRUE;
149: h = user->h;
150: c = user->kappa/(lambda-user->kappa);
151: d = n;
153: /*
154: Interior grid points
155: */
156: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
157: j[0] = i-1; j[1] = i; j[2] = i+1;
158: A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
159: MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES);
160: }
162: /*
163: Boundary points
164: */
165: if (FirstBlock) {
166: i = 0;
167: j[0] = 0; j[1] = 1;
168: A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
169: MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
170: }
172: if (LastBlock) {
173: i = n-1;
174: j[0] = n-2; j[1] = n-1;
175: A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
176: MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
177: }
179: /*
180: Assemble matrix
181: */
182: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
183: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
184: if (fun != B) {
185: MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
186: MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
187: }
188: return(0);
189: }
191: /* ------------------------------------------------------------------- */
192: /*
193: FormJacobian - Computes Jacobian matrix T'(lambda)
195: Input Parameters:
196: . nep - the NEP context
197: . lambda - the scalar argument
198: . ctx - optional user-defined context, as set by NEPSetJacobian()
200: Output Parameters:
201: . jac - Jacobian matrix
202: . B - optionally different preconditioning matrix
203: */
204: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
205: {
207: ApplicationCtx *user = (ApplicationCtx*)ctx;
208: PetscScalar A[3],c;
209: PetscReal h;
210: PetscInt i,n,j[3],Istart,Iend;
211: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
214: /*
215: Compute Jacobian entries and insert into matrix
216: */
217: MatGetSize(jac,&n,NULL);
218: MatGetOwnershipRange(jac,&Istart,&Iend);
219: if (Istart==0) FirstBlock=PETSC_TRUE;
220: if (Iend==n) LastBlock=PETSC_TRUE;
221: h = user->h;
222: c = user->kappa/(lambda-user->kappa);
224: /*
225: Interior grid points
226: */
227: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
228: j[0] = i-1; j[1] = i; j[2] = i+1;
229: A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
230: MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES);
231: }
233: /*
234: Boundary points
235: */
236: if (FirstBlock) {
237: i = 0;
238: j[0] = 0; j[1] = 1;
239: A[0] = -2.0*h/3.0; A[1] = -h/6.0;
240: MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
241: }
243: if (LastBlock) {
244: i = n-1;
245: j[0] = n-2; j[1] = n-1;
246: A[0] = -h/6.0; A[1] = -h/3.0-c*c;
247: MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
248: }
250: /*
251: Assemble matrix
252: */
253: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
254: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
255: return(0);
256: }
258: /*TEST
260: test:
261: suffix: 1
262: args: -nep_target 21 -terse
263: requires: !single
264: output_file: output/test1_1.out
266: TEST*/