Actual source code: davidson.c

slepc-3.13.3 2020-06-14
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:    Skeleton of Davidson solver. Actual solvers are GD and JD.
 12: */

 14: #include "davidson.h"

 16: static PetscBool  cited = PETSC_FALSE;
 17: static const char citation[] =
 18:   "@Article{slepc-davidson,\n"
 19:   "   author = \"E. Romero and J. E. Roman\",\n"
 20:   "   title = \"A parallel implementation of {Davidson} methods for large-scale eigenvalue problems in {SLEPc}\",\n"
 21:   "   journal = \"{ACM} Trans. Math. Software\",\n"
 22:   "   volume = \"40\",\n"
 23:   "   number = \"2\",\n"
 24:   "   pages = \"13:1--13:29\",\n"
 25:   "   year = \"2014,\"\n"
 26:   "   doi = \"https://doi.org/10.1145/2543696\"\n"
 27:   "}\n";

 29: PetscErrorCode EPSSetUp_XD(EPS eps)
 30: {
 32:   EPS_DAVIDSON   *data = (EPS_DAVIDSON*)eps->data;
 33:   dvdDashboard   *dvd = &data->ddb;
 34:   dvdBlackboard  b;
 35:   PetscInt       min_size_V,bs,initv,nmat;
 36:   Mat            A,B;
 37:   KSP            ksp;
 38:   PetscBool      ipB,ispositive;
 39:   HarmType_t     harm;
 40:   InitType_t     init;
 41:   PetscScalar    target;

 44:   /* Setup EPS options and get the problem specification */
 45:   bs = data->blocksize;
 46:   if (bs <= 0) bs = 1;
 47:   if (eps->ncv!=PETSC_DEFAULT) {
 48:     if (eps->ncv<eps->nev) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The value of ncv must be at least nev");
 49:   } else if (eps->mpd!=PETSC_DEFAULT) eps->ncv = eps->mpd + eps->nev + bs;
 50:   else if (eps->nev<500) eps->ncv = PetscMin(eps->n-bs,PetscMax(2*eps->nev,eps->nev+15))+bs;
 51:   else eps->ncv = PetscMin(eps->n-bs,eps->nev+500)+bs;
 52:   if (eps->mpd==PETSC_DEFAULT) eps->mpd = eps->ncv;
 53:   if (eps->mpd > eps->ncv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The mpd has to be less or equal than ncv");
 54:   if (eps->mpd < 2) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The mpd has to be greater than 2");
 55:   if (eps->max_it==PETSC_DEFAULT) eps->max_it = PetscMax(100*eps->ncv,2*eps->n);
 56:   if (!eps->which) eps->which = EPS_LARGEST_MAGNITUDE;
 57:   if (eps->ishermitian && (eps->which==EPS_LARGEST_IMAGINARY || eps->which==EPS_SMALLEST_IMAGINARY)) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Wrong value of eps->which");
 58:   if (!(eps->nev + bs <= eps->ncv)) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The ncv has to be greater than nev plus blocksize");
 59:   if (eps->trueres) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"-eps_true_residual is temporally disable in this solver.");

 61:   EPSXDSetRestart_XD(eps,data->minv,data->plusk);
 62:   min_size_V = data->minv;
 63:   if (!min_size_V) min_size_V = PetscMin(PetscMax(bs,5),eps->mpd/2);
 64:   if (!(min_size_V+bs <= eps->mpd)) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The value of minv must be less than mpd minus blocksize");
 65:   initv = data->initialsize;
 66:   if (eps->mpd < initv) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"The initv has to be less or equal than mpd");

 68:   /* Change the default sigma to inf if necessary */
 69:   if (eps->which == EPS_LARGEST_MAGNITUDE || eps->which == EPS_LARGEST_REAL || eps->which == EPS_LARGEST_IMAGINARY) {
 70:     STSetDefaultShift(eps->st,PETSC_MAX_REAL);
 71:   }

 73:   /* Set up preconditioner */
 74:   STSetUp(eps->st);

 76:   /* Setup problem specification in dvd */
 77:   STGetNumMatrices(eps->st,&nmat);
 78:   STGetMatrix(eps->st,0,&A);
 79:   if (nmat>1) { STGetMatrix(eps->st,1,&B); }
 80:   EPSReset_XD(eps);
 81:   PetscMemzero(dvd,sizeof(dvdDashboard));
 82:   dvd->A = A; dvd->B = eps->isgeneralized? B: NULL;
 83:   ispositive = eps->ispositive;
 84:   dvd->sA = DVD_MAT_IMPLICIT | (eps->ishermitian? DVD_MAT_HERMITIAN: 0) | ((ispositive && !eps->isgeneralized) ? DVD_MAT_POS_DEF: 0);
 85:   /* Asume -eps_hermitian means hermitian-definite in generalized problems */
 86:   if (!ispositive && !eps->isgeneralized && eps->ishermitian) ispositive = PETSC_TRUE;
 87:   if (!eps->isgeneralized) dvd->sB = DVD_MAT_IMPLICIT | DVD_MAT_HERMITIAN | DVD_MAT_IDENTITY | DVD_MAT_UNITARY | DVD_MAT_POS_DEF;
 88:   else dvd->sB = DVD_MAT_IMPLICIT | (eps->ishermitian? DVD_MAT_HERMITIAN: 0) | (ispositive? DVD_MAT_POS_DEF: 0);
 89:   ipB = (dvd->B && data->ipB && DVD_IS(dvd->sB,DVD_MAT_HERMITIAN))?PETSC_TRUE:PETSC_FALSE;
 90:   dvd->sEP = ((!eps->isgeneralized || (eps->isgeneralized && ipB))? DVD_EP_STD: 0) | (ispositive? DVD_EP_HERMITIAN: 0) | ((eps->problem_type == EPS_GHIEP && ipB) ? DVD_EP_INDEFINITE : 0);
 91:   if (data->ipB && !ipB) data->ipB = PETSC_FALSE;
 92:   dvd->correctXnorm = (dvd->B && (DVD_IS(dvd->sB,DVD_MAT_HERMITIAN)||DVD_IS(dvd->sEP,DVD_EP_INDEFINITE)))?PETSC_TRUE:PETSC_FALSE;
 93:   dvd->nev        = eps->nev;
 94:   dvd->which      = eps->which;
 95:   dvd->withTarget = PETSC_TRUE;
 96:   switch (eps->which) {
 97:     case EPS_TARGET_MAGNITUDE:
 98:     case EPS_TARGET_IMAGINARY:
 99:       dvd->target[0] = target = eps->target;
100:       dvd->target[1] = 1.0;
101:       break;
102:     case EPS_TARGET_REAL:
103:       dvd->target[0] = PetscRealPart(target = eps->target);
104:       dvd->target[1] = 1.0;
105:       break;
106:     case EPS_LARGEST_REAL:
107:     case EPS_LARGEST_MAGNITUDE:
108:     case EPS_LARGEST_IMAGINARY: /* TODO: think about this case */
109:       dvd->target[0] = 1.0;
110:       dvd->target[1] = target = 0.0;
111:       break;
112:     case EPS_SMALLEST_MAGNITUDE:
113:     case EPS_SMALLEST_REAL:
114:     case EPS_SMALLEST_IMAGINARY: /* TODO: think about this case */
115:       dvd->target[0] = target = 0.0;
116:       dvd->target[1] = 1.0;
117:       break;
118:     case EPS_WHICH_USER:
119:       STGetShift(eps->st,&target);
120:       dvd->target[0] = target;
121:       dvd->target[1] = 1.0;
122:       break;
123:     case EPS_ALL:
124:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported option: which == EPS_ALL");
125:     default:
126:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported value of option 'which'");
127:   }
128:   dvd->tol = (eps->tol==PETSC_DEFAULT)? SLEPC_DEFAULT_TOL: eps->tol;
129:   dvd->eps = eps;

