Actual source code: epssetup.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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>

 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: {
 24:   PetscFunctionBegin;
 25:   PetscTryTypeMethod(eps,setdefaultst);
 26:   if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STSHIFT));
 27:   PetscFunctionReturn(PETSC_SUCCESS);
 28: }

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

 38:   PetscFunctionBegin;
 39:   if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STPRECOND));
 40:   PetscCall(STGetKSP(eps->st,&ksp));
 41:   if (!((PetscObject)ksp)->type_name) PetscCall(KSPSetType(ksp,KSPPREONLY));
 42:   PetscFunctionReturn(PETSC_SUCCESS);
 43: }

 45: /*
 46:    This is done by preconditioned eigensolvers that can also use the KSP.
 47:    It sets STPRECOND with the default KSP (GMRES) and maxit=5.
 48: */
 49: PetscErrorCode EPSSetDefaultST_GMRES(EPS eps)
 50: {
 51:   KSP            ksp;

 53:   PetscFunctionBegin;
 54:   if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STPRECOND));
 55:   PetscCall(STPrecondSetKSPHasMat(eps->st,PETSC_TRUE));
 56:   PetscCall(STGetKSP(eps->st,&ksp));
 57:   if (!((PetscObject)ksp)->type_name) {
 58:     PetscCall(KSPSetType(ksp,KSPGMRES));
 59:     PetscCall(KSPSetTolerances(ksp,PETSC_CURRENT,PETSC_CURRENT,PETSC_CURRENT,5));
 60:   }
 61:   PetscFunctionReturn(PETSC_SUCCESS);
 62: }

 64: #if defined(SLEPC_HAVE_SCALAPACK) || defined(SLEPC_HAVE_ELPA) || defined(SLEPC_HAVE_ELEMENTAL) || defined(SLEPC_HAVE_EVSL)
 65: /*
 66:    This is for direct eigensolvers that work with A and B directly, so
 67:    no need to factorize B.
 68: */
 69: PetscErrorCode EPSSetDefaultST_NoFactor(EPS eps)
 70: {
 71:   KSP            ksp;
 72:   PC             pc;

 74:   PetscFunctionBegin;
 75:   if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STSHIFT));
 76:   PetscCall(STGetKSP(eps->st,&ksp));
 77:   if (!((PetscObject)ksp)->type_name) PetscCall(KSPSetType(ksp,KSPPREONLY));
 78:   PetscCall(KSPGetPC(ksp,&pc));
 79:   if (!((PetscObject)pc)->type_name) PetscCall(PCSetType(pc,PCNONE));
 80:   PetscFunctionReturn(PETSC_SUCCESS);
 81: }
 82: #endif

 84: /*
 85:    Check that the ST selected by the user is compatible with the EPS solver and options
 86: */
 87: static PetscErrorCode EPSCheckCompatibleST(EPS eps)
 88: {
 89:   PetscBool      precond,shift,sinvert,cayley,lyapii;
 90: #if defined(PETSC_USE_COMPLEX)
 91:   PetscScalar    sigma;
 92: #endif

 94:   PetscFunctionBegin;
 95:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&precond));
 96:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&shift));
 97:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&sinvert));
 98:   PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&cayley));

