Actual source code: test3.c
slepc-3.9.0 2018-04-12
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2018, 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: ST st;
44: KSP ksp;
45: PC pc;
46: Mat F,J;
47: ApplicationCtx ctx;
48: PetscInt n=128;
49: PetscBool terse;
52: SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
53: PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
54: PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%D\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: EPSSetTarget(eps,0.0);
64: EPSGetST(eps,&st);
65: STSetType(st,STSINVERT);
66: STGetKSP(st,&ksp);
67: KSPSetType(ksp,KSPBCGS);
68: KSPGetPC(ksp,&pc);
69: PCSetType(pc,PCBJACOBI);
70: EPSSetFromOptions(eps);
72: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
73: Prepare nonlinear eigensolver context
74: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
76: NEPCreate(PETSC_COMM_WORLD,&nep);
78: /* Create Function and Jacobian matrices; set evaluation routines */
79: MatCreate(PETSC_COMM_WORLD,&F);
80: MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n);
81: MatSetFromOptions(F);
82: MatSeqAIJSetPreallocation(F,3,NULL);
83: MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
84: MatSetUp(F);
85: NEPSetFunction(nep,F,F,FormFunction,&ctx);
87: MatCreate(PETSC_COMM_WORLD,&J);
88: MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);
89: MatSetFromOptions(J);
90: MatSeqAIJSetPreallocation(J,3,NULL);
91: MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
92: MatSetUp(J);
93: NEPSetJacobian(nep,J,FormJacobian,&ctx);
95: NEPSetType(nep,NEPSLP);
96: NEPSLPSetEPS(nep,eps);
97: NEPSetFromOptions(nep);
99: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100: Solve the eigensystem
101: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
103: NEPSolve(nep);
105: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
106: Display solution and clean up
107: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
109: /* show detailed info unless -terse option is given by user */
110: PetscOptionsHasName(NULL,NULL,"-terse",&terse);
111: if (terse) {
112: NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
113: } else {
114: PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
115: NEPReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
116: NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
117: PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
118: }
120: NEPDestroy(&nep);
121: EPSDestroy(&eps);
122: MatDestroy(&F);
123: MatDestroy(&J);
124: SlepcFinalize();
125: return ierr;
126: }
128: /* ------------------------------------------------------------------- */
129: /*
130: FormFunction - Computes Function matrix T(lambda)
132: Input Parameters:
133: . nep - the NEP context
134: . lambda - the scalar argument
135: . ctx - optional user-defined context, as set by NEPSetFunction()
137: Output Parameters:
138: . fun - Function matrix
139: . B - optionally different preconditioning matrix
140: */
141: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
142: {
144: ApplicationCtx *user = (ApplicationCtx*)ctx;
145: PetscScalar A[3],c,d;
146: PetscReal h;
147: PetscInt i,n,j[3],Istart,Iend;
148: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
151: /*
152: Compute Function entries and insert into matrix
153: */
154: MatGetSize(fun,&n,NULL);
155: MatGetOwnershipRange(fun,&Istart,&Iend);
156: if (Istart==0) FirstBlock=PETSC_TRUE;
157: if (Iend==n) LastBlock=PETSC_TRUE;
158: h = user->h;
159: c = user->kappa/(lambda-user->kappa);
160: d = n;
162: /*
163: Interior grid points
164: */
165: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
166: j[0] = i-1; j[1] = i; j[2] = i+1;
167: A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
168: MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES);
169: }
171: /*
172: Boundary points
173: */
174: if (FirstBlock) {
175: i = 0;
176: j[0] = 0; j[1] = 1;
177: A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
178: MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
179: }
181: if (LastBlock) {
182: i = n-1;
183: j[0] = n-2; j[1] = n-1;
184: A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
185: MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
186: }
188: /*
189: Assemble matrix
190: */
191: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
192: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
193: if (fun != B) {
194: MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
195: MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
196: }
197: return(0);
198: }
200: /* ------------------------------------------------------------------- */
201: /*
202: FormJacobian - Computes Jacobian matrix T'(lambda)
204: Input Parameters:
205: . nep - the NEP context
206: . lambda - the scalar argument
207: . ctx - optional user-defined context, as set by NEPSetJacobian()
209: Output Parameters:
210: . jac - Jacobian matrix
211: . B - optionally different preconditioning matrix
212: */
213: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
214: {
216: ApplicationCtx *user = (ApplicationCtx*)ctx;
217: PetscScalar A[3],c;
218: PetscReal h;
219: PetscInt i,n,j[3],Istart,Iend;
220: PetscBool FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;
223: /*
224: Compute Jacobian entries and insert into matrix
225: */
226: MatGetSize(jac,&n,NULL);
227: MatGetOwnershipRange(jac,&Istart,&Iend);
228: if (Istart==0) FirstBlock=PETSC_TRUE;
229: if (Iend==n) LastBlock=PETSC_TRUE;
230: h = user->h;
231: c = user->kappa/(lambda-user->kappa);
233: /*
234: Interior grid points
235: */
236: for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
237: j[0] = i-1; j[1] = i; j[2] = i+1;
238: A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
239: MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES);
240: }
242: /*
243: Boundary points
244: */
245: if (FirstBlock) {
246: i = 0;
247: j[0] = 0; j[1] = 1;
248: A[0] = -2.0*h/3.0; A[1] = -h/6.0;
249: MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
250: }
252: if (LastBlock) {
253: i = n-1;
254: j[0] = n-2; j[1] = n-1;
255: A[0] = -h/6.0; A[1] = -h/3.0-c*c;
256: MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
257: }
259: /*
260: Assemble matrix
261: */
262: MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
263: MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
264: return(0);
265: }