Actual source code: epssetup.c

slepc-3.11.0 2019-03-29
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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:    EPS routines related to problem setup
 12: */

 14: #include <slepc/private/epsimpl.h>       /*I "slepceps.h" I*/

 16: /*
 17:    Let the solver choose the ST type that should be used by default,
 18:    otherwise set it to SHIFT.
 19:    This is called at EPSSetFromOptions (before STSetFromOptions)
 20:    and also at EPSSetUp (in case EPSSetFromOptions was not called).
 21: */
 22: PetscErrorCode EPSSetDefaultST(EPS eps)
 23: {

 27:   if (eps->ops->setdefaultst) { (*eps->ops->setdefaultst)(eps); }
 28:   if (!((PetscObject)eps->st)->type_name) {
 29:     STSetType(eps->st,STSHIFT);
 30:   }
 31:   return(0);
 32: }

 34: /*
 35:    This is done by preconditioned eigensolvers that use the PC only.
 36:    It sets STPRECOND with KSPPREONLY.
 37: */
 38: PetscErrorCode EPSSetDefaultST_Precond(EPS eps)
 39: {
 41:   KSP            ksp;

 44:   if (!((PetscObject)eps->st)->type_name) {
 45:     STSetType(eps->st,STPRECOND);
 46:     STPrecondSetKSPHasMat(eps->st,PETSC_TRUE);
 47:   }
 48:   STGetKSP(eps->st,&ksp);
 49:   if (!((PetscObject)ksp)->type_name) {
 50:     KSPSetType(ksp,KSPPREONLY);
 51:   }
 52:   return(0);
 53: }

 55: /*
 56:    This is done by preconditioned eigensolvers that can also use the KSP.
 57:    It sets STPRECOND with the default KSP (GMRES).
 58: */
 59: PetscErrorCode EPSSetDefaultST_GMRES(EPS eps)
 60: {
 62:   KSP            ksp;

 65:   if (!((PetscObject)eps->st)->type_name) {
 66:     STSetType(eps->st,STPRECOND);
 67:     STPrecondSetKSPHasMat(eps->st,PETSC_TRUE);
 68:   }
 69:   STGetKSP(eps->st,&ksp);
 70:   KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,5);
 71:   return(0);
 72: }

 74: /*@
 75:    EPSSetUp - Sets up all the internal data structures necessary for the
 76:    execution of the eigensolver. Then calls STSetUp() for any set-up
 77:    operations associated to the ST object.

 79:    Collective on EPS

 81:    Input Parameter:
 82: .  eps   - eigenproblem solver context

 84:    Notes:
 85:    This function need not be called explicitly in most cases, since EPSSolve()
 86:    calls it. It can be useful when one wants to measure the set-up time
 87:    separately from the solve time.

 89:    Level: developer

 91: .seealso: EPSCreate(), EPSSolve(), EPSDestroy(), STSetUp(), EPSSetInitialSpace()
 92: @*/
 93: PetscErrorCode EPSSetUp(EPS eps)
 94: {
 96:   Mat            A,B;
 97:   SlepcSC        sc;
 98:   PetscInt       k,nmat;
 99:   PetscBool      flg,istrivial,precond;
100: #if defined(PETSC_USE_COMPLEX)
101:   PetscScalar    sigma;
102: #endif

106:   if (eps->state) return(0);
107:   PetscLogEventBegin(EPS_SetUp,eps,0,0,0);

109:   /* reset the convergence flag from the previous solves */
110:   eps->reason = EPS_CONVERGED_ITERATING;

112:   /* Set default solver type (EPSSetFromOptions was not called) */
113:   if (!((PetscObject)eps)->type_name) {
114:     EPSSetType(eps,EPSKRYLOVSCHUR);
115:   }
116:   if (!eps->st) { EPSGetST(eps,&eps->st); }
117:   EPSSetDefaultST(eps);

119:   PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&precond);
120:   if (eps->categ==EPS_CATEGORY_PRECOND && !precond) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires ST=PRECOND");
121:   if (eps->categ!=EPS_CATEGORY_PRECOND && precond) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"STPRECOND is intended for preconditioned eigensolvers only");

