Actual source code: test3.c

slepc-3.17.1 2022-04-11
Report Typos and Errors
  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[] = "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:   KSP            ksp;
 44:   PC             pc;
 45:   Mat            F,J;
 46:   ApplicationCtx ctx;
 47:   PetscInt       n=128;
 48:   PetscReal      deftol;
 49:   PetscBool      terse,flag,ts;

 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:         Create a standalone EPS with appropriate settings
 59:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 61:   EPSCreate(PETSC_COMM_WORLD,&eps);
 62:   EPSSetType(eps,EPSGD);
 63:   EPSSetFromOptions(eps);

 65:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 66:         Create a standalone KSP with appropriate settings
 67:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 69:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 70:   KSPSetType(ksp,KSPBCGS);
 71:   KSPGetPC(ksp,&pc);
 72:   PCSetType(pc,PCBJACOBI);
 73:   KSPSetFromOptions(ksp);

 75:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 76:                Prepare nonlinear eigensolver context
 77:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 79:   NEPCreate(PETSC_COMM_WORLD,&nep);

 81:   /* Create Function and Jacobian matrices; set evaluation routines */
 82:   MatCreate(PETSC_COMM_WORLD,&F);
 83:   MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n);
 84:   MatSetFromOptions(F);
 85:   MatSeqAIJSetPreallocation(F,3,NULL);
 86:   MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
 87:   MatSetUp(F);
 88:   NEPSetFunction(nep,F,F,FormFunction,&ctx);

 90:   MatCreate(PETSC_COMM_WORLD,&J);
 91:   MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);
 92:   MatSetFromOptions(J);
 93:   MatSeqAIJSetPreallocation(J,3,NULL);
 94:   MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
 95:   MatSetUp(J);
 96:   NEPSetJacobian(nep,J,FormJacobian,&ctx);

 98:   /* Set options */
 99:   NEPSetType(nep,NEPSLP);
100:   NEPSLPSetEPS(nep,eps);
101:   NEPSLPSetKSP(nep,ksp);
102:   NEPSetFromOptions(nep);

104:   /* Print some options */
105:   PetscObjectTypeCompare((PetscObject)nep,NEPSLP,&flag);
106:   if (flag) {
107:     NEPGetTwoSided(nep,&ts);
108:     if (ts) {
109:       NEPSLPGetDeflationThreshold(nep,&deftol);
110:       PetscPrintf(PETSC_COMM_WORLD," Two-sided solve with deflation threshold=%g\n",(double)deftol);
111:     }
112:   }

114:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
115:                       Solve the eigensystem
116:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

118:   NEPSolve(nep);

120:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121:                     Display solution and clean up
122:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

124:   /* show detailed info unless -terse option is given by user */
125:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
126:   if (terse) NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
127:   else {
128:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
129:     NEPConvergedReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
130:     NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
131:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
132:   }

134:   NEPDestroy(&nep);
135:   EPSDestroy(&eps);
136:   KSPDestroy(&ksp);
137:   MatDestroy(&F);
138:   MatDestroy(&J);
139:   SlepcFinalize();
140:   return 0;
141: }

143: /* ------------------------------------------------------------------- */
144: /*
145:    FormFunction - Computes Function matrix  T(lambda)

147:    Input Parameters:
148: .  nep    - the NEP context
149: .  lambda - the scalar argument
150: .  ctx    - optional user-defined context, as set by NEPSetFunction()

152:    Output Parameters:
153: .  fun - Function matrix
154: .  B   - optionally different preconditioning matrix
155: */
156: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
157: {
158:   ApplicationCtx *user = (ApplicationCtx*)ctx;
159:   PetscScalar    A[3],c,d;
160:   PetscReal      h;
161:   PetscInt       i,n,j[3],Istart,Iend;
162:   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;

165:   /*
166:      Compute Function entries and insert into matrix
167:   */
168:   MatGetSize(fun,&n,NULL);
169:   MatGetOwnershipRange(fun,&Istart,&Iend);
170:   if (Istart==0) FirstBlock=PETSC_TRUE;
171:   if (Iend==n) LastBlock=PETSC_TRUE;
172:   h = user->h;
173:   c = user->kappa/(lambda-user->kappa);
174:   d = n;

176:   /*
177:      Interior grid points
178:   */
179:   for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
180:     j[0] = i-1; j[1] = i; j[2] = i+1;
181:     A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
182:     MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES);
183:   }

185:   /*
186:      Boundary points
187:   */
188:   if (FirstBlock) {
189:     i = 0;
190:     j[0] = 0; j[1] = 1;
191:     A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
192:     MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
193:   }

195:   if (LastBlock) {
196:     i = n-1;
197:     j[0] = n-2; j[1] = n-1;
198:     A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
199:     MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
200:   }

202:   /*
203:      Assemble matrix
204:   */
205:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
206:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
207:   if (fun != B) {
208:     MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
209:     MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
210:   }
211:   PetscFunctionReturn(0);
212: }

214: /* ------------------------------------------------------------------- */
215: /*
216:    FormJacobian - Computes Jacobian matrix  T'(lambda)

218:    Input Parameters:
219: .  nep    - the NEP context
220: .  lambda - the scalar argument
221: .  ctx    - optional user-defined context, as set by NEPSetJacobian()

223:    Output Parameters:
224: .  jac - Jacobian matrix
225: .  B   - optionally different preconditioning matrix
226: */
227: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
228: {
229:   ApplicationCtx *user = (ApplicationCtx*)ctx;
230:   PetscScalar    A[3],c;
231:   PetscReal      h;
232:   PetscInt       i,n,j[3],Istart,Iend;
233:   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;

236:   /*
237:      Compute Jacobian entries and insert into matrix
238:   */
239:   MatGetSize(jac,&n,NULL);
240:   MatGetOwnershipRange(jac,&Istart,&Iend);
241:   if (Istart==0) FirstBlock=PETSC_TRUE;
242:   if (Iend==n) LastBlock=PETSC_TRUE;
243:   h = user->h;
244:   c = user->kappa/(lambda-user->kappa);

246:   /*
247:      Interior grid points
248:   */
249:   for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
250:     j[0] = i-1; j[1] = i; j[2] = i+1;
251:     A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
252:     MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES);
253:   }

255:   /*
256:      Boundary points
257:   */
258:   if (FirstBlock) {
259:     i = 0;
260:     j[0] = 0; j[1] = 1;
261:     A[0] = -2.0*h/3.0; A[1] = -h/6.0;
262:     MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
263:   }

265:   if (LastBlock) {
266:     i = n-1;
267:     j[0] = n-2; j[1] = n-1;
268:     A[0] = -h/6.0; A[1] = -h/3.0-c*c;
269:     MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
270:   }

272:   /*
273:      Assemble matrix
274:   */
275:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
276:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
277:   PetscFunctionReturn(0);
278: }

280: /*TEST

282:    test:
283:       args: -nep_target 21 -terse
284:       requires: !single
285:       test:
286:          suffix: 1
287:       test:
288:          suffix: 1_ts
289:          args: -nep_two_sided

291: TEST*/