100:   /* preconditioned eigensolvers */
101:   PetscCheck(eps->categ!=EPS_CATEGORY_PRECOND || precond,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires ST=PRECOND");
102:   PetscCheck(eps->categ==EPS_CATEGORY_PRECOND || !precond,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"STPRECOND is intended for preconditioned eigensolvers only");

104:   /* harmonic extraction */
105:   PetscCheck(precond || shift || !eps->extraction || eps->extraction==EPS_RITZ,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cannot use a spectral transformation combined with harmonic extraction");

107:   /* real shifts in Hermitian problems */
108: #if defined(PETSC_USE_COMPLEX)
109:   PetscCall(STGetShift(eps->st,&sigma));
110:   PetscCheck(!eps->ishermitian || PetscImaginaryPart(sigma)==0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Hermitian problems are not compatible with complex shifts");
111: #endif

113:   /* Cayley with PGNHEP */
114:   PetscCheck(!cayley || eps->problem_type!=EPS_PGNHEP,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cayley spectral transformation is not compatible with PGNHEP");

116:   /* make sure that the user does not specify smallest magnitude with shift-and-invert */
117:   if ((cayley || sinvert) && (eps->categ==EPS_CATEGORY_KRYLOV || eps->categ==EPS_CATEGORY_OTHER)) {
118:     PetscCall(PetscObjectTypeCompare((PetscObject)eps,EPSLYAPII,&lyapii));
119:     PetscCheck(lyapii || eps->which==EPS_TARGET_MAGNITUDE || eps->which==EPS_TARGET_REAL || eps->which==EPS_TARGET_IMAGINARY || eps->which==EPS_ALL || eps->which==EPS_WHICH_USER,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"Shift-and-invert requires a target 'which' (see EPSSetWhichEigenpairs), for instance -st_type sinvert -eps_target 0 -eps_target_magnitude");
120:   }
121:   PetscFunctionReturn(PETSC_SUCCESS);
122: }

124: /*
125:    MatEstimateSpectralRange_EPS: estimate the spectral range [left,right] of a
126:    symmetric/Hermitian matrix A using an auxiliary EPS object
127: */
128: PetscErrorCode MatEstimateSpectralRange_EPS(Mat A,PetscReal *left,PetscReal *right)
129: {
130:   PetscInt       nconv;
131:   PetscScalar    eig0;
132:   PetscReal      tol=1e-3,errest=tol;
133:   EPS            eps;

135:   PetscFunctionBegin;
136:   *left = 0.0; *right = 0.0;
137:   PetscCall(EPSCreate(PetscObjectComm((PetscObject)A),&eps));
138:   PetscCall(EPSSetOptionsPrefix(eps,"eps_filter_"));
139:   PetscCall(EPSSetOperators(eps,A,NULL));
140:   PetscCall(EPSSetProblemType(eps,EPS_HEP));
141:   PetscCall(EPSSetTolerances(eps,tol,50));
142:   PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_ABS));
143:   PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
144:   PetscCall(EPSSolve(eps));
145:   PetscCall(EPSGetConverged(eps,&nconv));
146:   if (nconv>0) {
147:     PetscCall(EPSGetEigenvalue(eps,0,&eig0,NULL));
148:     PetscCall(EPSGetErrorEstimate(eps,0,&errest));
149:   } else eig0 = eps->eigr[0];
150:   *left = PetscRealPart(eig0)-errest;
151:   PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL));
152:   PetscCall(EPSSolve(eps));
153:   PetscCall(EPSGetConverged(eps,&nconv));
154:   if (nconv>0) {
155:     PetscCall(EPSGetEigenvalue(eps,0,&eig0,NULL));
156:     PetscCall(EPSGetErrorEstimate(eps,0,&errest));
157:   } else eig0 = eps->eigr[0];
158:   *right = PetscRealPart(eig0)+errest;
159:   PetscCall(EPSDestroy(&eps));
160:   PetscFunctionReturn(PETSC_SUCCESS);
161: }

163: /*
164:    EPSSetUpSort_Basic: configure the EPS sorting criterion according to 'which'
165: */
166: PetscErrorCode EPSSetUpSort_Basic(EPS eps)
167: {
168:   PetscFunctionBegin;
169:   switch (eps->which) {
170:     case EPS_LARGEST_MAGNITUDE:
171:       eps->sc->comparison    = SlepcCompareLargestMagnitude;
172:       eps->sc->comparisonctx = NULL;
173:       break;
174:     case EPS_SMALLEST_MAGNITUDE:
175:       eps->sc->comparison    = SlepcCompareSmallestMagnitude;
176:       eps->sc->comparisonctx = NULL;
177:       break;
178:     case EPS_LARGEST_REAL:
179:       eps->sc->comparison    = SlepcCompareLargestReal;
180:       eps->sc->comparisonctx = NULL;
181:       break;
182:     case EPS_SMALLEST_REAL:
183:       eps->sc->comparison    = SlepcCompareSmallestReal;
184:       eps->sc->comparisonctx = NULL;
185:       break;
186:     case EPS_LARGEST_IMAGINARY:
187:       eps->sc->comparison    = SlepcCompareLargestImaginary;
188:       eps->sc->comparisonctx = NULL;
189:       break;
190:     case EPS_SMALLEST_IMAGINARY:
191:       eps->sc->comparison    = SlepcCompareSmallestImaginary;
192:       eps->sc->comparisonctx = NULL;
193:       break;
194:     case EPS_TARGET_MAGNITUDE:
195:       eps->sc->comparison    = SlepcCompareTargetMagnitude;
196:       eps->sc->comparisonctx = &eps->target;
197:       break;
198:     case EPS_TARGET_REAL:
199:       eps->sc->comparison    = SlepcCompareTargetReal;
200:       eps->sc->comparisonctx = &eps->target;
201:       break;
202:     case EPS_TARGET_IMAGINARY:
203:       eps->sc->comparison    = SlepcCompareTargetImaginary;
204:       eps->sc->comparisonctx = &eps->target;
205:       break;
206:     case EPS_ALL:
207:       eps->sc->comparison    = SlepcCompareSmallestReal;
208:       eps->sc->comparisonctx = NULL;
209:       break;
210:     case EPS_WHICH_USER:
211:       break;
212:   }
213:   eps->sc->map    = NULL;
214:   eps->sc->mapobj = NULL;
215:   PetscFunctionReturn(PETSC_SUCCESS);
216: }

