Actual source code: epssetup.c

slepc-3.7.4 2017-05-17
Report Typos and Errors
  1: /*
  2:       EPS routines related to problem setup.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

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

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

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

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

 28: /*@
 29:    EPSSetUp - Sets up all the internal data structures necessary for the
 30:    execution of the eigensolver. Then calls STSetUp() for any set-up
 31:    operations associated to the ST object.

 33:    Collective on EPS

 35:    Input Parameter:
 36: .  eps   - eigenproblem solver context

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

 43:    Level: developer

 45: .seealso: EPSCreate(), EPSSolve(), EPSDestroy(), STSetUp(), EPSSetInitialSpace()
 46: @*/
 47: PetscErrorCode EPSSetUp(EPS eps)
 48: {
 50:   Mat            A,B;
 51:   SlepcSC        sc;
 52:   PetscInt       k,nmat;
 53:   PetscBool      flg,istrivial;
 54: #if defined(PETSC_USE_COMPLEX)
 55:   PetscScalar    sigma;
 56: #endif

 60:   if (eps->state) return(0);
 61:   PetscLogEventBegin(EPS_SetUp,eps,0,0,0);

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

 66:   /* Set default solver type (EPSSetFromOptions was not called) */
 67:   if (!((PetscObject)eps)->type_name) {
 68:     EPSSetType(eps,EPSKRYLOVSCHUR);
 69:   }
 70:   if (!eps->st) { EPSGetST(eps,&eps->st); }
 71:   if (!((PetscObject)eps->st)->type_name) {
 72:     PetscObjectTypeCompareAny((PetscObject)eps,&flg,EPSGD,EPSJD,EPSRQCG,EPSBLOPEX,EPSPRIMME,"");
 73:     STSetType(eps->st,flg?STPRECOND:STSHIFT);
 74:   }
 75:   STSetTransform(eps->st,PETSC_TRUE);
 76:   if (!eps->ds) { EPSGetDS(eps,&eps->ds); }
 77:   DSReset(eps->ds);
 78:   if (!eps->rg) { EPSGetRG(eps,&eps->rg); }
 79:   if (!((PetscObject)eps->rg)->type_name) {
 80:     RGSetType(eps->rg,RGINTERVAL);
 81:   }

 83:   /* Set problem dimensions */
 84:   STGetNumMatrices(eps->st,&nmat);
 85:   if (!nmat) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"EPSSetOperators must be called first");
 86:   STMatGetSize(eps->st,&eps->n,NULL);
 87:   STMatGetLocalSize(eps->st,&eps->nloc,NULL);

 89:   /* Set default problem type */
 90:   if (!eps->problem_type) {
 91:     if (nmat==1) {
 92:       EPSSetProblemType(eps,EPS_NHEP);
 93:     } else {
 94:       EPSSetProblemType(eps,EPS_GNHEP);
 95:     }
 96:   } else if (nmat==1 && eps->isgeneralized) {
 97:     PetscInfo(eps,"Eigenproblem set as generalized but no matrix B was provided; reverting to a standard eigenproblem\n");
 98:     eps->isgeneralized = PETSC_FALSE;
 99:     eps->problem_type = eps->ishermitian? EPS_HEP: EPS_NHEP;
100:   } else if (nmat>1 && !eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Inconsistent EPS state");

102:   if (eps->nev > eps->n) eps->nev = eps->n;
103:   if (eps->ncv > eps->n) eps->ncv = eps->n;

105:   /* initialization of matrix norms */
106:   if (eps->conv==EPS_CONV_NORM) {
107:     if (!eps->nrma) {
108:       STGetOperators(eps->st,0,&A);
109:       MatNorm(A,NORM_INFINITY,&eps->nrma);
110:     }
111:     if (nmat>1 && !eps->nrmb) {
112:       STGetOperators(eps->st,1,&B);
113:       MatNorm(B,NORM_INFINITY,&eps->nrmb);
114:     }
115:   }

117:   /* call specific solver setup */
118:   (*eps->ops->setup)(eps);

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

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

127:   /* fill sorting criterion context */
128:   switch (eps->which) {
129:     case EPS_LARGEST_MAGNITUDE:
130:       eps->sc->comparison    = SlepcCompareLargestMagnitude;
131:       eps->sc->comparisonctx = NULL;
132:       break;
133:     case EPS_SMALLEST_MAGNITUDE:
134:       eps->sc->comparison    = SlepcCompareSmallestMagnitude;
135:       eps->sc->comparisonctx = NULL;
136:       break;
137:     case EPS_LARGEST_REAL:
138:       eps->sc->comparison    = SlepcCompareLargestReal;
139:       eps->sc->comparisonctx = NULL;
140:       break;
141:     case EPS_SMALLEST_REAL:
142:       eps->sc->comparison    = SlepcCompareSmallestReal;
143:       eps->sc->comparisonctx = NULL;
144:       break;
145:     case EPS_LARGEST_IMAGINARY:
146:       eps->sc->comparison    = SlepcCompareLargestImaginary;
147:       eps->sc->comparisonctx = NULL;
148:       break;
149:     case EPS_SMALLEST_IMAGINARY:
150:       eps->sc->comparison    = SlepcCompareSmallestImaginary;
151:       eps->sc->comparisonctx = NULL;
152:       break;
153:     case EPS_TARGET_MAGNITUDE:
154:       eps->sc->comparison    = SlepcCompareTargetMagnitude;
155:       eps->sc->comparisonctx = &eps->target;
156:       break;
157:     case EPS_TARGET_REAL:
158:       eps->sc->comparison    = SlepcCompareTargetReal;
159:       eps->sc->comparisonctx = &eps->target;
160:       break;
161:     case EPS_TARGET_IMAGINARY:
162:       eps->sc->comparison    = SlepcCompareTargetImaginary;
163:       eps->sc->comparisonctx = &eps->target;
164:       break;
165:     case EPS_ALL:
166:       eps->sc->comparison    = SlepcCompareSmallestReal;
167:       eps->sc->comparisonctx = NULL;
168:       break;
169:     case EPS_WHICH_USER:
170:       break;
171:   }
172:   eps->sc->map    = NULL;
173:   eps->sc->mapobj = NULL;

175:   /* fill sorting criterion for DS */
176:   DSGetSlepcSC(eps->ds,&sc);
177:   RGIsTrivial(eps->rg,&istrivial);
178:   if (eps->which==EPS_ALL) {
179:     sc->rg            = NULL;
180:     sc->comparison    = SlepcCompareLargestMagnitude;
181:     sc->comparisonctx = NULL;
182:     sc->map           = NULL;
183:     sc->mapobj        = NULL;
184:   } else {
185:     sc->rg            = istrivial? NULL: eps->rg;
186:     sc->comparison    = eps->sc->comparison;
187:     sc->comparisonctx = eps->sc->comparisonctx;
188:     sc->map           = SlepcMap_ST;
189:     sc->mapobj        = (PetscObject)eps->st;
190:   }

192:   /* Build balancing matrix if required */
193:   if (!eps->ishermitian && (eps->balance==EPS_BALANCE_ONESIDE || eps->balance==EPS_BALANCE_TWOSIDE)) {
194:     if (!eps->D) {
195:       BVCreateVec(eps->V,&eps->D);
196:       PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->D);
197:     } else {
198:       VecSet(eps->D,1.0);
199:     }
200:     EPSBuildBalance_Krylov(eps);
201:     STSetBalanceMatrix(eps->st,eps->D);
202:   }

