Actual source code: narnoldi.c

slepc-3.9.1 2018-05-02
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2018, 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>         /*I "slepcnep.h" I*/

 27: typedef struct {
 28:   KSP      ksp;              /* linear solver object */
 29: } NEP_NARNOLDI;

 31: PETSC_STATIC_INLINE PetscErrorCode NEPNArnoldi_KSPSolve(NEP nep,Vec b,Vec x)
 32: {
 34:   PetscInt       lits;
 35:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;

 38:   KSPSolve(ctx->ksp,b,x);
 39:   KSPGetIterationNumber(ctx->ksp,&lits);
 40:   PetscInfo2(nep,"iter=%D, linear solve iterations=%D\n",nep->its,lits);
 41:   return(0);
 42: }

 44: PetscErrorCode NEPSetUp_NArnoldi(NEP nep)
 45: {
 47:   PetscBool      istrivial;

 50:   NEPSetDimensions_Default(nep,nep->nev,&nep->ncv,&nep->mpd);
 51:   if (nep->ncv>nep->nev+nep->mpd) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must not be larger than nev+mpd");
 52:   if (nep->nev>1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Requested several eigenpairs but this solver can compute only one");
 53:   if (!nep->max_it) nep->max_it = nep->ncv;
 54:   if (nep->max_it < nep->ncv) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Current implementation is unrestarted, must set max_it >= ncv");
 55:   if (nep->which && nep->which!=NEP_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Wrong value of which");
 56:   if (nep->fui!=NEP_USER_INTERFACE_SPLIT) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NARNOLDI only available for split operator");

 58:   RGIsTrivial(nep->rg,&istrivial);
 59:   if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver does not support region filtering");

 61:   NEPAllocateSolution(nep,0);
 62:   NEPSetWorkVecs(nep,3);

 64:   /* set-up DS and transfer split operator functions */
 65:   DSSetType(nep->ds,DSNEP);
 66:   DSNEPSetFN(nep->ds,nep->nt,nep->f);
 67:   DSAllocate(nep->ds,nep->ncv);
 68:   return(0);
 69: }

 71: PetscErrorCode NEPSolve_NArnoldi(NEP nep)
 72: {
 73:   PetscErrorCode     ierr;
 74:   NEP_NARNOLDI       *ctx = (NEP_NARNOLDI*)nep->data;
 75:   Mat                T=nep->function,Tsigma;
 76:   Vec                f,r=nep->work[0],x=nep->work[1],w=nep->work[2];
 77:   PetscScalar        *X,lambda;
 78:   PetscReal          beta,resnorm=0.0,nrm;
 79:   PetscInt           n;
 80:   PetscBool          breakdown;
 81:   KSPConvergedReason kspreason;

 84:   /* get initial space and shift */
 85:   NEPGetDefaultShift(nep,&lambda);
 86:   if (!nep->nini) {
 87:     BVSetRandomColumn(nep->V,0);
 88:     BVNormColumn(nep->V,0,NORM_2,&nrm);
 89:     BVScaleColumn(nep->V,0,1.0/nrm);
 90:     n = 1;
 91:   } else n = nep->nini;

 93:   /* build projected matrices for initial space */
 94:   DSSetDimensions(nep->ds,n,0,0,0);
 95:   NEPProjectOperator(nep,0,n);

 97:   /* prepare linear solver */
 98:   if (!ctx->ksp) { NEPNArnoldiGetKSP(nep,&ctx->ksp); }
 99:   NEPComputeFunction(nep,lambda,T,T);
100:   MatDuplicate(T,MAT_COPY_VALUES,&Tsigma);
101:   KSPSetOperators(ctx->ksp,Tsigma,Tsigma);

103:   /* Restart loop */
104:   while (nep->reason == NEP_CONVERGED_ITERATING) {
105:     nep->its++;

107:     /* solve projected problem */
108:     DSSetDimensions(nep->ds,n,0,0,0);
109:     DSSetState(nep->ds,DS_STATE_RAW);
110:     DSSolve(nep->ds,nep->eigr,NULL);
111:     DSSynchronize(nep->ds,nep->eigr,NULL);
112:     lambda = nep->eigr[0];

114:     /* compute Ritz vector, x = V*s */
115:     DSGetArray(nep->ds,DS_MAT_X,&X);
116:     BVSetActiveColumns(nep->V,0,n);
117:     BVMultVec(nep->V,1.0,0.0,x,X);
118:     DSRestoreArray(nep->ds,DS_MAT_X,&X);

120:     /* compute the residual, r = T(lambda)*x */
121:     NEPApplyFunction(nep,lambda,x,w,r,NULL,NULL);

123:     /* convergence test */
124:     VecNorm(r,NORM_2,&resnorm);
125:     (*nep->converged)(nep,lambda,0,resnorm,&nep->errest[nep->nconv],nep->convergedctx);
126:     if (nep->errest[nep->nconv]<=nep->tol) {
127:       BVInsertVec(nep->V,nep->nconv,x);
128:       nep->nconv = nep->nconv + 1;
129:     }
130:     (*nep->stopping)(nep,nep->its,nep->max_it,nep->nconv,nep->nev,&nep->reason,nep->stoppingctx);
131:     NEPMonitor(nep,nep->its,nep->nconv,nep->eigr,nep->eigi,nep->errest,1);

133:     if (nep->reason == NEP_CONVERGED_ITERATING) {

135:       /* continuation vector: f = T(sigma)\r */
136:       BVGetColumn(nep->V,n,&f);
137:       NEPNArnoldi_KSPSolve(nep,r,f);
138:       BVRestoreColumn(nep->V,n,&f);
139:       KSPGetConvergedReason(ctx->ksp,&kspreason);
140:       if (kspreason<0) {
141:         PetscInfo1(nep,"iter=%D, linear solve failed, stopping solve\n",nep->its);
142:         nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
143:         break;
144:       }

146:       /* orthonormalize */
147:       BVOrthonormalizeColumn(nep->V,n,PETSC_FALSE,&beta,&breakdown);
148:       if (breakdown || beta==0.0) {
149:         PetscInfo1(nep,"iter=%D, orthogonalization failed, stopping solve\n",nep->its);
150:         nep->reason = NEP_DIVERGED_BREAKDOWN;
151:         break;
152:       }

154:       /* update projected matrices */
155:       DSSetDimensions(nep->ds,n+1,0,0,0);
156:       NEPProjectOperator(nep,n,n+1);
157:       n++;
158:     }
159:   }
160:   MatDestroy(&Tsigma);
161:   return(0);
162: }

164: PetscErrorCode NEPSetFromOptions_NArnoldi(PetscOptionItems *PetscOptionsObject,NEP nep)
165: {
167:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;

170:   PetscOptionsHead(PetscOptionsObject,"NEP N-Arnoldi Options");
171:   PetscOptionsTail();

173:   if (!ctx->ksp) { NEPNArnoldiGetKSP(nep,&ctx->ksp); }
174:   KSPSetFromOptions(ctx->ksp);
175:   return(0);
176: }

178: static PetscErrorCode NEPNArnoldiSetKSP_NArnoldi(NEP nep,KSP ksp)
179: {
181:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;

184:   PetscObjectReference((PetscObject)ksp);
185:   KSPDestroy(&ctx->ksp);
186:   ctx->ksp = ksp;
187:   PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp);
188:   nep->state = NEP_STATE_INITIAL;
189:   return(0);
190: }

192: /*@
193:    NEPNArnoldiSetKSP - Associate a linear solver object (KSP) to the nonlinear
194:    eigenvalue solver.

196:    Collective on NEP

198:    Input Parameters:
199: +  nep - eigenvalue solver
200: -  ksp - the linear solver object

202:    Level: advanced

204: .seealso: NEPNArnoldiGetKSP()
205: @*/
206: PetscErrorCode NEPNArnoldiSetKSP(NEP nep,KSP ksp)
207: {

214:   PetscTryMethod(nep,"NEPNArnoldiSetKSP_C",(NEP,KSP),(nep,ksp));
215:   return(0);
216: }

218: static PetscErrorCode NEPNArnoldiGetKSP_NArnoldi(NEP nep,KSP *ksp)
219: {
221:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;

224:   if (!ctx->ksp) {
225:     KSPCreate(PetscObjectComm((PetscObject)nep),&ctx->ksp);
226:     PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)nep,1);
227:     KSPSetOptionsPrefix(ctx->ksp,((PetscObject)nep)->prefix);
228:     KSPAppendOptionsPrefix(ctx->ksp,"nep_narnoldi_");
229:     PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp);
230:     KSPSetErrorIfNotConverged(ctx->ksp,PETSC_TRUE);
231:     KSPSetTolerances(ctx->ksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
232:   }
233:   *ksp = ctx->ksp;
234:   return(0);
235: }

