Actual source code: pepopts.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:    PEP routines related to options that can be set via the command-line
 12:    or procedurally
 13: */

 15: #include <slepc/private/pepimpl.h>
 16: #include <petscdraw.h>

 18: /*@C
 19:    PEPMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
 20:    indicated by the user.

 22:    Collective

 24:    Input Parameters:
 25: +  pep      - the polynomial eigensolver context
 26: .  opt      - the command line option for this monitor
 27: .  name     - the monitor type one is seeking
 28: .  ctx      - an optional user context for the monitor, or `NULL`
 29: -  trackall - whether this monitor tracks all eigenvalues or not

 31:    Level: developer

 33: .seealso: [](ch:pep), `PEPMonitorSet()`, `PEPSetTrackAll()`
 34: @*/
 35: PetscErrorCode PEPMonitorSetFromOptions(PEP pep,const char opt[],const char name[],void *ctx,PetscBool trackall)
 36: {
 37:   PetscErrorCode       (*mfunc)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
 38:   PetscErrorCode       (*cfunc)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
 39:   PetscErrorCode       (*dfunc)(PetscViewerAndFormat**);
 40:   PetscViewerAndFormat *vf;
 41:   PetscViewer          viewer;
 42:   PetscViewerFormat    format;
 43:   PetscViewerType      vtype;
 44:   char                 key[PETSC_MAX_PATH_LEN];
 45:   PetscBool            flg;

 47:   PetscFunctionBegin;
 48:   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->options,((PetscObject)pep)->prefix,opt,&viewer,&format,&flg));
 49:   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);

 51:   PetscCall(PetscViewerGetType(viewer,&vtype));
 52:   PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
 53:   PetscCall(PetscFunctionListFind(PEPMonitorList,key,&mfunc));
 54:   PetscCheck(mfunc,PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Specified viewer and format not supported");
 55:   PetscCall(PetscFunctionListFind(PEPMonitorCreateList,key,&cfunc));
 56:   PetscCall(PetscFunctionListFind(PEPMonitorDestroyList,key,&dfunc));
 57:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
 58:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

 60:   PetscCall((*cfunc)(viewer,format,ctx,&vf));
 61:   PetscCall(PetscViewerDestroy(&viewer));
 62:   PetscCall(PEPMonitorSet(pep,mfunc,vf,(PetscCtxDestroyFn*)dfunc));
 63:   if (trackall) PetscCall(PEPSetTrackAll(pep,PETSC_TRUE));
 64:   PetscFunctionReturn(PETSC_SUCCESS);
 65: }

 67: /*@
 68:    PEPSetFromOptions - Sets `PEP` options from the options database.
 69:    This routine must be called before `PEPSetUp()` if the user is to be
 70:    allowed to configure the solver.

 72:    Collective

 74:    Input Parameter:
 75: .  pep - the polynomial eigensolver context

 77:    Note:
 78:    To see all options, run your program with the `-help` option.

 80:    Level: beginner

 82: .seealso: [](ch:pep), `PEPSetOptionsPrefix()`
 83: @*/
 84: PetscErrorCode PEPSetFromOptions(PEP pep)
 85: {
 86:   char            type[256];
 87:   PetscBool       set,flg,flg1,flg2,flg3,flg4,flg5;
 88:   PetscReal       r,t,array[2]={0,0};
 89:   PetscScalar     s;
 90:   PetscInt        i,j,k;
 91:   PEPScale        scale;
 92:   PEPRefine       refine;
 93:   PEPRefineScheme scheme;

 95:   PetscFunctionBegin;
 97:   PetscCall(PEPRegisterAll());
 98:   PetscObjectOptionsBegin((PetscObject)pep);
 99:     PetscCall(PetscOptionsFList("-pep_type","Polynomial eigensolver method","PEPSetType",PEPList,(char*)(((PetscObject)pep)->type_name?((PetscObject)pep)->type_name:PEPTOAR),type,sizeof(type),&flg));
100:     if (flg) PetscCall(PEPSetType(pep,type));
101:     else if (!((PetscObject)pep)->type_name) PetscCall(PEPSetType(pep,PEPTOAR));

103:     PetscCall(PetscOptionsBoolGroupBegin("-pep_general","General polynomial eigenvalue problem","PEPSetProblemType",&flg));
104:     if (flg) PetscCall(PEPSetProblemType(pep,PEP_GENERAL));
105:     PetscCall(PetscOptionsBoolGroup("-pep_hermitian","Hermitian polynomial eigenvalue problem","PEPSetProblemType",&flg));
106:     if (flg) PetscCall(PEPSetProblemType(pep,PEP_HERMITIAN));
107:     PetscCall(PetscOptionsBoolGroup("-pep_hyperbolic","Hyperbolic polynomial eigenvalue problem","PEPSetProblemType",&flg));
108:     if (flg) PetscCall(PEPSetProblemType(pep,PEP_HYPERBOLIC));
109:     PetscCall(PetscOptionsBoolGroupEnd("-pep_gyroscopic","Gyroscopic polynomial eigenvalue problem","PEPSetProblemType",&flg));
110:     if (flg) PetscCall(PEPSetProblemType(pep,PEP_GYROSCOPIC));

112:     scale = pep->scale;
113:     PetscCall(PetscOptionsEnum("-pep_scale","Scaling strategy","PEPSetScale",PEPScaleTypes,(PetscEnum)scale,(PetscEnum*)&scale,&flg1));
114:     r = pep->sfactor;
115:     PetscCall(PetscOptionsReal("-pep_scale_factor","Scale factor","PEPSetScale",pep->sfactor,&r,&flg2));
116:     if (!flg2 && r==1.0) r = PETSC_DETERMINE;
117:     j = pep->sits;
118:     PetscCall(PetscOptionsInt("-pep_scale_its","Number of iterations in diagonal scaling","PEPSetScale",pep->sits,&j,&flg3));
119:     t = pep->slambda;
120:     PetscCall(PetscOptionsReal("-pep_scale_lambda","Estimate of eigenvalue (modulus) for diagonal scaling","PEPSetScale",pep->slambda,&t,&flg4));
121:     if (flg1 || flg2 || flg3 || flg4) PetscCall(PEPSetScale(pep,scale,r,NULL,NULL,j,t));

123:     PetscCall(PetscOptionsEnum("-pep_extract","Extraction method","PEPSetExtract",PEPExtractTypes,(PetscEnum)pep->extract,(PetscEnum*)&pep->extract,NULL));

125:     refine = pep->refine;
126:     PetscCall(PetscOptionsEnum("-pep_refine","Iterative refinement method","PEPSetRefine",PEPRefineTypes,(PetscEnum)refine,(PetscEnum*)&refine,&flg1));
127:     i = pep->npart;
128:     PetscCall(PetscOptionsInt("-pep_refine_partitions","Number of partitions of the communicator for iterative refinement","PEPSetRefine",pep->npart,&i,&flg2));
129:     r = pep->rtol;
130:     PetscCall(PetscOptionsReal("-pep_refine_tol","Tolerance for iterative refinement","PEPSetRefine",pep->rtol==(PetscReal)PETSC_DETERMINE?SLEPC_DEFAULT_TOL/1000:pep->rtol,&r,&flg3));
131:     j = pep->rits;
132:     PetscCall(PetscOptionsInt("-pep_refine_its","Maximum number of iterations for iterative refinement","PEPSetRefine",pep->rits,&j,&flg4));
133:     scheme = pep->scheme;
134:     PetscCall(PetscOptionsEnum("-pep_refine_scheme","Scheme used for linear systems within iterative refinement","PEPSetRefine",PEPRefineSchemes,(PetscEnum)scheme,(PetscEnum*)&scheme,&flg5));
135:     if (flg1 || flg2 || flg3 || flg4 || flg5) PetscCall(PEPSetRefine(pep,refine,i,r,j,scheme));

137:     i = pep->max_it;
138:     PetscCall(PetscOptionsInt("-pep_max_it","Maximum number of iterations","PEPSetTolerances",pep->max_it,&i,&flg1));
139:     r = pep->tol;
140:     PetscCall(PetscOptionsReal("-pep_tol","Tolerance","PEPSetTolerances",SlepcDefaultTol(pep->tol),&r,&flg2));
141:     if (flg1 || flg2) PetscCall(PEPSetTolerances(pep,r,i));

143:     PetscCall(PetscOptionsBoolGroupBegin("-pep_conv_rel","Relative error convergence test","PEPSetConvergenceTest",&flg));
144:     if (flg) PetscCall(PEPSetConvergenceTest(pep,PEP_CONV_REL));
145:     PetscCall(PetscOptionsBoolGroup("-pep_conv_norm","Convergence test relative to the matrix norms","PEPSetConvergenceTest",&flg));
146:     if (flg) PetscCall(PEPSetConvergenceTest(pep,PEP_CONV_NORM));
147:     PetscCall(PetscOptionsBoolGroup("-pep_conv_abs","Absolute error convergence test","PEPSetConvergenceTest",&flg));
148:     if (flg) PetscCall(PEPSetConvergenceTest(pep,PEP_CONV_ABS));
149:     PetscCall(PetscOptionsBoolGroupEnd("-pep_conv_user","User-defined convergence test","PEPSetConvergenceTest",&flg));
150:     if (flg) PetscCall(PEPSetConvergenceTest(pep,PEP_CONV_USER));

152:     PetscCall(PetscOptionsBoolGroupBegin("-pep_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","PEPSetStoppingTest",&flg));
153:     if (flg) PetscCall(PEPSetStoppingTest(pep,PEP_STOP_BASIC));
154:     PetscCall(PetscOptionsBoolGroupEnd("-pep_stop_user","User-defined stopping test","PEPSetStoppingTest",&flg));
155:     if (flg) PetscCall(PEPSetStoppingTest(pep,PEP_STOP_USER));

157:     i = pep->nev;
158:     PetscCall(PetscOptionsInt("-pep_nev","Number of eigenvalues to compute","PEPSetDimensions",pep->nev,&i,&flg1));
159:     j = pep->ncv;
160:     PetscCall(PetscOptionsInt("-pep_ncv","Number of basis vectors","PEPSetDimensions",pep->ncv,&j,&flg2));
161:     k = pep->mpd;
162:     PetscCall(PetscOptionsInt("-pep_mpd","Maximum dimension of projected problem","PEPSetDimensions",pep->mpd,&k,&flg3));
163:     if (flg1 || flg2 || flg3) PetscCall(PEPSetDimensions(pep,i,j,k));

165:     PetscCall(PetscOptionsEnum("-pep_basis","Polynomial basis","PEPSetBasis",PEPBasisTypes,(PetscEnum)pep->basis,(PetscEnum*)&pep->basis,NULL));

167:     PetscCall(PetscOptionsBoolGroupBegin("-pep_largest_magnitude","Compute largest eigenvalues in magnitude","PEPSetWhichEigenpairs",&flg));
168:     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_LARGEST_MAGNITUDE));
169:     PetscCall(PetscOptionsBoolGroup("-pep_smallest_magnitude","Compute smallest eigenvalues in magnitude","PEPSetWhichEigenpairs",&flg));
170:     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_SMALLEST_MAGNITUDE));
171:     PetscCall(PetscOptionsBoolGroup("-pep_largest_real","Compute eigenvalues with largest real parts","PEPSetWhichEigenpairs",&flg));
172:     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_LARGEST_REAL));
173:     PetscCall(PetscOptionsBoolGroup("-pep_smallest_real","Compute eigenvalues with smallest real parts","PEPSetWhichEigenpairs",&flg));
174:     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_SMALLEST_REAL));
175:     PetscCall(PetscOptionsBoolGroup("-pep_largest_imaginary","Compute eigenvalues with largest imaginary parts","PEPSetWhichEigenpairs",&flg));
176:     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_LARGEST_IMAGINARY));
177:     PetscCall(PetscOptionsBoolGroup("-pep_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","PEPSetWhichEigenpairs",&flg));
178:     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_SMALLEST_IMAGINARY));
179:     PetscCall(PetscOptionsBoolGroup("-pep_target_magnitude","Compute eigenvalues closest to target","PEPSetWhichEigenpairs",&flg));
180:     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE));
181:     PetscCall(PetscOptionsBoolGroup("-pep_target_real","Compute eigenvalues with real parts closest to target","PEPSetWhichEigenpairs",&flg));
182:     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_TARGET_REAL));
183:     PetscCall(PetscOptionsBoolGroup("-pep_target_imaginary","Compute eigenvalues with imaginary parts closest to target","PEPSetWhichEigenpairs",&flg));
184:     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_TARGET_IMAGINARY));
185:     PetscCall(PetscOptionsBoolGroupEnd("-pep_all","Compute all eigenvalues in an interval or a region","PEPSetWhichEigenpairs",&flg));
186:     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_ALL));