204:   /* Setup ST */
205:   STSetUp(eps->st);

207: #if defined(PETSC_USE_COMPLEX)
208:   STGetShift(eps->st,&sigma);
209:   if (eps->ishermitian && PetscImaginaryPart(sigma) != 0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Hermitian problems are not compatible with complex shifts");
210: #endif
211:   PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&flg);
212:   if (flg && eps->problem_type == EPS_PGNHEP) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cayley spectral transformation is not compatible with PGNHEP");

214:   /* process deflation and initial vectors */
215:   if (eps->nds<0) {
216:     k = -eps->nds;
217:     BVInsertConstraints(eps->V,&k,eps->defl);
218:     SlepcBasisDestroy_Private(&eps->nds,&eps->defl);
219:     eps->nds = k;
220:     STCheckNullSpace(eps->st,eps->V);
221:   }
222:   if (eps->nini<0) {
223:     k = -eps->nini;
224:     if (k>eps->ncv) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The number of initial vectors is larger than ncv");
225:     BVInsertVecs(eps->V,0,&k,eps->IS,PETSC_TRUE);
226:     SlepcBasisDestroy_Private(&eps->nini,&eps->IS);
227:     eps->nini = k;
228:   }

230:   PetscLogEventEnd(EPS_SetUp,eps,0,0,0);
231:   eps->state = EPS_STATE_SETUP;
232:   return(0);
233: }

