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:    Notes:
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:    Notes:
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_UMLIMITED 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 that
374:    (a) in cases where nev is small, the user sets ncv (a reasonable default is 2*nev); and
375:    (b) in cases where nev is large, the user sets mpd.

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

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

384:    Level: intermediate

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

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

419:    Logically Collective

421:    Input Parameters:
422: +  pep   - the polynomial eigensolver context
423: -  which - the portion of the spectrum to be sought

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

437:    Notes:
438:    The parameter 'which' can have one of these values

440: +     PEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
441: .     PEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
442: .     PEP_LARGEST_REAL - largest real parts
443: .     PEP_SMALLEST_REAL - smallest real parts
444: .     PEP_LARGEST_IMAGINARY - largest imaginary parts
445: .     PEP_SMALLEST_IMAGINARY - smallest imaginary parts
446: .     PEP_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
447: .     PEP_TARGET_REAL - eigenvalues with real part closest to target
448: .     PEP_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
449: .     PEP_ALL - all eigenvalues contained in a given interval or region
450: -     PEP_WHICH_USER - user defined ordering set with PEPSetEigenvalueComparison()

452:    Not all eigensolvers implemented in PEP account for all the possible values
453:    stated above. If SLEPc is compiled for real numbers PEP_LARGEST_IMAGINARY
454:    and PEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
455:    for eigenvalue selection.

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

459:    The criterion PEP_TARGET_IMAGINARY is available only in case PETSc and
460:    SLEPc have been built with complex scalars.

462:    PEP_ALL is intended for use in combination with an interval (see
463:    PEPSetInterval()), when all eigenvalues within the interval are requested,
464:    and also for computing all eigenvalues in a region with the CISS solver.
465:    In both cases, the number of eigenvalues is unknown, so the nev parameter
466:    has a different sense, see PEPSetDimensions().

468:    Level: intermediate

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

506: /*@
507:     PEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
508:     sought.

510:     Not Collective

512:     Input Parameter:
513: .   pep - the polynomial eigensolver context

515:     Output Parameter:
516: .   which - the portion of the spectrum to be sought

518:     Notes:
519:     See PEPSetWhichEigenpairs() for possible values of 'which'.

521:     Level: intermediate

523: .seealso: [](ch:pep), `PEPSetWhichEigenpairs()`, `PEPWhich`
524: @*/
525: PetscErrorCode PEPGetWhichEigenpairs(PEP pep,PEPWhich *which)
526: {
527:   PetscFunctionBegin;
529:   PetscAssertPointer(which,2);
530:   *which = pep->which;
531:   PetscFunctionReturn(PETSC_SUCCESS);
532: }

534: /*@C
535:    PEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
536:    when PEPSetWhichEigenpairs() is set to PEP_WHICH_USER.

538:    Logically Collective

540:    Input Parameters:
541: +  pep  - the polynomial eigensolver context
542: .  comp - a pointer to the comparison function
543: -  ctx  - a context pointer (the last parameter to the comparison function)

545:    Note:
546:    The returning parameter 'res' can be
547: +  negative - if the 1st eigenvalue is preferred to the 2st one
548: .  zero     - if both eigenvalues are equally preferred
549: -  positive - if the 2st eigenvalue is preferred to the 1st one

551:    Level: advanced

553: .seealso: [](ch:pep), `PEPSetWhichEigenpairs()`, `PEPWhich`
554: @*/
555: PetscErrorCode PEPSetEigenvalueComparison(PEP pep,SlepcEigenvalueComparisonFn *comp,void *ctx)
556: {
557:   PetscFunctionBegin;
559:   pep->sc->comparison    = comp;
560:   pep->sc->comparisonctx = ctx;
561:   pep->which             = PEP_WHICH_USER;
562:   PetscFunctionReturn(PETSC_SUCCESS);
563: }