123:   STSetTransform(eps->st,PETSC_TRUE);
124:   if (eps->useds && !eps->ds) { EPSGetDS(eps,&eps->ds); }
125:   if (eps->twosided) {
126:     if (!eps->hasts) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver does not support computing left eigenvectors (no two-sided variant)");
127:     if (eps->ishermitian) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Two-sided methods are not intended for symmetric problems");
128:     if (!eps->dsts) { DSDuplicate(eps->ds,&eps->dsts); }
129:   }
130:   if (!eps->rg) { EPSGetRG(eps,&eps->rg); }
131:   if (!((PetscObject)eps->rg)->type_name) {
132:     RGSetType(eps->rg,RGINTERVAL);
133:   }

135:   /* Set problem dimensions */
136:   STGetNumMatrices(eps->st,&nmat);
137:   if (!nmat) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"EPSSetOperators must be called first");
138:   STMatGetSize(eps->st,&eps->n,NULL);
139:   STMatGetLocalSize(eps->st,&eps->nloc,NULL);

141:   /* Set default problem type */
142:   if (!eps->problem_type) {
143:     if (nmat==1) {
144:       EPSSetProblemType(eps,EPS_NHEP);
145:     } else {
146:       EPSSetProblemType(eps,EPS_GNHEP);
147:     }
148:   } else if (nmat==1 && eps->isgeneralized) {
149:     PetscInfo(eps,"Eigenproblem set as generalized but no matrix B was provided; reverting to a standard eigenproblem\n");
150:     eps->isgeneralized = PETSC_FALSE;
151:     eps->problem_type = eps->ishermitian? EPS_HEP: EPS_NHEP;
152:   } else if (nmat>1 && !eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Inconsistent EPS state");

154:   if (eps->nev > eps->n) eps->nev = eps->n;
155:   if (eps->ncv > eps->n) eps->ncv = eps->n;

157:   /* initialization of matrix norms */
158:   if (eps->conv==EPS_CONV_NORM) {
159:     if (!eps->nrma) {
160:       STGetMatrix(eps->st,0,&A);
161:       MatNorm(A,NORM_INFINITY,&eps->nrma);
162:     }
163:     if (nmat>1 && !eps->nrmb) {
164:       STGetMatrix(eps->st,1,&B);
165:       MatNorm(B,NORM_INFINITY,&eps->nrmb);
166:     }
167:   }

169:   /* call specific solver setup */
170:   (*eps->ops->setup)(eps);

172:   /* if purification is set, check that it really makes sense */
173:   if (eps->purify) {
174:     if (eps->categ==EPS_CATEGORY_PRECOND || eps->categ==EPS_CATEGORY_CONTOUR) eps->purify = PETSC_FALSE;
175:     else {
176:       if (!eps->isgeneralized) eps->purify = PETSC_FALSE;
177:       else if (!eps->ishermitian && !eps->ispositive) eps->purify = PETSC_FALSE;
178:       else {
179:         PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&flg);
180:         if (flg) eps->purify = PETSC_FALSE;
181:       }
182:     }
183:   }

185:   /* check extraction */
186:   PetscObjectTypeCompareAny((PetscObject)eps->st,&flg,STPRECOND,STSHIFT,"");
187:   if (!flg && eps->extraction && eps->extraction!=EPS_RITZ) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cannot use a spectral transformation combined with harmonic extraction");

189:   /* set tolerance if not yet set */
190:   if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL;

192:   /* fill sorting criterion context */
193:   switch (eps->which) {
194:     case EPS_LARGEST_MAGNITUDE:
195:       eps->sc->comparison    = SlepcCompareLargestMagnitude;
196:       eps->sc->comparisonctx = NULL;
197:       break;
198:     case EPS_SMALLEST_MAGNITUDE:
199:       eps->sc->comparison    = SlepcCompareSmallestMagnitude;
200:       eps->sc->comparisonctx = NULL;
201:       break;
202:     case EPS_LARGEST_REAL:
203:       eps->sc->comparison    = SlepcCompareLargestReal;
204:       eps->sc->comparisonctx = NULL;
205:       break;
206:     case EPS_SMALLEST_REAL:
207:       eps->sc->comparison    = SlepcCompareSmallestReal;
208:       eps->sc->comparisonctx = NULL;
209:       break;
210:     case EPS_LARGEST_IMAGINARY:
211:       eps->sc->comparison    = SlepcCompareLargestImaginary;
212:       eps->sc->comparisonctx = NULL;
213:       break;
214:     case EPS_SMALLEST_IMAGINARY:
215:       eps->sc->comparison    = SlepcCompareSmallestImaginary;
216:       eps->sc->comparisonctx = NULL;
217:       break;
218:     case EPS_TARGET_MAGNITUDE:
219:       eps->sc->comparison    = SlepcCompareTargetMagnitude;
220:       eps->sc->comparisonctx = &eps->target;
221:       break;
222:     case EPS_TARGET_REAL:
223:       eps->sc->comparison    = SlepcCompareTargetReal;
224:       eps->sc->comparisonctx = &eps->target;
225:       break;
226:     case EPS_TARGET_IMAGINARY:
227: #if defined(PETSC_USE_COMPLEX)
228:       eps->sc->comparison    = SlepcCompareTargetImaginary;
229:       eps->sc->comparisonctx = &eps->target;
230: #endif
231:       break;
232:     case EPS_ALL:
233:       eps->sc->comparison    = SlepcCompareSmallestReal;
234:       eps->sc->comparisonctx = NULL;
235:       break;
236:     case EPS_WHICH_USER:
237:       break;
238:   }
239:   eps->sc->map    = NULL;
240:   eps->sc->mapobj = NULL;