131:   /* Setup the extraction technique */
132:   if (!eps->extraction) {
133:     if (ipB || ispositive) eps->extraction = EPS_RITZ;
134:     else {
135:       switch (eps->which) {
136:         case EPS_TARGET_REAL:
137:         case EPS_TARGET_MAGNITUDE:
138:         case EPS_TARGET_IMAGINARY:
139:         case EPS_SMALLEST_MAGNITUDE:
140:         case EPS_SMALLEST_REAL:
141:         case EPS_SMALLEST_IMAGINARY:
142:           eps->extraction = EPS_HARMONIC;
143:           break;
144:         case EPS_LARGEST_REAL:
145:         case EPS_LARGEST_MAGNITUDE:
146:         case EPS_LARGEST_IMAGINARY:
147:           eps->extraction = EPS_HARMONIC_LARGEST;
148:           break;
149:         default:
150:           eps->extraction = EPS_RITZ;
151:       }
152:     }
153:   }
154:   switch (eps->extraction) {
155:     case EPS_RITZ:              harm = DVD_HARM_NONE; break;
156:     case EPS_HARMONIC:          harm = DVD_HARM_RR; break;
157:     case EPS_HARMONIC_RELATIVE: harm = DVD_HARM_RRR; break;
158:     case EPS_HARMONIC_RIGHT:    harm = DVD_HARM_REIGS; break;
159:     case EPS_HARMONIC_LARGEST:  harm = DVD_HARM_LEIGS; break;
160:     default: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Unsupported extraction type");
161:   }

