Actual source code: epssetup.c

slepc-3.5.2 2014-10-10
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-2014, 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: advanced

 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:   if (!eps->ds) { EPSGetDS(eps,&eps->ds); }
 76:   DSReset(eps->ds);
 77:   if (!eps->rg) { EPSGetRG(eps,&eps->rg); }
 78:   if (!((PetscObject)eps->rg)->type_name) {
 79:     RGSetType(eps->rg,RGINTERVAL);
 80:   }
 81:   if (!((PetscObject)eps->rand)->type_name) {
 82:     PetscRandomSetFromOptions(eps->rand);
 83:   }

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

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

104:   if (eps->nev > eps->n) eps->nev = eps->n;
105:   if (eps->ncv > eps->n) eps->ncv = eps->n;

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

119:   /* call specific solver setup */
120:   (*eps->ops->setup)(eps);

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

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

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

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

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

206:   /* Setup ST */
207:   STSetUp(eps->st);

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

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

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

239: /*@
240:    EPSSetOperators - Sets the matrices associated with the eigenvalue problem.

242:    Collective on EPS and Mat

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

249:    Notes:
250:    To specify a standard eigenproblem, use NULL for parameter B.

252:    It must be called after EPSSetUp(). If it is called again after EPSSetUp() then
253:    the EPS object is reset.

255:    Level: beginner

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


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

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

298:    Collective on EPS and Mat

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

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

307:    Level: intermediate

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

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

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

337:    Collective on EPS and Vec

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

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

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

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

356:    Level: intermediate

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

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

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

379:    Collective on EPS and Vec

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

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

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

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

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

399:    Level: intermediate

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

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

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

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

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

454:    Collective on EPS

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

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

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

475:   requested = eps->ncv + extra;

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

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

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

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