242:   /* fill sorting criterion for DS */
243:   if (eps->useds && eps->which!=EPS_ALL) {
244:     DSGetSlepcSC(eps->ds,&sc);
245:     RGIsTrivial(eps->rg,&istrivial);
246:     sc->rg            = istrivial? NULL: eps->rg;
247:     sc->comparison    = eps->sc->comparison;
248:     sc->comparisonctx = eps->sc->comparisonctx;
249:     sc->map           = SlepcMap_ST;
250:     sc->mapobj        = (PetscObject)eps->st;
251:     if (eps->twosided) {
252:       DSSetSlepcSC(eps->dsts,sc);
253:     }
254:   }

256:   /* Build balancing matrix if required */
257:   if (!eps->ishermitian && (eps->balance==EPS_BALANCE_ONESIDE || eps->balance==EPS_BALANCE_TWOSIDE)) {
258:     if (!eps->D) {
259:       BVCreateVec(eps->V,&eps->D);
260:       PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->D);
261:     } else {
262:       VecSet(eps->D,1.0);
263:     }
264:     EPSBuildBalance_Krylov(eps);
265:     STSetBalanceMatrix(eps->st,eps->D);
266:   }

268:   /* Setup ST */
269:   STSetUp(eps->st);

271: #if defined(PETSC_USE_COMPLEX)
272:   STGetShift(eps->st,&sigma);
273:   if (eps->ishermitian && PetscImaginaryPart(sigma) != 0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Hermitian problems are not compatible with complex shifts");
274: #endif
275:   PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&flg);
276:   if (flg && eps->problem_type == EPS_PGNHEP) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cayley spectral transformation is not compatible with PGNHEP");
277:   PetscObjectTypeCompareAny((PetscObject)eps->st,&flg,STSINVERT,STCAYLEY,"");
278:   if (flg && (eps->categ==EPS_CATEGORY_KRYLOV || eps->categ==EPS_CATEGORY_OTHER) && (eps->which!=EPS_TARGET_MAGNITUDE && eps->which!=EPS_TARGET_REAL && eps->which!=EPS_TARGET_IMAGINARY && eps->which!=EPS_ALL && eps->which!=EPS_WHICH_USER)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Shift-and-invert requires a target 'which' (see EPSSetWhichEigenpairs), for instance -st_type sinvert -eps_target 0 -eps_target_magnitude");

280:   /* process deflation and initial vectors */
281:   if (eps->nds<0) {
282:     k = -eps->nds;
283:     BVInsertConstraints(eps->V,&k,eps->defl);
284:     SlepcBasisDestroy_Private(&eps->nds,&eps->defl);
285:     eps->nds = k;
286:     STCheckNullSpace(eps->st,eps->V);
287:   }
288:   if (eps->nini<0) {
289:     k = -eps->nini;
290:     if (k>eps->ncv) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The number of initial vectors is larger than ncv");
291:     BVInsertVecs(eps->V,0,&k,eps->IS,PETSC_TRUE);
292:     SlepcBasisDestroy_Private(&eps->nini,&eps->IS);
293:     eps->nini = k;
294:   }

296:   PetscLogEventEnd(EPS_SetUp,eps,0,0,0);
297:   eps->state = EPS_STATE_SETUP;
298:   return(0);
299: }

