Actual source code: narnoldi.c

slepc-3.15.1 2021-05-28
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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: */
 10: /*
 11:    SLEPc nonlinear eigensolver: "narnoldi"

 13:    Method: Nonlinear Arnoldi

 15:    Algorithm:

 17:        Arnoldi for nonlinear eigenproblems.

 19:    References:

 21:        [1] H. Voss, "An Arnoldi method for nonlinear eigenvalue problems",
 22:            BIT 44:387-401, 2004.
 23: */

 25: #include <slepc/private/nepimpl.h>
 26: #include <../src/nep/impls/nepdefl.h>

 28: typedef struct {
 29:   PetscInt lag;             /* interval to rebuild preconditioner */
 30:   KSP      ksp;             /* linear solver object */
 31: } NEP_NARNOLDI;

 33: PetscErrorCode NEPSetUp_NArnoldi(NEP nep)
 34: {

 38:   NEPSetDimensions_Default(nep,nep->nev,&nep->ncv,&nep->mpd);
 39:   if (nep->ncv>nep->nev+nep->mpd) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must not be larger than nev+mpd");
 40:   if (nep->max_it==PETSC_DEFAULT) nep->max_it = nep->nev*nep->ncv;
 41:   if (!nep->which) nep->which = NEP_TARGET_MAGNITUDE;
 42:   if (nep->which!=NEP_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver supports only target magnitude eigenvalues");
 43:   NEPCheckUnsupported(nep,NEP_FEATURE_CALLBACK | NEP_FEATURE_REGION | NEP_FEATURE_TWOSIDED);
 44:   NEPAllocateSolution(nep,0);
 45:   NEPSetWorkVecs(nep,3);
 46:   return(0);
 47: }

 49: PetscErrorCode NEPSolve_NArnoldi(NEP nep)
 50: {
 51:   PetscErrorCode     ierr;
 52:   NEP_NARNOLDI       *ctx = (NEP_NARNOLDI*)nep->data;
 53:   Mat                T,H;
 54:   Vec                f,r,u,uu;
 55:   PetscScalar        *X,lambda=0.0,lambda2=0.0,*eigr,*Ap,sigma;
 56:   const PetscScalar  *Hp;
 57:   PetscReal          beta,resnorm=0.0,nrm,perr=0.0;
 58:   PetscInt           n,i,j,ldds,ldh;
 59:   PetscBool          breakdown,skip=PETSC_FALSE;
 60:   BV                 Vext;
 61:   DS                 ds;
 62:   NEP_EXT_OP         extop=NULL;
 63:   SlepcSC            sc;
 64:   KSPConvergedReason kspreason;

 67:   /* get initial space and shift */
 68:   NEPGetDefaultShift(nep,&sigma);
 69:   if (!nep->nini) {
 70:     BVSetRandomColumn(nep->V,0);
 71:     BVNormColumn(nep->V,0,NORM_2,&nrm);
 72:     BVScaleColumn(nep->V,0,1.0/nrm);
 73:     n = 1;
 74:   } else n = nep->nini;

 76:   if (!ctx->ksp) { NEPNArnoldiGetKSP(nep,&ctx->ksp); }
 77:   NEPDeflationInitialize(nep,nep->V,ctx->ksp,PETSC_FALSE,nep->nev,&extop);
 78:   NEPDeflationCreateBV(extop,nep->ncv,&Vext);

 80:   /* prepare linear solver */
 81:   NEPDeflationSolveSetUp(extop,sigma);

 83:   BVGetColumn(Vext,0,&f);
 84:   VecDuplicate(f,&r);
 85:   VecDuplicate(f,&u);
 86:   BVGetColumn(nep->V,0,&uu);
 87:   NEPDeflationCopyToExtendedVec(extop,uu,NULL,f,PETSC_FALSE);
 88:   BVRestoreColumn(nep->V,0,&uu);
 89:   VecCopy(f,r);
 90:   NEPDeflationFunctionSolve(extop,r,f);
 91:   VecNorm(f,NORM_2,&nrm);
 92:   VecScale(f,1.0/nrm);
 93:   BVRestoreColumn(Vext,0,&f);

 95:   DSCreate(PetscObjectComm((PetscObject)nep),&ds);
 96:   PetscLogObjectParent((PetscObject)nep,(PetscObject)ds);
 97:   DSSetType(ds,DSNEP);
 98:   DSNEPSetFN(ds,nep->nt,nep->f);
 99:   DSAllocate(ds,nep->ncv);
100:   DSGetSlepcSC(ds,&sc);
101:   sc->comparison    = nep->sc->comparison;
102:   sc->comparisonctx = nep->sc->comparisonctx;
103:   DSSetFromOptions(ds);

105:   /* build projected matrices for initial space */
106:   DSSetDimensions(ds,n,0,0,0);
107:   NEPDeflationProjectOperator(extop,Vext,ds,0,n);

109:   PetscMalloc1(nep->ncv,&eigr);

111:   /* Restart loop */
112:   while (nep->reason == NEP_CONVERGED_ITERATING) {
113:     nep->its++;

115:     /* solve projected problem */
116:     DSSetDimensions(ds,n,0,0,0);
117:     DSSetState(ds,DS_STATE_RAW);
118:     DSSolve(ds,eigr,NULL);
119:     DSSynchronize(ds,eigr,NULL);
120:     if (nep->its>1) lambda2 = lambda;
121:     lambda = eigr[0];
122:     nep->eigr[nep->nconv] = lambda;

124:     /* compute Ritz vector, x = V*s */
125:     DSGetArray(ds,DS_MAT_X,&X);
126:     BVSetActiveColumns(Vext,0,n);
127:     BVMultVec(Vext,1.0,0.0,u,X);
128:     DSRestoreArray(ds,DS_MAT_X,&X);

130:     /* compute the residual, r = T(lambda)*x */
131:     NEPDeflationComputeFunction(extop,lambda,&T);
132:     MatMult(T,u,r);

134:     /* convergence test */
135:     VecNorm(r,NORM_2,&resnorm);
136:     if (nep->its>1) perr=nep->errest[nep->nconv];
137:     (*nep->converged)(nep,lambda,0,resnorm,&nep->errest[nep->nconv],nep->convergedctx);
138:     if (nep->errest[nep->nconv]<=nep->tol) {
139:       nep->nconv = nep->nconv + 1;
140:       NEPDeflationLocking(extop,u,lambda);
141:       skip = PETSC_TRUE;
142:     }
143:     (*nep->stopping)(nep,nep->its,nep->max_it,nep->nconv,nep->nev,&nep->reason,nep->stoppingctx);
144:     if (!skip || nep->reason>0) {
145:       NEPMonitor(nep,nep->its,nep->nconv,nep->eigr,nep->eigi,nep->errest,(nep->reason>0)?nep->nconv:nep->nconv+1);
146:     }

148:     if (nep->reason == NEP_CONVERGED_ITERATING) {
149:       if (!skip) {
150:         if (n>=nep->ncv) {
151:           nep->reason = NEP_DIVERGED_SUBSPACE_EXHAUSTED;
152:           break;
153:         }
154:         if (ctx->lag && !(nep->its%ctx->lag) && nep->its>=2*ctx->lag && perr && nep->errest[nep->nconv]>.5*perr) {
155:           NEPDeflationSolveSetUp(extop,lambda2);
156:         }

158:         /* continuation vector: f = T(sigma)\r */
159:         BVGetColumn(Vext,n,&f);
160:         NEPDeflationFunctionSolve(extop,r,f);
161:         BVRestoreColumn(Vext,n,&f);
162:         KSPGetConvergedReason(ctx->ksp,&kspreason);
163:         if (kspreason<0) {
164:           PetscInfo1(nep,"iter=%D, linear solve failed, stopping solve\n",nep->its);
165:           nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
166:           break;
167:         }

169:         /* orthonormalize */
170:         BVOrthonormalizeColumn(Vext,n,PETSC_FALSE,&beta,&breakdown);
171:         if (breakdown || beta==0.0) {
172:           PetscInfo1(nep,"iter=%D, orthogonalization failed, stopping solve\n",nep->its);
173:           nep->reason = NEP_DIVERGED_BREAKDOWN;
174:           break;
175:         }

177:         /* update projected matrices */
178:         DSSetDimensions(ds,n+1,0,0,0);
179:         NEPDeflationProjectOperator(extop,Vext,ds,n,n+1);
180:         n++;
181:       } else {
182:         nep->its--;  /* do not count this as a full iteration */
183:         BVGetColumn(Vext,0,&f);
184:         NEPDeflationSetRandomVec(extop,f);
185:         NEPDeflationSolveSetUp(extop,sigma);
186:         VecCopy(f,r);
187:         NEPDeflationFunctionSolve(extop,r,f);
188:         VecNorm(f,NORM_2,&nrm);
189:         VecScale(f,1.0/nrm);
190:         BVRestoreColumn(Vext,0,&f);
191:         n = 1;
192:         DSSetDimensions(ds,n,0,0,0);
193:         NEPDeflationProjectOperator(extop,Vext,ds,n-1,n);
194:         skip = PETSC_FALSE;
195:       }
196:     }
197:   }

199:   NEPDeflationGetInvariantPair(extop,NULL,&H);
200:   MatGetSize(H,NULL,&ldh);
201:   DSSetType(nep->ds,DSNHEP);
202:   DSAllocate(nep->ds,PetscMax(nep->nconv,1));
203:   DSGetLeadingDimension(nep->ds,&ldds);
204:   MatDenseGetArrayRead(H,&Hp);
205:   DSGetArray(nep->ds,DS_MAT_A,&Ap);
206:   for (j=0;j<nep->nconv;j++)
207:     for (i=0;i<nep->nconv;i++) Ap[j*ldds+i] = Hp[j*ldh+i];
208:   DSRestoreArray(nep->ds,DS_MAT_A,&Ap);
209:   MatDenseRestoreArrayRead(H,&Hp);
210:   MatDestroy(&H);
211:   DSSetDimensions(nep->ds,nep->nconv,0,0,nep->nconv);
212:   DSSolve(nep->ds,nep->eigr,nep->eigi);
213:   NEPDeflationReset(extop);
214:   VecDestroy(&u);
215:   VecDestroy(&r);
216:   BVDestroy(&Vext);
217:   DSDestroy(&ds);
218:   PetscFree(eigr);
219:   return(0);
220: }

222: static PetscErrorCode NEPNArnoldiSetLagPreconditioner_NArnoldi(NEP nep,PetscInt lag)
223: {
224:   NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;

227:   if (lag<0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Lag must be non-negative");
228:   ctx->lag = lag;
229:   return(0);
230: }

232: /*@
233:    NEPNArnoldiSetLagPreconditioner - Determines when the preconditioner is rebuilt in the
234:    nonlinear solve.

236:    Logically Collective on nep

238:    Input Parameters:
239: +  nep - nonlinear eigenvalue solver
240: -  lag - 0 indicates NEVER rebuild, 1 means rebuild every time the Jacobian is
241:           computed within the nonlinear iteration, 2 means every second time
242:           the Jacobian is built, etc.

244:    Options Database Keys:
245: .  -nep_narnoldi_lag_preconditioner <lag>

247:    Notes:
248:    The default is 1.
249:    The preconditioner is ALWAYS built in the first iteration of a nonlinear solve.

251:    Level: intermediate

253: .seealso: NEPNArnoldiGetLagPreconditioner()
254: @*/
255: PetscErrorCode NEPNArnoldiSetLagPreconditioner(NEP nep,PetscInt lag)
256: {

262:   PetscTryMethod(nep,"NEPNArnoldiSetLagPreconditioner_C",(NEP,PetscInt),(nep,lag));
263:   return(0);
264: }

266: static PetscErrorCode NEPNArnoldiGetLagPreconditioner_NArnoldi(NEP nep,PetscInt *lag)
267: {
268:   NEP_NARNOLDI *ctx = (NEP_NARNOLDI*)nep->data;

271:   *lag = ctx->lag;
272:   return(0);
273: }

275: /*@
276:    NEPNArnoldiGetLagPreconditioner - Indicates how often the preconditioner is rebuilt.

278:    Not Collective

280:    Input Parameter:
281: .  nep - nonlinear eigenvalue solver

283:    Output Parameter:
284: .  lag - the lag parameter

286:    Level: intermediate

288: .seealso: NEPNArnoldiSetLagPreconditioner()
289: @*/
290: PetscErrorCode NEPNArnoldiGetLagPreconditioner(NEP nep,PetscInt *lag)
291: {

297:   PetscUseMethod(nep,"NEPNArnoldiGetLagPreconditioner_C",(NEP,PetscInt*),(nep,lag));
298:   return(0);
299: }

301: PetscErrorCode NEPSetFromOptions_NArnoldi(PetscOptionItems *PetscOptionsObject,NEP nep)
302: {
304:   PetscInt       i;
305:   PetscBool      flg;
306:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;

309:   PetscOptionsHead(PetscOptionsObject,"NEP N-Arnoldi Options");
310:     i = 0;
311:     PetscOptionsInt("-nep_narnoldi_lag_preconditioner","Interval to rebuild preconditioner","NEPNArnoldiSetLagPreconditioner",ctx->lag,&i,&flg);
312:     if (flg) { NEPNArnoldiSetLagPreconditioner(nep,i); }

314:   PetscOptionsTail();

316:   if (!ctx->ksp) { NEPNArnoldiGetKSP(nep,&ctx->ksp); }
317:   KSPSetFromOptions(ctx->ksp);
318:   return(0);
319: }

321: static PetscErrorCode NEPNArnoldiSetKSP_NArnoldi(NEP nep,KSP ksp)
322: {
324:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;

327:   PetscObjectReference((PetscObject)ksp);
328:   KSPDestroy(&ctx->ksp);
329:   ctx->ksp = ksp;
330:   PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp);
331:   nep->state = NEP_STATE_INITIAL;
332:   return(0);
333: }

335: /*@
336:    NEPNArnoldiSetKSP - Associate a linear solver object (KSP) to the nonlinear
337:    eigenvalue solver.

339:    Collective on nep

341:    Input Parameters:
342: +  nep - eigenvalue solver
343: -  ksp - the linear solver object

345:    Level: advanced

347: .seealso: NEPNArnoldiGetKSP()
348: @*/
349: PetscErrorCode NEPNArnoldiSetKSP(NEP nep,KSP ksp)
350: {

357:   PetscTryMethod(nep,"NEPNArnoldiSetKSP_C",(NEP,KSP),(nep,ksp));
358:   return(0);
359: }

361: static PetscErrorCode NEPNArnoldiGetKSP_NArnoldi(NEP nep,KSP *ksp)
362: {
364:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;

367:   if (!ctx->ksp) {
368:     KSPCreate(PetscObjectComm((PetscObject)nep),&ctx->ksp);
369:     PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)nep,1);
370:     KSPSetOptionsPrefix(ctx->ksp,((PetscObject)nep)->prefix);
371:     KSPAppendOptionsPrefix(ctx->ksp,"nep_narnoldi_");
372:     PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp);
373:     PetscObjectSetOptions((PetscObject)ctx->ksp,((PetscObject)nep)->options);
374:     KSPSetErrorIfNotConverged(ctx->ksp,PETSC_TRUE);
375:     KSPSetTolerances(ctx->ksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
376:   }
377:   *ksp = ctx->ksp;
378:   return(0);
379: }

