Actual source code: pepopts.c

slepc-3.16.1 2021-11-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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 on pep

 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;
 46:   PetscErrorCode       ierr;

 49:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->options,((PetscObject)pep)->prefix,opt,&viewer,&format,&flg);
 50:   if (!flg) return(0);

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

 61:   (*cfunc)(viewer,format,ctx,&vf);
 62:   PetscObjectDereference((PetscObject)viewer);
 63:   PEPMonitorSet(pep,mfunc,vf,(PetscErrorCode(*)(void **))dfunc);
 64:   if (trackall) {
 65:     PEPSetTrackAll(pep,PETSC_TRUE);
 66:   }
 67:   return(0);
 68: }

 70: /*@
 71:    PEPSetFromOptions - Sets PEP options from the options database.
 72:    This routine must be called before PEPSetUp() if the user is to be
 73:    allowed to set the solver type.

 75:    Collective on pep

 77:    Input Parameters:
 78: .  pep - the polynomial eigensolver context

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

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

 99:   PEPRegisterAll();
100:   PetscObjectOptionsBegin((PetscObject)pep);
101:     PetscOptionsFList("-pep_type","Polynomial eigensolver method","PEPSetType",PEPList,(char*)(((PetscObject)pep)->type_name?((PetscObject)pep)->type_name:PEPTOAR),type,sizeof(type),&flg);
102:     if (flg) {
103:       PEPSetType(pep,type);
104:     } else if (!((PetscObject)pep)->type_name) {
105:       PEPSetType(pep,PEPTOAR);
106:     }

108:     PetscOptionsBoolGroupBegin("-pep_general","General polynomial eigenvalue problem","PEPSetProblemType",&flg);
109:     if (flg) { PEPSetProblemType(pep,PEP_GENERAL); }
110:     PetscOptionsBoolGroup("-pep_hermitian","Hermitian polynomial eigenvalue problem","PEPSetProblemType",&flg);
111:     if (flg) { PEPSetProblemType(pep,PEP_HERMITIAN); }
112:     PetscOptionsBoolGroup("-pep_hyperbolic","Hyperbolic polynomial eigenvalue problem","PEPSetProblemType",&flg);
113:     if (flg) { PEPSetProblemType(pep,PEP_HYPERBOLIC); }
114:     PetscOptionsBoolGroupEnd("-pep_gyroscopic","Gyroscopic polynomial eigenvalue problem","PEPSetProblemType",&flg);
115:     if (flg) { PEPSetProblemType(pep,PEP_GYROSCOPIC); }

117:     scale = pep->scale;
118:     PetscOptionsEnum("-pep_scale","Scaling strategy","PEPSetScale",PEPScaleTypes,(PetscEnum)scale,(PetscEnum*)&scale,&flg1);
119:     r = pep->sfactor;
120:     PetscOptionsReal("-pep_scale_factor","Scale factor","PEPSetScale",pep->sfactor,&r,&flg2);
121:     if (!flg2 && r==1.0) r = PETSC_DEFAULT;
122:     j = pep->sits;
123:     PetscOptionsInt("-pep_scale_its","Number of iterations in diagonal scaling","PEPSetScale",pep->sits,&j,&flg3);
124:     t = pep->slambda;
125:     PetscOptionsReal("-pep_scale_lambda","Estimate of eigenvalue (modulus) for diagonal scaling","PEPSetScale",pep->slambda,&t,&flg4);
126:     if (flg1 || flg2 || flg3 || flg4) { PEPSetScale(pep,scale,r,NULL,NULL,j,t); }

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

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

142:     i = pep->max_it;
143:     PetscOptionsInt("-pep_max_it","Maximum number of iterations","PEPSetTolerances",pep->max_it,&i,&flg1);
144:     r = pep->tol;
145:     PetscOptionsReal("-pep_tol","Tolerance","PEPSetTolerances",SlepcDefaultTol(pep->tol),&r,&flg2);
146:     if (flg1 || flg2) { PEPSetTolerances(pep,r,i); }

148:     PetscOptionsBoolGroupBegin("-pep_conv_rel","Relative error convergence test","PEPSetConvergenceTest",&flg);
149:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_REL); }
150:     PetscOptionsBoolGroup("-pep_conv_norm","Convergence test relative to the matrix norms","PEPSetConvergenceTest",&flg);
151:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_NORM); }
152:     PetscOptionsBoolGroup("-pep_conv_abs","Absolute error convergence test","PEPSetConvergenceTest",&flg);
153:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_ABS); }
154:     PetscOptionsBoolGroupEnd("-pep_conv_user","User-defined convergence test","PEPSetConvergenceTest",&flg);
155:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_USER); }

157:     PetscOptionsBoolGroupBegin("-pep_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","PEPSetStoppingTest",&flg);
158:     if (flg) { PEPSetStoppingTest(pep,PEP_STOP_BASIC); }
159:     PetscOptionsBoolGroupEnd("-pep_stop_user","User-defined stopping test","PEPSetStoppingTest",&flg);
160:     if (flg) { PEPSetStoppingTest(pep,PEP_STOP_USER); }

162:     i = pep->nev;
163:     PetscOptionsInt("-pep_nev","Number of eigenvalues to compute","PEPSetDimensions",pep->nev,&i,&flg1);
164:     j = pep->ncv;
165:     PetscOptionsInt("-pep_ncv","Number of basis vectors","PEPSetDimensions",pep->ncv,&j,&flg2);
166:     k = pep->mpd;
167:     PetscOptionsInt("-pep_mpd","Maximum dimension of projected problem","PEPSetDimensions",pep->mpd,&k,&flg3);
168:     if (flg1 || flg2 || flg3) { PEPSetDimensions(pep,i,j,k); }

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

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

193:     PetscOptionsScalar("-pep_target","Value of the target","PEPSetTarget",pep->target,&s,&flg);
194:     if (flg) {
195:       if (pep->which!=PEP_TARGET_REAL && pep->which!=PEP_TARGET_IMAGINARY) {
196:         PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE);
197:       }
198:       PEPSetTarget(pep,s);
199:     }

201:     k = 2;
202:     PetscOptionsRealArray("-pep_interval","Computational interval (two real values separated with a comma without spaces)","PEPSetInterval",array,&k,&flg);
203:     if (flg) {
204:       if (k<2) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_SIZ,"Must pass two values in -pep_interval (comma-separated without spaces)");
205:       PEPSetWhichEigenpairs(pep,PEP_ALL);
206:       PEPSetInterval(pep,array[0],array[1]);
207:     }

209:     /* -----------------------------------------------------------------------*/
210:     /*
211:       Cancels all monitors hardwired into code before call to PEPSetFromOptions()
212:     */
213:     PetscOptionsBool("-pep_monitor_cancel","Remove any hardwired monitor routines","PEPMonitorCancel",PETSC_FALSE,&flg,&set);
214:     if (set && flg) { PEPMonitorCancel(pep); }
215:     PEPMonitorSetFromOptions(pep,"-pep_monitor","first_approximation",NULL,PETSC_FALSE);
216:     PEPMonitorSetFromOptions(pep,"-pep_monitor_all","all_approximations",NULL,PETSC_TRUE);
217:     PEPMonitorSetFromOptions(pep,"-pep_monitor_conv","convergence_history",NULL,PETSC_FALSE);

