Actual source code: test4.c

slepc-3.11.2 2019-07-30
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 RII solver with a user-provided KSP.\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:   KSP            ksp;
 43:   PC             pc;
 44:   Mat            F,J;
 45:   ApplicationCtx ctx;
 46:   PetscInt       n=128,lag,its;
 47:   PetscBool      terse,flg,cct;

 50:   SlepcInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
 51:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 52:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%D\n\n",n);
 53:   ctx.h = 1.0/(PetscReal)n;
 54:   ctx.kappa = 1.0;

 56:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 57:         Create a standalone KSP with appropriate settings
 58:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 60:   KSPCreate(PETSC_COMM_WORLD,&ksp);
 61:   KSPSetType(ksp,KSPBCGS);
 62:   KSPGetPC(ksp,&pc);
 63:   PCSetType(pc,PCBJACOBI);
 64:   KSPSetFromOptions(ksp);

 66:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 67:                Prepare nonlinear eigensolver context
 68:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 70:   NEPCreate(PETSC_COMM_WORLD,&nep);

 72:   /* Create Function and Jacobian matrices; set evaluation routines */
 73:   MatCreate(PETSC_COMM_WORLD,&F);
 74:   MatSetSizes(F,PETSC_DECIDE,PETSC_DECIDE,n,n);
 75:   MatSetFromOptions(F);
 76:   MatSeqAIJSetPreallocation(F,3,NULL);
 77:   MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
 78:   MatSetUp(F);
 79:   NEPSetFunction(nep,F,F,FormFunction,&ctx);

 81:   MatCreate(PETSC_COMM_WORLD,&J);
 82:   MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,n,n);
 83:   MatSetFromOptions(J);
 84:   MatSeqAIJSetPreallocation(J,3,NULL);
 85:   MatMPIAIJSetPreallocation(F,3,NULL,1,NULL);
 86:   MatSetUp(J);
 87:   NEPSetJacobian(nep,J,FormJacobian,&ctx);

 89:   NEPSetType(nep,NEPRII);
 90:   NEPRIISetKSP(nep,ksp);
 91:   NEPRIISetMaximumIterations(nep,6);
 92:   NEPRIISetConstCorrectionTol(nep,PETSC_TRUE);
 93:   NEPSetFromOptions(nep);

 95:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 96:                       Solve the eigensystem
 97:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

 99:   NEPSolve(nep);
100:   PetscObjectTypeCompare((PetscObject)nep,NEPRII,&flg);
101:   if (flg) {
102:     NEPRIIGetMaximumIterations(nep,&its);
103:     NEPRIIGetLagPreconditioner(nep,&lag);
104:     NEPRIIGetConstCorrectionTol(nep,&cct);
105:     PetscPrintf(PETSC_COMM_WORLD," Maximum inner iterations of RII is %D\n",its);
106:     PetscPrintf(PETSC_COMM_WORLD," Preconditioner rebuilt every %D iterations\n",lag);
107:     if (cct) { PetscPrintf(PETSC_COMM_WORLD," Using a constant correction tolerance\n"); }
108:     PetscPrintf(PETSC_COMM_WORLD,"\n");
109:   }

111:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
112:                     Display solution and clean up
113:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

273: /*TEST

275:    test:
276:       suffix: 1
277:       args: -nep_target 21 -nep_rii_lag_preconditioner 2 -terse
278:       requires: !single

280: TEST*/