Actual source code: test3.c

slepc-3.12.2 2020-01-13
Report Typos and Errors
  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:   KSP            ksp;
 44:   PC             pc;
 45:   Mat            F,J;
 46:   ApplicationCtx ctx;
 47:   PetscInt       n=128;
 48:   PetscBool      terse;

 51:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 52:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 53:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%D\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,PCSOR);
 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:   NEPSetType(nep,NEPSLP);
 99:   NEPSLPSetEPS(nep,eps);
100:   NEPSLPSetKSP(nep,ksp);
101:   NEPSetFromOptions(nep);

103:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
104:                       Solve the eigensystem
105:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

107:   NEPSolve(nep);

109:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
110:                     Display solution and clean up
111:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

113:   /* show detailed info unless -terse option is given by user */
114:   PetscOptionsHasName(NULL,NULL,"-terse",&terse);
115:   if (terse) {
116:     NEPErrorView(nep,NEP_ERROR_RELATIVE,NULL);
117:   } else {
118:     PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD,PETSC_VIEWER_ASCII_INFO_DETAIL);
119:     NEPReasonView(nep,PETSC_VIEWER_STDOUT_WORLD);
120:     NEPErrorView(nep,NEP_ERROR_RELATIVE,PETSC_VIEWER_STDOUT_WORLD);
121:     PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD);
122:   }

124:   NEPDestroy(&nep);
125:   EPSDestroy(&eps);
126:   KSPDestroy(&ksp);
127:   MatDestroy(&F);
128:   MatDestroy(&J);
129:   SlepcFinalize();
130:   return ierr;
131: }

133: /* ------------------------------------------------------------------- */
134: /*
135:    FormFunction - Computes Function matrix  T(lambda)

137:    Input Parameters:
138: .  nep    - the NEP context
139: .  lambda - the scalar argument
140: .  ctx    - optional user-defined context, as set by NEPSetFunction()

142:    Output Parameters:
143: .  fun - Function matrix
144: .  B   - optionally different preconditioning matrix
145: */
146: PetscErrorCode FormFunction(NEP nep,PetscScalar lambda,Mat fun,Mat B,void *ctx)
147: {
149:   ApplicationCtx *user = (ApplicationCtx*)ctx;
150:   PetscScalar    A[3],c,d;
151:   PetscReal      h;
152:   PetscInt       i,n,j[3],Istart,Iend;
153:   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;

156:   /*
157:      Compute Function entries and insert into matrix
158:   */
159:   MatGetSize(fun,&n,NULL);
160:   MatGetOwnershipRange(fun,&Istart,&Iend);
161:   if (Istart==0) FirstBlock=PETSC_TRUE;
162:   if (Iend==n) LastBlock=PETSC_TRUE;
163:   h = user->h;
164:   c = user->kappa/(lambda-user->kappa);
165:   d = n;

167:   /*
168:      Interior grid points
169:   */
170:   for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
171:     j[0] = i-1; j[1] = i; j[2] = i+1;
172:     A[0] = A[2] = -d-lambda*h/6.0; A[1] = 2.0*(d-lambda*h/3.0);
173:     MatSetValues(fun,1,&i,3,j,A,INSERT_VALUES);
174:   }

176:   /*
177:      Boundary points
178:   */
179:   if (FirstBlock) {
180:     i = 0;
181:     j[0] = 0; j[1] = 1;
182:     A[0] = 2.0*(d-lambda*h/3.0); A[1] = -d-lambda*h/6.0;
183:     MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
184:   }

186:   if (LastBlock) {
187:     i = n-1;
188:     j[0] = n-2; j[1] = n-1;
189:     A[0] = -d-lambda*h/6.0; A[1] = d-lambda*h/3.0+lambda*c;
190:     MatSetValues(fun,1,&i,2,j,A,INSERT_VALUES);
191:   }

193:   /*
194:      Assemble matrix
195:   */
196:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
197:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
198:   if (fun != B) {
199:     MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
200:     MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
201:   }
202:   return(0);
203: }

205: /* ------------------------------------------------------------------- */
206: /*
207:    FormJacobian - Computes Jacobian matrix  T'(lambda)

209:    Input Parameters:
210: .  nep    - the NEP context
211: .  lambda - the scalar argument
212: .  ctx    - optional user-defined context, as set by NEPSetJacobian()

214:    Output Parameters:
215: .  jac - Jacobian matrix
216: .  B   - optionally different preconditioning matrix
217: */
218: PetscErrorCode FormJacobian(NEP nep,PetscScalar lambda,Mat jac,void *ctx)
219: {
221:   ApplicationCtx *user = (ApplicationCtx*)ctx;
222:   PetscScalar    A[3],c;
223:   PetscReal      h;
224:   PetscInt       i,n,j[3],Istart,Iend;
225:   PetscBool      FirstBlock=PETSC_FALSE,LastBlock=PETSC_FALSE;

228:   /*
229:      Compute Jacobian entries and insert into matrix
230:   */
231:   MatGetSize(jac,&n,NULL);
232:   MatGetOwnershipRange(jac,&Istart,&Iend);
233:   if (Istart==0) FirstBlock=PETSC_TRUE;
234:   if (Iend==n) LastBlock=PETSC_TRUE;
235:   h = user->h;
236:   c = user->kappa/(lambda-user->kappa);

238:   /*
239:      Interior grid points
240:   */
241:   for (i=(FirstBlock? Istart+1: Istart);i<(LastBlock? Iend-1: Iend);i++) {
242:     j[0] = i-1; j[1] = i; j[2] = i+1;
243:     A[0] = A[2] = -h/6.0; A[1] = -2.0*h/3.0;
244:     MatSetValues(jac,1,&i,3,j,A,INSERT_VALUES);
245:   }

247:   /*
248:      Boundary points
249:   */
250:   if (FirstBlock) {
251:     i = 0;
252:     j[0] = 0; j[1] = 1;
253:     A[0] = -2.0*h/3.0; A[1] = -h/6.0;
254:     MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
255:   }

257:   if (LastBlock) {
258:     i = n-1;
259:     j[0] = n-2; j[1] = n-1;
260:     A[0] = -h/6.0; A[1] = -h/3.0-c*c;
261:     MatSetValues(jac,1,&i,2,j,A,INSERT_VALUES);
262:   }

264:   /*
265:      Assemble matrix
266:   */
267:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
268:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
269:   return(0);
270: }

272: /*TEST

274:    test:
275:       suffix: 1
276:       args: -nep_target 21 -terse
277:       requires: !single

279: TEST*/