218: /*
219:    EPSSetUpSort_Default: configure both EPS and DS sorting criterion
220: */
221: PetscErrorCode EPSSetUpSort_Default(EPS eps)
222: {
223:   SlepcSC        sc;
224:   PetscBool      istrivial;

226:   PetscFunctionBegin;
227:   /* fill sorting criterion context */
228:   PetscCall(EPSSetUpSort_Basic(eps));
229:   /* fill sorting criterion for DS */
230:   PetscCall(DSGetSlepcSC(eps->ds,&sc));
231:   PetscCall(RGIsTrivial(eps->rg,&istrivial));
232:   sc->rg            = istrivial? NULL: eps->rg;
233:   sc->comparison    = eps->sc->comparison;
234:   sc->comparisonctx = eps->sc->comparisonctx;
235:   sc->map           = SlepcMap_ST;
236:   sc->mapobj        = (PetscObject)eps->st;
237:   PetscFunctionReturn(PETSC_SUCCESS);
238: }

240: /*@
241:    EPSSetDSType - Sets the type of the internal DS object based on the current
242:    settings of the eigenvalue solver.

244:    Collective

246:    Input Parameter:
247: .  eps - eigenproblem solver context

249:    Note:
250:    This function need not be called explicitly, since it will be called at
251:    both EPSSetFromOptions() and EPSSetUp().

253:    Level: developer

255: .seealso: `EPSSetFromOptions()`, `EPSSetUp()`
256: @*/
257: PetscErrorCode EPSSetDSType(EPS eps)
258: {
259:   PetscFunctionBegin;
261:   PetscTryTypeMethod(eps,setdstype);
262:   PetscFunctionReturn(PETSC_SUCCESS);
263: }