219:     /* -----------------------------------------------------------------------*/
220:     PetscOptionsName("-pep_view","Print detailed information on solver used","PEPView",NULL);
221:     PetscOptionsName("-pep_view_vectors","View computed eigenvectors","PEPVectorsView",NULL);
222:     PetscOptionsName("-pep_view_values","View computed eigenvalues","PEPValuesView",NULL);
223:     PetscOptionsName("-pep_converged_reason","Print reason for convergence, and number of iterations","PEPConvergedReasonView",NULL);
224:     PetscOptionsName("-pep_error_absolute","Print absolute errors of each eigenpair","PEPErrorView",NULL);
225:     PetscOptionsName("-pep_error_relative","Print relative errors of each eigenpair","PEPErrorView",NULL);
226:     PetscOptionsName("-pep_error_backward","Print backward errors of each eigenpair","PEPErrorView",NULL);

228:     if (pep->ops->setfromoptions) {
229:       (*pep->ops->setfromoptions)(PetscOptionsObject,pep);
230:     }
231:     PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)pep);
232:   PetscOptionsEnd();

234:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
235:   BVSetFromOptions(pep->V);
236:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
237:   RGSetFromOptions(pep->rg);
238:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
239:   DSSetFromOptions(pep->ds);
240:   if (!pep->st) { PEPGetST(pep,&pep->st); }
241:   PEPSetDefaultST(pep);
242:   STSetFromOptions(pep->st);
243:   if (!pep->refineksp) { PEPRefineGetKSP(pep,&pep->refineksp); }
244:   KSPSetFromOptions(pep->refineksp);
245:   return(0);
246: }

248: /*@C
249:    PEPGetTolerances - Gets the tolerance and maximum iteration count used
250:    by the PEP convergence tests.

252:    Not Collective

254:    Input Parameter:
255: .  pep - the polynomial eigensolver context

257:    Output Parameters:
258: +  tol - the convergence tolerance
259: -  maxits - maximum number of iterations

261:    Notes:
262:    The user can specify NULL for any parameter that is not needed.

264:    Level: intermediate

266: .seealso: PEPSetTolerances()
267: @*/
268: PetscErrorCode PEPGetTolerances(PEP pep,PetscReal *tol,PetscInt *maxits)
269: {
272:   if (tol)    *tol    = pep->tol;
273:   if (maxits) *maxits = pep->max_it;
274:   return(0);
275: }