381: /*@
382:    NEPNArnoldiGetKSP - Retrieve the linear solver object (KSP) associated with
383:    the nonlinear eigenvalue solver.

385:    Not Collective

387:    Input Parameter:
388: .  nep - nonlinear eigenvalue solver

390:    Output Parameter:
391: .  ksp - the linear solver object

393:    Level: advanced

395: .seealso: NEPNArnoldiSetKSP()
396: @*/
397: PetscErrorCode NEPNArnoldiGetKSP(NEP nep,KSP *ksp)
398: {

404:   PetscUseMethod(nep,"NEPNArnoldiGetKSP_C",(NEP,KSP*),(nep,ksp));
405:   return(0);
406: }

408: PetscErrorCode NEPView_NArnoldi(NEP nep,PetscViewer viewer)
409: {
411:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;
412:   PetscBool      isascii;

415:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
416:   if (isascii) {
417:     if (ctx->lag) {
418:       PetscViewerASCIIPrintf(viewer,"  updating the preconditioner every %D iterations\n",ctx->lag);
419:     }
420:     if (!ctx->ksp) { NEPNArnoldiGetKSP(nep,&ctx->ksp); }
421:     PetscViewerASCIIPushTab(viewer);
422:     KSPView(ctx->ksp,viewer);
423:     PetscViewerASCIIPopTab(viewer);
424:   }
425:   return(0);
426: }