237: /*@
238:    NEPNArnoldiGetKSP - Retrieve the linear solver object (KSP) associated with
239:    the nonlinear eigenvalue solver.

241:    Not Collective

243:    Input Parameter:
244: .  nep - nonlinear eigenvalue solver

246:    Output Parameter:
247: .  ksp - the linear solver object

249:    Level: advanced

251: .seealso: NEPNArnoldiSetKSP()
252: @*/
253: PetscErrorCode NEPNArnoldiGetKSP(NEP nep,KSP *ksp)
254: {

260:   PetscUseMethod(nep,"NEPNArnoldiGetKSP_C",(NEP,KSP*),(nep,ksp));
261:   return(0);
262: }

264: PetscErrorCode NEPView_NArnoldi(NEP nep,PetscViewer viewer)
265: {
267:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;
268:   PetscBool      isascii;

271:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
272:   if (isascii) {
273:     if (!ctx->ksp) { NEPNArnoldiGetKSP(nep,&ctx->ksp); }
274:     KSPView(ctx->ksp,viewer);
275:   }
276:   return(0);
277: }

279: PetscErrorCode NEPReset_NArnoldi(NEP nep)
280: {
282:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;

285:   KSPReset(ctx->ksp);
286:   return(0);
287: }

289: PetscErrorCode NEPDestroy_NArnoldi(NEP nep)
290: {
292:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;

295:   KSPDestroy(&ctx->ksp);
296:   PetscFree(nep->data);
297:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetKSP_C",NULL);
298:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetKSP_C",NULL);
299:   return(0);
300: }

302: PETSC_EXTERN PetscErrorCode NEPCreate_NArnoldi(NEP nep)
303: {
305:   NEP_NARNOLDI   *ctx;

308:   PetscNewLog(nep,&ctx);
309:   nep->data = (void*)ctx;

311:   nep->ops->solve          = NEPSolve_NArnoldi;
312:   nep->ops->setup          = NEPSetUp_NArnoldi;
313:   nep->ops->setfromoptions = NEPSetFromOptions_NArnoldi;
314:   nep->ops->reset          = NEPReset_NArnoldi;
315:   nep->ops->destroy        = NEPDestroy_NArnoldi;
316:   nep->ops->view           = NEPView_NArnoldi;

318:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetKSP_C",NEPNArnoldiSetKSP_NArnoldi);
319:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetKSP_C",NEPNArnoldiGetKSP_NArnoldi);
320:   return(0);
321: }