Actual source code: test4.c

slepc-3.17.2 2022-08-09
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 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,herm;
 48:   PetscReal      thres;

 50:   SlepcInitialize(&argc,&argv,(char*)0,help);
 51:   PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL);
 52:   PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%" PetscInt_FMT "\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:   NEPRIISetHermitian(nep,PETSC_TRUE);
 94:   NEPSetFromOptions(nep);

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

100:   NEPSolve(nep);
101:   PetscObjectTypeCompare((PetscObject)nep,NEPRII,&flg);
102:   if (flg) {
103:     NEPRIIGetMaximumIterations(nep,&its);
104:     NEPRIIGetLagPreconditioner(nep,&lag);
105:     NEPRIIGetDeflationThreshold(nep,&thres);
106:     NEPRIIGetConstCorrectionTol(nep,&cct);
107:     NEPRIIGetHermitian(nep,&herm);
108:     PetscPrintf(PETSC_COMM_WORLD," Maximum inner iterations of RII is %" PetscInt_FMT "\n",its);
109:     PetscPrintf(PETSC_COMM_WORLD," Preconditioner rebuilt every %" PetscInt_FMT " iterations\n",lag);
110:     if (thres>0.0) PetscPrintf(PETSC_COMM_WORLD," Using deflation threshold=%g\n",(double)thres);
111:     if (cct) PetscPrintf(PETSC_COMM_WORLD," Using a constant correction tolerance\n");
112:     if (herm) PetscPrintf(PETSC_COMM_WORLD," Hermitian version of scalar equation\n");
113:     PetscPrintf(PETSC_COMM_WORLD,"\n");
114:   }

116:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
117:                     Display solution and clean up
118:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

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

130:   NEPDestroy(&nep);
131:   KSPDestroy(&ksp);
132:   MatDestroy(&F);
133:   MatDestroy(&J);
134:   SlepcFinalize();
135:   return 0;
136: }

138: /* ------------------------------------------------------------------- */
139: /*
140:    FormFunction - Computes Function matrix  T(lambda)

142:    Input Parameters:
143: .  nep    - the NEP context
144: .  lambda - the scalar argument
145: .  ctx    - optional user-defined context, as set by NEPSetFunction()

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

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

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

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

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

197:   /*
198:      Assemble matrix
199:   */
200:   MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
201:   MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
202:   if (fun != B) {
203:     MatAssemblyBegin(fun,MAT_FINAL_ASSEMBLY);
204:     MatAssemblyEnd(fun,MAT_FINAL_ASSEMBLY);
205:   }
206:   PetscFunctionReturn(0);
207: }

209: /* ------------------------------------------------------------------- */
210: /*
211:    FormJacobian - Computes Jacobian matrix  T'(lambda)

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

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

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

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

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

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

267:   /*
268:      Assemble matrix
269:   */
270:   MatAssemblyBegin(jac,MAT_FINAL_ASSEMBLY);
271:   MatAssemblyEnd(jac,MAT_FINAL_ASSEMBLY);
272:   PetscFunctionReturn(0);
273: }

275: /*TEST

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

282: TEST*/