428: PetscErrorCode NEPReset_NArnoldi(NEP nep)
429: {
431:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;

434:   KSPReset(ctx->ksp);
435:   return(0);
436: }

438: PetscErrorCode NEPDestroy_NArnoldi(NEP nep)
439: {
441:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;

444:   KSPDestroy(&ctx->ksp);
445:   PetscFree(nep->data);
446:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetLagPreconditioner_C",NULL);
447:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetLagPreconditioner_C",NULL);
448:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetKSP_C",NULL);
449:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetKSP_C",NULL);
450:   return(0);
451: }

453: SLEPC_EXTERN PetscErrorCode NEPCreate_NArnoldi(NEP nep)
454: {
456:   NEP_NARNOLDI   *ctx;

459:   PetscNewLog(nep,&ctx);
460:   nep->data = (void*)ctx;
461:   ctx->lag  = 1;

463:   nep->useds = PETSC_TRUE;

465:   nep->ops->solve          = NEPSolve_NArnoldi;
466:   nep->ops->setup          = NEPSetUp_NArnoldi;
467:   nep->ops->setfromoptions = NEPSetFromOptions_NArnoldi;
468:   nep->ops->reset          = NEPReset_NArnoldi;
469:   nep->ops->destroy        = NEPDestroy_NArnoldi;
470:   nep->ops->view           = NEPView_NArnoldi;
471:   nep->ops->computevectors = NEPComputeVectors_Schur;

473:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetLagPreconditioner_C",NEPNArnoldiSetLagPreconditioner_NArnoldi);
474:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetLagPreconditioner_C",NEPNArnoldiGetLagPreconditioner_NArnoldi);
475:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetKSP_C",NEPNArnoldiSetKSP_NArnoldi);
476:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetKSP_C",NEPNArnoldiGetKSP_NArnoldi);
477:   return(0);
478: }