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: 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 set the solver type.

 72:    Collective

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

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

 80:    Level: beginner

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

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

511:     Not Collective

513:     Input Parameter:
514: .   pep - eigensolver context obtained from PEPCreate()

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

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

522:     Level: intermediate

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

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

539:    Logically Collective

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

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

552:    Level: advanced

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

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

569:    Logically Collective

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

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

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

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

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

594:    Level: intermediate

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

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

614:    Not Collective

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

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

622:    Level: intermediate

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

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

639:    Logically Collective

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

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

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

655:    Level: intermediate

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

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

671:    Not Collective

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

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

679:    Level: intermediate

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

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

696:    Logically Collective

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

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

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

710:    Level: developer

712: .seealso: PEPGetTrackAll()
713: @*/
714: PetscErrorCode PEPSetTrackAll(PEP pep,PetscBool trackall)
715: {
716:   PetscFunctionBegin;
719:   pep->trackall = trackall;
720:   PetscFunctionReturn(PETSC_SUCCESS);
721: }

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

727:    Not Collective

729:    Input Parameter:
730: .  pep - the eigensolver context

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

735:    Level: developer

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

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

752:    Logically Collective

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

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

764:    Level: advanced

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

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

790:    Logically Collective

792:    Input Parameters:
793: +  pep  - eigensolver context obtained from PEPCreate()
794: -  conv - the type of convergence test

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

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

809:    Level: intermediate

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

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

837:    Not Collective

839:    Input Parameters:
840: .  pep   - eigensolver context obtained from PEPCreate()

842:    Output Parameters:
843: .  conv  - the type of convergence test

845:    Level: intermediate

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

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

862:    Logically Collective

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

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

876:    Level: advanced

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

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

900:    Logically Collective

902:    Input Parameters:
903: +  pep  - eigensolver context obtained from PEPCreate()
904: -  stop - the type of stopping test

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

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

915:    Level: advanced

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

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

941:    Not Collective

943:    Input Parameters:
944: .  pep   - eigensolver context obtained from PEPCreate()

946:    Output Parameters:
947: .  stop  - the type of stopping test

949:    Level: advanced

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

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

965:    Collective

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

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

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

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

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

1003:    Level: intermediate

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

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

1060:    Not Collective

1062:    Input Parameter:
1063: .  pep - the eigensolver context

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

1073:    Level: intermediate

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

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

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

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

1099:    Logically Collective

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

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

1108:    Level: intermediate

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

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

1124:    Not Collective

1126:    Input Parameter:
1127: .  pep - the eigensolver context

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

1132:    Level: intermediate

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

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

1149:    Logically Collective

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

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

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

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

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

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

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

1192:    Level: intermediate

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

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

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

1242:    Not Collective

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

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

1254:    Level: intermediate

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

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

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

1277:    Logically Collective

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

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

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

1295:    Level: advanced

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

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

1319:    Logically Collective

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

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

1329:    Level: advanced

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

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

1353:    Not Collective

1355:    Input Parameters:
1356: .  pep - the polynomial eigensolver context

1358:    Output Parameters:
1359: .  prefix - pointer to the prefix string used is returned

1361:    Note:
1362:    On the Fortran side, the user should pass in a string 'prefix' of
1363:    sufficient length to hold the prefix.

1365:    Level: advanced

1367: .seealso: PEPSetOptionsPrefix(), PEPAppendOptionsPrefix()
1368: @*/
1369: PetscErrorCode PEPGetOptionsPrefix(PEP pep,const char *prefix[])
1370: {
1371:   PetscFunctionBegin;
1373:   PetscAssertPointer(prefix,2);
1374:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)pep,prefix));
1375:   PetscFunctionReturn(PETSC_SUCCESS);
1376: }