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[],PetscCtx 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(PetscOptionsBoolGroup("-pep_all","Compute all eigenvalues in an interval or a region","PEPSetWhichEigenpairs",&flg));
186:     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_ALL));
187:     PetscCall(PetscOptionsBoolGroupEnd("-pep_which_user","Select the user-defined selection criterion","PEPSetWhichEigenpairs",&flg));
188:     if (flg) PetscCall(PEPSetWhichEigenpairs(pep,PEP_WHICH_USER));

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

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

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

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

223:     PetscTryTypeMethod(pep,setfromoptions,PetscOptionsObject);
224:     PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)pep,PetscOptionsObject));
225:   PetscOptionsEnd();

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

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

246:    Not Collective

248:    Input Parameter:
249: .  pep - the polynomial eigensolver context

251:    Output Parameters:
252: +  tol - the convergence tolerance
253: -  maxits - maximum number of iterations

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

258:    Level: intermediate

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

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

275:    Logically Collective

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

282:    Options Database Keys:
283: +  -pep_tol tol       - sets the convergence tolerance
284: -  -pep_max_it maxits - sets the maximum number of iterations allowed

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

292:    Level: intermediate

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

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

325:    Not Collective

327:    Input Parameter:
328: .  pep - the polynomial eigensolver context

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

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

338:    Level: intermediate

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

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

356:    Logically Collective

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

364:    Options Database Keys:
365: +  -pep_nev nev - sets the number of eigenvalues
366: .  -pep_ncv ncv - sets the dimension of the subspace
367: -  -pep_mpd mpd - sets the maximum projected dimension

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

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

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

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

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

388:    Level: intermediate

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

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

423:    Logically Collective

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

429:    Options Database Keys:
430: +  -pep_largest_magnitude  - sets largest eigenvalues in magnitude
431: .  -pep_smallest_magnitude - sets smallest eigenvalues in magnitude
432: .  -pep_largest_real       - sets largest real parts
433: .  -pep_smallest_real      - sets smallest real parts
434: .  -pep_largest_imaginary  - sets largest imaginary parts
435: .  -pep_smallest_imaginary - sets smallest imaginary parts
436: .  -pep_target_magnitude   - sets eigenvalues closest to target
437: .  -pep_target_real        - sets real parts closest to target
438: .  -pep_target_imaginary   - sets imaginary parts closest to target
439: .  -pep_all                - sets all eigenvalues in an interval or region
440: -  -pep_which_user         - select the user-defined selection criterion

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

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

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

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

458:    Level: intermediate

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

496: /*@
497:     PEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
498:     sought.

500:     Not Collective

502:     Input Parameter:
503: .   pep - the polynomial eigensolver context

505:     Output Parameter:
506: .   which - the portion of the spectrum to be sought

508:     Level: intermediate

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

521: /*@C
522:    PEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
523:    when `PEPSetWhichEigenpairs()` is set to `PEP_WHICH_USER`.

525:    Logically Collective

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

532:    Level: advanced

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

546: /*@
547:    PEPSetProblemType - Specifies the type of the polynomial eigenvalue problem.

549:    Logically Collective

551:    Input Parameters:
552: +  pep  - the polynomial eigensolver context
553: -  type - a known type of polynomial eigenvalue problem

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

561:    Notes:
562:    See `PEPProblemType` for possible problem types.

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

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

573:    Level: intermediate

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

590: /*@
591:    PEPGetProblemType - Gets the problem type from the `PEP` object.

593:    Not Collective

595:    Input Parameter:
596: .  pep - the polynomial eigensolver context

598:    Output Parameter:
599: .  type - the problem type

601:    Level: intermediate

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

614: /*@
615:    PEPSetBasis - Specifies the type of polynomial basis used to describe the
616:    polynomial eigenvalue problem.

618:    Logically Collective

620:    Input Parameters:
621: +  pep   - the polynomial eigensolver context
622: -  basis - the type of polynomial basis, see `PEPBasis` for possible values

624:    Options Database Key:
625: .  -pep_basis (monomial|chebyshev1|chebyshev2|legendre|laguerre|hermite) - select the basis type

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

634:    Level: intermediate

636: .seealso: [](ch:pep), `PEPSetOperators()`, `PEPGetBasis()`, `PEPBasis`
637: @*/
638: PetscErrorCode PEPSetBasis(PEP pep,PEPBasis basis)
639: {
640:   PetscFunctionBegin;
643:   pep->basis = basis;
644:   PetscFunctionReturn(PETSC_SUCCESS);
645: }

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

650:    Not Collective

652:    Input Parameter:
653: .  pep - the polynomial eigensolver context

655:    Output Parameter:
656: .  basis - the polynomial basis

658:    Level: intermediate

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

671: /*@
672:    PEPSetTrackAll - Specifies if the solver must compute the residual of all
673:    approximate eigenpairs or not.

675:    Logically Collective

677:    Input Parameters:
678: +  pep      - the polynomial eigensolver context
679: -  trackall - whether compute all residuals or not

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

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

689:    Level: developer

691: .seealso: [](ch:pep), `PEPGetTrackAll()`
692: @*/
693: PetscErrorCode PEPSetTrackAll(PEP pep,PetscBool trackall)
694: {
695:   PetscFunctionBegin;
698:   pep->trackall = trackall;
699:   PetscFunctionReturn(PETSC_SUCCESS);
700: }

