Actual source code: arnoldi.c

slepc-3.14.1 2020-12-08
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2020, 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 eigensolver: "arnoldi"

 13:    Method: Explicitly Restarted Arnoldi

 15:    Algorithm:

 17:        Arnoldi method with explicit restart and deflation.

 19:    References:

 21:        [1] "Arnoldi Methods in SLEPc", SLEPc Technical Report STR-4,
 22:            available at https://slepc.upv.es.
 23: */

 25: #include <slepc/private/epsimpl.h>

 27: typedef struct {
 28:   PetscBool delayed;
 29: } EPS_ARNOLDI;

 31: PetscErrorCode EPSSetUp_Arnoldi(EPS eps)
 32: {

 36:   EPSCheckDefinite(eps);
 37:   EPSSetDimensions_Default(eps,eps->nev,&eps->ncv,&eps->mpd);
 38:   if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd");
 39:   if (eps->max_it==PETSC_DEFAULT) eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
 40:   if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
 41:   if (eps->which==EPS_ALL) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support computing all eigenvalues");
 42:   EPSCheckUnsupported(eps,EPS_FEATURE_ARBITRARY | EPS_FEATURE_TWOSIDED);

 44:   EPSAllocateSolution(eps,1);
 45:   EPS_SetInnerProduct(eps);
 46:   DSSetType(eps->ds,DSNHEP);
 47:   if (eps->extraction==EPS_REFINED || eps->extraction==EPS_REFINED_HARMONIC) {
 48:     DSSetRefined(eps->ds,PETSC_TRUE);
 49:   }
 50:   DSSetExtraRow(eps->ds,PETSC_TRUE);
 51:   DSAllocate(eps->ds,eps->ncv+1);
 52:   return(0);
 53: }

 55: PetscErrorCode EPSSolve_Arnoldi(EPS eps)
 56: {
 57:   PetscErrorCode     ierr;
 58:   PetscInt           k,nv,ld;
 59:   Mat                U,Op;
 60:   PetscScalar        *H;
 61:   PetscReal          beta,gamma=1.0;
 62:   PetscBool          breakdown,harmonic,refined;
 63:   BVOrthogRefineType orthog_ref;
 64:   EPS_ARNOLDI        *arnoldi = (EPS_ARNOLDI*)eps->data;

 67:   DSGetLeadingDimension(eps->ds,&ld);
 68:   DSGetRefined(eps->ds,&refined);
 69:   harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
 70:   BVGetOrthogonalization(eps->V,NULL,&orthog_ref,NULL,NULL);

 72:   /* Get the starting Arnoldi vector */
 73:   EPSGetStartVector(eps,0,NULL);

 75:   /* Restart loop */
 76:   while (eps->reason == EPS_CONVERGED_ITERATING) {
 77:     eps->its++;

 79:     /* Compute an nv-step Arnoldi factorization */
 80:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
 81:     DSSetDimensions(eps->ds,nv,0,eps->nconv,0);
 82:     DSGetArray(eps->ds,DS_MAT_A,&H);
 83:     if (!arnoldi->delayed) {
 84:       STGetOperator(eps->st,&Op);
 85:       BVMatArnoldi(eps->V,Op,H,ld,eps->nconv,&nv,&beta,&breakdown);
 86:       STRestoreOperator(eps->st,&Op);
 87:     } else if (orthog_ref == BV_ORTHOG_REFINE_NEVER) {
 88:       EPSDelayedArnoldi1(eps,H,ld,eps->nconv,&nv,&beta,&breakdown);
 89:     } else {
 90:       EPSDelayedArnoldi(eps,H,ld,eps->nconv,&nv,&beta,&breakdown);
 91:     }
 92:     DSRestoreArray(eps->ds,DS_MAT_A,&H);
 93:     DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
 94:     BVSetActiveColumns(eps->V,eps->nconv,nv);

 96:     /* Compute translation of Krylov decomposition if harmonic extraction used */
 97:     if (harmonic) {
 98:       DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,NULL,&gamma);
 99:     }

101:     /* Solve projected problem */
102:     DSSolve(eps->ds,eps->eigr,eps->eigi);
103:     DSSort(eps->ds,eps->eigr,eps->eigi,NULL,NULL,NULL);
104:     DSUpdateExtraRow(eps->ds);
105:     DSSynchronize(eps->ds,eps->eigr,eps->eigi);

107:     /* Check convergence */
108:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,beta,0.0,gamma,&k);
109:     if (refined) {
110:       DSGetMat(eps->ds,DS_MAT_X,&U);
111:       BVMultInPlace(eps->V,U,eps->nconv,k+1);
112:       MatDestroy(&U);
113:       BVOrthonormalizeColumn(eps->V,k,PETSC_FALSE,NULL,NULL);
114:     } else {
115:       DSGetMat(eps->ds,DS_MAT_Q,&U);
116:       BVMultInPlace(eps->V,U,eps->nconv,PetscMin(k+1,nv));
117:       MatDestroy(&U);
118:     }
119:     (*eps->stopping)(eps,eps->its,eps->max_it,k,eps->nev,&eps->reason,eps->stoppingctx);
120:     if (eps->reason == EPS_CONVERGED_ITERATING && breakdown) {
121:       PetscInfo2(eps,"Breakdown in Arnoldi method (it=%D norm=%g)\n",eps->its,(double)beta);
122:       EPSGetStartVector(eps,k,&breakdown);
123:       if (breakdown) {
124:         eps->reason = EPS_DIVERGED_BREAKDOWN;
125:         PetscInfo(eps,"Unable to generate more start vectors\n");
126:       }
127:     }
128:     eps->nconv = k;
129:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
130:   }