301: /*@
302:    EPSSetOperators - Sets the matrices associated with the eigenvalue problem.

304:    Collective on EPS and Mat

306:    Input Parameters:
307: +  eps - the eigenproblem solver context
308: .  A  - the matrix associated with the eigensystem
309: -  B  - the second matrix in the case of generalized eigenproblems

311:    Notes:
312:    To specify a standard eigenproblem, use NULL for parameter B.

314:    It must be called before EPSSetUp(). If it is called again after EPSSetUp() then
315:    the EPS object is reset.

317:    Level: beginner

319: .seealso: EPSSolve(), EPSSetUp(), EPSReset(), EPSGetST(), STGetMatrix()
320: @*/
321: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
322: {
324:   PetscInt       m,n,m0,nmat;
325:   Mat            mat[2];


334:   /* Check for square matrices */
335:   MatGetSize(A,&m,&n);
336:   if (m!=n) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"A is a non-square matrix");
337:   if (B) {
338:     MatGetSize(B,&m0,&n);
339:     if (m0!=n) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"B is a non-square matrix");
340:     if (m!=m0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Dimensions of A and B do not match");
341:   }
342:   if (eps->state && n!=eps->n) { EPSReset(eps); }
343:   eps->nrma = 0.0;
344:   eps->nrmb = 0.0;
345:   if (!eps->st) { EPSGetST(eps,&eps->st); }
346:   mat[0] = A;
347:   if (B) {
348:     mat[1] = B;
349:     nmat = 2;
350:   } else nmat = 1;
351:   STSetMatrices(eps->st,nmat,mat);
352:   eps->state = EPS_STATE_INITIAL;
353:   return(0);
354: }

356: /*@
357:    EPSGetOperators - Gets the matrices associated with the eigensystem.

359:    Collective on EPS and Mat

361:    Input Parameter:
362: .  eps - the EPS context

364:    Output Parameters:
365: +  A  - the matrix associated with the eigensystem
366: -  B  - the second matrix in the case of generalized eigenproblems

368:    Level: intermediate

370: .seealso: EPSSolve(), EPSGetST(), STGetMatrix(), STSetMatrices()
371: @*/
372: PetscErrorCode EPSGetOperators(EPS eps,Mat *A,Mat *B)
373: {
375:   ST             st;
376:   PetscInt       k;

380:   EPSGetST(eps,&st);
381:   if (A) { STGetMatrix(st,0,A); }
382:   if (B) {
383:     STGetNumMatrices(st,&k);
384:     if (k==1) B = NULL;
385:     else {
386:       STGetMatrix(st,1,B);
387:     }
388:   }
389:   return(0);
390: }

392: /*@
393:    EPSSetDeflationSpace - Specify a basis of vectors that constitute the deflation
394:    space.

396:    Collective on EPS and Vec

398:    Input Parameter:
399: +  eps - the eigenproblem solver context
400: .  n   - number of vectors
401: -  v   - set of basis vectors of the deflation space

403:    Notes:
404:    When a deflation space is given, the eigensolver seeks the eigensolution
405:    in the restriction of the problem to the orthogonal complement of this
406:    space. This can be used for instance in the case that an invariant
407:    subspace is known beforehand (such as the nullspace of the matrix).

409:    These vectors do not persist from one EPSSolve() call to the other, so the
410:    deflation space should be set every time.

412:    The vectors do not need to be mutually orthonormal, since they are explicitly
413:    orthonormalized internally.

415:    Level: intermediate

417: .seealso: EPSSetInitialSpace()
418: @*/
419: PetscErrorCode EPSSetDeflationSpace(EPS eps,PetscInt n,Vec *v)
420: {

426:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n out of range");
427:   SlepcBasisReference_Private(n,v,&eps->nds,&eps->defl);
428:   if (n>0) eps->state = EPS_STATE_INITIAL;
429:   return(0);
430: }

432: /*@C
433:    EPSSetInitialSpace - Specify a basis of vectors that constitute the initial
434:    space, that is, the subspace from which the solver starts to iterate.

436:    Collective on EPS and Vec

438:    Input Parameter:
439: +  eps - the eigenproblem solver context
440: .  n   - number of vectors
441: -  is  - set of basis vectors of the initial space

443:    Notes:
444:    Some solvers start to iterate on a single vector (initial vector). In that case,
445:    the other vectors are ignored.

447:    These vectors do not persist from one EPSSolve() call to the other, so the
448:    initial space should be set every time.

450:    The vectors do not need to be mutually orthonormal, since they are explicitly
451:    orthonormalized internally.

453:    Common usage of this function is when the user can provide a rough approximation
454:    of the wanted eigenspace. Then, convergence may be faster.

456:    Level: intermediate

458: .seealso: EPSSetDeflationSpace()
459: @*/
460: PetscErrorCode EPSSetInitialSpace(EPS eps,PetscInt n,Vec *is)
461: {

467:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
468:   SlepcBasisReference_Private(n,is,&eps->nini,&eps->IS);
469:   if (n>0) eps->state = EPS_STATE_INITIAL;
470:   return(0);
471: }

