Actual source code: krylovschur.c

  1: /*

  3:    SLEPc eigensolver: "krylovschur"

  5:    Method: Krylov-Schur

  7:    Algorithm:

  9:        Single-vector Krylov-Schur method for non-symmetric problems,
 10:        including harmonic extraction.

 12:    References:

 14:        [1] "Krylov-Schur Methods in SLEPc", SLEPc Technical Report STR-7,
 15:            available at http://www.grycap.upv.es/slepc.

 17:        [2] G.W. Stewart, "A Krylov-Schur Algorithm for Large Eigenproblems",
 18:            SIAM J. Matrix Anal. App. 23(3):601-614, 2001.

 20:        [3] "Practical Implementation of Harmonic Krylov-Schur", SLEPc Technical
 21:             Report STR-9, available at http://www.grycap.upv.es/slepc.

 23:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 24:    SLEPc - Scalable Library for Eigenvalue Problem Computations
 25:    Copyright (c) 2002-2013, Universitat Politecnica de Valencia, Spain

 27:    This file is part of SLEPc.

 29:    SLEPc is free software: you can redistribute it and/or modify it under  the
 30:    terms of version 3 of the GNU Lesser General Public License as published by
 31:    the Free Software Foundation.

 33:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 34:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 35:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 36:    more details.

 38:    You  should have received a copy of the GNU Lesser General  Public  License
 39:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 40:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 41: */

 43: #include <slepc-private/epsimpl.h>                /*I "slepceps.h" I*/
 44: #include <slepcblaslapack.h>
 45:  #include krylovschur.h

 49: PetscErrorCode EPSGetArbitraryValues(EPS eps,PetscScalar *rr,PetscScalar *ri)
 50: {
 52:   PetscInt       i,newi,ld,n,l;
 53:   Vec            xr=eps->work[1],xi=eps->work[2];
 54:   PetscScalar    re,im,*Zr,*Zi,*X;

 57:   DSGetLeadingDimension(eps->ds,&ld);
 58:   DSGetDimensions(eps->ds,&n,NULL,&l,NULL,NULL);
 59:   for (i=l;i<n;i++) {
 60:     re = eps->eigr[i];
 61:     im = eps->eigi[i];
 62:     STBackTransform(eps->st,1,&re,&im);
 63:     newi = i;
 64:     DSVectors(eps->ds,DS_MAT_X,&newi,NULL);
 65:     DSGetArray(eps->ds,DS_MAT_X,&X);
 66:     Zr = X+i*ld;
 67:     if (newi==i+1) Zi = X+newi*ld;
 68:     else Zi = NULL;
 69:     EPSComputeRitzVector(eps,Zr,Zi,eps->V,n,xr,xi);
 70:     DSRestoreArray(eps->ds,DS_MAT_X,&X);
 71:     (*eps->arbitrary)(re,im,xr,xi,rr+i,ri+i,eps->arbitraryctx);
 72:   }
 73:   return(0);
 74: }

 78: PetscErrorCode EPSSetUp_KrylovSchur(EPS eps)
 79: {
 80:   PetscErrorCode  ierr;
 81:   PetscBool       issinv;
 82:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
 83:   enum { EPS_KS_DEFAULT,EPS_KS_SYMM,EPS_KS_SLICE,EPS_KS_INDEF } variant;

 86:   /* spectrum slicing requires special treatment of default values */
 87:   if (eps->which==EPS_ALL) {
 88:     if (eps->inta==0.0 && eps->intb==0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Must define a computational interval when using EPS_ALL");
 89:     if (!eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing only available for symmetric/Hermitian eigenproblems");
 90:     if (eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs cannot be used with spectrum slicing");
 91:     if (!((PetscObject)(eps->st))->type_name) { /* default to shift-and-invert */
 92:       STSetType(eps->st,STSINVERT);
 93:     }
 94:     PetscObjectTypeCompareAny((PetscObject)eps->st,&issinv,STSINVERT,STCAYLEY,"");
 95:     if (!issinv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Shift-and-invert or Cayley ST is needed for spectrum slicing");
 96: #if defined(PETSC_USE_REAL_DOUBLE)
 97:     if (eps->tol==PETSC_DEFAULT) eps->tol = 1e-10;  /* use tighter tolerance */
 98: #endif
 99:     if (eps->intb >= PETSC_MAX_REAL) { /* right-open interval */
100:       if (eps->inta <= PETSC_MIN_REAL) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"The defined computational interval should have at least one of their sides bounded");
101:       STSetDefaultShift(eps->st,eps->inta);
102:     } else {
103:       STSetDefaultShift(eps->st,eps->intb);
104:     }

106:     if (eps->nev==1) eps->nev = 40;  /* nev not set, use default value */
107:     if (eps->nev<10) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"nev cannot be less than 10 in spectrum slicing runs");
108:     eps->ops->backtransform = NULL;
109:   }

111:   if (eps->isgeneralized && eps->ishermitian && !eps->ispositive && eps->arbitrary) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Arbitrary selection of eigenpairs not implemented for indefinite problems");

113:   /* proceed with the general case */
114:   if (eps->ncv) { /* ncv set */
115:     if (eps->ncv<eps->nev+1) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev+1");
116:   } else if (eps->mpd) { /* mpd set */
117:     eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
118:   } else { /* neither set: defaults depend on nev being small or large */
119:     if (eps->nev<500) eps->ncv = PetscMin(eps->n,PetscMax(2*eps->nev,eps->nev+15));
120:     else {
121:       eps->mpd = 500;
122:       eps->ncv = PetscMin(eps->n,eps->nev+eps->mpd);
123:     }
124:   }
125:   if (!eps->mpd) eps->mpd = eps->ncv;
126:   if (eps->ncv>eps->nev+eps->mpd) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must not be larger than nev+mpd");
127:   if (!eps->max_it) {
128:     if (eps->which==EPS_ALL) eps->max_it = 100;  /* special case for spectrum slicing */
129:     else eps->max_it = PetscMax(100,2*eps->n/eps->ncv);
130:   }
131:   if (!eps->which) { EPSSetWhichEigenpairs_Default(eps); }
132:   if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Wrong value of eps->which");

134:   if (!eps->extraction) {
135:     EPSSetExtraction(eps,EPS_RITZ);
136:   } else if (eps->extraction!=EPS_RITZ && eps->extraction!=EPS_HARMONIC)
137:     SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");

139:   if (!ctx->keep) ctx->keep = 0.5;

141:   EPSAllocateSolution(eps);
142:   if (eps->arbitrary) {
143:     EPSSetWorkVecs(eps,3);
144:   } else if (eps->ishermitian && !eps->ispositive){
145:     EPSSetWorkVecs(eps,2);
146:   } else {
147:     EPSSetWorkVecs(eps,1);
148:   }

150:   /* dispatch solve method */
151:   if (eps->leftvecs) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Left vectors not supported in this solver");
152:   if (eps->ishermitian) {
153:     if (eps->which==EPS_ALL) {
154:       if (eps->isgeneralized && !eps->ispositive) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Spectrum slicing not implemented for indefinite problems");
155:       else variant = EPS_KS_SLICE;
156:     } else if (eps->isgeneralized && !eps->ispositive) {
157:       variant = EPS_KS_INDEF;
158:     } else {
159:       switch (eps->extraction) {
160:         case EPS_RITZ:     variant = EPS_KS_SYMM; break;
161:         case EPS_HARMONIC: variant = EPS_KS_DEFAULT; break;
162:         default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
163:       }
164:     }
165:   } else {
166:     switch (eps->extraction) {
167:       case EPS_RITZ:     variant = EPS_KS_DEFAULT; break;
168:       case EPS_HARMONIC: variant = EPS_KS_DEFAULT; break;
169:       default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
170:     }
171:   }
172:   switch (variant) {
173:     case EPS_KS_DEFAULT:
174:       eps->ops->solve = EPSSolve_KrylovSchur_Default;
175:       DSSetType(eps->ds,DSNHEP);
176:       break;
177:     case EPS_KS_SYMM:
178:       eps->ops->solve = EPSSolve_KrylovSchur_Symm;
179:       DSSetType(eps->ds,DSHEP);
180:       DSSetCompact(eps->ds,PETSC_TRUE);
181:       DSSetExtraRow(eps->ds,PETSC_TRUE);
182:       break;
183:     case EPS_KS_SLICE:
184:       eps->ops->solve = EPSSolve_KrylovSchur_Slice;
185:       DSSetType(eps->ds,DSHEP);
186:       DSSetCompact(eps->ds,PETSC_TRUE);
187:       break;
188:     case EPS_KS_INDEF:
189:       eps->ops->solve = EPSSolve_KrylovSchur_Indefinite;
190:       DSSetType(eps->ds,DSGHIEP);
191:       DSSetCompact(eps->ds,PETSC_TRUE);
192:       break;
193:     default: SETERRQ(PetscObjectComm((PetscObject)eps),1,"Unexpected error");
194:   }
195:   DSAllocate(eps->ds,eps->ncv+1);
196:   return(0);
197: }

201: PetscErrorCode EPSSolve_KrylovSchur_Default(EPS eps)
202: {
203:   PetscErrorCode  ierr;
204:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
205:   PetscInt        i,j,*pj,k,l,nv,ld;
206:   Vec             u=eps->work[0];
207:   PetscScalar     *S,*Q,*g;
208:   PetscReal       beta,gamma=1.0;
209:   PetscBool       breakdown,harmonic;

212:   DSGetLeadingDimension(eps->ds,&ld);
213:   harmonic = (eps->extraction==EPS_HARMONIC || eps->extraction==EPS_REFINED_HARMONIC)?PETSC_TRUE:PETSC_FALSE;
214:   if (harmonic) { PetscMalloc(ld*sizeof(PetscScalar),&g); }
215:   if (eps->arbitrary) pj = &j;
216:   else pj = NULL;

218:   /* Get the starting Arnoldi vector */
219:   EPSGetStartVector(eps,0,eps->V[0],NULL);
220:   l = 0;

222:   /* Restart loop */
223:   while (eps->reason == EPS_CONVERGED_ITERATING) {
224:     eps->its++;

226:     /* Compute an nv-step Arnoldi factorization */
227:     nv = PetscMin(eps->nconv+eps->mpd,eps->ncv);
228:     DSGetArray(eps->ds,DS_MAT_A,&S);
229:     EPSBasicArnoldi(eps,PETSC_FALSE,S,ld,eps->V,eps->nconv+l,&nv,u,&beta,&breakdown);
230:     VecScale(u,1.0/beta);
231:     DSRestoreArray(eps->ds,DS_MAT_A,&S);
232:     DSSetDimensions(eps->ds,nv,0,eps->nconv,eps->nconv+l);
233:     if (l==0) {
234:       DSSetState(eps->ds,DS_STATE_INTERMEDIATE);
235:     } else {
236:       DSSetState(eps->ds,DS_STATE_RAW);
237:     }

239:     /* Compute translation of Krylov decomposition if harmonic extraction used */
240:     if (harmonic) {
241:       DSTranslateHarmonic(eps->ds,eps->target,beta,PETSC_FALSE,g,&gamma);
242:     }

244:     /* Solve projected problem */
245:     DSSolve(eps->ds,eps->eigr,eps->eigi);
246:     if (eps->arbitrary) {
247:       EPSGetArbitraryValues(eps,eps->rr,eps->ri);
248:       j=1;
249:     }
250:     DSSort(eps->ds,eps->eigr,eps->eigi,eps->rr,eps->ri,pj);

252:     /* Check convergence */
253:     EPSKrylovConvergence(eps,PETSC_FALSE,eps->nconv,nv-eps->nconv,eps->V,nv,beta,gamma,&k);
254:     if (eps->its >= eps->max_it) eps->reason = EPS_DIVERGED_ITS;
255:     if (k >= eps->nev) eps->reason = EPS_CONVERGED_TOL;

257:     /* Update l */
258:     if (eps->reason != EPS_CONVERGED_ITERATING || breakdown) l = 0;
259:     else {
260:       l = PetscMax(1,(PetscInt)((nv-k)*ctx->keep));
261: #if !defined(PETSC_USE_COMPLEX)
262:       DSGetArray(eps->ds,DS_MAT_A,&S);
263:       if (S[k+l+(k+l-1)*ld] != 0.0) {
264:         if (k+l<nv-1) l = l+1;
265:         else l = l-1;
266:       }
267:       DSRestoreArray(eps->ds,DS_MAT_A,&S);
268: #endif
269:     }

271:     if (eps->reason == EPS_CONVERGED_ITERATING) {
272:       if (breakdown) {
273:         /* Start a new Arnoldi factorization */
274:         PetscInfo2(eps,"Breakdown in Krylov-Schur method (it=%D norm=%G)\n",eps->its,beta);
275:         EPSGetStartVector(eps,k,eps->V[k],&breakdown);
276:         if (breakdown) {
277:           eps->reason = EPS_DIVERGED_BREAKDOWN;
278:           PetscInfo(eps,"Unable to generate more start vectors\n");
279:         }
280:       } else {
281:         /* Undo translation of Krylov decomposition */
282:         if (harmonic) {
283:           DSSetDimensions(eps->ds,nv,0,k,l);
284:           DSTranslateHarmonic(eps->ds,0.0,beta,PETSC_TRUE,g,&gamma);
285:           /* gamma u^ = u - U*g~ */
286:           SlepcVecMAXPBY(u,1.0,-1.0,nv,g,eps->V);
287:           VecScale(u,1.0/gamma);
288:         }
289:         /* Prepare the Rayleigh quotient for restart */
290:         DSGetArray(eps->ds,DS_MAT_A,&S);
291:         DSGetArray(eps->ds,DS_MAT_Q,&Q);
292:         for (i=k;i<k+l;i++) {
293:           S[k+l+i*ld] = Q[nv-1+i*ld]*beta*gamma;
294:         }
295:         DSRestoreArray(eps->ds,DS_MAT_A,&S);
296:         DSRestoreArray(eps->ds,DS_MAT_Q,&Q);
297:       }
298:     }
299:     /* Update the corresponding vectors V(:,idx) = V*Q(:,idx) */
300:     DSGetArray(eps->ds,DS_MAT_Q,&Q);
301:     SlepcUpdateVectors(nv,eps->V,eps->nconv,k+l,Q,ld,PETSC_FALSE);
302:     DSRestoreArray(eps->ds,DS_MAT_Q,&Q);

304:     if (eps->reason == EPS_CONVERGED_ITERATING && !breakdown) {
305:       VecCopy(u,eps->V[k+l]);
306:     }
307:     eps->nconv = k;
308:     EPSMonitor(eps,eps->its,eps->nconv,eps->eigr,eps->eigi,eps->errest,nv);
309:   }

311:   if (harmonic) { PetscFree(g); }
312:   /* truncate Schur decomposition and change the state to raw so that
313:      PSVectors() computes eigenvectors from scratch */
314:   DSSetDimensions(eps->ds,eps->nconv,0,0,0);
315:   DSSetState(eps->ds,DS_STATE_RAW);
316:   return(0);
317: }

321: static PetscErrorCode EPSKrylovSchurSetRestart_KrylovSchur(EPS eps,PetscReal keep)
322: {
323:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

326:   if (keep==PETSC_DEFAULT) ctx->keep = 0.5;
327:   else {
328:     if (keep<0.1 || keep>0.9) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"The keep argument must be in the range [0.1,0.9]");
329:     ctx->keep = keep;
330:   }
331:   return(0);
332: }

336: /*@
337:    EPSKrylovSchurSetRestart - Sets the restart parameter for the Krylov-Schur
338:    method, in particular the proportion of basis vectors that must be kept
339:    after restart.

341:    Logically Collective on EPS

343:    Input Parameters:
344: +  eps - the eigenproblem solver context
345: -  keep - the number of vectors to be kept at restart

347:    Options Database Key:
348: .  -eps_krylovschur_restart - Sets the restart parameter

350:    Notes:
351:    Allowed values are in the range [0.1,0.9]. The default is 0.5.

353:    Level: advanced

355: .seealso: EPSKrylovSchurGetRestart()
356: @*/
357: PetscErrorCode EPSKrylovSchurSetRestart(EPS eps,PetscReal keep)
358: {

364:   PetscTryMethod(eps,"EPSKrylovSchurSetRestart_C",(EPS,PetscReal),(eps,keep));
365:   return(0);
366: }

370: static PetscErrorCode EPSKrylovSchurGetRestart_KrylovSchur(EPS eps,PetscReal *keep)
371: {
372:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

375:   *keep = ctx->keep;
376:   return(0);
377: }

381: /*@
382:    EPSKrylovSchurGetRestart - Gets the restart parameter used in the
383:    Krylov-Schur method.

385:    Not Collective

387:    Input Parameter:
388: .  eps - the eigenproblem solver context

390:    Output Parameter:
391: .  keep - the restart parameter

393:    Level: advanced

395: .seealso: EPSKrylovSchurSetRestart()
396: @*/
397: PetscErrorCode EPSKrylovSchurGetRestart(EPS eps,PetscReal *keep)
398: {

404:   PetscTryMethod(eps,"EPSKrylovSchurGetRestart_C",(EPS,PetscReal*),(eps,keep));
405:   return(0);
406: }

410: PetscErrorCode EPSSetFromOptions_KrylovSchur(EPS eps)
411: {
413:   PetscBool      flg;
414:   PetscReal      keep;

417:   PetscOptionsHead("EPS Krylov-Schur Options");
418:   PetscOptionsReal("-eps_krylovschur_restart","Proportion of vectors kept after restart","EPSKrylovSchurSetRestart",0.5,&keep,&flg);
419:   if (flg) {
420:     EPSKrylovSchurSetRestart(eps,keep);
421:   }
422:   PetscOptionsTail();
423:   return(0);
424: }

428: PetscErrorCode EPSView_KrylovSchur(EPS eps,PetscViewer viewer)
429: {
430:   PetscErrorCode  ierr;
431:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;
432:   PetscBool       isascii;

435:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
436:   if (isascii) {
437:     PetscViewerASCIIPrintf(viewer,"  Krylov-Schur: %d%% of basis vectors kept after restart\n",(int)(100*ctx->keep));
438:   }
439:   return(0);
440: }

444: PetscErrorCode EPSReset_KrylovSchur(EPS eps)
445: {
446:   PetscErrorCode  ierr;
447:   EPS_KRYLOVSCHUR *ctx = (EPS_KRYLOVSCHUR*)eps->data;

450:   ctx->keep = 0.0;
451:   EPSReset_Default(eps);
452:   return(0);
453: }

457: PetscErrorCode EPSDestroy_KrylovSchur(EPS eps)
458: {

462:   PetscFree(eps->data);
463:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",NULL);
464:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",NULL);
465:   return(0);
466: }

470: PETSC_EXTERN PetscErrorCode EPSCreate_KrylovSchur(EPS eps)
471: {

475:   PetscNewLog(eps,EPS_KRYLOVSCHUR,&eps->data);
476:   eps->ops->setup          = EPSSetUp_KrylovSchur;
477:   eps->ops->setfromoptions = EPSSetFromOptions_KrylovSchur;
478:   eps->ops->destroy        = EPSDestroy_KrylovSchur;
479:   eps->ops->reset          = EPSReset_KrylovSchur;
480:   eps->ops->view           = EPSView_KrylovSchur;
481:   eps->ops->backtransform  = EPSBackTransform_Default;
482:   eps->ops->computevectors = EPSComputeVectors_Schur;
483:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurSetRestart_C",EPSKrylovSchurSetRestart_KrylovSchur);
484:   PetscObjectComposeFunction((PetscObject)eps,"EPSKrylovSchurGetRestart_C",EPSKrylovSchurGetRestart_KrylovSchur);
485:   return(0);
486: }