Actual source code: test4.c

slepc-3.20.1 2023-11-27
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:   PetscFunctionBeginUser;
 51:   PetscCall(SlepcInitialize(&argc,&argv,(char*)0,help));
 52:   PetscCall(PetscOptionsGetInt(NULL,NULL,"-n",&n,NULL));
 53:   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n1-D Nonlinear Eigenproblem, n=%" PetscInt_FMT "\n\n",n));
 54:   ctx.h = 1.0/(PetscReal)n;
 55:   ctx.kappa = 1.0;

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

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

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

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

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

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

 90:   PetscCall(NEPSetType(nep,NEPRII));
 91:   PetscCall(NEPRIISetKSP(nep,ksp));
 92:   PetscCall(NEPRIISetMaximumIterations(nep,6));
 93:   PetscCall(NEPRIISetConstCorrectionTol(nep,PETSC_TRUE));
 94:   PetscCall(NEPRIISetHermitian(nep,PETSC_TRUE));
 95:   PetscCall(NEPSetFromOptions(nep));

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

276: /*TEST

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

283: TEST*/