188:     PetscCall(PetscOptionsScalar("-pep_target","Value of the target","PEPSetTarget",pep->target,&s,&flg));
189:     if (flg) {
190:       if (pep->which!=PEP_TARGET_REAL && pep->which!=PEP_TARGET_IMAGINARY) PetscCall(PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE));
191:       PetscCall(PEPSetTarget(pep,s));
192:     }

194:     k = 2;
195:     PetscCall(PetscOptionsRealArray("-pep_interval","Computational interval (two real values separated with a comma without spaces)","PEPSetInterval",array,&k,&flg));
196:     if (flg) {
197:       PetscCheck(k>1,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_SIZ,"Must pass two values in -pep_interval (comma-separated without spaces)");
198:       PetscCall(PEPSetWhichEigenpairs(pep,PEP_ALL));
199:       PetscCall(PEPSetInterval(pep,array[0],array[1]));
200:     }

202:     /* -----------------------------------------------------------------------*/
203:     /*
204:       Cancels all monitors hardwired into code before call to PEPSetFromOptions()
205:     */
206:     PetscCall(PetscOptionsBool("-pep_monitor_cancel","Remove any hardwired monitor routines","PEPMonitorCancel",PETSC_FALSE,&flg,&set));
207:     if (set && flg) PetscCall(PEPMonitorCancel(pep));
208:     PetscCall(PEPMonitorSetFromOptions(pep,"-pep_monitor","first_approximation",NULL,PETSC_FALSE));
209:     PetscCall(PEPMonitorSetFromOptions(pep,"-pep_monitor_all","all_approximations",NULL,PETSC_TRUE));
210:     PetscCall(PEPMonitorSetFromOptions(pep,"-pep_monitor_conv","convergence_history",NULL,PETSC_FALSE));