565: /*@
566:    PEPSetProblemType - Specifies the type of the polynomial eigenvalue problem.

568:    Logically Collective

570:    Input Parameters:
571: +  pep  - the polynomial eigensolver context
572: -  type - a known type of polynomial eigenvalue problem

574:    Options Database Keys:
575: +  -pep_general - general problem with no particular structure
576: .  -pep_hermitian - problem whose coefficient matrices are Hermitian
577: .  -pep_hyperbolic - Hermitian problem that satisfies the definition of hyperbolic
578: -  -pep_gyroscopic - problem with Hamiltonian structure

580:    Notes:
581:    Allowed values for the problem type are general (PEP_GENERAL), Hermitian
582:    (PEP_HERMITIAN), hyperbolic (PEP_HYPERBOLIC), and gyroscopic (PEP_GYROSCOPIC).

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

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

593:    Level: intermediate

595: .seealso: [](ch:pep), `PEPSetOperators()`, `PEPSetType()`, `PEPGetProblemType()`, `PEPProblemType`
596: @*/
597: PetscErrorCode PEPSetProblemType(PEP pep,PEPProblemType type)
598: {
599:   PetscFunctionBegin;
602:   PetscCheck(type==PEP_GENERAL || type==PEP_HERMITIAN || type==PEP_HYPERBOLIC || type==PEP_GYROSCOPIC,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
603:   if (type != pep->problem_type) {
604:     pep->problem_type = type;
605:     pep->state = PEP_STATE_INITIAL;
606:   }
607:   PetscFunctionReturn(PETSC_SUCCESS);
608: }

610: /*@
611:    PEPGetProblemType - Gets the problem type from the PEP object.

613:    Not Collective

615:    Input Parameter:
616: .  pep - the polynomial eigensolver context

618:    Output Parameter:
619: .  type - the problem type

621:    Level: intermediate

623: .seealso: [](ch:pep), `PEPSetProblemType()`, `PEPProblemType`
624: @*/
625: PetscErrorCode PEPGetProblemType(PEP pep,PEPProblemType *type)
626: {
627:   PetscFunctionBegin;
629:   PetscAssertPointer(type,2);
630:   *type = pep->problem_type;
631:   PetscFunctionReturn(PETSC_SUCCESS);
632: }

634: /*@
635:    PEPSetBasis - Specifies the type of polynomial basis used to describe the
636:    polynomial eigenvalue problem.

638:    Logically Collective

640:    Input Parameters:
641: +  pep   - the polynomial eigensolver context
642: -  basis - the type of polynomial basis

644:    Options Database Key:
645: .  -pep_basis <basis> - Select the basis type

647:    Notes:
648:    By default, the coefficient matrices passed via PEPSetOperators() are
649:    expressed in the monomial basis, i.e.
650:    P(lambda) = A_0 + lambda*A_1 + lambda^2*A_2 + ... + lambda^d*A_d.
651:    Other polynomial bases may have better numerical behavior, but the user
652:    must then pass the coefficient matrices accordingly.

654:    Level: intermediate

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

667: /*@
668:    PEPGetBasis - Gets the type of polynomial basis from the PEP object.

670:    Not Collective

672:    Input Parameter:
673: .  pep - the polynomial eigensolver context

675:    Output Parameter:
676: .  basis - the polynomial basis

678:    Level: intermediate

680: .seealso: [](ch:pep), `PEPSetBasis()`, `PEPBasis`
681: @*/
682: PetscErrorCode PEPGetBasis(PEP pep,PEPBasis *basis)
683: {
684:   PetscFunctionBegin;
686:   PetscAssertPointer(basis,2);
687:   *basis = pep->basis;
688:   PetscFunctionReturn(PETSC_SUCCESS);
689: }

691: /*@
692:    PEPSetTrackAll - Specifies if the solver must compute the residual of all
693:    approximate eigenpairs or not.

695:    Logically Collective

697:    Input Parameters:
698: +  pep      - the polynomial eigensolver context
699: -  trackall - whether compute all residuals or not

701:    Notes:
702:    If the user sets trackall=PETSC_TRUE then the solver explicitly computes
703:    the residual for each eigenpair approximation. Computing the residual is
704:    usually an expensive operation and solvers commonly compute the associated
705:    residual to the first unconverged eigenpair.

707:    The option '-pep_monitor_all' automatically activates this option.

709:    Level: developer

711: .seealso: [](ch:pep), `PEPGetTrackAll()`
712: @*/
713: PetscErrorCode PEPSetTrackAll(PEP pep,PetscBool trackall)
714: {
715:   PetscFunctionBegin;
718:   pep->trackall = trackall;
719:   PetscFunctionReturn(PETSC_SUCCESS);
720: }

722: /*@
723:    PEPGetTrackAll - Returns the flag indicating whether all residual norms must
724:    be computed or not.

726:    Not Collective

728:    Input Parameter:
729: .  pep - the polynomial eigensolver context

731:    Output Parameter:
732: .  trackall - the returned flag

734:    Level: developer

736: .seealso: [](ch:pep), `PEPSetTrackAll()`
737: @*/
738: PetscErrorCode PEPGetTrackAll(PEP pep,PetscBool *trackall)
739: {
740:   PetscFunctionBegin;
742:   PetscAssertPointer(trackall,2);
743:   *trackall = pep->trackall;
744:   PetscFunctionReturn(PETSC_SUCCESS);
745: }

747: /*@C
748:    PEPSetConvergenceTestFunction - Sets a function to compute the error estimate
749:    used in the convergence test.

751:    Logically Collective

753:    Input Parameters:
754: +  pep     - the polynomial eigensolver context
755: .  conv    - convergence test function, see PEPConvergenceTestFn for the calling sequence
756: .  ctx     - context for private data for the convergence routine (may be NULL)
757: -  destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence

759:    Note:
760:    If the error estimate returned by the convergence test function is less than
761:    the tolerance, then the eigenvalue is accepted as converged.

763:    Level: advanced

765: .seealso: [](ch:pep), `PEPSetConvergenceTest()`, `PEPSetTolerances()`
766: @*/
767: PetscErrorCode PEPSetConvergenceTestFunction(PEP pep,PEPConvergenceTestFn *conv,void *ctx,PetscCtxDestroyFn *destroy)
768: {
769:   PetscFunctionBegin;
771:   if (pep->convergeddestroy) PetscCall((*pep->convergeddestroy)(&pep->convergedctx));
772:   pep->convergeduser    = conv;
773:   pep->convergeddestroy = destroy;
774:   pep->convergedctx     = ctx;
775:   if (conv == PEPConvergedRelative) pep->conv = PEP_CONV_REL;
776:   else if (conv == PEPConvergedNorm) pep->conv = PEP_CONV_NORM;
777:   else if (conv == PEPConvergedAbsolute) pep->conv = PEP_CONV_ABS;
778:   else {
779:     pep->conv      = PEP_CONV_USER;
780:     pep->converged = pep->convergeduser;
781:   }
782:   PetscFunctionReturn(PETSC_SUCCESS);
783: }

785: /*@
786:    PEPSetConvergenceTest - Specifies how to compute the error estimate
787:    used in the convergence test.

789:    Logically Collective

791:    Input Parameters:
792: +  pep  - the polynomial eigensolver context
793: -  conv - the type of convergence test

795:    Options Database Keys:
796: +  -pep_conv_abs    - Sets the absolute convergence test
797: .  -pep_conv_rel    - Sets the convergence test relative to the eigenvalue
798: .  -pep_conv_norm   - Sets the convergence test relative to the matrix norms
799: -  -pep_conv_user   - Selects the user-defined convergence test

801:    Note:
802:    The parameter 'conv' can have one of these values
803: +     PEP_CONV_ABS    - absolute error ||r||
804: .     PEP_CONV_REL    - error relative to the eigenvalue l, ||r||/|l|
805: .     PEP_CONV_NORM   - error relative matrix norms, ||r||/sum_i(l^i*||A_i||)
806: -     PEP_CONV_USER   - function set by PEPSetConvergenceTestFunction()

808:    Level: intermediate

810: .seealso: [](ch:pep), `PEPGetConvergenceTest()`, `PEPSetConvergenceTestFunction()`, `PEPSetStoppingTest()`, `PEPConv`
811: @*/
812: PetscErrorCode PEPSetConvergenceTest(PEP pep,PEPConv conv)
813: {
814:   PetscFunctionBegin;
817:   switch (conv) {
818:     case PEP_CONV_ABS:  pep->converged = PEPConvergedAbsolute; break;
819:     case PEP_CONV_REL:  pep->converged = PEPConvergedRelative; break;
820:     case PEP_CONV_NORM: pep->converged = PEPConvergedNorm; break;
821:     case PEP_CONV_USER:
822:       PetscCheck(pep->convergeduser,PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetConvergenceTestFunction() first");
823:       pep->converged = pep->convergeduser;
824:       break;
825:     default:
826:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
827:   }
828:   pep->conv = conv;
829:   PetscFunctionReturn(PETSC_SUCCESS);
830: }

832: /*@
833:    PEPGetConvergenceTest - Gets the method used to compute the error estimate
834:    used in the convergence test.

836:    Not Collective

838:    Input Parameter:
839: .  pep   - the polynomial eigensolver context

841:    Output Parameter:
842: .  conv  - the type of convergence test

844:    Level: intermediate

846: .seealso: [](ch:pep), `PEPSetConvergenceTest()`, `PEPConv`
847: @*/
848: PetscErrorCode PEPGetConvergenceTest(PEP pep,PEPConv *conv)
849: {
850:   PetscFunctionBegin;
852:   PetscAssertPointer(conv,2);
853:   *conv = pep->conv;
854:   PetscFunctionReturn(PETSC_SUCCESS);
855: }

857: /*@C
858:    PEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
859:    iteration of the eigensolver.

861:    Logically Collective

863:    Input Parameters:
864: +  pep     - the polynomial eigensolver context
865: .  stop    - stopping test function, see PEPStoppingTestFn for the calling sequence
866: .  ctx     - context for private data for the stopping routine (may be NULL)
867: -  destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence

869:    Note:
870:    Normal usage is to first call the default routine PEPStoppingBasic() and then
871:    set reason to PEP_CONVERGED_USER if some user-defined conditions have been
872:    met. To let the eigensolver continue iterating, the result must be left as
873:    PEP_CONVERGED_ITERATING.

875:    Level: advanced

877: .seealso: [](ch:pep), `PEPSetStoppingTest()`, `PEPStoppingBasic()`
878: @*/
879: PetscErrorCode PEPSetStoppingTestFunction(PEP pep,PEPStoppingTestFn *stop,void *ctx,PetscCtxDestroyFn *destroy)
880: {
881:   PetscFunctionBegin;
883:   if (pep->stoppingdestroy) PetscCall((*pep->stoppingdestroy)(&pep->stoppingctx));
884:   pep->stoppinguser    = stop;
885:   pep->stoppingdestroy = destroy;
886:   pep->stoppingctx     = ctx;
887:   if (stop == PEPStoppingBasic) pep->stop = PEP_STOP_BASIC;
888:   else {
889:     pep->stop     = PEP_STOP_USER;
890:     pep->stopping = pep->stoppinguser;
891:   }
892:   PetscFunctionReturn(PETSC_SUCCESS);
893: }

895: /*@
896:    PEPSetStoppingTest - Specifies how to decide the termination of the outer
897:    loop of the eigensolver.

899:    Logically Collective

901:    Input Parameters:
902: +  pep  - the polynomial eigensolver context
903: -  stop - the type of stopping test

905:    Options Database Keys:
906: +  -pep_stop_basic - Sets the default stopping test
907: -  -pep_stop_user  - Selects the user-defined stopping test

909:    Note:
910:    The parameter 'stop' can have one of these values
911: +     PEP_STOP_BASIC - default stopping test
912: -     PEP_STOP_USER  - function set by PEPSetStoppingTestFunction()

914:    Level: advanced

916: .seealso: [](ch:pep), `PEPGetStoppingTest()`, `PEPSetStoppingTestFunction()`, `PEPSetConvergenceTest()`, `PEPStop`
917: @*/
918: PetscErrorCode PEPSetStoppingTest(PEP pep,PEPStop stop)
919: {
920:   PetscFunctionBegin;
923:   switch (stop) {
924:     case PEP_STOP_BASIC: pep->stopping = PEPStoppingBasic; break;
925:     case PEP_STOP_USER:
926:       PetscCheck(pep->stoppinguser,PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetStoppingTestFunction() first");
927:       pep->stopping = pep->stoppinguser;
928:       break;
929:     default:
930:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
931:   }
932:   pep->stop = stop;
933:   PetscFunctionReturn(PETSC_SUCCESS);
934: }

936: /*@
937:    PEPGetStoppingTest - Gets the method used to decide the termination of the outer
938:    loop of the eigensolver.

940:    Not Collective

942:    Input Parameter:
943: .  pep   - the polynomial eigensolver context

945:    Output Parameter:
946: .  stop  - the type of stopping test

948:    Level: advanced

950: .seealso: [](ch:pep), `PEPSetStoppingTest()`, `PEPStop`
951: @*/
952: PetscErrorCode PEPGetStoppingTest(PEP pep,PEPStop *stop)
953: {
954:   PetscFunctionBegin;
956:   PetscAssertPointer(stop,2);
957:   *stop = pep->stop;
958:   PetscFunctionReturn(PETSC_SUCCESS);
959: }

961: /*@
962:    PEPSetScale - Specifies the scaling strategy to be used.

964:    Collective

966:    Input Parameters:
967: +  pep    - the polynomial eigensolver context
968: .  scale  - scaling strategy
969: .  alpha  - the scaling factor used in the scalar strategy
970: .  Dl     - the left diagonal matrix of the diagonal scaling algorithm
971: .  Dr     - the right diagonal matrix of the diagonal scaling algorithm
972: .  its    - number of iterations of the diagonal scaling algorithm
973: -  lambda - approximation to wanted eigenvalues (modulus)

975:    Options Database Keys:
976: +  -pep_scale <type> - scaling type, one of <none,scalar,diagonal,both>
977: .  -pep_scale_factor <alpha> - the scaling factor
978: .  -pep_scale_its <its> - number of iterations
979: -  -pep_scale_lambda <lambda> - approximation to eigenvalues

981:    Notes:
982:    There are two non-exclusive scaling strategies, scalar and diagonal.

984:    In the scalar strategy, scaling is applied to the eigenvalue, that is,
985:    mu = lambda/alpha is the new eigenvalue and all matrices are scaled
986:    accordingly. After solving the scaled problem, the original lambda is
987:    recovered. Parameter 'alpha' must be positive. Use PETSC_DETERMINE to let
988:    the solver compute a reasonable scaling factor, and PETSC_CURRENT to
989:    retain a previously set value.

991:    In the diagonal strategy, the solver works implicitly with matrix Dl*A*Dr,
992:    where Dl and Dr are appropriate diagonal matrices. This improves the accuracy
993:    of the computed results in some cases. The user may provide the Dr and Dl
994:    matrices represented as Vec objects storing diagonal elements. If not
995:    provided, these matrices are computed internally. This option requires
996:    that the polynomial coefficient matrices are of MATAIJ type.
997:    The parameter 'its' is the number of iterations performed by the method.
998:    Parameter 'lambda' must be positive. Use PETSC_DETERMINE or set lambda = 1.0
999:    if no information about eigenvalues is available. PETSC_CURRENT can also
1000:    be used to leave its and lambda unchanged.

1002:    Level: intermediate

1004: .seealso: [](ch:pep), `PEPGetScale()`
1005: @*/
1006: PetscErrorCode PEPSetScale(PEP pep,PEPScale scale,PetscReal alpha,Vec Dl,Vec Dr,PetscInt its,PetscReal lambda)
1007: {
1008:   PetscFunctionBegin;
1011:   pep->scale = scale;
1012:   if (scale==PEP_SCALE_SCALAR || scale==PEP_SCALE_BOTH) {
1014:     if (alpha == (PetscReal)PETSC_DETERMINE) {
1015:       pep->sfactor = 0.0;
1016:       pep->sfactor_set = PETSC_FALSE;
1017:     } else if (alpha != (PetscReal)PETSC_CURRENT) {
1018:       PetscCheck(alpha>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of alpha. Must be > 0");
1019:       pep->sfactor = alpha;
1020:       pep->sfactor_set = PETSC_TRUE;
1021:     }
1022:   }
1023:   if (scale==PEP_SCALE_DIAGONAL || scale==PEP_SCALE_BOTH) {
1024:     if (Dl) {
1026:       PetscCheckSameComm(pep,1,Dl,4);
1027:       PetscCall(PetscObjectReference((PetscObject)Dl));
1028:       PetscCall(VecDestroy(&pep->Dl));
1029:       pep->Dl = Dl;
1030:     }
1031:     if (Dr) {
1033:       PetscCheckSameComm(pep,1,Dr,5);
1034:       PetscCall(PetscObjectReference((PetscObject)Dr));
1035:       PetscCall(VecDestroy(&pep->Dr));
1036:       pep->Dr = Dr;
1037:     }
1040:     if (its==PETSC_DETERMINE) pep->sits = 5;
1041:     else if (its!=PETSC_CURRENT) {
1042:       PetscCheck(its>0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be > 0");
1043:       pep->sits = its;
1044:     }
1045:     if (lambda == (PetscReal)PETSC_DETERMINE) pep->slambda = 1.0;
1046:     else if (lambda != (PetscReal)PETSC_CURRENT) {
1047:       PetscCheck(lambda>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of lambda. Must be > 0");
1048:       pep->slambda = lambda;
1049:     }
1050:   }
1051:   pep->state = PEP_STATE_INITIAL;
1052:   PetscFunctionReturn(PETSC_SUCCESS);
1053: }

1055: /*@
1056:    PEPGetScale - Gets the scaling strategy used by the PEP object, and the
1057:    associated parameters.

1059:    Not Collective

1061:    Input Parameter:
1062: .  pep - the polynomial eigensolver context

1064:    Output Parameters:
1065: +  scale  - scaling strategy
1066: .  alpha  - the scaling factor used in the scalar strategy
1067: .  Dl     - the left diagonal matrix of the diagonal scaling algorithm
1068: .  Dr     - the right diagonal matrix of the diagonal scaling algorithm
1069: .  its    - number of iterations of the diagonal scaling algorithm
1070: -  lambda - approximation to wanted eigenvalues (modulus)

1072:    Level: intermediate

1074:    Note:
1075:    The user can specify NULL for any parameter that is not needed.

1077:    If Dl or Dr were not set by the user, then the ones computed internally are
1078:    returned (or a null pointer if called before PEPSetUp).

1080: .seealso: [](ch:pep), `PEPSetScale()`, `PEPSetUp()`
1081: @*/
1082: PetscErrorCode PEPGetScale(PEP pep,PEPScale *scale,PetscReal *alpha,Vec *Dl,Vec *Dr,PetscInt *its,PetscReal *lambda)
1083: {
1084:   PetscFunctionBegin;
1086:   if (scale)  *scale  = pep->scale;
1087:   if (alpha)  *alpha  = pep->sfactor;
1088:   if (Dl)     *Dl     = pep->Dl;
1089:   if (Dr)     *Dr     = pep->Dr;
1090:   if (its)    *its    = pep->sits;
1091:   if (lambda) *lambda = pep->slambda;
1092:   PetscFunctionReturn(PETSC_SUCCESS);
1093: }

1095: /*@
1096:    PEPSetExtract - Specifies the extraction strategy to be used.

1098:    Logically Collective

1100:    Input Parameters:
1101: +  pep     - the polynomial eigensolver context
1102: -  extract - extraction strategy

1104:    Options Database Key:
1105: .  -pep_extract <type> - extraction type, one of <none,norm,residual,structured>

1107:    Level: intermediate

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

1120: /*@
1121:    PEPGetExtract - Gets the extraction strategy used by the PEP object.

1123:    Not Collective

1125:    Input Parameter:
1126: .  pep - the polynomial eigensolver context

1128:    Output Parameter:
1129: .  extract - extraction strategy

1131:    Level: intermediate

1133: .seealso: [](ch:pep), `PEPSetExtract()`, `PEPExtract`
1134: @*/
1135: PetscErrorCode PEPGetExtract(PEP pep,PEPExtract *extract)
1136: {
1137:   PetscFunctionBegin;
1139:   PetscAssertPointer(extract,2);
1140:   *extract = pep->extract;
1141:   PetscFunctionReturn(PETSC_SUCCESS);
1142: }

1144: /*@
1145:    PEPSetRefine - Specifies the refinement type (and options) to be used
1146:    after the solve.

1148:    Logically Collective

1150:    Input Parameters:
1151: +  pep    - the polynomial eigensolver context
1152: .  refine - refinement type
1153: .  npart  - number of partitions of the communicator
1154: .  tol    - the convergence tolerance
1155: .  its    - maximum number of refinement iterations
1156: -  scheme - which scheme to be used for solving the involved linear systems

1158:    Options Database Keys:
1159: +  -pep_refine <type> - refinement type, one of <none,simple,multiple>
1160: .  -pep_refine_partitions <n> - the number of partitions
1161: .  -pep_refine_tol <tol> - the tolerance
1162: .  -pep_refine_its <its> - number of iterations
1163: -  -pep_refine_scheme - to set the scheme for the linear solves

1165:    Notes:
1166:    By default, iterative refinement is disabled, since it may be very
1167:    costly. There are two possible refinement strategies, simple and multiple.
1168:    The simple approach performs iterative refinement on each of the
1169:    converged eigenpairs individually, whereas the multiple strategy works
1170:    with the invariant pair as a whole, refining all eigenpairs simultaneously.
1171:    The latter may be required for the case of multiple eigenvalues.

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

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

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

1188:    Use PETSC_CURRENT to retain the current value of npart, tol or its. Use
1189:    PETSC_DETERMINE to assign a default value.

1191:    Level: intermediate

1193: .seealso: [](ch:pep), `PEPGetRefine()`
1194: @*/
1195: PetscErrorCode PEPSetRefine(PEP pep,PEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,PEPRefineScheme scheme)
1196: {
1197:   PetscMPIInt    size;

1199:   PetscFunctionBegin;
1206:   pep->refine = refine;
1207:   if (refine) {  /* process parameters only if not REFINE_NONE */
1208:     if (npart!=pep->npart) {
1209:       PetscCall(PetscSubcommDestroy(&pep->refinesubc));
1210:       PetscCall(KSPDestroy(&pep->refineksp));
1211:     }
1212:     if (npart == PETSC_DETERMINE) {
1213:       pep->npart = 1;
1214:     } else if (npart != PETSC_CURRENT) {
1215:       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)pep),&size));
1216:       PetscCheck(npart>0 && npart<=size,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
1217:       pep->npart = npart;
1218:     }
1219:     if (tol == (PetscReal)PETSC_DETERMINE) {
1220:       pep->rtol = PETSC_DETERMINE;
1221:     } else if (tol != (PetscReal)PETSC_CURRENT) {
1222:       PetscCheck(tol>0.0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1223:       pep->rtol = tol;
1224:     }
1225:     if (its==PETSC_DETERMINE) {
1226:       pep->rits = PETSC_DETERMINE;
1227:     } else if (its != PETSC_CURRENT) {
1228:       PetscCheck(its>=0,PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be >= 0");
1229:       pep->rits = its;
1230:     }
1231:     pep->scheme = scheme;
1232:   }
1233:   pep->state = PEP_STATE_INITIAL;
1234:   PetscFunctionReturn(PETSC_SUCCESS);
1235: }

1237: /*@
1238:    PEPGetRefine - Gets the refinement strategy used by the PEP object, and the
1239:    associated parameters.

1241:    Not Collective

1243:    Input Parameter:
1244: .  pep - the polynomial eigensolver context

1246:    Output Parameters:
1247: +  refine - refinement type
1248: .  npart  - number of partitions of the communicator
1249: .  tol    - the convergence tolerance
1250: .  its    - maximum number of refinement iterations
1251: -  scheme - the scheme used for solving linear systems

1253:    Level: intermediate

1255:    Note:
1256:    The user can specify NULL for any parameter that is not needed.

1258: .seealso: [](ch:pep), `PEPSetRefine()`
1259: @*/
1260: PetscErrorCode PEPGetRefine(PEP pep,PEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,PEPRefineScheme *scheme)
1261: {
1262:   PetscFunctionBegin;
1264:   if (refine) *refine = pep->refine;
1265:   if (npart)  *npart  = pep->npart;
1266:   if (tol)    *tol    = pep->rtol;
1267:   if (its)    *its    = pep->rits;
1268:   if (scheme) *scheme = pep->scheme;
1269:   PetscFunctionReturn(PETSC_SUCCESS);
1270: }

1272: /*@
1273:    PEPSetOptionsPrefix - Sets the prefix used for searching for all
1274:    PEP options in the database.

1276:    Logically Collective

1278:    Input Parameters:
1279: +  pep - the polynomial eigensolver context
1280: -  prefix - the prefix string to prepend to all PEP option requests

1282:    Notes:
1283:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1284:    The first character of all runtime options is AUTOMATICALLY the
1285:    hyphen.

1287:    For example, to distinguish between the runtime options for two
1288:    different PEP contexts, one could call
1289: .vb
1290:       PEPSetOptionsPrefix(pep1,"qeig1_")
1291:       PEPSetOptionsPrefix(pep2,"qeig2_")
1292: .ve

1294:    Level: advanced

1296: .seealso: [](ch:pep), `PEPAppendOptionsPrefix()`, `PEPGetOptionsPrefix()`
1297: @*/
1298: PetscErrorCode PEPSetOptionsPrefix(PEP pep,const char prefix[])
1299: {
1300:   PetscFunctionBegin;
1302:   if (!pep->st) PetscCall(PEPGetST(pep,&pep->st));
1303:   PetscCall(STSetOptionsPrefix(pep->st,prefix));
1304:   if (!pep->V) PetscCall(PEPGetBV(pep,&pep->V));
1305:   PetscCall(BVSetOptionsPrefix(pep->V,prefix));
1306:   if (!pep->ds) PetscCall(PEPGetDS(pep,&pep->ds));
1307:   PetscCall(DSSetOptionsPrefix(pep->ds,prefix));
1308:   if (!pep->rg) PetscCall(PEPGetRG(pep,&pep->rg));
1309:   PetscCall(RGSetOptionsPrefix(pep->rg,prefix));
1310:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)pep,prefix));
1311:   PetscFunctionReturn(PETSC_SUCCESS);
1312: }

1314: /*@
1315:    PEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1316:    PEP options in the database.

1318:    Logically Collective

1320:    Input Parameters:
1321: +  pep - the polynomial eigensolver context
1322: -  prefix - the prefix string to prepend to all PEP option requests

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

1328:    Level: advanced

1330: .seealso: [](ch:pep), `PEPSetOptionsPrefix()`, `PEPGetOptionsPrefix()`
1331: @*/
1332: PetscErrorCode PEPAppendOptionsPrefix(PEP pep,const char prefix[])
1333: {
1334:   PetscFunctionBegin;
1336:   if (!pep->st) PetscCall(PEPGetST(pep,&pep->st));
1337:   PetscCall(STAppendOptionsPrefix(pep->st,prefix));
1338:   if (!pep->V) PetscCall(PEPGetBV(pep,&pep->V));
1339:   PetscCall(BVAppendOptionsPrefix(pep->V,prefix));
1340:   if (!pep->ds) PetscCall(PEPGetDS(pep,&pep->ds));
1341:   PetscCall(DSAppendOptionsPrefix(pep->ds,prefix));
1342:   if (!pep->rg) PetscCall(PEPGetRG(pep,&pep->rg));
1343:   PetscCall(RGAppendOptionsPrefix(pep->rg,prefix));
1344:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)pep,prefix));
1345:   PetscFunctionReturn(PETSC_SUCCESS);
1346: }

1348: /*@
1349:    PEPGetOptionsPrefix - Gets the prefix used for searching for all
1350:    `PEP` options in the database.

1352:    Not Collective

1354:    Input Parameter:
1355: .  pep - the polynomial eigensolver context

1357:    Output Parameter:
1358: .  prefix - pointer to the prefix string used is returned

1360:    Level: advanced

1362: .seealso: [](ch:pep), `PEPSetOptionsPrefix()`, `PEPAppendOptionsPrefix()`
1363: @*/
1364: PetscErrorCode PEPGetOptionsPrefix(PEP pep,const char *prefix[])
1365: {
1366:   PetscFunctionBegin;
1368:   PetscAssertPointer(prefix,2);
1369:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pep,prefix));
1370:   PetscFunctionReturn(PETSC_SUCCESS);
1371: }