237: /*@
238:    EPSSetOperators - Sets the matrices associated with the eigenvalue problem.

240:    Collective on EPS and Mat

242:    Input Parameters:
243: +  eps - the eigenproblem solver context
244: .  A  - the matrix associated with the eigensystem
245: -  B  - the second matrix in the case of generalized eigenproblems

247:    Notes:
248:    To specify a standard eigenproblem, use NULL for parameter B.

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

253:    Level: beginner

255: .seealso: EPSSolve(), EPSSetUp(), EPSReset(), EPSGetST(), STGetOperators()
256: @*/
257: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
258: {
260:   PetscInt       m,n,m0,nmat;
261:   Mat            mat[2];


270:   /* Check for square matrices */
271:   MatGetSize(A,&m,&n);
272:   if (m!=n) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"A is a non-square matrix");
273:   if (B) {
274:     MatGetSize(B,&m0,&n);
275:     if (m0!=n) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"B is a non-square matrix");
276:     if (m!=m0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Dimensions of A and B do not match");
277:   }
278:   if (eps->state && n!=eps->n) { EPSReset(eps); }
279:   eps->nrma = 0.0;
280:   eps->nrmb = 0.0;
281:   if (!eps->st) { EPSGetST(eps,&eps->st); }
282:   mat[0] = A;
283:   if (B) {
284:     mat[1] = B;
285:     nmat = 2;
286:   } else nmat = 1;
287:   STSetOperators(eps->st,nmat,mat);
288:   eps->state = EPS_STATE_INITIAL;
289:   return(0);
290: }

294: /*@
295:    EPSGetOperators - Gets the matrices associated with the eigensystem.

297:    Collective on EPS and Mat

299:    Input Parameter:
300: .  eps - the EPS context

302:    Output Parameters:
303: +  A  - the matrix associated with the eigensystem
304: -  B  - the second matrix in the case of generalized eigenproblems

306:    Level: intermediate

308: .seealso: EPSSolve(), EPSGetST(), STGetOperators(), STSetOperators()
309: @*/
310: PetscErrorCode EPSGetOperators(EPS eps,Mat *A,Mat *B)
311: {
313:   ST             st;
314:   PetscInt       k;

318:   EPSGetST(eps,&st);
319:   if (A) { STGetOperators(st,0,A); }
320:   if (B) {
321:     STGetNumMatrices(st,&k);
322:     if (k==1) B = NULL;
323:     else {
324:       STGetOperators(st,1,B);
325:     }
326:   }
327:   return(0);
328: }

332: /*@
333:    EPSSetDeflationSpace - Specify a basis of vectors that constitute the deflation
334:    space.

336:    Collective on EPS and Vec

338:    Input Parameter:
339: +  eps - the eigenproblem solver context
340: .  n   - number of vectors
341: -  v   - set of basis vectors of the deflation space

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

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

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

355:    Level: intermediate

357: .seealso: EPSSetInitialSpace()
358: @*/
359: PetscErrorCode EPSSetDeflationSpace(EPS eps,PetscInt n,Vec *v)
360: {

366:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n out of range");
367:   SlepcBasisReference_Private(n,v,&eps->nds,&eps->defl);
368:   if (n>0) eps->state = EPS_STATE_INITIAL;
369:   return(0);
370: }