212:     /* -----------------------------------------------------------------------*/
213:     PetscCall(PetscOptionsName("-pep_view","Print detailed information on solver used","PEPView",&set));
214:     PetscCall(PetscOptionsName("-pep_view_vectors","View computed eigenvectors","PEPVectorsView",&set));
215:     PetscCall(PetscOptionsName("-pep_view_values","View computed eigenvalues","PEPValuesView",&set));
216:     PetscCall(PetscOptionsName("-pep_converged_reason","Print reason for convergence, and number of iterations","PEPConvergedReasonView",&set));
217:     PetscCall(PetscOptionsName("-pep_error_absolute","Print absolute errors of each eigenpair","PEPErrorView",&set));
218:     PetscCall(PetscOptionsName("-pep_error_relative","Print relative errors of each eigenpair","PEPErrorView",&set));
219:     PetscCall(PetscOptionsName("-pep_error_backward","Print backward errors of each eigenpair","PEPErrorView",&set));

221:     PetscTryTypeMethod(pep,setfromoptions,PetscOptionsObject);
222:     PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)pep,PetscOptionsObject));
223:   PetscOptionsEnd();

225:   if (!pep->V) PetscCall(PEPGetBV(pep,&pep->V));
226:   PetscCall(BVSetFromOptions(pep->V));
227:   if (!pep->rg) PetscCall(PEPGetRG(pep,&pep->rg));
228:   PetscCall(RGSetFromOptions(pep->rg));
229:   if (!pep->ds) PetscCall(PEPGetDS(pep,&pep->ds));
230:   PetscCall(PEPSetDSType(pep));
231:   PetscCall(DSSetFromOptions(pep->ds));
232:   if (!pep->st) PetscCall(PEPGetST(pep,&pep->st));
233:   PetscCall(PEPSetDefaultST(pep));
234:   PetscCall(STSetFromOptions(pep->st));
235:   if (!pep->refineksp) PetscCall(PEPRefineGetKSP(pep,&pep->refineksp));
236:   PetscCall(KSPSetFromOptions(pep->refineksp));
237:   PetscFunctionReturn(PETSC_SUCCESS);
238: }

240: /*@
241:    PEPGetTolerances - Gets the tolerance and maximum iteration count used
242:    by the `PEP` convergence tests.

244:    Not Collective

246:    Input Parameter:
247: .  pep - the polynomial eigensolver context

249:    Output Parameters:
250: +  tol - the convergence tolerance
251: -  maxits - maximum number of iterations

253:    Note:
254:    The user can specify `NULL` for any parameter that is not needed.

256:    Level: intermediate

258: .seealso: [](ch:pep), `PEPSetTolerances()`
259: @*/
260: PetscErrorCode PEPGetTolerances(PEP pep,PetscReal *tol,PetscInt *maxits)
261: {
262:   PetscFunctionBegin;
264:   if (tol)    *tol    = pep->tol;
265:   if (maxits) *maxits = pep->max_it;
266:   PetscFunctionReturn(PETSC_SUCCESS);
267: }