163:   /* Setup the type of starting subspace */
164:   init = data->krylovstart? DVD_INITV_KRYLOV: DVD_INITV_CLASSIC;

166:   /* Preconfigure dvd */
167:   STGetKSP(eps->st,&ksp);
168:   dvd_schm_basic_preconf(dvd,&b,eps->mpd,min_size_V,bs,initv,PetscAbs(eps->nini),data->plusk,harm,ksp,init,eps->trackall,data->ipB,data->doubleexp);

170:   /* Allocate memory */
171:   EPSAllocateSolution(eps,0);

173:   /* Setup orthogonalization */
174:   EPS_SetInnerProduct(eps);
175:   if (!(ipB && dvd->B)) {
176:     BVSetMatrix(eps->V,NULL,PETSC_FALSE);
177:   }

179:   /* Configure dvd for a basic GD */
180:   dvd_schm_basic_conf(dvd,&b,eps->mpd,min_size_V,bs,initv,PetscAbs(eps->nini),data->plusk,harm,dvd->withTarget,target,ksp,data->fix,init,eps->trackall,data->ipB,data->dynamic,data->doubleexp);
181:   return(0);
182: }

184: PetscErrorCode EPSSolve_XD(EPS eps)
185: {
186:   EPS_DAVIDSON   *data = (EPS_DAVIDSON*)eps->data;
187:   dvdDashboard   *d = &data->ddb;
188:   PetscInt       l,k;

192:   PetscCitationsRegister(citation,&cited);
193:   /* Call the starting routines */
194:   EPSDavidsonFLCall(d->startList,d);

196:   while (eps->reason == EPS_CONVERGED_ITERATING) {

198:     /* Initialize V, if it is needed */
199:     BVGetActiveColumns(d->eps->V,&l,&k);
200:     if (l == k) { d->initV(d); }

202:     /* Find the best approximated eigenpairs in V, X */
203:     d->calcPairs(d);

205:     /* Test for convergence */
206:     (*eps->stopping)(eps,eps->its,eps->max_it,eps->nconv,eps->nev,&eps->reason,eps->stoppingctx);
207:     if (eps->reason != EPS_CONVERGED_ITERATING) break;

209:     /* Expand the subspace */
210:     d->updateV(d);

212:     /* Monitor progress */
213:     eps->nconv = d->nconv;
214:     eps->its++;
215:     BVGetActiveColumns(d->eps->V,NULL,&k);
216:     EPSMonitor(eps,eps->its,eps->nconv+d->npreconv,eps->eigr,eps->eigi,eps->errest,PetscMin(k,eps->nev));
217:   }

219:   /* Call the ending routines */
220:   EPSDavidsonFLCall(d->endList,d);
221:   return(0);
222: }

224: PetscErrorCode EPSReset_XD(EPS eps)
225: {
226:   EPS_DAVIDSON   *data = (EPS_DAVIDSON*)eps->data;
227:   dvdDashboard   *dvd = &data->ddb;

231:   /* Call step destructors and destroys the list */
232:   EPSDavidsonFLCall(dvd->destroyList,dvd);
233:   EPSDavidsonFLDestroy(&dvd->destroyList);
234:   EPSDavidsonFLDestroy(&dvd->startList);
235:   EPSDavidsonFLDestroy(&dvd->endList);
236:   return(0);
237: }