374: /*@
375:    EPSSetInitialSpace - Specify a basis of vectors that constitute the initial
376:    space, that is, the subspace from which the solver starts to iterate.

378:    Collective on EPS and Vec

380:    Input Parameter:
381: +  eps - the eigenproblem solver context
382: .  n   - number of vectors
383: -  is  - set of basis vectors of the initial space

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

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

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

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

398:    Level: intermediate

400: .seealso: EPSSetDeflationSpace()
401: @*/
402: PetscErrorCode EPSSetInitialSpace(EPS eps,PetscInt n,Vec *is)
403: {

409:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
410:   SlepcBasisReference_Private(n,is,&eps->nini,&eps->IS);
411:   if (n>0) eps->state = EPS_STATE_INITIAL;
412:   return(0);
413: }

417: /*
418:   EPSSetDimensions_Default - Set reasonable values for ncv, mpd if not set
419:   by the user. This is called at setup.
420:  */
421: PetscErrorCode EPSSetDimensions_Default(EPS eps,PetscInt nev,PetscInt *ncv,PetscInt *mpd)
422: {
424:   PetscBool      krylov;

427:   if (*ncv) { /* ncv set */
428:     PetscObjectTypeCompareAny((PetscObject)eps,&krylov,EPSKRYLOVSCHUR,EPSARNOLDI,EPSLANCZOS,"");
429:     if (krylov) {
430:       if (*ncv<nev+1 && !(*ncv==nev && *ncv==eps->n)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev+1");
431:     } else {
432:       if (*ncv<nev) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev");
433:     }
434:   } else if (*mpd) { /* mpd set */
435:     *ncv = PetscMin(eps->n,nev+(*mpd));
436:   } else { /* neither set: defaults depend on nev being small or large */
437:     if (nev<500) *ncv = PetscMin(eps->n,PetscMax(2*nev,nev+15));
438:     else {
439:       *mpd = 500;
440:       *ncv = PetscMin(eps->n,nev+(*mpd));
441:     }
442:   }
443:   if (!*mpd) *mpd = *ncv;
444:   return(0);
445: }

449: /*@
450:    EPSAllocateSolution - Allocate memory storage for common variables such
451:    as eigenvalues and eigenvectors.

453:    Collective on EPS

455:    Input Parameters:
456: +  eps   - eigensolver context
457: -  extra - number of additional positions, used for methods that require a
458:            working basis slightly larger than ncv

460:    Developers Note:
461:    This is PETSC_EXTERN because it may be required by user plugin EPS
462:    implementations.

464:    Level: developer
465: @*/
466: PetscErrorCode EPSAllocateSolution(EPS eps,PetscInt extra)
467: {
469:   PetscInt       oldsize,newc,requested;
470:   PetscLogDouble cnt;
471:   Vec            t;

474:   requested = eps->ncv + extra;

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

480:   /* allocate space for eigenvalues and friends */
481:   if (requested != oldsize || !eps->eigr) {
482:     if (oldsize) {
483:       PetscFree4(eps->eigr,eps->eigi,eps->errest,eps->perm);
484:     }
485:     PetscMalloc4(requested,&eps->eigr,requested,&eps->eigi,requested,&eps->errest,requested,&eps->perm);
486:     cnt = 2*newc*sizeof(PetscScalar) + 2*newc*sizeof(PetscReal) + newc*sizeof(PetscInt);
487:     PetscLogObjectMemory((PetscObject)eps,cnt);
488:   }

490:   /* workspace for the case of arbitrary selection */
491:   if (eps->arbitrary) {
492:     if (eps->rr) {
493:       PetscFree2(eps->rr,eps->ri);
494:     }
495:     PetscMalloc2(requested,&eps->rr,requested,&eps->ri);
496:     PetscLogObjectMemory((PetscObject)eps,2*newc*sizeof(PetscScalar));
497:   }

499:   /* allocate V */
500:   if (!eps->V) { EPSGetBV(eps,&eps->V); }
501:   if (!oldsize) {
502:     if (!((PetscObject)(eps->V))->type_name) {
503:       BVSetType(eps->V,BVSVEC);
504:     }
505:     STMatCreateVecs(eps->st,&t,NULL);
506:     BVSetSizesFromVec(eps->V,t,requested);
507:     VecDestroy(&t);
508:   } else {
509:     BVResize(eps->V,requested,PETSC_FALSE);
510:   }
511:   return(0);
512: }