269: /*@
270:    PEPSetTolerances - Sets the tolerance and maximum iteration count used
271:    by the `PEP` convergence tests.

273:    Logically Collective

275:    Input Parameters:
276: +  pep    - the polynomial eigensolver context
277: .  tol    - the convergence tolerance
278: -  maxits - maximum number of iterations to use

280:    Options Database Keys:
281: +  -pep_tol \<tol\>       - sets the convergence tolerance
282: -  -pep_max_it \<maxits\> - sets the maximum number of iterations allowed

284:    Note:
285:    Use `PETSC_CURRENT` to retain the current value of any of the parameters.
286:    Use `PETSC_DETERMINE` for either argument to assign a default value computed
287:    internally (may be different in each solver).
288:    For `maxits` use `PETSC_UNLIMITED` to indicate there is no upper bound on this value.

290:    Level: intermediate

292: .seealso: [](ch:pep), `PEPGetTolerances()`
293: @*/
294: PetscErrorCode PEPSetTolerances(PEP pep,PetscReal tol,PetscInt maxits)
295: {
296:   PetscFunctionBegin;
300:   if (tol == (PetscReal)PETSC_DETERMINE) {
301:     pep->tol   = PETSC_DETERMINE;
302:     pep->state = PEP_STATE_INITIAL;
303:   } else if (tol != (PetscReal)PETSC_CURRENT) {
304:     PetscCheck(tol>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
305:     pep->tol = tol;
306:   }
307:   if (maxits == PETSC_DETERMINE) {
308:     pep->max_it = PETSC_DETERMINE;
309:     pep->state  = PEP_STATE_INITIAL;
310:   } else if (maxits == PETSC_UNLIMITED) {
311:     pep->max_it = PETSC_INT_MAX;
312:   } else if (maxits != PETSC_CURRENT) {
313:     PetscCheck(maxits>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
314:     pep->max_it = maxits;
315:   }
316:   PetscFunctionReturn(PETSC_SUCCESS);
317: }

319: /*@
320:    PEPGetDimensions - Gets the number of eigenvalues to compute
321:    and the dimension of the subspace.

323:    Not Collective

325:    Input Parameter:
326: .  pep - the polynomial eigensolver context

328:    Output Parameters:
329: +  nev - number of eigenvalues to compute
330: .  ncv - the maximum dimension of the subspace to be used by the solver
331: -  mpd - the maximum dimension allowed for the projected problem

333:    Notes:
334:    The user can specify `NULL` for any parameter that is not needed.

336:    Level: intermediate

338: .seealso: [](ch:pep), `PEPSetDimensions()`
339: @*/
340: PetscErrorCode PEPGetDimensions(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
341: {
342:   PetscFunctionBegin;
344:   if (nev) *nev = pep->nev;
345:   if (ncv) *ncv = pep->ncv;
346:   if (mpd) *mpd = pep->mpd;
347:   PetscFunctionReturn(PETSC_SUCCESS);
348: }

350: /*@
351:    PEPSetDimensions - Sets the number of eigenvalues to compute
352:    and the dimension of the subspace.

354:    Logically Collective

356:    Input Parameters:
357: +  pep - the polynomial eigensolver context
358: .  nev - number of eigenvalues to compute
359: .  ncv - the maximum dimension of the subspace to be used by the solver
360: -  mpd - the maximum dimension allowed for the projected problem

362:    Options Database Keys:
363: +  -pep_nev \<nev\> - sets the number of eigenvalues
364: .  -pep_ncv \<ncv\> - sets the dimension of the subspace
365: -  -pep_mpd \<mpd\> - sets the maximum projected dimension

367:    Notes:
368:    Use `PETSC_DETERMINE` for `ncv` and `mpd` to assign a reasonably good value, which is
369:    dependent on the solution method. For any of the arguments, use `PETSC_CURRENT`
370:    to preserve the current value.

372:    The parameters `ncv` and `mpd` are intimately related, so that the user is advised
373:    to set one of them at most. Normal usage is\:

375:     1. in cases where `nev` is small, the user sets `ncv` (a reasonable default is `2*nev`).
376:     2. in cases where `nev` is large, the user sets `mpd`.

378:    The value of `ncv` should always be between `nev` and `(nev+mpd)`, typically
379:    `ncv=nev+mpd`. If `nev` is not too large, `mpd=nev` is a reasonable choice, otherwise
380:    a smaller value should be used.

382:    When computing all eigenvalues in an interval, see `PEPSetInterval()`, these
383:    parameters lose relevance, and tuning must be done with
384:    `PEPSTOARSetDimensions()`.

386:    Level: intermediate

388: .seealso: [](ch:pep), `PEPGetDimensions()`, `PEPSetInterval()`, `PEPSTOARSetDimensions()`
389: @*/
390: PetscErrorCode PEPSetDimensions(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)
391: {
392:   PetscFunctionBegin;
397:   if (nev != PETSC_CURRENT) {
398:     PetscCheck(nev>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
399:     pep->nev = nev;
400:   }
401:   if (ncv == PETSC_DETERMINE) {
402:     pep->ncv = PETSC_DETERMINE;
403:   } else if (ncv != PETSC_CURRENT) {
404:     PetscCheck(ncv>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
405:     pep->ncv = ncv;
406:   }
407:   if (mpd == PETSC_DETERMINE) {
408:     pep->mpd = PETSC_DETERMINE;
409:   } else if (mpd != PETSC_CURRENT) {
410:     PetscCheck(mpd>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
411:     pep->mpd = mpd;
412:   }
413:   pep->state = PEP_STATE_INITIAL;
414:   PetscFunctionReturn(PETSC_SUCCESS);
415: }

417: /*@
418:    PEPSetWhichEigenpairs - Specifies which portion of the spectrum is
419:    to be sought.

421:    Logically Collective

423:    Input Parameters:
424: +  pep   - the polynomial eigensolver context
425: -  which - the portion of the spectrum to be sought, see `PEPWhich` for possible values

427:    Options Database Keys:
428: +  -pep_largest_magnitude  - sets largest eigenvalues in magnitude
429: .  -pep_smallest_magnitude - sets smallest eigenvalues in magnitude
430: .  -pep_largest_real       - sets largest real parts
431: .  -pep_smallest_real      - sets smallest real parts
432: .  -pep_largest_imaginary  - sets largest imaginary parts
433: .  -pep_smallest_imaginary - sets smallest imaginary parts
434: .  -pep_target_magnitude   - sets eigenvalues closest to target
435: .  -pep_target_real        - sets real parts closest to target
436: .  -pep_target_imaginary   - sets imaginary parts closest to target
437: -  -pep_all                - sets all eigenvalues in an interval or region

439:    Notes:
440:    Not all eigensolvers implemented in `PEP` account for all the possible values
441:    of `which`. Also, some values make sense only for certain types of
442:    problems. If SLEPc is compiled for real numbers `PEP_LARGEST_IMAGINARY`
443:    and `PEP_SMALLEST_IMAGINARY` use the absolute value of the imaginary part
444:    for eigenvalue selection.

446:    The target is a scalar value provided with `PEPSetTarget()`.

448:    The criterion `PEP_TARGET_IMAGINARY` is available only in case PETSc and
449:    SLEPc have been built with complex scalars.

451:    `PEP_ALL` is intended for use in combination with an interval (see
452:    `PEPSetInterval()`), when all eigenvalues within the interval are requested,
453:    and also for computing all eigenvalues in a region with the `PEPCISS` solver.

455:    Level: intermediate

457: .seealso: [](ch:pep), `PEPGetWhichEigenpairs()`, `PEPSetTarget()`, `PEPSetInterval()`, `PEPSetDimensions()`, `PEPSetEigenvalueComparison()`, `PEPWhich`
458: @*/
459: PetscErrorCode PEPSetWhichEigenpairs(PEP pep,PEPWhich which)
460: {
461:   PetscFunctionBegin;
464:   switch (which) {
465:     case PEP_LARGEST_MAGNITUDE:
466:     case PEP_SMALLEST_MAGNITUDE:
467:     case PEP_LARGEST_REAL:
468:     case PEP_SMALLEST_REAL:
469:     case PEP_LARGEST_IMAGINARY:
470:     case PEP_SMALLEST_IMAGINARY:
471:     case PEP_TARGET_MAGNITUDE:
472:     case PEP_TARGET_REAL:
473: #if defined(PETSC_USE_COMPLEX)
474:     case PEP_TARGET_IMAGINARY:
475: #endif
476:     case PEP_ALL:
477:     case PEP_WHICH_USER:
478:       if (pep->which != which) {
479:         pep->state = PEP_STATE_INITIAL;
480:         pep->which = which;
481:       }
482:       break;
483: #if !defined(PETSC_USE_COMPLEX)
484:     case PEP_TARGET_IMAGINARY:
485:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"PEP_TARGET_IMAGINARY can be used only with complex scalars");
486: #endif
487:     default:
488:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
489:   }
490:   PetscFunctionReturn(PETSC_SUCCESS);
491: }

493: /*@
494:     PEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
495:     sought.

497:     Not Collective

499:     Input Parameter:
500: .   pep - the polynomial eigensolver context

502:     Output Parameter:
503: .   which - the portion of the spectrum to be sought

505:     Level: intermediate

507: .seealso: [](ch:pep), `PEPSetWhichEigenpairs()`, `PEPWhich`
508: @*/
509: PetscErrorCode PEPGetWhichEigenpairs(PEP pep,PEPWhich *which)
510: {
511:   PetscFunctionBegin;
513:   PetscAssertPointer(which,2);
514:   *which = pep->which;
515:   PetscFunctionReturn(PETSC_SUCCESS);
516: }

518: /*@C
519:    PEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
520:    when `PEPSetWhichEigenpairs()` is set to `PEP_WHICH_USER`.

522:    Logically Collective

524:    Input Parameters:
525: +  pep  - the polynomial eigensolver context
526: .  comp - the comparison function, see `SlepcEigenvalueComparisonFn` for the calling sequence
527: -  ctx  - a context pointer (the last parameter to the comparison function)

529:    Level: advanced

531: .seealso: [](ch:pep), `PEPSetWhichEigenpairs()`, `PEPWhich`
532: @*/
533: PetscErrorCode PEPSetEigenvalueComparison(PEP pep,SlepcEigenvalueComparisonFn *comp,void *ctx)
534: {
535:   PetscFunctionBegin;
537:   pep->sc->comparison    = comp;
538:   pep->sc->comparisonctx = ctx;
539:   pep->which             = PEP_WHICH_USER;
540:   PetscFunctionReturn(PETSC_SUCCESS);
541: }

543: /*@
544:    PEPSetProblemType - Specifies the type of the polynomial eigenvalue problem.

546:    Logically Collective

548:    Input Parameters:
549: +  pep  - the polynomial eigensolver context
550: -  type - a known type of polynomial eigenvalue problem

552:    Options Database Keys:
553: +  -pep_general    - general problem with no particular structure
554: .  -pep_hermitian  - problem whose coefficient matrices are Hermitian
555: .  -pep_hyperbolic - Hermitian problem that satisfies the definition of hyperbolic
556: -  -pep_gyroscopic - problem with gyroscopic structure

558:    Notes:
559:    See `PEPProblemType` for possible problem types.

561:    This function is used to instruct SLEPc to exploit certain structure in
562:    the polynomial eigenproblem. By default, no particular structure is assumed.

564:    If the problem matrices are Hermitian (symmetric in the real case) or
565:    Hermitian/skew-Hermitian then the solver can exploit this fact to perform
566:    less operations or provide better stability. Hyperbolic problems are a
567:    particular case of Hermitian problems, some solvers may treat them simply as
568:    Hermitian.

570:    Level: intermediate

572: .seealso: [](ch:pep), `PEPSetOperators()`, `PEPSetType()`, `PEPGetProblemType()`, `PEPProblemType`
573: @*/
574: PetscErrorCode PEPSetProblemType(PEP pep,PEPProblemType type)
575: {
576:   PetscFunctionBegin;
579:   PetscCheck(type==PEP_GENERAL || type==PEP_HERMITIAN || type==PEP_HYPERBOLIC || type==PEP_GYROSCOPIC,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
580:   if (type != pep->problem_type) {
581:     pep->problem_type = type;
582:     pep->state = PEP_STATE_INITIAL;
583:   }
584:   PetscFunctionReturn(PETSC_SUCCESS);
585: }

587: /*@
588:    PEPGetProblemType - Gets the problem type from the `PEP` object.

590:    Not Collective

592:    Input Parameter:
593: .  pep - the polynomial eigensolver context

595:    Output Parameter:
596: .  type - the problem type

598:    Level: intermediate

600: .seealso: [](ch:pep), `PEPSetProblemType()`, `PEPProblemType`
601: @*/
602: PetscErrorCode PEPGetProblemType(PEP pep,PEPProblemType *type)
603: {
604:   PetscFunctionBegin;
606:   PetscAssertPointer(type,2);
607:   *type = pep->problem_type;
608:   PetscFunctionReturn(PETSC_SUCCESS);
609: }

611: /*@
612:    PEPSetBasis - Specifies the type of polynomial basis used to describe the
613:    polynomial eigenvalue problem.

615:    Logically Collective

617:    Input Parameters:
618: +  pep   - the polynomial eigensolver context
619: -  basis - the type of polynomial basis, see `PEPBasis` for possible values

621:    Options Database Key:
622: .  -pep_basis \<basis\> - Select the basis type

624:    Note:
625:    By default, the coefficient matrices passed via `PEPSetOperators()` are
626:    expressed in the monomial basis, i.e.
627:    $P(\lambda) = A_0 + \lambda A_1 + \lambda^2 A_2 + \dots + \lambda^d A_d$.
628:    Other polynomial bases may have better numerical behavior, but the user
629:    must then pass the coefficient matrices accordingly.

631:    Level: intermediate

633: .seealso: [](ch:pep), `PEPSetOperators()`, `PEPGetBasis()`, `PEPBasis`
634: @*/
635: PetscErrorCode PEPSetBasis(PEP pep,PEPBasis basis)
636: {
637:   PetscFunctionBegin;
640:   pep->basis = basis;
641:   PetscFunctionReturn(PETSC_SUCCESS);
642: }

644: /*@
645:    PEPGetBasis - Gets the type of polynomial basis from the `PEP` object.

647:    Not Collective

649:    Input Parameter:
650: .  pep - the polynomial eigensolver context

652:    Output Parameter:
653: .  basis - the polynomial basis

655:    Level: intermediate

657: .seealso: [](ch:pep), `PEPSetBasis()`, `PEPBasis`
658: @*/
659: PetscErrorCode PEPGetBasis(PEP pep,PEPBasis *basis)
660: {
661:   PetscFunctionBegin;
663:   PetscAssertPointer(basis,2);
664:   *basis = pep->basis;
665:   PetscFunctionReturn(PETSC_SUCCESS);
666: }

668: /*@
669:    PEPSetTrackAll - Specifies if the solver must compute the residual of all
670:    approximate eigenpairs or not.

672:    Logically Collective

674:    Input Parameters:
675: +  pep      - the polynomial eigensolver context
676: -  trackall - whether compute all residuals or not

678:    Notes:
679:    If the user sets `trackall`=`PETSC_TRUE` then the solver explicitly computes
680:    the residual for each eigenpair approximation. Computing the residual is
681:    usually an expensive operation and solvers commonly compute the associated
682:    residual to the first unconverged eigenpair.

684:    The option `-pep_monitor_all` automatically activates this option.

686:    Level: developer

688: .seealso: [](ch:pep), `PEPGetTrackAll()`
689: @*/
690: PetscErrorCode PEPSetTrackAll(PEP pep,PetscBool trackall)
691: {
692:   PetscFunctionBegin;
695:   pep->trackall = trackall;
696:   PetscFunctionReturn(PETSC_SUCCESS);
697: }

699: /*@
700:    PEPGetTrackAll - Returns the flag indicating whether all residual norms must
701:    be computed or not.

703:    Not Collective

705:    Input Parameter:
706: .  pep - the polynomial eigensolver context

708:    Output Parameter:
709: .  trackall - the returned flag

711:    Level: developer

713: .seealso: [](ch:pep), `PEPSetTrackAll()`
714: @*/
715: PetscErrorCode PEPGetTrackAll(PEP pep,PetscBool *trackall)
716: {
717:   PetscFunctionBegin;
719:   PetscAssertPointer(trackall,2);
720:   *trackall = pep->trackall;
721:   PetscFunctionReturn(PETSC_SUCCESS);
722: }

724: /*@C
725:    PEPSetConvergenceTestFunction - Sets a function to compute the error estimate
726:    used in the convergence test.

728:    Logically Collective

730:    Input Parameters:
731: +  pep     - the polynomial eigensolver context
732: .  conv    - convergence test function, see `PEPConvergenceTestFn` for the calling sequence
733: .  ctx     - context for private data for the convergence routine (may be `NULL`)
734: -  destroy - a routine for destroying the context (may be `NULL`), see `PetscCtxDestroyFn`
735:              for the calling sequence

737:    Notes:
738:    When this is called with a user-defined function, then the convergence
739:    criterion is set to `PEP_CONV_USER`, see `PEPSetConvergenceTest()`.

741:    If the error estimate returned by the convergence test function is less than
742:    the tolerance, then the eigenvalue is accepted as converged.

744:    Level: advanced

746: .seealso: [](ch:pep), `PEPSetConvergenceTest()`, `PEPSetTolerances()`
747: @*/
748: PetscErrorCode PEPSetConvergenceTestFunction(PEP pep,PEPConvergenceTestFn *conv,void *ctx,PetscCtxDestroyFn *destroy)
749: {
750:   PetscFunctionBegin;
752:   if (pep->convergeddestroy) PetscCall((*pep->convergeddestroy)(&pep->convergedctx));
753:   pep->convergeduser    = conv;
754:   pep->convergeddestroy = destroy;
755:   pep->convergedctx     = ctx;
756:   if (conv == PEPConvergedRelative) pep->conv = PEP_CONV_REL;
757:   else if (conv == PEPConvergedNorm) pep->conv = PEP_CONV_NORM;
758:   else if (conv == PEPConvergedAbsolute) pep->conv = PEP_CONV_ABS;
759:   else {
760:     pep->conv      = PEP_CONV_USER;
761:     pep->converged = pep->convergeduser;
762:   }
763:   PetscFunctionReturn(PETSC_SUCCESS);
764: }

766: /*@
767:    PEPSetConvergenceTest - Specifies how to compute the error estimate
768:    used in the convergence test.

770:    Logically Collective

772:    Input Parameters:
773: +  pep  - the polynomial eigensolver context
774: -  conv - the type of convergence test, see `PEPConv` for possible values

776:    Options Database Keys:
777: +  -pep_conv_abs  - sets the absolute convergence test
778: .  -pep_conv_rel  - sets the convergence test relative to the eigenvalue
779: .  -pep_conv_norm - sets the convergence test relative to the matrix norms
780: -  -pep_conv_user - selects the user-defined convergence test

782:    Level: intermediate

784: .seealso: [](ch:pep), `PEPGetConvergenceTest()`, `PEPSetConvergenceTestFunction()`, `PEPSetStoppingTest()`, `PEPConv`
785: @*/
786: PetscErrorCode PEPSetConvergenceTest(PEP pep,PEPConv conv)
787: {
788:   PetscFunctionBegin;
791:   switch (conv) {
792:     case PEP_CONV_ABS:  pep->converged = PEPConvergedAbsolute; break;
793:     case PEP_CONV_REL:  pep->converged = PEPConvergedRelative; break;
794:     case PEP_CONV_NORM: pep->converged = PEPConvergedNorm; break;
795:     case PEP_CONV_USER:
796:       PetscCheck(pep->convergeduser,PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetConvergenceTestFunction() first");
797:       pep->converged = pep->convergeduser;
798:       break;
799:     default:
800:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
801:   }
802:   pep->conv = conv;
803:   PetscFunctionReturn(PETSC_SUCCESS);
804: }

806: /*@
807:    PEPGetConvergenceTest - Gets the method used to compute the error estimate
808:    used in the convergence test.

810:    Not Collective

812:    Input Parameter:
813: .  pep   - the polynomial eigensolver context

815:    Output Parameter:
816: .  conv  - the type of convergence test

818:    Level: intermediate

820: .seealso: [](ch:pep), `PEPSetConvergenceTest()`, `PEPConv`
821: @*/
822: PetscErrorCode PEPGetConvergenceTest(PEP pep,PEPConv *conv)
823: {
824:   PetscFunctionBegin;
826:   PetscAssertPointer(conv,2);
827:   *conv = pep->conv;
828:   PetscFunctionReturn(PETSC_SUCCESS);
829: }

831: /*@C
832:    PEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
833:    iteration of the eigensolver.

835:    Logically Collective

837:    Input Parameters:
838: +  pep     - the polynomial eigensolver context
839: .  stop    - stopping test function, see `PEPStoppingTestFn` for the calling sequence
840: .  ctx     - context for private data for the stopping routine (may be `NULL`)
841: -  destroy - a routine for destroying the context (may be `NULL`), see `PetscCtxDestroyFn`
842:              for the calling sequence

844:    Note:
845:    When implementing a function for this, normal usage is to first call the
846:    default routine `PEPStoppingBasic()` and then set `reason` to `PEP_CONVERGED_USER`
847:    if some user-defined conditions have been met. To let the eigensolver continue
848:    iterating, the result must be left as `PEP_CONVERGED_ITERATING`.

850:    Level: advanced

852: .seealso: [](ch:pep), `PEPSetStoppingTest()`, `PEPStoppingBasic()`
853: @*/
854: PetscErrorCode PEPSetStoppingTestFunction(PEP pep,PEPStoppingTestFn *stop,void *ctx,PetscCtxDestroyFn *destroy)
855: {
856:   PetscFunctionBegin;
858:   if (pep->stoppingdestroy) PetscCall((*pep->stoppingdestroy)(&pep->stoppingctx));
859:   pep->stoppinguser    = stop;
860:   pep->stoppingdestroy = destroy;
861:   pep->stoppingctx     = ctx;
862:   if (stop == PEPStoppingBasic) pep->stop = PEP_STOP_BASIC;
863:   else {
864:     pep->stop     = PEP_STOP_USER;
865:     pep->stopping = pep->stoppinguser;
866:   }
867:   PetscFunctionReturn(PETSC_SUCCESS);
868: }

870: /*@
871:    PEPSetStoppingTest - Specifies how to decide the termination of the outer
872:    loop of the eigensolver.

874:    Logically Collective

876:    Input Parameters:
877: +  pep  - the polynomial eigensolver context
878: -  stop - the type of stopping test, see `PEPStop`

880:    Options Database Keys:
881: +  -pep_stop_basic - sets the default stopping test
882: -  -pep_stop_user  - selects the user-defined stopping test

884:    Level: advanced

886: .seealso: [](ch:pep), `PEPGetStoppingTest()`, `PEPSetStoppingTestFunction()`, `PEPSetConvergenceTest()`, `PEPStop`
887: @*/
888: PetscErrorCode PEPSetStoppingTest(PEP pep,PEPStop stop)
889: {
890:   PetscFunctionBegin;
893:   switch (stop) {
894:     case PEP_STOP_BASIC: pep->stopping = PEPStoppingBasic; break;
895:     case PEP_STOP_USER:
896:       PetscCheck(pep->stoppinguser,PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetStoppingTestFunction() first");
897:       pep->stopping = pep->stoppinguser;
898:       break;
899:     default:
900:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
901:   }
902:   pep->stop = stop;
903:   PetscFunctionReturn(PETSC_SUCCESS);
904: }

906: /*@
907:    PEPGetStoppingTest - Gets the method used to decide the termination of the outer
908:    loop of the eigensolver.

910:    Not Collective

912:    Input Parameter:
913: .  pep   - the polynomial eigensolver context

915:    Output Parameter:
916: .  stop  - the type of stopping test

918:    Level: advanced

920: .seealso: [](ch:pep), `PEPSetStoppingTest()`, `PEPStop`
921: @*/
922: PetscErrorCode PEPGetStoppingTest(PEP pep,PEPStop *stop)
923: {
924:   PetscFunctionBegin;
926:   PetscAssertPointer(stop,2);
927:   *stop = pep->stop;
928:   PetscFunctionReturn(PETSC_SUCCESS);
929: }

931: /*@
932:    PEPSetScale - Specifies the scaling strategy to be used.

934:    Collective

936:    Input Parameters:
937: +  pep    - the polynomial eigensolver context
938: .  scale  - scaling strategy, see `PEPScale` for possible values
939: .  alpha  - the scaling factor used in the scalar strategy
940: .  Dl     - the left diagonal matrix of the diagonal scaling algorithm
941: .  Dr     - the right diagonal matrix of the diagonal scaling algorithm
942: .  its    - number of iterations of the diagonal scaling algorithm
943: -  lambda - approximation to wanted eigenvalues (modulus)

945:    Options Database Keys:
946: +  -pep_scale \<scale\>         - set the scaling type, one of `none`,`scalar`,`diagonal`,`both`
947: .  -pep_scale_factor \<alpha\>  - set the scaling factor
948: .  -pep_scale_its \<its\>       - set the number of iterations
949: -  -pep_scale_lambda \<lambda\> - set the approximation to eigenvalues

951:    Notes:
952:    There are two non-exclusive scaling strategies, scalar and diagonal.
953:    See discussion in section [](#sec:scaling).

955:    In the scalar strategy, scaling is applied to the eigenvalue, that is,
956:    $\mu = \lambda/\alpha$ is the new eigenvalue and all matrices are scaled
957:    accordingly. After solving the scaled problem, the original $\lambda$ is
958:    recovered. Parameter `alpha` must be positive. Use `PETSC_DETERMINE` to let
959:    the solver compute a reasonable scaling factor, and `PETSC_CURRENT` to
960:    retain a previously set value.

962:    In the diagonal strategy, the solver works implicitly with matrix $D_\ell P(\lambda)D_r$,
963:    where $D_\ell$ and $D_r$ are appropriate diagonal matrices. This improves the accuracy
964:    of the computed results in some cases. The user may provide the `Dl` and `Dr`
965:    matrices represented as `Vec` objects storing diagonal elements. If not
966:    provided, these matrices are computed internally. This option requires
967:    that the polynomial coefficient matrices are of `MATAIJ` type.
968:    The parameter `its` is the number of iterations performed by the method.
969:    Parameter `lambda` must be positive. Use `PETSC_DETERMINE` or set `lambda` = 1.0
970:    if no information about eigenvalues is available. `PETSC_CURRENT` can also
971:    be used to leave `its` and `lambda` unchanged.

973:    Level: intermediate

975: .seealso: [](ch:pep), [](#sec:scaling), `PEPGetScale()`
976: @*/
977: PetscErrorCode PEPSetScale(PEP pep,PEPScale scale,PetscReal alpha,Vec Dl,Vec Dr,PetscInt its,PetscReal lambda)
978: {
979:   PetscFunctionBegin;
982:   pep->scale = scale;
983:   if (scale==PEP_SCALE_SCALAR || scale==PEP_SCALE_BOTH) {
985:     if (alpha == (PetscReal)PETSC_DETERMINE) {
986:       pep->sfactor = 0.0;
987:       pep->sfactor_set = PETSC_FALSE;
988:     } else if (alpha != (PetscReal)PETSC_CURRENT) {
989:       PetscCheck(alpha>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of alpha. Must be > 0");
990:       pep->sfactor = alpha;
991:       pep->sfactor_set = PETSC_TRUE;
992:     }
993:   }
994:   if (scale==PEP_SCALE_DIAGONAL || scale==PEP_SCALE_BOTH) {
995:     if (Dl) {
997:       PetscCheckSameComm(pep,1,Dl,4);
998:       PetscCall(PetscObjectReference((PetscObject)Dl));
999:       PetscCall(VecDestroy(&pep->Dl));
1000:       pep->Dl = Dl;
1001:     }
1002:     if (Dr) {
1004:       PetscCheckSameComm(pep,1,Dr,5);
1005:       PetscCall(PetscObjectReference((PetscObject)Dr));
1006:       PetscCall(VecDestroy(&pep->Dr));
1007:       pep->Dr = Dr;
1008:     }
1011:     if (its==PETSC_DETERMINE) pep->sits = 5;
1012:     else if (its!=PETSC_CURRENT) {
1013:       PetscCheck(its>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be > 0");
1014:       pep->sits = its;
1015:     }
1016:     if (lambda == (PetscReal)PETSC_DETERMINE) pep->slambda = 1.0;
1017:     else if (lambda != (PetscReal)PETSC_CURRENT) {
1018:       PetscCheck(lambda>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of lambda. Must be > 0");
1019:       pep->slambda = lambda;
1020:     }
1021:   }
1022:   pep->state = PEP_STATE_INITIAL;
1023:   PetscFunctionReturn(PETSC_SUCCESS);
1024: }

1026: /*@
1027:    PEPGetScale - Gets the scaling strategy used by the `PEP` object, and the
1028:    associated parameters.

1030:    Not Collective

1032:    Input Parameter:
1033: .  pep - the polynomial eigensolver context

1035:    Output Parameters:
1036: +  scale  - scaling strategy
1037: .  alpha  - the scaling factor used in the scalar strategy
1038: .  Dl     - the left diagonal matrix of the diagonal scaling algorithm
1039: .  Dr     - the right diagonal matrix of the diagonal scaling algorithm
1040: .  its    - number of iterations of the diagonal scaling algorithm
1041: -  lambda - approximation to wanted eigenvalues (modulus)

1043:    Notes:
1044:    The user can specify `NULL` for any parameter that is not needed.

1046:    If `Dl` or `Dr` were not set by the user, then the ones computed internally are
1047:    returned (or a `NULL` pointer if called before `PEPSetUp()`).

1049:    Level: intermediate

1051: .seealso: [](ch:pep), `PEPSetScale()`, `PEPSetUp()`
1052: @*/
1053: PetscErrorCode PEPGetScale(PEP pep,PEPScale *scale,PetscReal *alpha,Vec *Dl,Vec *Dr,PetscInt *its,PetscReal *lambda)
1054: {
1055:   PetscFunctionBegin;
1057:   if (scale)  *scale  = pep->scale;
1058:   if (alpha)  *alpha  = pep->sfactor;
1059:   if (Dl)     *Dl     = pep->Dl;
1060:   if (Dr)     *Dr     = pep->Dr;
1061:   if (its)    *its    = pep->sits;
1062:   if (lambda) *lambda = pep->slambda;
1063:   PetscFunctionReturn(PETSC_SUCCESS);
1064: }

1066: /*@
1067:    PEPSetExtract - Specifies the extraction strategy to be used.

1069:    Logically Collective

1071:    Input Parameters:
1072: +  pep     - the polynomial eigensolver context
1073: -  extract - extraction strategy, see `PEPExtract` for possible values

1075:    Options Database Key:
1076: .  -pep_extract \<extract\> - extraction type, one of `none`,`norm`,`residual`,`structured`

1078:    Note:
1079:    This is relevant for solvers based on linearization. Once the solver has
1080:    converged, the polynomial eigenvectors can be extracted from the
1081:    eigenvectors of the linearized problem in different ways. See the
1082:    discussion in section [](#sec:pepextr).

1084:    Level: intermediate

1086: .seealso: [](ch:pep), [](#sec:pepextr), `PEPExtract`, `PEPGetExtract()`
1087: @*/
1088: PetscErrorCode PEPSetExtract(PEP pep,PEPExtract extract)
1089: {
1090:   PetscFunctionBegin;
1093:   pep->extract = extract;
1094:   PetscFunctionReturn(PETSC_SUCCESS);
1095: }

1097: /*@
1098:    PEPGetExtract - Gets the extraction strategy used by the `PEP` object.

1100:    Not Collective

1102:    Input Parameter:
1103: .  pep - the polynomial eigensolver context

1105:    Output Parameter:
1106: .  extract - extraction strategy

1108:    Level: intermediate

1110: .seealso: [](ch:pep), `PEPSetExtract()`, `PEPExtract`
1111: @*/
1112: PetscErrorCode PEPGetExtract(PEP pep,PEPExtract *extract)
1113: {
1114:   PetscFunctionBegin;
1116:   PetscAssertPointer(extract,2);
1117:   *extract = pep->extract;
1118:   PetscFunctionReturn(PETSC_SUCCESS);
1119: }

1121: /*@
1122:    PEPSetRefine - Specifies the refinement type (and options) to be used
1123:    after the solve.

1125:    Logically Collective

1127:    Input Parameters:
1128: +  pep    - the polynomial eigensolver context
1129: .  refine - refinement type, see `PEPRefine` for possible values
1130: .  npart  - number of partitions of the communicator
1131: .  tol    - the convergence tolerance
1132: .  its    - maximum number of refinement iterations
1133: -  scheme - the scheme for solving the involved linear systems, see `PEPRefineScheme`
1134:             for possible values

1136:    Options Database Keys:
1137: +  -pep_refine \<refine\>           - set the refinement type, one of `none`,`simple`,`multiple`
1138: .  -pep_refine_partitions \<npart\> - set the number of partitions
1139: .  -pep_refine_tol \<tol\>         - set the tolerance
1140: .  -pep_refine_its \<its\>         - set the number of iterations
1141: -  -pep_refine_scheme \<scheme\>   - set the scheme for the linear solves, `schur`,`mbe`, or `explicit`

1143:    Notes:
1144:    This function configures the parameters of Newton iterative refinement,
1145:    see section [](#sec:refine) for a discussion of the different strategies.

1147:    By default, iterative refinement is disabled, since it may be very
1148:    costly. There are two possible refinement strategies, simple and multiple.
1149:    The simple approach performs iterative refinement on each of the
1150:    converged eigenpairs individually, whereas the multiple strategy works
1151:    with the invariant pair as a whole, refining all eigenpairs simultaneously.
1152:    The latter may be required for the case of multiple eigenvalues.

1154:    In some cases, especially when using direct solvers within the
1155:    iterative refinement method, it may be helpful for improved scalability
1156:    to split the communicator in several partitions. The `npart` parameter
1157:    indicates how many partitions to use (defaults to 1).

1159:    The `tol` and `its` parameters specify the stopping criterion. In the simple
1160:    method, refinement continues until the residual of each eigenpair is
1161:    below the tolerance (`tol` defaults to the `PEP` tolerance, but may be set to a
1162:    different value). In contrast, the multiple method simply performs its
1163:    refinement iterations (just one by default).

1165:    The `scheme` argument is used to change the way in which linear systems are
1166:    solved. Possible choices are explicit, mixed block elimination (MBE),
1167:    and Schur complement.

1169:    Use `PETSC_CURRENT` to retain the current value of `npart`, `tol` or `its`. Use
1170:    `PETSC_DETERMINE` to assign a default value.

1172:    Level: intermediate

1174: .seealso: [](ch:pep), [](#sec:refine), `PEPGetRefine()`
1175: @*/
1176: PetscErrorCode PEPSetRefine(PEP pep,PEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,PEPRefineScheme scheme)
1177: {
1178:   PetscMPIInt    size;

1180:   PetscFunctionBegin;
1187:   pep->refine = refine;
1188:   if (refine) {  /* process parameters only if not REFINE_NONE */
1189:     if (npart!=pep->npart) {
1190:       PetscCall(PetscSubcommDestroy(&pep->refinesubc));
1191:       PetscCall(KSPDestroy(&pep->refineksp));
1192:     }
1193:     if (npart == PETSC_DETERMINE) {
1194:       pep->npart = 1;
1195:     } else if (npart != PETSC_CURRENT) {
1196:       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pep),&size));
1197:       PetscCheck(npart>0 && npart<=size,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
1198:       pep->npart = npart;
1199:     }
1200:     if (tol == (PetscReal)PETSC_DETERMINE) {
1201:       pep->rtol = PETSC_DETERMINE;
1202:     } else if (tol != (PetscReal)PETSC_CURRENT) {
1203:       PetscCheck(tol>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1204:       pep->rtol = tol;
1205:     }
1206:     if (its==PETSC_DETERMINE) {
1207:       pep->rits = PETSC_DETERMINE;
1208:     } else if (its != PETSC_CURRENT) {
1209:       PetscCheck(its>=0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be >= 0");
1210:       pep->rits = its;
1211:     }
1212:     pep->scheme = scheme;
1213:   }
1214:   pep->state = PEP_STATE_INITIAL;
1215:   PetscFunctionReturn(PETSC_SUCCESS);
1216: }

1218: /*@
1219:    PEPGetRefine - Gets the refinement strategy used by the `PEP` object, and the
1220:    associated parameters.

1222:    Not Collective

1224:    Input Parameter:
1225: .  pep - the polynomial eigensolver context

1227:    Output Parameters:
1228: +  refine - refinement type
1229: .  npart  - number of partitions of the communicator
1230: .  tol    - the convergence tolerance
1231: .  its    - maximum number of refinement iterations
1232: -  scheme - the scheme used for solving linear systems

1234:    Note:
1235:    The user can specify `NULL` for any parameter that is not needed.

1237:    Level: intermediate

1239: .seealso: [](ch:pep), `PEPSetRefine()`
1240: @*/
1241: PetscErrorCode PEPGetRefine(PEP pep,PEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,PEPRefineScheme *scheme)
1242: {
1243:   PetscFunctionBegin;
1245:   if (refine) *refine = pep->refine;
1246:   if (npart)  *npart  = pep->npart;
1247:   if (tol)    *tol    = pep->rtol;
1248:   if (its)    *its    = pep->rits;
1249:   if (scheme) *scheme = pep->scheme;
1250:   PetscFunctionReturn(PETSC_SUCCESS);
1251: }

1253: /*@
1254:    PEPSetOptionsPrefix - Sets the prefix used for searching for all
1255:    `PEP` options in the database.

1257:    Logically Collective

1259:    Input Parameters:
1260: +  pep    - the polynomial eigensolver context
1261: -  prefix - the prefix string to prepend to all `PEP` option requests

1263:    Notes:
1264:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1265:    The first character of all runtime options is AUTOMATICALLY the
1266:    hyphen.

1268:    For example, to distinguish between the runtime options for two
1269:    different `PEP` contexts, one could call
1270: .vb
1271:    PEPSetOptionsPrefix(pep1,"qeig1_")
1272:    PEPSetOptionsPrefix(pep2,"qeig2_")
1273: .ve

1275:    Level: advanced

1277: .seealso: [](ch:pep), `PEPAppendOptionsPrefix()`, `PEPGetOptionsPrefix()`
1278: @*/
1279: PetscErrorCode PEPSetOptionsPrefix(PEP pep,const char prefix[])
1280: {
1281:   PetscFunctionBegin;
1283:   if (!pep->st) PetscCall(PEPGetST(pep,&pep->st));
1284:   PetscCall(STSetOptionsPrefix(pep->st,prefix));
1285:   if (!pep->V) PetscCall(PEPGetBV(pep,&pep->V));
1286:   PetscCall(BVSetOptionsPrefix(pep->V,prefix));
1287:   if (!pep->ds) PetscCall(PEPGetDS(pep,&pep->ds));
1288:   PetscCall(DSSetOptionsPrefix(pep->ds,prefix));
1289:   if (!pep->rg) PetscCall(PEPGetRG(pep,&pep->rg));
1290:   PetscCall(RGSetOptionsPrefix(pep->rg,prefix));
1291:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pep,prefix));
1292:   PetscFunctionReturn(PETSC_SUCCESS);
1293: }

1295: /*@
1296:    PEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1297:    `PEP` options in the database.

1299:    Logically Collective

1301:    Input Parameters:
1302: +  pep    - the polynomial eigensolver context
1303: -  prefix - the prefix string to prepend to all `PEP` option requests

1305:    Notes:
1306:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1307:    The first character of all runtime options is AUTOMATICALLY the hyphen.

1309:    Level: advanced

1311: .seealso: [](ch:pep), `PEPSetOptionsPrefix()`, `PEPGetOptionsPrefix()`
1312: @*/
1313: PetscErrorCode PEPAppendOptionsPrefix(PEP pep,const char prefix[])
1314: {
1315:   PetscFunctionBegin;
1317:   if (!pep->st) PetscCall(PEPGetST(pep,&pep->st));
1318:   PetscCall(STAppendOptionsPrefix(pep->st,prefix));
1319:   if (!pep->V) PetscCall(PEPGetBV(pep,&pep->V));
1320:   PetscCall(BVAppendOptionsPrefix(pep->V,prefix));
1321:   if (!pep->ds) PetscCall(PEPGetDS(pep,&pep->ds));
1322:   PetscCall(DSAppendOptionsPrefix(pep->ds,prefix));
1323:   if (!pep->rg) PetscCall(PEPGetRG(pep,&pep->rg));
1324:   PetscCall(RGAppendOptionsPrefix(pep->rg,prefix));
1325:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pep,prefix));
1326:   PetscFunctionReturn(PETSC_SUCCESS);
1327: }

1329: /*@
1330:    PEPGetOptionsPrefix - Gets the prefix used for searching for all
1331:    `PEP` options in the database.

1333:    Not Collective

1335:    Input Parameter:
1336: .  pep - the polynomial eigensolver context

1338:    Output Parameter:
1339: .  prefix - pointer to the prefix string used is returned

1341:    Level: advanced

1343: .seealso: [](ch:pep), `PEPSetOptionsPrefix()`, `PEPAppendOptionsPrefix()`
1344: @*/
1345: PetscErrorCode PEPGetOptionsPrefix(PEP pep,const char *prefix[])
1346: {
1347:   PetscFunctionBegin;
1349:   PetscAssertPointer(prefix,2);
1350:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pep,prefix));
1351:   PetscFunctionReturn(PETSC_SUCCESS);
1352: }