277: /*@
278:    PEPSetTolerances - Sets the tolerance and maximum iteration count used
279:    by the PEP convergence tests.

281:    Logically Collective on pep

283:    Input Parameters:
284: +  pep - the polynomial eigensolver context
285: .  tol - the convergence tolerance
286: -  maxits - maximum number of iterations to use

288:    Options Database Keys:
289: +  -pep_tol <tol> - Sets the convergence tolerance
290: -  -pep_max_it <maxits> - Sets the maximum number of iterations allowed

292:    Notes:
293:    Use PETSC_DEFAULT for either argument to assign a reasonably good value.

295:    Level: intermediate

297: .seealso: PEPGetTolerances()
298: @*/
299: PetscErrorCode PEPSetTolerances(PEP pep,PetscReal tol,PetscInt maxits)
300: {
305:   if (tol == PETSC_DEFAULT) {
306:     pep->tol   = PETSC_DEFAULT;
307:     pep->state = PEP_STATE_INITIAL;
308:   } else {
309:     if (tol <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
310:     pep->tol = tol;
311:   }
312:   if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
313:     pep->max_it = PETSC_DEFAULT;
314:     pep->state  = PEP_STATE_INITIAL;
315:   } else {
316:     if (maxits <= 0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
317:     pep->max_it = maxits;
318:   }
319:   return(0);
320: }

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

326:    Not Collective

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

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

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

339:    Level: intermediate

341: .seealso: PEPSetDimensions()
342: @*/
343: PetscErrorCode PEPGetDimensions(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
344: {
347:   if (nev) *nev = pep->nev;
348:   if (ncv) *ncv = pep->ncv;
349:   if (mpd) *mpd = pep->mpd;
350:   return(0);
351: }

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

357:    Logically Collective on pep

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

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

370:    Notes:
371:    Use PETSC_DEFAULT for ncv and mpd to assign a reasonably good value, which is
372:    dependent on the solution method.

374:    The parameters ncv and mpd are intimately related, so that the user is advised
375:    to set one of them at most. Normal usage is that
376:    (a) in cases where nev is small, the user sets ncv (a reasonable default is 2*nev); and
377:    (b) in cases where nev is large, the user sets mpd.

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

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

386:    Level: intermediate

388: .seealso: PEPGetDimensions(), PEPSetInterval(), PEPSTOARSetDimensions()
389: @*/
390: PetscErrorCode PEPSetDimensions(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)
391: {
397:   if (nev<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
398:   pep->nev = nev;
399:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
400:     pep->ncv = PETSC_DEFAULT;
401:   } else {
402:     if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
403:     pep->ncv = ncv;
404:   }
405:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
406:     pep->mpd = PETSC_DEFAULT;
407:   } else {
408:     if (mpd<1) SETERRQ(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:   return(0);
413: }

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

419:    Logically Collective on pep

421:    Input Parameters:
422: +  pep   - eigensolver context obtained from PEPCreate()
423: -  which - the portion of the spectrum to be sought

425:    Possible values:
426:    The parameter 'which' can have one of these values

428: +     PEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
429: .     PEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
430: .     PEP_LARGEST_REAL - largest real parts
431: .     PEP_SMALLEST_REAL - smallest real parts
432: .     PEP_LARGEST_IMAGINARY - largest imaginary parts
433: .     PEP_SMALLEST_IMAGINARY - smallest imaginary parts
434: .     PEP_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
435: .     PEP_TARGET_REAL - eigenvalues with real part closest to target
436: .     PEP_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
437: .     PEP_ALL - all eigenvalues contained in a given interval or region
438: -     PEP_WHICH_USER - user defined ordering set with PEPSetEigenvalueComparison()

440:    Options Database Keys:
441: +   -pep_largest_magnitude - Sets largest eigenvalues in magnitude
442: .   -pep_smallest_magnitude - Sets smallest eigenvalues in magnitude
443: .   -pep_largest_real - Sets largest real parts
444: .   -pep_smallest_real - Sets smallest real parts
445: .   -pep_largest_imaginary - Sets largest imaginary parts
446: .   -pep_smallest_imaginary - Sets smallest imaginary parts
447: .   -pep_target_magnitude - Sets eigenvalues closest to target
448: .   -pep_target_real - Sets real parts closest to target
449: .   -pep_target_imaginary - Sets imaginary parts closest to target
450: -   -pep_all - Sets all eigenvalues in an interval or region

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

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

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

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

469:    Level: intermediate

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

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

512:     Not Collective

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

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

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

523:     Level: intermediate

525: .seealso: PEPSetWhichEigenpairs(), PEPWhich
526: @*/
527: PetscErrorCode PEPGetWhichEigenpairs(PEP pep,PEPWhich *which)
528: {
532:   *which = pep->which;
533:   return(0);
534: }

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

540:    Logically Collective on pep

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

547:    Calling Sequence of func:
548: $   func(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)

550: +   ar     - real part of the 1st eigenvalue
551: .   ai     - imaginary part of the 1st eigenvalue
552: .   br     - real part of the 2nd eigenvalue
553: .   bi     - imaginary part of the 2nd eigenvalue
554: .   res    - result of comparison
555: -   ctx    - optional context, as set by PEPSetEigenvalueComparison()

557:    Note:
558:    The returning parameter 'res' can be
559: +  negative - if the 1st eigenvalue is preferred to the 2st one
560: .  zero     - if both eigenvalues are equally preferred
561: -  positive - if the 2st eigenvalue is preferred to the 1st one

563:    Level: advanced

565: .seealso: PEPSetWhichEigenpairs(), PEPWhich
566: @*/
567: PetscErrorCode PEPSetEigenvalueComparison(PEP pep,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)
568: {
571:   pep->sc->comparison    = func;
572:   pep->sc->comparisonctx = ctx;
573:   pep->which             = PEP_WHICH_USER;
574:   return(0);
575: }

577: /*@
578:    PEPSetProblemType - Specifies the type of the polynomial eigenvalue problem.

580:    Logically Collective on pep

582:    Input Parameters:
583: +  pep  - the polynomial eigensolver context
584: -  type - a known type of polynomial eigenvalue problem

586:    Options Database Keys:
587: +  -pep_general - general problem with no particular structure
588: .  -pep_hermitian - problem whose coefficient matrices are Hermitian
589: .  -pep_hyperbolic - Hermitian problem that satisfies the definition of hyperbolic
590: -  -pep_gyroscopic - problem with Hamiltonian structure

592:    Notes:
593:    Allowed values for the problem type are: general (PEP_GENERAL), Hermitian
594:    (PEP_HERMITIAN), hyperbolic (PEP_HYPERBOLIC), and gyroscopic (PEP_GYROSCOPIC).

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

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

605:    Level: intermediate

607: .seealso: PEPSetOperators(), PEPSetType(), PEPGetProblemType(), PEPProblemType
608: @*/
609: PetscErrorCode PEPSetProblemType(PEP pep,PEPProblemType type)
610: {
614:   if (type!=PEP_GENERAL && type!=PEP_HERMITIAN && type!=PEP_HYPERBOLIC && type!=PEP_GYROSCOPIC) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
615:   if (type != pep->problem_type) {
616:     pep->problem_type = type;
617:     pep->state = PEP_STATE_INITIAL;
618:   }
619:   return(0);
620: }

622: /*@
623:    PEPGetProblemType - Gets the problem type from the PEP object.

625:    Not Collective

627:    Input Parameter:
628: .  pep - the polynomial eigensolver context

630:    Output Parameter:
631: .  type - the problem type

633:    Level: intermediate

635: .seealso: PEPSetProblemType(), PEPProblemType
636: @*/
637: PetscErrorCode PEPGetProblemType(PEP pep,PEPProblemType *type)
638: {
642:   *type = pep->problem_type;
643:   return(0);
644: }

646: /*@
647:    PEPSetBasis - Specifies the type of polynomial basis used to describe the
648:    polynomial eigenvalue problem.

650:    Logically Collective on pep

652:    Input Parameters:
653: +  pep   - the polynomial eigensolver context
654: -  basis - the type of polynomial basis

656:    Options Database Key:
657: .  -pep_basis <basis> - Select the basis type

659:    Notes:
660:    By default, the coefficient matrices passed via PEPSetOperators() are
661:    expressed in the monomial basis, i.e.
662:    P(lambda) = A_0 + lambda*A_1 + lambda^2*A_2 + ... + lambda^d*A_d.
663:    Other polynomial bases may have better numerical behaviour, but the user
664:    must then pass the coefficient matrices accordingly.

666:    Level: intermediate

668: .seealso: PEPSetOperators(), PEPGetBasis(), PEPBasis
669: @*/
670: PetscErrorCode PEPSetBasis(PEP pep,PEPBasis basis)
671: {
675:   pep->basis = basis;
676:   return(0);
677: }

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

682:    Not Collective

684:    Input Parameter:
685: .  pep - the polynomial eigensolver context

687:    Output Parameter:
688: .  basis - the polynomial basis

690:    Level: intermediate

692: .seealso: PEPSetBasis(), PEPBasis
693: @*/
694: PetscErrorCode PEPGetBasis(PEP pep,PEPBasis *basis)
695: {
699:   *basis = pep->basis;
700:   return(0);
701: }

703: /*@
704:    PEPSetTrackAll - Specifies if the solver must compute the residual of all
705:    approximate eigenpairs or not.

707:    Logically Collective on pep

709:    Input Parameters:
710: +  pep      - the eigensolver context
711: -  trackall - whether compute all residuals or not

713:    Notes:
714:    If the user sets trackall=PETSC_TRUE then the solver explicitly computes
715:    the residual for each eigenpair approximation. Computing the residual is
716:    usually an expensive operation and solvers commonly compute the associated
717:    residual to the first unconverged eigenpair.

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

721:    Level: developer

723: .seealso: PEPGetTrackAll()
724: @*/
725: PetscErrorCode PEPSetTrackAll(PEP pep,PetscBool trackall)
726: {
730:   pep->trackall = trackall;
731:   return(0);
732: }

734: /*@
735:    PEPGetTrackAll - Returns the flag indicating whether all residual norms must
736:    be computed or not.

738:    Not Collective

740:    Input Parameter:
741: .  pep - the eigensolver context

743:    Output Parameter:
744: .  trackall - the returned flag

746:    Level: developer

748: .seealso: PEPSetTrackAll()
749: @*/
750: PetscErrorCode PEPGetTrackAll(PEP pep,PetscBool *trackall)
751: {
755:   *trackall = pep->trackall;
756:   return(0);
757: }

759: /*@C
760:    PEPSetConvergenceTestFunction - Sets a function to compute the error estimate
761:    used in the convergence test.

763:    Logically Collective on pep

765:    Input Parameters:
766: +  pep     - eigensolver context obtained from PEPCreate()
767: .  func    - a pointer to the convergence test function
768: .  ctx     - context for private data for the convergence routine (may be null)
769: -  destroy - a routine for destroying the context (may be null)

771:    Calling Sequence of func:
772: $   func(PEP pep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)

774: +   pep    - eigensolver context obtained from PEPCreate()
775: .   eigr   - real part of the eigenvalue
776: .   eigi   - imaginary part of the eigenvalue
777: .   res    - residual norm associated to the eigenpair
778: .   errest - (output) computed error estimate
779: -   ctx    - optional context, as set by PEPSetConvergenceTestFunction()

781:    Note:
782:    If the error estimate returned by the convergence test function is less than
783:    the tolerance, then the eigenvalue is accepted as converged.

785:    Level: advanced

787: .seealso: PEPSetConvergenceTest(), PEPSetTolerances()
788: @*/
789: PetscErrorCode PEPSetConvergenceTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
790: {

795:   if (pep->convergeddestroy) {
796:     (*pep->convergeddestroy)(pep->convergedctx);
797:   }
798:   pep->convergeduser    = func;
799:   pep->convergeddestroy = destroy;
800:   pep->convergedctx     = ctx;
801:   if (func == PEPConvergedRelative) pep->conv = PEP_CONV_REL;
802:   else if (func == PEPConvergedNorm) pep->conv = PEP_CONV_NORM;
803:   else if (func == PEPConvergedAbsolute) pep->conv = PEP_CONV_ABS;
804:   else {
805:     pep->conv      = PEP_CONV_USER;
806:     pep->converged = pep->convergeduser;
807:   }
808:   return(0);
809: }

811: /*@
812:    PEPSetConvergenceTest - Specifies how to compute the error estimate
813:    used in the convergence test.

815:    Logically Collective on pep

817:    Input Parameters:
818: +  pep  - eigensolver context obtained from PEPCreate()
819: -  conv - the type of convergence test

821:    Options Database Keys:
822: +  -pep_conv_abs    - Sets the absolute convergence test
823: .  -pep_conv_rel    - Sets the convergence test relative to the eigenvalue
824: .  -pep_conv_norm   - Sets the convergence test relative to the matrix norms
825: -  -pep_conv_user   - Selects the user-defined convergence test

827:    Note:
828:    The parameter 'conv' can have one of these values
829: +     PEP_CONV_ABS    - absolute error ||r||
830: .     PEP_CONV_REL    - error relative to the eigenvalue l, ||r||/|l|
831: .     PEP_CONV_NORM   - error relative matrix norms, ||r||/sum_i(l^i*||A_i||)
832: -     PEP_CONV_USER   - function set by PEPSetConvergenceTestFunction()

834:    Level: intermediate

836: .seealso: PEPGetConvergenceTest(), PEPSetConvergenceTestFunction(), PEPSetStoppingTest(), PEPConv
837: @*/
838: PetscErrorCode PEPSetConvergenceTest(PEP pep,PEPConv conv)
839: {
843:   switch (conv) {
844:     case PEP_CONV_ABS:  pep->converged = PEPConvergedAbsolute; break;
845:     case PEP_CONV_REL:  pep->converged = PEPConvergedRelative; break;
846:     case PEP_CONV_NORM: pep->converged = PEPConvergedNorm; break;
847:     case PEP_CONV_USER:
848:       if (!pep->convergeduser) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetConvergenceTestFunction() first");
849:       pep->converged = pep->convergeduser;
850:       break;
851:     default:
852:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
853:   }
854:   pep->conv = conv;
855:   return(0);
856: }

858: /*@
859:    PEPGetConvergenceTest - Gets the method used to compute the error estimate
860:    used in the convergence test.

862:    Not Collective

864:    Input Parameters:
865: .  pep   - eigensolver context obtained from PEPCreate()

867:    Output Parameters:
868: .  conv  - the type of convergence test

870:    Level: intermediate

872: .seealso: PEPSetConvergenceTest(), PEPConv
873: @*/
874: PetscErrorCode PEPGetConvergenceTest(PEP pep,PEPConv *conv)
875: {
879:   *conv = pep->conv;
880:   return(0);
881: }

883: /*@C
884:    PEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
885:    iteration of the eigensolver.

887:    Logically Collective on pep

889:    Input Parameters:
890: +  pep     - eigensolver context obtained from PEPCreate()
891: .  func    - pointer to the stopping test function
892: .  ctx     - context for private data for the stopping routine (may be null)
893: -  destroy - a routine for destroying the context (may be null)

895:    Calling Sequence of func:
896: $   func(PEP pep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,PEPConvergedReason *reason,void *ctx)

898: +   pep    - eigensolver context obtained from PEPCreate()
899: .   its    - current number of iterations
900: .   max_it - maximum number of iterations
901: .   nconv  - number of currently converged eigenpairs
902: .   nev    - number of requested eigenpairs
903: .   reason - (output) result of the stopping test
904: -   ctx    - optional context, as set by PEPSetStoppingTestFunction()

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

912:    Level: advanced

914: .seealso: PEPSetStoppingTest(), PEPStoppingBasic()
915: @*/
916: PetscErrorCode PEPSetStoppingTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
917: {

922:   if (pep->stoppingdestroy) {
923:     (*pep->stoppingdestroy)(pep->stoppingctx);
924:   }
925:   pep->stoppinguser    = func;
926:   pep->stoppingdestroy = destroy;
927:   pep->stoppingctx     = ctx;
928:   if (func == PEPStoppingBasic) pep->stop = PEP_STOP_BASIC;
929:   else {
930:     pep->stop     = PEP_STOP_USER;
931:     pep->stopping = pep->stoppinguser;
932:   }
933:   return(0);
934: }

936: /*@
937:    PEPSetStoppingTest - Specifies how to decide the termination of the outer
938:    loop of the eigensolver.

940:    Logically Collective on pep

942:    Input Parameters:
943: +  pep  - eigensolver context obtained from PEPCreate()
944: -  stop - the type of stopping test

946:    Options Database Keys:
947: +  -pep_stop_basic - Sets the default stopping test
948: -  -pep_stop_user  - Selects the user-defined stopping test

950:    Note:
951:    The parameter 'stop' can have one of these values
952: +     PEP_STOP_BASIC - default stopping test
953: -     PEP_STOP_USER  - function set by PEPSetStoppingTestFunction()

955:    Level: advanced

957: .seealso: PEPGetStoppingTest(), PEPSetStoppingTestFunction(), PEPSetConvergenceTest(), PEPStop
958: @*/
959: PetscErrorCode PEPSetStoppingTest(PEP pep,PEPStop stop)
960: {
964:   switch (stop) {
965:     case PEP_STOP_BASIC: pep->stopping = PEPStoppingBasic; break;
966:     case PEP_STOP_USER:
967:       if (!pep->stoppinguser) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetStoppingTestFunction() first");
968:       pep->stopping = pep->stoppinguser;
969:       break;
970:     default:
971:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
972:   }
973:   pep->stop = stop;
974:   return(0);
975: }

977: /*@
978:    PEPGetStoppingTest - Gets the method used to decide the termination of the outer
979:    loop of the eigensolver.

981:    Not Collective

983:    Input Parameters:
984: .  pep   - eigensolver context obtained from PEPCreate()

986:    Output Parameters:
987: .  stop  - the type of stopping test

989:    Level: advanced

991: .seealso: PEPSetStoppingTest(), PEPStop
992: @*/
993: PetscErrorCode PEPGetStoppingTest(PEP pep,PEPStop *stop)
994: {
998:   *stop = pep->stop;
999:   return(0);
1000: }

1002: /*@
1003:    PEPSetScale - Specifies the scaling strategy to be used.

1005:    Logically Collective on pep

1007:    Input Parameters:
1008: +  pep    - the eigensolver context
1009: .  scale  - scaling strategy
1010: .  alpha  - the scaling factor used in the scalar strategy
1011: .  Dl     - the left diagonal matrix of the diagonal scaling algorithm
1012: .  Dr     - the right diagonal matrix of the diagonal scaling algorithm
1013: .  its    - number of iterations of the diagonal scaling algorithm
1014: -  lambda - approximation to wanted eigenvalues (modulus)

1016:    Options Database Keys:
1017: +  -pep_scale <type> - scaling type, one of <none,scalar,diagonal,both>
1018: .  -pep_scale_factor <alpha> - the scaling factor
1019: .  -pep_scale_its <its> - number of iterations
1020: -  -pep_scale_lambda <lambda> - approximation to eigenvalues

1022:    Notes:
1023:    There are two non-exclusive scaling strategies: scalar and diagonal.

1025:    In the scalar strategy, scaling is applied to the eigenvalue, that is,
1026:    mu = lambda/alpha is the new eigenvalue and all matrices are scaled
1027:    accordingly. After solving the scaled problem, the original lambda is
1028:    recovered. Parameter 'alpha' must be positive. Use PETSC_DEFAULT to let
1029:    the solver compute a reasonable scaling factor.

1031:    In the diagonal strategy, the solver works implicitly with matrix Dl*A*Dr,
1032:    where Dl and Dr are appropriate diagonal matrices. This improves the accuracy
1033:    of the computed results in some cases. The user may provide the Dr and Dl
1034:    matrices represented as Vec objects storing diagonal elements. If not
1035:    provided, these matrices are computed internally. This option requires
1036:    that the polynomial coefficient matrices are of MATAIJ type.
1037:    The parameter 'its' is the number of iterations performed by the method.
1038:    Parameter 'lambda' must be positive. Use PETSC_DEFAULT or set lambda = 1.0 if
1039:    no information about eigenvalues is available.

1041:    Level: intermediate

1043: .seealso: PEPGetScale()
1044: @*/
1045: PetscErrorCode PEPSetScale(PEP pep,PEPScale scale,PetscReal alpha,Vec Dl,Vec Dr,PetscInt its,PetscReal lambda)
1046: {

1052:   pep->scale = scale;
1053:   if (scale==PEP_SCALE_SCALAR || scale==PEP_SCALE_BOTH) {
1055:     if (alpha == PETSC_DEFAULT || alpha == PETSC_DECIDE) {
1056:       pep->sfactor = 0.0;
1057:       pep->sfactor_set = PETSC_FALSE;
1058:     } else {
1059:       if (alpha<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of alpha. Must be > 0");
1060:       pep->sfactor = alpha;
1061:       pep->sfactor_set = PETSC_TRUE;
1062:     }
1063:   }
1064:   if (scale==PEP_SCALE_DIAGONAL || scale==PEP_SCALE_BOTH) {
1065:     if (Dl) {
1068:       PetscObjectReference((PetscObject)Dl);
1069:       VecDestroy(&pep->Dl);
1070:       pep->Dl = Dl;
1071:     }
1072:     if (Dr) {
1075:       PetscObjectReference((PetscObject)Dr);
1076:       VecDestroy(&pep->Dr);
1077:       pep->Dr = Dr;
1078:     }
1081:     if (its==PETSC_DECIDE || its==PETSC_DEFAULT) pep->sits = 5;
1082:     else pep->sits = its;
1083:     if (lambda==PETSC_DECIDE || lambda==PETSC_DEFAULT) pep->slambda = 1.0;
1084:     else if (lambda<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of lambda. Must be > 0");
1085:     else pep->slambda = lambda;
1086:   }
1087:   pep->state = PEP_STATE_INITIAL;
1088:   return(0);
1089: }

1091: /*@C
1092:    PEPGetScale - Gets the scaling strategy used by the PEP object, and the
1093:    associated parameters.

1095:    Not Collectiv, but vectors are shared by all processors that share the PEP

1097:    Input Parameter:
1098: .  pep - the eigensolver context

1100:    Output Parameters:
1101: +  scale  - scaling strategy
1102: .  alpha  - the scaling factor used in the scalar strategy
1103: .  Dl     - the left diagonal matrix of the diagonal scaling algorithm
1104: .  Dr     - the right diagonal matrix of the diagonal scaling algorithm
1105: .  its    - number of iterations of the diagonal scaling algorithm
1106: -  lambda - approximation to wanted eigenvalues (modulus)

1108:    Level: intermediate

1110:    Note:
1111:    The user can specify NULL for any parameter that is not needed.

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

1116: .seealso: PEPSetScale(), PEPSetUp()
1117: @*/
1118: PetscErrorCode PEPGetScale(PEP pep,PEPScale *scale,PetscReal *alpha,Vec *Dl,Vec *Dr,PetscInt *its,PetscReal *lambda)
1119: {
1122:   if (scale)  *scale  = pep->scale;
1123:   if (alpha)  *alpha  = pep->sfactor;
1124:   if (Dl)     *Dl     = pep->Dl;
1125:   if (Dr)     *Dr     = pep->Dr;
1126:   if (its)    *its    = pep->sits;
1127:   if (lambda) *lambda = pep->slambda;
1128:   return(0);
1129: }

1131: /*@
1132:    PEPSetExtract - Specifies the extraction strategy to be used.

1134:    Logically Collective on pep

1136:    Input Parameters:
1137: +  pep     - the eigensolver context
1138: -  extract - extraction strategy

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

1143:    Level: intermediate

1145: .seealso: PEPGetExtract()
1146: @*/
1147: PetscErrorCode PEPSetExtract(PEP pep,PEPExtract extract)
1148: {
1152:   pep->extract = extract;
1153:   return(0);
1154: }

1156: /*@
1157:    PEPGetExtract - Gets the extraction strategy used by the PEP object.

1159:    Not Collective

1161:    Input Parameter:
1162: .  pep - the eigensolver context

1164:    Output Parameter:
1165: .  extract - extraction strategy

1167:    Level: intermediate

1169: .seealso: PEPSetExtract(), PEPExtract
1170: @*/
1171: PetscErrorCode PEPGetExtract(PEP pep,PEPExtract *extract)
1172: {
1176:   *extract = pep->extract;
1177:   return(0);
1178: }

1180: /*@
1181:    PEPSetRefine - Specifies the refinement type (and options) to be used
1182:    after the solve.

1184:    Logically Collective on pep

1186:    Input Parameters:
1187: +  pep    - the polynomial eigensolver context
1188: .  refine - refinement type
1189: .  npart  - number of partitions of the communicator
1190: .  tol    - the convergence tolerance
1191: .  its    - maximum number of refinement iterations
1192: -  scheme - which scheme to be used for solving the involved linear systems

1194:    Options Database Keys:
1195: +  -pep_refine <type> - refinement type, one of <none,simple,multiple>
1196: .  -pep_refine_partitions <n> - the number of partitions
1197: .  -pep_refine_tol <tol> - the tolerance
1198: .  -pep_refine_its <its> - number of iterations
1199: -  -pep_refine_scheme - to set the scheme for the linear solves

1201:    Notes:
1202:    By default, iterative refinement is disabled, since it may be very
1203:    costly. There are two possible refinement strategies: simple and multiple.
1204:    The simple approach performs iterative refinement on each of the
1205:    converged eigenpairs individually, whereas the multiple strategy works
1206:    with the invariant pair as a whole, refining all eigenpairs simultaneously.
1207:    The latter may be required for the case of multiple eigenvalues.

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

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

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

1224:    Level: intermediate

1226: .seealso: PEPGetRefine()
1227: @*/
1228: PetscErrorCode PEPSetRefine(PEP pep,PEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,PEPRefineScheme scheme)
1229: {
1231:   PetscMPIInt    size;

1240:   pep->refine = refine;
1241:   if (refine) {  /* process parameters only if not REFINE_NONE */
1242:     if (npart!=pep->npart) {
1243:       PetscSubcommDestroy(&pep->refinesubc);
1244:       KSPDestroy(&pep->refineksp);
1245:     }
1246:     if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
1247:       pep->npart = 1;
1248:     } else {
1249:       MPI_Comm_size(PetscObjectComm((PetscObject)pep),&size);
1250:       if (npart<1 || npart>size) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
1251:       pep->npart = npart;
1252:     }
1253:     if (tol == PETSC_DEFAULT || tol == PETSC_DECIDE) {
1254:       pep->rtol = PETSC_DEFAULT;
1255:     } else {
1256:       if (tol<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1257:       pep->rtol = tol;
1258:     }
1259:     if (its==PETSC_DECIDE || its==PETSC_DEFAULT) {
1260:       pep->rits = PETSC_DEFAULT;
1261:     } else {
1262:       if (its<0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be >= 0");
1263:       pep->rits = its;
1264:     }
1265:     pep->scheme = scheme;
1266:   }
1267:   pep->state = PEP_STATE_INITIAL;
1268:   return(0);
1269: }

1271: /*@C
1272:    PEPGetRefine - Gets the refinement strategy used by the PEP object, and the
1273:    associated parameters.

1275:    Not Collective

1277:    Input Parameter:
1278: .  pep - the polynomial eigensolver context

1280:    Output Parameters:
1281: +  refine - refinement type
1282: .  npart  - number of partitions of the communicator
1283: .  tol    - the convergence tolerance
1284: .  its    - maximum number of refinement iterations
1285: -  scheme - the scheme used for solving linear systems

1287:    Level: intermediate

1289:    Note:
1290:    The user can specify NULL for any parameter that is not needed.

1292: .seealso: PEPSetRefine()
1293: @*/
1294: PetscErrorCode PEPGetRefine(PEP pep,PEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,PEPRefineScheme *scheme)
1295: {
1298:   if (refine) *refine = pep->refine;
1299:   if (npart)  *npart  = pep->npart;
1300:   if (tol)    *tol    = pep->rtol;
1301:   if (its)    *its    = pep->rits;
1302:   if (scheme) *scheme = pep->scheme;
1303:   return(0);
1304: }

1306: /*@C
1307:    PEPSetOptionsPrefix - Sets the prefix used for searching for all
1308:    PEP options in the database.

1310:    Logically Collective on pep

1312:    Input Parameters:
1313: +  pep - the polynomial eigensolver context
1314: -  prefix - the prefix string to prepend to all PEP option requests

1316:    Notes:
1317:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1318:    The first character of all runtime options is AUTOMATICALLY the
1319:    hyphen.

1321:    For example, to distinguish between the runtime options for two
1322:    different PEP contexts, one could call
1323: .vb
1324:       PEPSetOptionsPrefix(pep1,"qeig1_")
1325:       PEPSetOptionsPrefix(pep2,"qeig2_")
1326: .ve

1328:    Level: advanced

1330: .seealso: PEPAppendOptionsPrefix(), PEPGetOptionsPrefix()
1331: @*/
1332: PetscErrorCode PEPSetOptionsPrefix(PEP pep,const char *prefix)
1333: {

1338:   if (!pep->st) { PEPGetST(pep,&pep->st); }
1339:   STSetOptionsPrefix(pep->st,prefix);
1340:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
1341:   BVSetOptionsPrefix(pep->V,prefix);
1342:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
1343:   DSSetOptionsPrefix(pep->ds,prefix);
1344:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
1345:   RGSetOptionsPrefix(pep->rg,prefix);
1346:   PetscObjectSetOptionsPrefix((PetscObject)pep,prefix);
1347:   return(0);
1348: }

1350: /*@C
1351:    PEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1352:    PEP options in the database.

1354:    Logically Collective on pep

1356:    Input Parameters:
1357: +  pep - the polynomial eigensolver context
1358: -  prefix - the prefix string to prepend to all PEP option requests

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

1364:    Level: advanced

1366: .seealso: PEPSetOptionsPrefix(), PEPGetOptionsPrefix()
1367: @*/
1368: PetscErrorCode PEPAppendOptionsPrefix(PEP pep,const char *prefix)
1369: {

1374:   if (!pep->st) { PEPGetST(pep,&pep->st); }
1375:   STAppendOptionsPrefix(pep->st,prefix);
1376:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
1377:   BVAppendOptionsPrefix(pep->V,prefix);
1378:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
1379:   DSAppendOptionsPrefix(pep->ds,prefix);
1380:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
1381:   RGAppendOptionsPrefix(pep->rg,prefix);
1382:   PetscObjectAppendOptionsPrefix((PetscObject)pep,prefix);
1383:   return(0);
1384: }

1386: /*@C
1387:    PEPGetOptionsPrefix - Gets the prefix used for searching for all
1388:    PEP options in the database.

1390:    Not Collective

1392:    Input Parameters:
1393: .  pep - the polynomial eigensolver context

1395:    Output Parameters:
1396: .  prefix - pointer to the prefix string used is returned

1398:    Note:
1399:    On the Fortran side, the user should pass in a string 'prefix' of
1400:    sufficient length to hold the prefix.

1402:    Level: advanced

1404: .seealso: PEPSetOptionsPrefix(), PEPAppendOptionsPrefix()
1405: @*/
1406: PetscErrorCode PEPGetOptionsPrefix(PEP pep,const char *prefix[])
1407: {

1413:   PetscObjectGetOptionsPrefix((PetscObject)pep,prefix);
1414:   return(0);
1415: }