132:   /* truncate Schur decomposition and change the state to raw so that
133:      DSVectors() computes eigenvectors from scratch */
134:   DSSetDimensions(eps->ds,eps->nconv,0,0,0);
135:   DSSetState(eps->ds,DS_STATE_RAW);
136:   return(0);
137: }

139: PetscErrorCode EPSSetFromOptions_Arnoldi(PetscOptionItems *PetscOptionsObject,EPS eps)
140: {
142:   PetscBool      set,val;
143:   EPS_ARNOLDI    *arnoldi = (EPS_ARNOLDI*)eps->data;

146:   PetscOptionsHead(PetscOptionsObject,"EPS Arnoldi Options");

148:     PetscOptionsBool("-eps_arnoldi_delayed","Use delayed reorthogonalization","EPSArnoldiSetDelayed",arnoldi->delayed,&val,&set);
149:     if (set) { EPSArnoldiSetDelayed(eps,val); }

151:   PetscOptionsTail();
152:   return(0);
153: }

155: static PetscErrorCode EPSArnoldiSetDelayed_Arnoldi(EPS eps,PetscBool delayed)
156: {
157:   EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data;

160:   arnoldi->delayed = delayed;
161:   return(0);
162: }

164: /*@
165:    EPSArnoldiSetDelayed - Activates or deactivates delayed reorthogonalization
166:    in the Arnoldi iteration.

168:    Logically Collective on eps

170:    Input Parameters:
171: +  eps - the eigenproblem solver context
172: -  delayed - boolean flag

174:    Options Database Key:
175: .  -eps_arnoldi_delayed - Activates delayed reorthogonalization in Arnoldi

177:    Note:
178:    Delayed reorthogonalization is an aggressive optimization for the Arnoldi
179:    eigensolver than may provide better scalability, but sometimes makes the
180:    solver converge less than the default algorithm.

182:    Level: advanced

184: .seealso: EPSArnoldiGetDelayed()
185: @*/
186: PetscErrorCode EPSArnoldiSetDelayed(EPS eps,PetscBool delayed)
187: {

193:   PetscTryMethod(eps,"EPSArnoldiSetDelayed_C",(EPS,PetscBool),(eps,delayed));
194:   return(0);
195: }

197: static PetscErrorCode EPSArnoldiGetDelayed_Arnoldi(EPS eps,PetscBool *delayed)
198: {
199:   EPS_ARNOLDI *arnoldi = (EPS_ARNOLDI*)eps->data;

202:   *delayed = arnoldi->delayed;
203:   return(0);
204: }

206: /*@
207:    EPSArnoldiGetDelayed - Gets the type of reorthogonalization used during the Arnoldi
208:    iteration.

210:    Not Collective

212:    Input Parameter:
213: .  eps - the eigenproblem solver context

215:    Input Parameter:
216: .  delayed - boolean flag indicating if delayed reorthogonalization has been enabled

218:    Level: advanced

220: .seealso: EPSArnoldiSetDelayed()
221: @*/
222: PetscErrorCode EPSArnoldiGetDelayed(EPS eps,PetscBool *delayed)
223: {

229:   PetscUseMethod(eps,"EPSArnoldiGetDelayed_C",(EPS,PetscBool*),(eps,delayed));
230:   return(0);
231: }

233: PetscErrorCode EPSDestroy_Arnoldi(EPS eps)
234: {

238:   PetscFree(eps->data);
239:   PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiSetDelayed_C",NULL);
240:   PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiGetDelayed_C",NULL);
241:   return(0);
242: }

244: PetscErrorCode EPSView_Arnoldi(EPS eps,PetscViewer viewer)
245: {
247:   PetscBool      isascii;
248:   EPS_ARNOLDI    *arnoldi = (EPS_ARNOLDI*)eps->data;

251:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
252:   if (isascii && arnoldi->delayed) {
253:     PetscViewerASCIIPrintf(viewer,"  using delayed reorthogonalization\n");
254:   }
255:   return(0);
256: }

258: SLEPC_EXTERN PetscErrorCode EPSCreate_Arnoldi(EPS eps)
259: {
260:   EPS_ARNOLDI    *ctx;

264:   PetscNewLog(eps,&ctx);
265:   eps->data = (void*)ctx;

267:   eps->useds = PETSC_TRUE;

269:   eps->ops->solve          = EPSSolve_Arnoldi;
270:   eps->ops->setup          = EPSSetUp_Arnoldi;
271:   eps->ops->setupsort      = EPSSetUpSort_Default;
272:   eps->ops->setfromoptions = EPSSetFromOptions_Arnoldi;
273:   eps->ops->destroy        = EPSDestroy_Arnoldi;
274:   eps->ops->view           = EPSView_Arnoldi;
275:   eps->ops->backtransform  = EPSBackTransform_Default;
276:   eps->ops->computevectors = EPSComputeVectors_Schur;

278:   PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiSetDelayed_C",EPSArnoldiSetDelayed_Arnoldi);
279:   PetscObjectComposeFunction((PetscObject)eps,"EPSArnoldiGetDelayed_C",EPSArnoldiGetDelayed_Arnoldi);
280:   return(0);
281: }