Actual source code: narnoldi.c
 
   slepc-3.7.4 2017-05-17
   
  1: /*
  3:    SLEPc nonlinear eigensolver: "narnoldi"
  5:    Method: Nonlinear Arnoldi
  7:    Algorithm:
  9:        Arnoldi for nonlinear eigenproblems.
 11:    References:
 13:        [1] H. Voss, "An Arnoldi method for nonlinear eigenvalue problems",
 14:            BIT 44:387-401, 2004.
 16:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 17:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 18:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain
 20:    This file is part of SLEPc.
 22:    SLEPc is free software: you can redistribute it and/or modify it under  the
 23:    terms of version 3 of the GNU Lesser General Public License as published by
 24:    the Free Software Foundation.
 26:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 27:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 28:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 29:    more details.
 31:    You  should have received a copy of the GNU Lesser General  Public  License
 32:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 33:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 34: */
 36: #include <slepc/private/nepimpl.h>         /*I "slepcnep.h" I*/
 38: typedef struct {
 39:   KSP      ksp;              /* linear solver object */
 40: } NEP_NARNOLDI;
 44: PETSC_STATIC_INLINE PetscErrorCode NEPNArnoldi_KSPSolve(NEP nep,Vec b,Vec x)
 45: {
 47:   PetscInt       lits;
 48:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;
 51:   KSPSolve(ctx->ksp,b,x);
 52:   KSPGetIterationNumber(ctx->ksp,&lits);
 53:   PetscInfo2(nep,"iter=%D, linear solve iterations=%D\n",nep->its,lits);
 54:   return(0);
 55: }
 59: PetscErrorCode NEPSetUp_NArnoldi(NEP nep)
 60: {
 62:   PetscBool      istrivial;
 65:   NEPSetDimensions_Default(nep,nep->nev,&nep->ncv,&nep->mpd);
 66:   if (nep->ncv>nep->nev+nep->mpd) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must not be larger than nev+mpd");
 67:   if (nep->nev>1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Requested several eigenpairs but this solver can compute only one");
 68:   if (!nep->max_it) nep->max_it = nep->ncv;
 69:   if (nep->max_it < nep->ncv) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Current implementation is unrestarted, must set max_it >= ncv");
 70:   if (nep->which && nep->which!=NEP_TARGET_MAGNITUDE) SETERRQ(PetscObjectComm((PetscObject)nep),1,"Wrong value of which");
 71:   if (nep->fui!=NEP_USER_INTERFACE_SPLIT) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NARNOLDI only available for split operator");
 73:   RGIsTrivial(nep->rg,&istrivial);
 74:   if (!istrivial) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver does not support region filtering");
 76:   NEPAllocateSolution(nep,0);
 77:   NEPSetWorkVecs(nep,3);
 79:   /* set-up DS and transfer split operator functions */
 80:   DSSetType(nep->ds,DSNEP);
 81:   DSNEPSetFN(nep->ds,nep->nt,nep->f);
 82:   DSAllocate(nep->ds,nep->ncv);
 83:   return(0);
 84: }
 88: PetscErrorCode NEPSolve_NArnoldi(NEP nep)
 89: {
 90:   PetscErrorCode     ierr;
 91:   NEP_NARNOLDI       *ctx = (NEP_NARNOLDI*)nep->data;
 92:   Mat                T=nep->function,Tsigma;
 93:   Vec                f,r=nep->work[0],x=nep->work[1],w=nep->work[2];
 94:   PetscScalar        *X,lambda;
 95:   PetscReal          beta,resnorm=0.0,nrm;
 96:   PetscInt           n;
 97:   PetscBool          breakdown;
 98:   KSPConvergedReason kspreason;
101:   /* get initial space and shift */
102:   NEPGetDefaultShift(nep,&lambda);
103:   if (!nep->nini) {
104:     BVSetRandomColumn(nep->V,0);
105:     BVNormColumn(nep->V,0,NORM_2,&nrm);
106:     BVScaleColumn(nep->V,0,1.0/nrm);
107:     n = 1;
108:   } else n = nep->nini;
110:   /* build projected matrices for initial space */
111:   DSSetDimensions(nep->ds,n,0,0,0);
112:   NEPProjectOperator(nep,0,n);
114:   /* prepare linear solver */
115:   if (!ctx->ksp) { NEPNArnoldiGetKSP(nep,&ctx->ksp); }
116:   NEPComputeFunction(nep,lambda,T,T);
117:   MatDuplicate(T,MAT_COPY_VALUES,&Tsigma);
118:   KSPSetOperators(ctx->ksp,Tsigma,Tsigma);
120:   /* Restart loop */
121:   while (nep->reason == NEP_CONVERGED_ITERATING) {
122:     nep->its++;
124:     /* solve projected problem */
125:     DSSetDimensions(nep->ds,n,0,0,0);
126:     DSSetState(nep->ds,DS_STATE_RAW);
127:     DSSolve(nep->ds,nep->eigr,NULL);
128:     lambda = nep->eigr[0];
130:     /* compute Ritz vector, x = V*s */
131:     DSGetArray(nep->ds,DS_MAT_X,&X);
132:     BVSetActiveColumns(nep->V,0,n);
133:     BVMultVec(nep->V,1.0,0.0,x,X);
134:     DSRestoreArray(nep->ds,DS_MAT_X,&X);
136:     /* compute the residual, r = T(lambda)*x */
137:     NEPApplyFunction(nep,lambda,x,w,r,NULL,NULL);
139:     /* convergence test */
140:     VecNorm(r,NORM_2,&resnorm);
141:     (*nep->converged)(nep,lambda,0,resnorm,&nep->errest[nep->nconv],nep->convergedctx);
142:     if (nep->errest[nep->nconv]<=nep->tol) {
143:       BVInsertVec(nep->V,nep->nconv,x);
144:       nep->nconv = nep->nconv + 1;
145:     }
146:     (*nep->stopping)(nep,nep->its,nep->max_it,nep->nconv,nep->nev,&nep->reason,nep->stoppingctx);
147:     NEPMonitor(nep,nep->its,nep->nconv,nep->eigr,nep->eigi,nep->errest,1);
149:     if (nep->reason == NEP_CONVERGED_ITERATING) {
151:       /* continuation vector: f = T(sigma)\r */
152:       BVGetColumn(nep->V,n,&f);
153:       NEPNArnoldi_KSPSolve(nep,r,f);
154:       BVRestoreColumn(nep->V,n,&f);
155:       KSPGetConvergedReason(ctx->ksp,&kspreason);
156:       if (kspreason<0) {
157:         PetscInfo1(nep,"iter=%D, linear solve failed, stopping solve\n",nep->its);
158:         nep->reason = NEP_DIVERGED_LINEAR_SOLVE;
159:         break;
160:       }
162:       /* orthonormalize */
163:       BVOrthogonalizeColumn(nep->V,n,NULL,&beta,&breakdown);
164:       if (breakdown || beta==0.0) {
165:         PetscInfo1(nep,"iter=%D, orthogonalization failed, stopping solve\n",nep->its);
166:         nep->reason = NEP_DIVERGED_BREAKDOWN;
167:         break;
168:       }
169:       BVScaleColumn(nep->V,n,1.0/beta);
171:       /* update projected matrices */
172:       DSSetDimensions(nep->ds,n+1,0,0,0);
173:       NEPProjectOperator(nep,n,n+1);
174:       n++;
175:     }
176:   }
177:   MatDestroy(&Tsigma);
178:   return(0);
179: }
183: PetscErrorCode NEPSetFromOptions_NArnoldi(PetscOptionItems *PetscOptionsObject,NEP nep)
184: {
186:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;
189:   if (!ctx->ksp) { NEPNArnoldiGetKSP(nep,&ctx->ksp); }
190:   KSPSetOperators(ctx->ksp,nep->function,nep->function_pre);
191:   KSPSetFromOptions(ctx->ksp);
192:   return(0);
193: }
197: static PetscErrorCode NEPNArnoldiSetKSP_NArnoldi(NEP nep,KSP ksp)
198: {
200:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;
203:   PetscObjectReference((PetscObject)ksp);
204:   KSPDestroy(&ctx->ksp);
205:   ctx->ksp = ksp;
206:   PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp);
207:   nep->state = NEP_STATE_INITIAL;
208:   return(0);
209: }
213: /*@
214:    NEPNArnoldiSetKSP - Associate a linear solver object (KSP) to the nonlinear
215:    eigenvalue solver.
217:    Collective on NEP
219:    Input Parameters:
220: +  nep - eigenvalue solver
221: -  ksp - the linear solver object
223:    Level: advanced
225: .seealso: NEPNArnoldiGetKSP()
226: @*/
227: PetscErrorCode NEPNArnoldiSetKSP(NEP nep,KSP ksp)
228: {
235:   PetscTryMethod(nep,"NEPNArnoldiSetKSP_C",(NEP,KSP),(nep,ksp));
236:   return(0);
237: }
241: static PetscErrorCode NEPNArnoldiGetKSP_NArnoldi(NEP nep,KSP *ksp)
242: {
244:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;
247:   if (!ctx->ksp) {
248:     KSPCreate(PetscObjectComm((PetscObject)nep),&ctx->ksp);
249:     KSPSetOptionsPrefix(ctx->ksp,((PetscObject)nep)->prefix);
250:     KSPAppendOptionsPrefix(ctx->ksp,"nep_narnoldi_");
251:     PetscObjectIncrementTabLevel((PetscObject)ctx->ksp,(PetscObject)nep,1);
252:     PetscLogObjectParent((PetscObject)nep,(PetscObject)ctx->ksp);
253:     KSPSetErrorIfNotConverged(ctx->ksp,PETSC_TRUE);
254:   }
255:   *ksp = ctx->ksp;
256:   return(0);
257: }
261: /*@
262:    NEPNArnoldiGetKSP - Retrieve the linear solver object (KSP) associated with
263:    the nonlinear eigenvalue solver.
265:    Not Collective
267:    Input Parameter:
268: .  nep - nonlinear eigenvalue solver
270:    Output Parameter:
271: .  ksp - the linear solver object
273:    Level: advanced
275: .seealso: NEPNArnoldiSetKSP()
276: @*/
277: PetscErrorCode NEPNArnoldiGetKSP(NEP nep,KSP *ksp)
278: {
284:   PetscUseMethod(nep,"NEPNArnoldiGetKSP_C",(NEP,KSP*),(nep,ksp));
285:   return(0);
286: }
290: PetscErrorCode NEPView_NArnoldi(NEP nep,PetscViewer viewer)
291: {
293:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;
294:   PetscBool      isascii;
297:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
298:   if (isascii) {
299:     if (!ctx->ksp) { NEPNArnoldiGetKSP(nep,&ctx->ksp); }
300:     PetscViewerASCIIPushTab(viewer);
301:     KSPView(ctx->ksp,viewer);
302:     PetscViewerASCIIPopTab(viewer);
303:   }
304:   return(0);
305: }
309: PetscErrorCode NEPDestroy_NArnoldi(NEP nep)
310: {
312:   NEP_NARNOLDI   *ctx = (NEP_NARNOLDI*)nep->data;
315:   KSPDestroy(&ctx->ksp);
316:   PetscFree(nep->data);
317:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetKSP_C",NULL);
318:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetKSP_C",NULL);
319:   return(0);
320: }
324: PETSC_EXTERN PetscErrorCode NEPCreate_NArnoldi(NEP nep)
325: {
327:   NEP_NARNOLDI   *ctx;
330:   PetscNewLog(nep,&ctx);
331:   nep->data = (void*)ctx;
333:   nep->ops->solve          = NEPSolve_NArnoldi;
334:   nep->ops->setup          = NEPSetUp_NArnoldi;
335:   nep->ops->setfromoptions = NEPSetFromOptions_NArnoldi;
336:   nep->ops->destroy        = NEPDestroy_NArnoldi;
337:   nep->ops->view           = NEPView_NArnoldi;
338:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiSetKSP_C",NEPNArnoldiSetKSP_NArnoldi);
339:   PetscObjectComposeFunction((PetscObject)nep,"NEPNArnoldiGetKSP_C",NEPNArnoldiGetKSP_NArnoldi);
340:   return(0);
341: }