473: /*
474:   EPSSetDimensions_Default - Set reasonable values for ncv, mpd if not set
475:   by the user. This is called at setup.
476:  */
477: PetscErrorCode EPSSetDimensions_Default(EPS eps,PetscInt nev,PetscInt *ncv,PetscInt *mpd)
478: {
480:   PetscBool      krylov;

483:   if (*ncv) { /* ncv set */
484:     PetscObjectTypeCompareAny((PetscObject)eps,&krylov,EPSKRYLOVSCHUR,EPSARNOLDI,EPSLANCZOS,"");
485:     if (krylov) {
486:       if (*ncv<nev+1 && !(*ncv==nev && *ncv==eps->n)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev+1");
487:     } else {
488:       if (*ncv<nev) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev");
489:     }
490:   } else if (*mpd) { /* mpd set */
491:     *ncv = PetscMin(eps->n,nev+(*mpd));
492:   } else { /* neither set: defaults depend on nev being small or large */
493:     if (nev<500) *ncv = PetscMin(eps->n,PetscMax(2*nev,nev+15));
494:     else {
495:       *mpd = 500;
496:       *ncv = PetscMin(eps->n,nev+(*mpd));
497:     }
498:   }
499:   if (!*mpd) *mpd = *ncv;
500:   return(0);
501: }

503: /*@
504:    EPSAllocateSolution - Allocate memory storage for common variables such
505:    as eigenvalues and eigenvectors.

507:    Collective on EPS

509:    Input Parameters:
510: +  eps   - eigensolver context
511: -  extra - number of additional positions, used for methods that require a
512:            working basis slightly larger than ncv

514:    Developers Note:
515:    This is SLEPC_EXTERN because it may be required by user plugin EPS
516:    implementations.

518:    Level: developer
519: @*/
520: PetscErrorCode EPSAllocateSolution(EPS eps,PetscInt extra)
521: {
523:   PetscInt       oldsize,newc,requested;
524:   PetscLogDouble cnt;
525:   Vec            t;

528:   requested = eps->ncv + extra;

530:   /* oldsize is zero if this is the first time setup is called */
531:   BVGetSizes(eps->V,NULL,NULL,&oldsize);
532:   newc = PetscMax(0,requested-oldsize);

534:   /* allocate space for eigenvalues and friends */
535:   if (requested != oldsize || !eps->eigr) {
536:     PetscFree4(eps->eigr,eps->eigi,eps->errest,eps->perm);
537:     PetscMalloc4(requested,&eps->eigr,requested,&eps->eigi,requested,&eps->errest,requested,&eps->perm);
538:     cnt = 2*newc*sizeof(PetscScalar) + 2*newc*sizeof(PetscReal) + newc*sizeof(PetscInt);
539:     PetscLogObjectMemory((PetscObject)eps,cnt);
540:   }

542:   /* workspace for the case of arbitrary selection */
543:   if (eps->arbitrary) {
544:     if (eps->rr) {
545:       PetscFree2(eps->rr,eps->ri);
546:     }
547:     PetscMalloc2(requested,&eps->rr,requested,&eps->ri);
548:     PetscLogObjectMemory((PetscObject)eps,2*newc*sizeof(PetscScalar));
549:   }

551:   /* allocate V */
552:   if (!eps->V) { EPSGetBV(eps,&eps->V); }
553:   if (!oldsize) {
554:     if (!((PetscObject)(eps->V))->type_name) {
555:       BVSetType(eps->V,BVSVEC);
556:     }
557:     STMatCreateVecsEmpty(eps->st,&t,NULL);
558:     BVSetSizesFromVec(eps->V,t,requested);
559:     VecDestroy(&t);
560:   } else {
561:     BVResize(eps->V,requested,PETSC_FALSE);
562:   }

564:   /* allocate W */
565:   if (eps->twosided) {
566:     BVDestroy(&eps->W);
567:     BVDuplicate(eps->V,&eps->W);
568:   }
569:   return(0);
570: }