Actual source code: test3.c

slepc-3.15.2 2021-09-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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;

 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:   EPSSetType(eps,EPSGD);
 64:   EPSSetFromOptions(eps);

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

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

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

 80:   NEPCreate(PETSC_COMM_WORLD,&nep);

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

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

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

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

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

119:   NEPSolve(nep);

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

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

136:   NEPDestroy(&nep);
137:   EPSDestroy(&eps);
138:   KSPDestroy(&ksp);
139:   MatDestroy(&F);
140:   MatDestroy(&J);
141:   SlepcFinalize();
142:   return ierr;
143: }

145: /* ------------------------------------------------------------------- */
146: /*
147:    FormFunction - Computes Function matrix  T(lambda)

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

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

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

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

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

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

205:   /*
206:      Assemble matrix
207:   */
208:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
209:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
210:   if (fun != B) {
211:     MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
212:     MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
213:   }
214:   return(0);
215: }

217: /* ------------------------------------------------------------------- */
218: /*
219:    FormJacobian - Computes Jacobian matrix  T'(lambda)

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

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

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

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

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

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

276:   /*
277:      Assemble matrix
278:   */
279:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
280:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
281:   return(0);
282: }

284: /*TEST

286:    test:
287:       args: -nep_target 21 -terse
288:       requires: !single
289:       test:
290:          suffix: 1
291:       test:
292:          suffix: 1_ts
293:          args: -nep_two_sided

295: TEST*/