702: /*@
703:    PEPGetTrackAll - Returns the flag indicating whether all residual norms must
704:    be computed or not.

706:    Not Collective

708:    Input Parameter:
709: .  pep - the polynomial eigensolver context

711:    Output Parameter:
712: .  trackall - the returned flag

714:    Level: developer

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

727: /*@C
728:    PEPSetConvergenceTestFunction - Sets a function to compute the error estimate
729:    used in the convergence test.

731:    Logically Collective

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

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

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

747:    Level: advanced

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

769: /*@
770:    PEPSetConvergenceTest - Specifies how to compute the error estimate
771:    used in the convergence test.

773:    Logically Collective

775:    Input Parameters:
776: +  pep  - the polynomial eigensolver context
777: -  conv - the type of convergence test, see `PEPConv` for possible values

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

785:    Level: intermediate

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

809: /*@
810:    PEPGetConvergenceTest - Gets the method used to compute the error estimate
811:    used in the convergence test.

813:    Not Collective

815:    Input Parameter:
816: .  pep   - the polynomial eigensolver context

818:    Output Parameter:
819: .  conv  - the type of convergence test

821:    Level: intermediate

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

834: /*@C
835:    PEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
836:    iteration of the eigensolver.

838:    Logically Collective

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

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

853:    Level: advanced

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

873: /*@
874:    PEPSetStoppingTest - Specifies how to decide the termination of the outer
875:    loop of the eigensolver.

877:    Logically Collective

879:    Input Parameters:
880: +  pep  - the polynomial eigensolver context
881: -  stop - the type of stopping test, see `PEPStop`

883:    Options Database Keys:
884: +  -pep_stop_basic - sets the default stopping test
885: -  -pep_stop_user  - selects the user-defined stopping test

887:    Level: advanced

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

909: /*@
910:    PEPGetStoppingTest - Gets the method used to decide the termination of the outer
911:    loop of the eigensolver.

913:    Not Collective

915:    Input Parameter:
916: .  pep   - the polynomial eigensolver context

918:    Output Parameter:
919: .  stop  - the type of stopping test

921:    Level: advanced

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

934: /*@
935:    PEPSetScale - Specifies the scaling strategy to be used.

937:    Collective

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

948:    Options Database Keys:
949: +  -pep_scale (none|scalar|diagonal|both) - set the scaling type
950: .  -pep_scale_factor alpha                - set the scaling factor
951: .  -pep_scale_its its                     - set the number of iterations
952: -  -pep_scale_lambda lambda               - set the approximation to eigenvalues

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

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

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

976:    Level: intermediate

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

1029: /*@
1030:    PEPGetScale - Gets the scaling strategy used by the `PEP` object, and the
1031:    associated parameters.

1033:    Not Collective

1035:    Input Parameter:
1036: .  pep - the polynomial eigensolver context

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

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

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

1052:    Level: intermediate

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

1069: /*@
1070:    PEPSetExtract - Specifies the extraction strategy to be used.

1072:    Logically Collective

1074:    Input Parameters:
1075: +  pep     - the polynomial eigensolver context
1076: -  extract - extraction strategy, see `PEPExtract` for possible values

1078:    Options Database Key:
1079: .  -pep_extract (none|norm|residual|structured) - set the extraction type

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

1087:    Level: intermediate

1089: .seealso: [](ch:pep), [](#sec:pepextr), `PEPExtract`, `PEPGetExtract()`
1090: @*/
1091: PetscErrorCode PEPSetExtract(PEP pep,PEPExtract extract)
1092: {
1093:   PetscFunctionBegin;
1096:   pep->extract = extract;
1097:   PetscFunctionReturn(PETSC_SUCCESS);
1098: }

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

1103:    Not Collective

1105:    Input Parameter:
1106: .  pep - the polynomial eigensolver context

1108:    Output Parameter:
1109: .  extract - extraction strategy

1111:    Level: intermediate

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

1124: /*@
1125:    PEPSetRefine - Specifies the refinement type (and options) to be used
1126:    after the solve.

1128:    Logically Collective

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

1139:    Options Database Keys:
1140: +  -pep_refine (none|simple|multiple)      - set the refinement type
1141: .  -pep_refine_partitions npart            - set the number of partitions
1142: .  -pep_refine_tol tol                     - set the tolerance
1143: .  -pep_refine_its its                     - set the number of iterations
1144: -  -pep_refine_scheme (schur|mbe|explicit) - set the scheme for the linear solves

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

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

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

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

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

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

1175:    Level: intermediate

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

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

1221: /*@
1222:    PEPGetRefine - Gets the refinement strategy used by the `PEP` object, and the
1223:    associated parameters.

1225:    Not Collective

1227:    Input Parameter:
1228: .  pep - the polynomial eigensolver context

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

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

1240:    Level: intermediate

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

1256: /*@
1257:    PEPSetOptionsPrefix - Sets the prefix used for searching for all
1258:    `PEP` options in the database.

1260:    Logically Collective

1262:    Input Parameters:
1263: +  pep    - the polynomial eigensolver context
1264: -  prefix - the prefix string to prepend to all `PEP` option requests

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

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

1278:    Level: advanced

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

1298: /*@
1299:    PEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1300:    `PEP` options in the database.

1302:    Logically Collective

1304:    Input Parameters:
1305: +  pep    - the polynomial eigensolver context
1306: -  prefix - the prefix string to prepend to all `PEP` option requests

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

1312:    Level: advanced

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

1332: /*@
1333:    PEPGetOptionsPrefix - Gets the prefix used for searching for all
1334:    `PEP` options in the database.

1336:    Not Collective

1338:    Input Parameter:
1339: .  pep - the polynomial eigensolver context

1341:    Output Parameter:
1342: .  prefix - pointer to the prefix string used is returned

1344:    Level: advanced

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