239: PetscErrorCode EPSXDSetKrylovStart_XD(EPS eps,PetscBool krylovstart)
240: {
241:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

244:   data->krylovstart = krylovstart;
245:   return(0);
246: }

248: PetscErrorCode EPSXDGetKrylovStart_XD(EPS eps,PetscBool *krylovstart)
249: {
250:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

253:   *krylovstart = data->krylovstart;
254:   return(0);
255: }

257: PetscErrorCode EPSXDSetBlockSize_XD(EPS eps,PetscInt blocksize)
258: {
259:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

262:   if (blocksize == PETSC_DEFAULT || blocksize == PETSC_DECIDE) blocksize = 1;
263:   if (blocksize <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid blocksize value");
264:   data->blocksize = blocksize;
265:   return(0);
266: }

268: PetscErrorCode EPSXDGetBlockSize_XD(EPS eps,PetscInt *blocksize)
269: {
270:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

273:   *blocksize = data->blocksize;
274:   return(0);
275: }

277: PetscErrorCode EPSXDSetRestart_XD(EPS eps,PetscInt minv,PetscInt plusk)
278: {
279:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

282:   if (minv == PETSC_DEFAULT || minv == PETSC_DECIDE) minv = 5;
283:   if (minv <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid minv value");
284:   if (plusk == PETSC_DEFAULT || plusk == PETSC_DECIDE) plusk = eps->problem_type == EPS_GHIEP?0:1;
285:   if (plusk < 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid plusk value");
286:   data->minv = minv;
287:   data->plusk = plusk;
288:   return(0);
289: }

291: PetscErrorCode EPSXDGetRestart_XD(EPS eps,PetscInt *minv,PetscInt *plusk)
292: {
293:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

296:   if (minv) *minv = data->minv;
297:   if (plusk) *plusk = data->plusk;
298:   return(0);
299: }

301: PetscErrorCode EPSXDGetInitialSize_XD(EPS eps,PetscInt *initialsize)
302: {
303:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

306:   *initialsize = data->initialsize;
307:   return(0);
308: }

310: PetscErrorCode EPSXDSetInitialSize_XD(EPS eps,PetscInt initialsize)
311: {
312:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

315:   if (initialsize == PETSC_DEFAULT || initialsize == PETSC_DECIDE) initialsize = 5;
316:   if (initialsize <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid initial size value");
317:   data->initialsize = initialsize;
318:   return(0);
319: }

321: PetscErrorCode EPSXDSetBOrth_XD(EPS eps,PetscBool borth)
322: {
323:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

326:   data->ipB = borth;
327:   return(0);
328: }

330: PetscErrorCode EPSXDGetBOrth_XD(EPS eps,PetscBool *borth)
331: {
332:   EPS_DAVIDSON *data = (EPS_DAVIDSON*)eps->data;

335:   *borth = data->ipB;
336:   return(0);
337: }

339: /*
340:   EPSComputeVectors_XD - Compute eigenvectors from the vectors
341:   provided by the eigensolver. This version is intended for solvers
342:   that provide Schur vectors from the QZ decomposition. Given the partial
343:   Schur decomposition OP*V=V*T, the following steps are performed:
344:       1) compute eigenvectors of (S,T): S*Z=T*Z*D
345:       2) compute eigenvectors of OP: X=V*Z
346:  */
347: PetscErrorCode EPSComputeVectors_XD(EPS eps)
348: {
350:   Mat            X;
351:   PetscBool      symm;

354:   PetscObjectTypeCompare((PetscObject)eps->ds,DSHEP,&symm);
355:   if (symm) return(0);
356:   DSVectors(eps->ds,DS_MAT_X,NULL,NULL);

358:   /* V <- V * X */
359:   DSGetMat(eps->ds,DS_MAT_X,&X);
360:   BVSetActiveColumns(eps->V,0,eps->nconv);
361:   BVMultInPlace(eps->V,X,0,eps->nconv);
362:   MatDestroy(&X);
363:   return(0);
364: }