265: /*@
266:    EPSSetUp - Sets up all the internal data structures necessary for the
267:    execution of the eigensolver. Then calls STSetUp() for any set-up
268:    operations associated to the ST object.

270:    Collective

272:    Input Parameter:
273: .  eps   - eigenproblem solver context

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

280:    Level: developer

282: .seealso: `EPSCreate()`, `EPSSolve()`, `EPSDestroy()`, `STSetUp()`, `EPSSetInitialSpace()`
283: @*/
284: PetscErrorCode EPSSetUp(EPS eps)
285: {
286:   Mat            A,B;
287:   PetscInt       k,nmat;
288:   PetscBool      flg;
289:   EPSStoppingCtx ctx;

291:   PetscFunctionBegin;
293:   if (eps->state) PetscFunctionReturn(PETSC_SUCCESS);
294:   PetscCall(PetscLogEventBegin(EPS_SetUp,eps,0,0,0));

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

299:   /* Set default solver type (EPSSetFromOptions was not called) */
300:   if (!((PetscObject)eps)->type_name) PetscCall(EPSSetType(eps,EPSKRYLOVSCHUR));
301:   if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
302:   PetscCall(EPSSetDefaultST(eps));

304:   PetscCall(STSetTransform(eps->st,PETSC_TRUE));
305:   PetscCall(STSetStructured(eps->st,PETSC_FALSE));
306:   if (eps->useds && !eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
307:   if (eps->useds) PetscCall(EPSSetDSType(eps));
308:   if (eps->twosided) {
309:     PetscCheck(!eps->ishermitian || (eps->isgeneralized && !eps->ispositive),PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Two-sided methods are not intended for %s problems",SLEPC_STRING_HERMITIAN);
310:   }
311:   if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
312:   if (!((PetscObject)eps->rg)->type_name) PetscCall(RGSetType(eps->rg,RGINTERVAL));

314:   /* Set problem dimensions */
315:   PetscCall(STGetNumMatrices(eps->st,&nmat));
316:   PetscCheck(nmat,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"EPSSetOperators must be called first");
317:   PetscCall(STMatGetSize(eps->st,&eps->n,NULL));
318:   PetscCall(STMatGetLocalSize(eps->st,&eps->nloc,NULL));

320:   /* Set default problem type */
321:   if (!eps->problem_type) {
322:     if (nmat==1) PetscCall(EPSSetProblemType(eps,EPS_NHEP));
323:     else PetscCall(EPSSetProblemType(eps,EPS_GNHEP));
324:   } else if (nmat==1 && eps->isgeneralized) {
325:     PetscCall(PetscInfo(eps,"Eigenproblem set as generalized but no matrix B was provided; reverting to a standard eigenproblem\n"));
326:     eps->isgeneralized = PETSC_FALSE;
327:     eps->problem_type = eps->ishermitian? EPS_HEP: EPS_NHEP;
328:   } else PetscCheck(nmat==1 || eps->isgeneralized,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Inconsistent EPS state: the problem type does not match the number of matrices");

330:   if (eps->isstructured) {
331:     /* make sure the user has set the appropriate matrix */
332:     PetscCall(STGetMatrix(eps->st,0,&A));
333:     if (eps->problem_type==EPS_BSE) PetscCall(SlepcCheckMatStruct(A,SLEPC_MAT_STRUCT_BSE,NULL));
334:     if (eps->problem_type==EPS_HAMILT) PetscCall(SlepcCheckMatStruct(A,SLEPC_MAT_STRUCT_HAMILT,NULL));
335:   }

337:   /* safeguard for small problems */
338:   if (eps->n == 0) PetscFunctionReturn(PETSC_SUCCESS);
339:   if (eps->nev > eps->n) eps->nev = eps->n;
340:   if (eps->problem_type == EPS_BSE) {
341:     if (2*eps->ncv > eps->n) eps->ncv = eps->n/2;
342:   } else {
343:     if (eps->ncv > eps->n) eps->ncv = eps->n;
344:   }

346:   /* check some combinations of eps->which */
347:   PetscCheck(!eps->ishermitian || (eps->isgeneralized && !eps->ispositive) || (eps->which!=EPS_LARGEST_IMAGINARY && eps->which!=EPS_SMALLEST_IMAGINARY && eps->which!=EPS_TARGET_IMAGINARY),PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Sorting the eigenvalues along the imaginary axis is not allowed when all eigenvalues are real");

349:   /* initialization of matrix norms */
350:   if (eps->conv==EPS_CONV_NORM) {
351:     if (!eps->nrma) {
352:       PetscCall(STGetMatrix(eps->st,0,&A));
353:       PetscCall(MatNorm(A,NORM_INFINITY,&eps->nrma));
354:     }
355:     if (nmat>1 && !eps->nrmb) {
356:       PetscCall(STGetMatrix(eps->st,1,&B));
357:       PetscCall(MatNorm(B,NORM_INFINITY,&eps->nrmb));
358:     }
359:   }

361:   /* call specific solver setup */
362:   PetscUseTypeMethod(eps,setup);

364:   /* threshold stopping test */
365:   if (eps->stop==EPS_STOP_THRESHOLD) {
366:     PetscCheck(eps->thres!=PETSC_MIN_REAL,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"Must give a threshold value with EPSSetThreshold()");
367:     PetscCheck(eps->which==EPS_LARGEST_MAGNITUDE || eps->which==EPS_SMALLEST_MAGNITUDE || eps->which==EPS_LARGEST_REAL || eps->which==EPS_SMALLEST_REAL || eps->which==EPS_TARGET_MAGNITUDE,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Threshold stopping test can only be used with largest/smallest/target magnitude or largest/smallest real selection of eigenvalues");
368:     if (eps->which==EPS_LARGEST_REAL || eps->which==EPS_SMALLEST_REAL) PetscCheck(eps->problem_type==EPS_HEP || eps->problem_type==EPS_GHEP || eps->problem_type==EPS_BSE,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Threshold stopping test with largest/smallest real can only be used in problems that have all eigenvaues real");
369:     else PetscCheck(eps->thres>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"In case of largest/smallest/target magnitude the threshold value must be positive");
370:     PetscCheck(eps->which==EPS_LARGEST_MAGNITUDE || eps->which==EPS_TARGET_MAGNITUDE || !eps->threlative,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"Can only use a relative threshold with largest/target magnitude selection of eigenvalues");
371:     PetscCall(PetscNew(&ctx));
372:     ctx->thres      = eps->thres;
373:     ctx->threlative = eps->threlative;
374:     ctx->which      = eps->which;
375:     PetscCall(EPSSetStoppingTestFunction(eps,EPSStoppingThreshold,ctx,PetscCtxDestroyDefault));
376:   }

378:   /* if purification is set, check that it really makes sense */
379:   if (eps->purify) {
380:     if (eps->categ==EPS_CATEGORY_PRECOND || eps->categ==EPS_CATEGORY_CONTOUR) eps->purify = PETSC_FALSE;
381:     else {
382:       if (!eps->isgeneralized) eps->purify = PETSC_FALSE;
383:       else if (!eps->ishermitian && !eps->ispositive) eps->purify = PETSC_FALSE;
384:       else {
385:         PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&flg));
386:         if (flg) eps->purify = PETSC_FALSE;
387:       }
388:     }
389:   }

391:   /* set tolerance if not yet set */
392:   if (eps->tol==(PetscReal)PETSC_DETERMINE) eps->tol = SLEPC_DEFAULT_TOL;

394:   /* set up sorting criterion */
395:   PetscTryTypeMethod(eps,setupsort);

397:   /* Build balancing matrix if required */
398:   if (eps->balance!=EPS_BALANCE_USER) {
399:     PetscCall(STSetBalanceMatrix(eps->st,NULL));
400:     if (!eps->ishermitian && (eps->balance==EPS_BALANCE_ONESIDE || eps->balance==EPS_BALANCE_TWOSIDE)) {
401:       if (!eps->D) PetscCall(BVCreateVec(eps->V,&eps->D));
402:       PetscCall(EPSBuildBalance_Krylov(eps));
403:       PetscCall(STSetBalanceMatrix(eps->st,eps->D));
404:     }
405:   }

407:   /* Setup ST */
408:   PetscCall(STSetUp(eps->st));
409:   PetscCall(EPSCheckCompatibleST(eps));

411:   /* process deflation and initial vectors */
412:   if (eps->nds<0) {
413:     PetscCheck(!eps->isstructured,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Deflation space is not supported in structured eigenproblems");
414:     k = -eps->nds;
415:     PetscCall(BVInsertConstraints(eps->V,&k,eps->defl));
416:     PetscCall(SlepcBasisDestroy_Private(&eps->nds,&eps->defl));
417:     eps->nds = k;
418:     PetscCall(STCheckNullSpace(eps->st,eps->V));
419:   }
420:   if (eps->nini<0) {
421:     k = -eps->nini;
422:     PetscCheck(k<=eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The number of initial vectors is larger than ncv");
423:     PetscCall(BVInsertVecs(eps->V,0,&k,eps->IS,PETSC_TRUE));
424:     PetscCall(SlepcBasisDestroy_Private(&eps->nini,&eps->IS));
425:     eps->nini = k;
426:   }
427:   if (eps->twosided && eps->ninil<0) {
428:     k = -eps->ninil;
429:     PetscCheck(k<=eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The number of left initial vectors is larger than ncv");
430:     PetscCall(BVInsertVecs(eps->W,0,&k,eps->ISL,PETSC_TRUE));
431:     PetscCall(SlepcBasisDestroy_Private(&eps->ninil,&eps->ISL));
432:     eps->ninil = k;
433:   }

435:   PetscCall(PetscLogEventEnd(EPS_SetUp,eps,0,0,0));
436:   eps->state = EPS_STATE_SETUP;
437:   PetscFunctionReturn(PETSC_SUCCESS);
438: }

440: /*@
441:    EPSSetOperators - Sets the matrices associated with the eigenvalue problem.

443:    Collective

445:    Input Parameters:
446: +  eps - the eigenproblem solver context
447: .  A  - the matrix associated with the eigensystem
448: -  B  - the second matrix in the case of generalized eigenproblems

450:    Notes:
451:    To specify a standard eigenproblem, use NULL for parameter B.

453:    It must be called before EPSSetUp(). If it is called again after EPSSetUp() and
454:    the matrix sizes have changed then the EPS object is reset.

456:    For structured eigenproblem types such as EPS_BSE (see EPSSetProblemType()), the
457:    provided matrices must have been created with the corresponding helper function,
458:    i.e., MatCreateBSE().

460:    Level: beginner

462: .seealso: `EPSSolve()`, `EPSSetUp()`, `EPSReset()`, `EPSGetST()`, `STGetMatrix()`, `EPSSetProblemType()`
463: @*/
464: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
465: {
466:   PetscInt       m,n,m0,mloc,nloc,mloc0,nmat;
467:   Mat            mat[2];

469:   PetscFunctionBegin;
473:   PetscCheckSameComm(eps,1,A,2);
474:   if (B) PetscCheckSameComm(eps,1,B,3);

476:   /* Check matrix sizes */
477:   PetscCall(MatGetSize(A,&m,&n));
478:   PetscCall(MatGetLocalSize(A,&mloc,&nloc));
479:   PetscCheck(m==n,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"A is a non-square matrix (%" PetscInt_FMT " rows, %" PetscInt_FMT " cols)",m,n);
480:   PetscCheck(mloc==nloc,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"A does not have equal row and column sizes (%" PetscInt_FMT ", %" PetscInt_FMT ")",mloc,nloc);
481:   if (B) {
482:     PetscCall(MatGetSize(B,&m0,&n));
483:     PetscCall(MatGetLocalSize(B,&mloc0,&nloc));
484:     PetscCheck(m0==n,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"B is a non-square matrix (%" PetscInt_FMT " rows, %" PetscInt_FMT " cols)",m0,n);
485:     PetscCheck(mloc0==nloc,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"B does not have equal row and column local sizes (%" PetscInt_FMT ", %" PetscInt_FMT ")",mloc0,nloc);
486:     PetscCheck(m==m0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Dimensions of A and B do not match (%" PetscInt_FMT ", %" PetscInt_FMT ")",m,m0);
487:     PetscCheck(mloc==mloc0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Local dimensions of A and B do not match (%" PetscInt_FMT ", %" PetscInt_FMT ")",mloc,mloc0);
488:   }
489:   if (eps->state && (n!=eps->n || nloc!=eps->nloc)) PetscCall(EPSReset(eps));
490:   eps->nrma = 0.0;
491:   eps->nrmb = 0.0;
492:   if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
493:   mat[0] = A;
494:   if (B) {
495:     mat[1] = B;
496:     nmat = 2;
497:   } else nmat = 1;
498:   PetscCall(STSetMatrices(eps->st,nmat,mat));
499:   eps->state = EPS_STATE_INITIAL;
500:   PetscFunctionReturn(PETSC_SUCCESS);
501: }

503: /*@
504:    EPSGetOperators - Gets the matrices associated with the eigensystem.

506:    Collective

508:    Input Parameter:
509: .  eps - the EPS context

511:    Output Parameters:
512: +  A  - the matrix associated with the eigensystem
513: -  B  - the second matrix in the case of generalized eigenproblems

515:    Note:
516:    Does not increase the reference count of the matrices, so you should not destroy them.

518:    Level: intermediate

520: .seealso: `EPSSolve()`, `EPSGetST()`, `STGetMatrix()`, `STSetMatrices()`
521: @*/
522: PetscErrorCode EPSGetOperators(EPS eps,Mat *A,Mat *B)
523: {
524:   ST             st;
525:   PetscInt       k;

527:   PetscFunctionBegin;
529:   PetscCall(EPSGetST(eps,&st));
530:   PetscCall(STGetNumMatrices(st,&k));
531:   if (A) {
532:     if (k<1) *A = NULL;
533:     else PetscCall(STGetMatrix(st,0,A));
534:   }
535:   if (B) {
536:     if (k<2) *B = NULL;
537:     else PetscCall(STGetMatrix(st,1,B));
538:   }
539:   PetscFunctionReturn(PETSC_SUCCESS);
540: }

542: /*@
543:    EPSSetDeflationSpace - Specify a basis of vectors that constitute the deflation
544:    space.

546:    Collective

548:    Input Parameters:
549: +  eps - the eigenproblem solver context
550: .  n   - number of vectors
551: -  v   - set of basis vectors of the deflation space

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

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

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

565:    Level: intermediate

567: .seealso: `EPSSetInitialSpace()`
568: @*/
569: PetscErrorCode EPSSetDeflationSpace(EPS eps,PetscInt n,Vec v[])
570: {
571:   PetscFunctionBegin;
574:   PetscCheck(n>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
575:   if (n>0) {
576:     PetscAssertPointer(v,3);
578:   }
579:   PetscCall(SlepcBasisReference_Private(n,v,&eps->nds,&eps->defl));
580:   if (n>0) eps->state = EPS_STATE_INITIAL;
581:   PetscFunctionReturn(PETSC_SUCCESS);
582: }

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

588:    Collective

590:    Input Parameters:
591: +  eps - the eigenproblem solver context
592: .  n   - number of vectors
593: -  is  - set of basis vectors of the initial space

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

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

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

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

608:    Level: intermediate

610: .seealso: `EPSSetLeftInitialSpace()`, `EPSSetDeflationSpace()`
611: @*/
612: PetscErrorCode EPSSetInitialSpace(EPS eps,PetscInt n,Vec is[])
613: {
614:   PetscFunctionBegin;
617:   PetscCheck(n>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
618:   if (n>0) {
619:     PetscAssertPointer(is,3);
621:   }
622:   PetscCall(SlepcBasisReference_Private(n,is,&eps->nini,&eps->IS));
623:   if (n>0) eps->state = EPS_STATE_INITIAL;
624:   PetscFunctionReturn(PETSC_SUCCESS);
625: }

627: /*@
628:    EPSSetLeftInitialSpace - Specify a basis of vectors that constitute the left
629:    initial space, used by two-sided solvers to start the left subspace.

631:    Collective

633:    Input Parameters:
634: +  eps - the eigenproblem solver context
635: .  n   - number of vectors
636: -  isl - set of basis vectors of the left initial space

638:    Notes:
639:    Left initial vectors are used to initiate the left search space in two-sided
640:    eigensolvers. Users should pass here an approximation of the left eigenspace,
641:    if available.

643:    The same comments in EPSSetInitialSpace() are applicable here.

645:    Level: intermediate

647: .seealso: `EPSSetInitialSpace()`, `EPSSetTwoSided()`
648: @*/
649: PetscErrorCode EPSSetLeftInitialSpace(EPS eps,PetscInt n,Vec isl[])
650: {
651:   PetscFunctionBegin;
654:   PetscCheck(n>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
655:   if (n>0) {
656:     PetscAssertPointer(isl,3);
658:   }
659:   PetscCall(SlepcBasisReference_Private(n,isl,&eps->ninil,&eps->ISL));
660:   if (n>0) eps->state = EPS_STATE_INITIAL;
661:   PetscFunctionReturn(PETSC_SUCCESS);
662: }

664: /*
665:   EPSSetDimensions_Default - Set reasonable values for ncv, mpd if not set
666:   by the user. This is called at setup.
667:  */
668: PetscErrorCode EPSSetDimensions_Default(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
669: {
670:   PetscBool      krylov;
671:   PetscInt       nev2, n = eps->problem_type==EPS_BSE? eps->n/2: eps->n;

673:   PetscFunctionBegin;
674:   if (*nev==0 && eps->stop!=EPS_STOP_THRESHOLD) *nev = 1;
675:   nev2 = eps->problem_type==EPS_BSE? (*nev+1)/2: *nev;
676:   if (*ncv!=PETSC_DETERMINE) { /* ncv set */
677:     PetscCall(PetscObjectTypeCompareAny((PetscObject)eps,&krylov,EPSKRYLOVSCHUR,EPSARNOLDI,EPSLANCZOS,""));
678:     if (krylov) {
679:       PetscCheck(*ncv>=nev2+1 || (*ncv==nev2 && *ncv==n),PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev+1");
680:     } else {
681:       PetscCheck(*ncv>=nev2,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev");
682:     }
683:   } else if (*mpd!=PETSC_DETERMINE) { /* mpd set */
684:     *ncv = PetscMin(n,nev2+(*mpd));
685:   } else { /* neither set: defaults depend on nev being small or large */
686:     if (nev2<500) *ncv = PetscMin(n,PetscMax(2*(nev2),nev2+15));
687:     else {
688:       *mpd = 500;
689:       *ncv = PetscMin(n,nev2+(*mpd));
690:     }
691:   }
692:   if (*mpd==PETSC_DETERMINE) *mpd = *ncv;
693:   PetscFunctionReturn(PETSC_SUCCESS);
694: }

696: /*@
697:    EPSAllocateSolution - Allocate memory storage for common variables such
698:    as eigenvalues and eigenvectors.

700:    Collective

702:    Input Parameters:
703: +  eps   - eigensolver context
704: -  extra - number of additional positions, used for methods that require a
705:            working basis slightly larger than ncv

707:    Developer Notes:
708:    This is SLEPC_EXTERN because it may be required by user plugin EPS
709:    implementations.

711:    Level: developer

713: .seealso: `EPSSetUp()`
714: @*/
715: PetscErrorCode EPSAllocateSolution(EPS eps,PetscInt extra)
716: {
717:   PetscInt       oldsize,requested;
718:   PetscRandom    rand;
719:   Vec            t;

721:   PetscFunctionBegin;
722:   requested = eps->ncv + extra;

724:   /* oldsize is zero if this is the first time setup is called */
725:   PetscCall(BVGetSizes(eps->V,NULL,NULL,&oldsize));

727:   /* allocate space for eigenvalues and friends */
728:   if (requested != oldsize || !eps->eigr) {
729:     PetscCall(PetscFree4(eps->eigr,eps->eigi,eps->errest,eps->perm));
730:     PetscCall(PetscMalloc4(requested,&eps->eigr,requested,&eps->eigi,requested,&eps->errest,requested,&eps->perm));
731:   }

733:   /* workspace for the case of arbitrary selection */
734:   if (eps->arbitrary) {
735:     if (eps->rr) PetscCall(PetscFree2(eps->rr,eps->ri));
736:     PetscCall(PetscMalloc2(requested,&eps->rr,requested,&eps->ri));
737:   }

739:   /* allocate V */
740:   if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
741:   if (!oldsize) {
742:     if (!((PetscObject)eps->V)->type_name) PetscCall(BVSetType(eps->V,BVMAT));
743:     PetscCall(STMatCreateVecsEmpty(eps->st,&t,NULL));
744:     PetscCall(BVSetSizesFromVec(eps->V,t,requested));
745:     PetscCall(VecDestroy(&t));
746:   } else PetscCall(BVResize(eps->V,requested,PETSC_FALSE));

748:   /* allocate W */
749:   if (eps->twosided) {
750:     PetscCall(BVGetRandomContext(eps->V,&rand));  /* make sure the random context is available when duplicating */
751:     PetscCall(BVDestroy(&eps->W));
752:     PetscCall(BVDuplicate(eps->V,&eps->W));
753:   }
754:   PetscFunctionReturn(PETSC_SUCCESS);
755: }

757: /*@
758:    EPSReallocateSolution - Reallocate memory storage for common variables such
759:    as the eigenvalues and the basis vectors.

761:    Collective

763:    Input Parameters:
764: +  eps     - eigensolver context
765: -  newsize - new size

767:    Developer Notes:
768:    This is SLEPC_EXTERN because it may be required by user plugin EPS
769:    implementations.

771:    This is called during the iteration in case the threshold stopping test has
772:    been selected.

774:    Level: developer

776: .seealso: `EPSAllocateSolution()`, `EPSSetThreshold()`
777: @*/
778: PetscErrorCode EPSReallocateSolution(EPS eps,PetscInt newsize)
779: {
780:   PetscInt    oldsize,*nperm;
781:   PetscReal   *nerrest;
782:   PetscScalar *neigr,*neigi;

784:   PetscFunctionBegin;
785:   PetscCall(BVGetSizes(eps->V,NULL,NULL,&oldsize));
786:   if (oldsize>=newsize) PetscFunctionReturn(PETSC_SUCCESS);
787:   PetscCall(PetscInfo(eps,"Reallocating basis vectors to %" PetscInt_FMT " columns\n",newsize));

789:   /* reallocate eigenvalues */
790:   PetscCall(PetscMalloc4(newsize,&neigr,newsize,&neigi,newsize,&nerrest,newsize,&nperm));
791:   PetscCall(PetscArraycpy(neigr,eps->eigr,oldsize));
792:   PetscCall(PetscArraycpy(neigi,eps->eigi,oldsize));
793:   PetscCall(PetscArraycpy(nerrest,eps->errest,oldsize));
794:   PetscCall(PetscArraycpy(nperm,eps->perm,oldsize));
795:   PetscCall(PetscFree4(eps->eigr,eps->eigi,eps->errest,eps->perm));
796:   eps->eigr   = neigr;
797:   eps->eigi   = neigi;
798:   eps->errest = nerrest;
799:   eps->perm   = nperm;
800:   /* reallocate V,W */
801:   PetscCall(BVResize(eps->V,newsize,PETSC_TRUE));
802:   if (eps->twosided) PetscCall(BVResize(eps->W,newsize,PETSC_TRUE));
803:   PetscFunctionReturn(PETSC_SUCCESS);
804: }