Actual source code: test3.c

slepc-3.21.1 2024-04-26
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:   PetscFunctionBeginUser;
 52:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 53:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 54:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%" PetscInt_FMT "\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:   PetscCall(EPSCreate(PETSC_COMM_WORLD,&eps));
 63:   PetscCall(EPSSetType(eps,EPSGD));
 64:   PetscCall(EPSSetFromOptions(eps));

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

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

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

 80:   PetscCall(NEPCreate(PETSC_COMM_WORLD,&nep));

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

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

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

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

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

117:   PetscCall(NEPSolve(nep));

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

279: /*TEST

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

290: TEST*/