Actual source code: pepopts.c

slepc-3.18.1 2022-11-02
Report Typos and Errors
  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 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;

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

 50:   PetscViewerGetType(viewer,&vtype);
 51:   SlepcMonitorMakeKey_Internal(name,vtype,format,key);
 52:   PetscFunctionListFind(PEPMonitorList,key,&mfunc);
 54:   PetscFunctionListFind(PEPMonitorCreateList,key,&cfunc);
 55:   PetscFunctionListFind(PEPMonitorDestroyList,key,&dfunc);
 56:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
 57:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

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

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

 71:    Collective on pep

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

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

 79:    Level: beginner

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

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

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

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

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

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

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

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

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

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

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

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

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

192:     k = 2;
193:     PetscOptionsRealArray("-pep_interval","Computational interval (two real values separated with a comma without spaces)","PEPSetInterval",array,&k,&flg);
194:     if (flg) {
196:       PEPSetWhichEigenpairs(pep,PEP_ALL);
197:       PEPSetInterval(pep,array[0],array[1]);
198:     }

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

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

219:     PetscTryTypeMethod(pep,setfromoptions,PetscOptionsObject);
220:     PetscObjectProcessOptionsHandlers((PetscObject)pep,PetscOptionsObject);
221:   PetscOptionsEnd();

223:   if (!pep->V) PEPGetBV(pep,&pep->V);
224:   BVSetFromOptions(pep->V);
225:   if (!pep->rg) PEPGetRG(pep,&pep->rg);
226:   RGSetFromOptions(pep->rg);
227:   if (!pep->ds) PEPGetDS(pep,&pep->ds);
228:   DSSetFromOptions(pep->ds);
229:   if (!pep->st) PEPGetST(pep,&pep->st);
230:   PEPSetDefaultST(pep);
231:   STSetFromOptions(pep->st);
232:   if (!pep->refineksp) PEPRefineGetKSP(pep,&pep->refineksp);
233:   KSPSetFromOptions(pep->refineksp);
234:   return 0;
235: }

237: /*@C
238:    PEPGetTolerances - Gets the tolerance and maximum iteration count used
239:    by the PEP convergence tests.

241:    Not Collective

243:    Input Parameter:
244: .  pep - the polynomial eigensolver context

246:    Output Parameters:
247: +  tol - the convergence tolerance
248: -  maxits - maximum number of iterations

250:    Notes:
251:    The user can specify NULL for any parameter that is not needed.

253:    Level: intermediate

255: .seealso: PEPSetTolerances()
256: @*/
257: PetscErrorCode PEPGetTolerances(PEP pep,PetscReal *tol,PetscInt *maxits)
258: {
260:   if (tol)    *tol    = pep->tol;
261:   if (maxits) *maxits = pep->max_it;
262:   return 0;
263: }

265: /*@
266:    PEPSetTolerances - Sets the tolerance and maximum iteration count used
267:    by the PEP convergence tests.

269:    Logically Collective on pep

271:    Input Parameters:
272: +  pep - the polynomial eigensolver context
273: .  tol - the convergence tolerance
274: -  maxits - maximum number of iterations to use

276:    Options Database Keys:
277: +  -pep_tol <tol> - Sets the convergence tolerance
278: -  -pep_max_it <maxits> - Sets the maximum number of iterations allowed

280:    Notes:
281:    Use PETSC_DEFAULT for either argument to assign a reasonably good value.

283:    Level: intermediate

285: .seealso: PEPGetTolerances()
286: @*/
287: PetscErrorCode PEPSetTolerances(PEP pep,PetscReal tol,PetscInt maxits)
288: {
292:   if (tol == PETSC_DEFAULT) {
293:     pep->tol   = PETSC_DEFAULT;
294:     pep->state = PEP_STATE_INITIAL;
295:   } else {
297:     pep->tol = tol;
298:   }
299:   if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
300:     pep->max_it = PETSC_DEFAULT;
301:     pep->state  = PEP_STATE_INITIAL;
302:   } else {
304:     pep->max_it = maxits;
305:   }
306:   return 0;
307: }

309: /*@C
310:    PEPGetDimensions - Gets the number of eigenvalues to compute
311:    and the dimension of the subspace.

313:    Not Collective

315:    Input Parameter:
316: .  pep - the polynomial eigensolver context

318:    Output Parameters:
319: +  nev - number of eigenvalues to compute
320: .  ncv - the maximum dimension of the subspace to be used by the solver
321: -  mpd - the maximum dimension allowed for the projected problem

323:    Notes:
324:    The user can specify NULL for any parameter that is not needed.

326:    Level: intermediate

328: .seealso: PEPSetDimensions()
329: @*/
330: PetscErrorCode PEPGetDimensions(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
331: {
333:   if (nev) *nev = pep->nev;
334:   if (ncv) *ncv = pep->ncv;
335:   if (mpd) *mpd = pep->mpd;
336:   return 0;
337: }

339: /*@
340:    PEPSetDimensions - Sets the number of eigenvalues to compute
341:    and the dimension of the subspace.

343:    Logically Collective on pep

345:    Input Parameters:
346: +  pep - the polynomial eigensolver context
347: .  nev - number of eigenvalues to compute
348: .  ncv - the maximum dimension of the subspace to be used by the solver
349: -  mpd - the maximum dimension allowed for the projected problem

351:    Options Database Keys:
352: +  -pep_nev <nev> - Sets the number of eigenvalues
353: .  -pep_ncv <ncv> - Sets the dimension of the subspace
354: -  -pep_mpd <mpd> - Sets the maximum projected dimension

356:    Notes:
357:    Use PETSC_DEFAULT for ncv and mpd to assign a reasonably good value, which is
358:    dependent on the solution method.

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

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

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

372:    Level: intermediate

374: .seealso: PEPGetDimensions(), PEPSetInterval(), PEPSTOARSetDimensions()
375: @*/
376: PetscErrorCode PEPSetDimensions(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)
377: {
383:   pep->nev = nev;
384:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
385:     pep->ncv = PETSC_DEFAULT;
386:   } else {
388:     pep->ncv = ncv;
389:   }
390:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
391:     pep->mpd = PETSC_DEFAULT;
392:   } else {
394:     pep->mpd = mpd;
395:   }
396:   pep->state = PEP_STATE_INITIAL;
397:   return 0;
398: }

400: /*@
401:    PEPSetWhichEigenpairs - Specifies which portion of the spectrum is
402:    to be sought.

404:    Logically Collective on pep

406:    Input Parameters:
407: +  pep   - eigensolver context obtained from PEPCreate()
408: -  which - the portion of the spectrum to be sought

410:    Possible values:
411:    The parameter 'which' can have one of these values

413: +     PEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
414: .     PEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
415: .     PEP_LARGEST_REAL - largest real parts
416: .     PEP_SMALLEST_REAL - smallest real parts
417: .     PEP_LARGEST_IMAGINARY - largest imaginary parts
418: .     PEP_SMALLEST_IMAGINARY - smallest imaginary parts
419: .     PEP_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
420: .     PEP_TARGET_REAL - eigenvalues with real part closest to target
421: .     PEP_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
422: .     PEP_ALL - all eigenvalues contained in a given interval or region
423: -     PEP_WHICH_USER - user defined ordering set with PEPSetEigenvalueComparison()

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:    Not all eigensolvers implemented in PEP account for all the possible values
439:    stated above. If SLEPc is compiled for real numbers PEP_LARGEST_IMAGINARY
440:    and PEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
441:    for eigenvalue selection.

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

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

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

454:    Level: intermediate

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

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

496:     Not Collective

498:     Input Parameter:
499: .   pep - eigensolver context obtained from PEPCreate()

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

504:     Notes:
505:     See PEPSetWhichEigenpairs() for possible values of 'which'.

507:     Level: intermediate

509: .seealso: PEPSetWhichEigenpairs(), PEPWhich
510: @*/
511: PetscErrorCode PEPGetWhichEigenpairs(PEP pep,PEPWhich *which)
512: {
515:   *which = pep->which;
516:   return 0;
517: }

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

523:    Logically Collective on pep

525:    Input Parameters:
526: +  pep  - eigensolver context obtained from PEPCreate()
527: .  func - a pointer to the comparison function
528: -  ctx  - a context pointer (the last parameter to the comparison function)

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

533: +   ar     - real part of the 1st eigenvalue
534: .   ai     - imaginary part of the 1st eigenvalue
535: .   br     - real part of the 2nd eigenvalue
536: .   bi     - imaginary part of the 2nd eigenvalue
537: .   res    - result of comparison
538: -   ctx    - optional context, as set by PEPSetEigenvalueComparison()

540:    Note:
541:    The returning parameter 'res' can be
542: +  negative - if the 1st eigenvalue is preferred to the 2st one
543: .  zero     - if both eigenvalues are equally preferred
544: -  positive - if the 2st eigenvalue is preferred to the 1st one

546:    Level: advanced

548: .seealso: PEPSetWhichEigenpairs(), PEPWhich
549: @*/
550: PetscErrorCode PEPSetEigenvalueComparison(PEP pep,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)
551: {
553:   pep->sc->comparison    = func;
554:   pep->sc->comparisonctx = ctx;
555:   pep->which             = PEP_WHICH_USER;
556:   return 0;
557: }

559: /*@
560:    PEPSetProblemType - Specifies the type of the polynomial eigenvalue problem.

562:    Logically Collective on pep

564:    Input Parameters:
565: +  pep  - the polynomial eigensolver context
566: -  type - a known type of polynomial eigenvalue problem

568:    Options Database Keys:
569: +  -pep_general - general problem with no particular structure
570: .  -pep_hermitian - problem whose coefficient matrices are Hermitian
571: .  -pep_hyperbolic - Hermitian problem that satisfies the definition of hyperbolic
572: -  -pep_gyroscopic - problem with Hamiltonian structure

574:    Notes:
575:    Allowed values for the problem type are general (PEP_GENERAL), Hermitian
576:    (PEP_HERMITIAN), hyperbolic (PEP_HYPERBOLIC), and gyroscopic (PEP_GYROSCOPIC).

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

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

587:    Level: intermediate

589: .seealso: PEPSetOperators(), PEPSetType(), PEPGetProblemType(), PEPProblemType
590: @*/
591: PetscErrorCode PEPSetProblemType(PEP pep,PEPProblemType type)
592: {
596:   if (type != pep->problem_type) {
597:     pep->problem_type = type;
598:     pep->state = PEP_STATE_INITIAL;
599:   }
600:   return 0;
601: }

603: /*@
604:    PEPGetProblemType - Gets the problem type from the PEP object.

606:    Not Collective

608:    Input Parameter:
609: .  pep - the polynomial eigensolver context

611:    Output Parameter:
612: .  type - the problem type

614:    Level: intermediate

616: .seealso: PEPSetProblemType(), PEPProblemType
617: @*/
618: PetscErrorCode PEPGetProblemType(PEP pep,PEPProblemType *type)
619: {
622:   *type = pep->problem_type;
623:   return 0;
624: }

626: /*@
627:    PEPSetBasis - Specifies the type of polynomial basis used to describe the
628:    polynomial eigenvalue problem.

630:    Logically Collective on pep

632:    Input Parameters:
633: +  pep   - the polynomial eigensolver context
634: -  basis - the type of polynomial basis

636:    Options Database Key:
637: .  -pep_basis <basis> - Select the basis type

639:    Notes:
640:    By default, the coefficient matrices passed via PEPSetOperators() are
641:    expressed in the monomial basis, i.e.
642:    P(lambda) = A_0 + lambda*A_1 + lambda^2*A_2 + ... + lambda^d*A_d.
643:    Other polynomial bases may have better numerical behaviour, but the user
644:    must then pass the coefficient matrices accordingly.

646:    Level: intermediate

648: .seealso: PEPSetOperators(), PEPGetBasis(), PEPBasis
649: @*/
650: PetscErrorCode PEPSetBasis(PEP pep,PEPBasis basis)
651: {
654:   pep->basis = basis;
655:   return 0;
656: }

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

661:    Not Collective

663:    Input Parameter:
664: .  pep - the polynomial eigensolver context

666:    Output Parameter:
667: .  basis - the polynomial basis

669:    Level: intermediate

671: .seealso: PEPSetBasis(), PEPBasis
672: @*/
673: PetscErrorCode PEPGetBasis(PEP pep,PEPBasis *basis)
674: {
677:   *basis = pep->basis;
678:   return 0;
679: }

681: /*@
682:    PEPSetTrackAll - Specifies if the solver must compute the residual of all
683:    approximate eigenpairs or not.

685:    Logically Collective on pep

687:    Input Parameters:
688: +  pep      - the eigensolver context
689: -  trackall - whether compute all residuals or not

691:    Notes:
692:    If the user sets trackall=PETSC_TRUE then the solver explicitly computes
693:    the residual for each eigenpair approximation. Computing the residual is
694:    usually an expensive operation and solvers commonly compute the associated
695:    residual to the first unconverged eigenpair.

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

699:    Level: developer

701: .seealso: PEPGetTrackAll()
702: @*/
703: PetscErrorCode PEPSetTrackAll(PEP pep,PetscBool trackall)
704: {
707:   pep->trackall = trackall;
708:   return 0;
709: }

711: /*@
712:    PEPGetTrackAll - Returns the flag indicating whether all residual norms must
713:    be computed or not.

715:    Not Collective

717:    Input Parameter:
718: .  pep - the eigensolver context

720:    Output Parameter:
721: .  trackall - the returned flag

723:    Level: developer

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

735: /*@C
736:    PEPSetConvergenceTestFunction - Sets a function to compute the error estimate
737:    used in the convergence test.

739:    Logically Collective on pep

741:    Input Parameters:
742: +  pep     - eigensolver context obtained from PEPCreate()
743: .  func    - a pointer to the convergence test function
744: .  ctx     - context for private data for the convergence routine (may be null)
745: -  destroy - a routine for destroying the context (may be null)

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

750: +   pep    - eigensolver context obtained from PEPCreate()
751: .   eigr   - real part of the eigenvalue
752: .   eigi   - imaginary part of the eigenvalue
753: .   res    - residual norm associated to the eigenpair
754: .   errest - (output) computed error estimate
755: -   ctx    - optional context, as set by PEPSetConvergenceTestFunction()

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

761:    Level: advanced

763: .seealso: PEPSetConvergenceTest(), PEPSetTolerances()
764: @*/
765: PetscErrorCode PEPSetConvergenceTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
766: {
768:   if (pep->convergeddestroy) (*pep->convergeddestroy)(pep->convergedctx);
769:   pep->convergeduser    = func;
770:   pep->convergeddestroy = destroy;
771:   pep->convergedctx     = ctx;
772:   if (func == PEPConvergedRelative) pep->conv = PEP_CONV_REL;
773:   else if (func == PEPConvergedNorm) pep->conv = PEP_CONV_NORM;
774:   else if (func == PEPConvergedAbsolute) pep->conv = PEP_CONV_ABS;
775:   else {
776:     pep->conv      = PEP_CONV_USER;
777:     pep->converged = pep->convergeduser;
778:   }
779:   return 0;
780: }

782: /*@
783:    PEPSetConvergenceTest - Specifies how to compute the error estimate
784:    used in the convergence test.

786:    Logically Collective on pep

788:    Input Parameters:
789: +  pep  - eigensolver context obtained from PEPCreate()
790: -  conv - the type of convergence test

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

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

805:    Level: intermediate

807: .seealso: PEPGetConvergenceTest(), PEPSetConvergenceTestFunction(), PEPSetStoppingTest(), PEPConv
808: @*/
809: PetscErrorCode PEPSetConvergenceTest(PEP pep,PEPConv conv)
810: {
813:   switch (conv) {
814:     case PEP_CONV_ABS:  pep->converged = PEPConvergedAbsolute; break;
815:     case PEP_CONV_REL:  pep->converged = PEPConvergedRelative; break;
816:     case PEP_CONV_NORM: pep->converged = PEPConvergedNorm; break;
817:     case PEP_CONV_USER:
819:       pep->converged = pep->convergeduser;
820:       break;
821:     default:
822:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
823:   }
824:   pep->conv = conv;
825:   return 0;
826: }

828: /*@
829:    PEPGetConvergenceTest - Gets the method used to compute the error estimate
830:    used in the convergence test.

832:    Not Collective

834:    Input Parameters:
835: .  pep   - eigensolver context obtained from PEPCreate()

837:    Output Parameters:
838: .  conv  - the type of convergence test

840:    Level: intermediate

842: .seealso: PEPSetConvergenceTest(), PEPConv
843: @*/
844: PetscErrorCode PEPGetConvergenceTest(PEP pep,PEPConv *conv)
845: {
848:   *conv = pep->conv;
849:   return 0;
850: }

852: /*@C
853:    PEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
854:    iteration of the eigensolver.

856:    Logically Collective on pep

858:    Input Parameters:
859: +  pep     - eigensolver context obtained from PEPCreate()
860: .  func    - pointer to the stopping test function
861: .  ctx     - context for private data for the stopping routine (may be null)
862: -  destroy - a routine for destroying the context (may be null)

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

867: +   pep    - eigensolver context obtained from PEPCreate()
868: .   its    - current number of iterations
869: .   max_it - maximum number of iterations
870: .   nconv  - number of currently converged eigenpairs
871: .   nev    - number of requested eigenpairs
872: .   reason - (output) result of the stopping test
873: -   ctx    - optional context, as set by PEPSetStoppingTestFunction()

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

881:    Level: advanced

883: .seealso: PEPSetStoppingTest(), PEPStoppingBasic()
884: @*/
885: PetscErrorCode PEPSetStoppingTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
886: {
888:   if (pep->stoppingdestroy) (*pep->stoppingdestroy)(pep->stoppingctx);
889:   pep->stoppinguser    = func;
890:   pep->stoppingdestroy = destroy;
891:   pep->stoppingctx     = ctx;
892:   if (func == PEPStoppingBasic) pep->stop = PEP_STOP_BASIC;
893:   else {
894:     pep->stop     = PEP_STOP_USER;
895:     pep->stopping = pep->stoppinguser;
896:   }
897:   return 0;
898: }

900: /*@
901:    PEPSetStoppingTest - Specifies how to decide the termination of the outer
902:    loop of the eigensolver.

904:    Logically Collective on pep

906:    Input Parameters:
907: +  pep  - eigensolver context obtained from PEPCreate()
908: -  stop - the type of stopping test

910:    Options Database Keys:
911: +  -pep_stop_basic - Sets the default stopping test
912: -  -pep_stop_user  - Selects the user-defined stopping test

914:    Note:
915:    The parameter 'stop' can have one of these values
916: +     PEP_STOP_BASIC - default stopping test
917: -     PEP_STOP_USER  - function set by PEPSetStoppingTestFunction()

919:    Level: advanced

921: .seealso: PEPGetStoppingTest(), PEPSetStoppingTestFunction(), PEPSetConvergenceTest(), PEPStop
922: @*/
923: PetscErrorCode PEPSetStoppingTest(PEP pep,PEPStop stop)
924: {
927:   switch (stop) {
928:     case PEP_STOP_BASIC: pep->stopping = PEPStoppingBasic; break;
929:     case PEP_STOP_USER:
931:       pep->stopping = pep->stoppinguser;
932:       break;
933:     default:
934:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
935:   }
936:   pep->stop = stop;
937:   return 0;
938: }

940: /*@
941:    PEPGetStoppingTest - Gets the method used to decide the termination of the outer
942:    loop of the eigensolver.

944:    Not Collective

946:    Input Parameters:
947: .  pep   - eigensolver context obtained from PEPCreate()

949:    Output Parameters:
950: .  stop  - the type of stopping test

952:    Level: advanced

954: .seealso: PEPSetStoppingTest(), PEPStop
955: @*/
956: PetscErrorCode PEPGetStoppingTest(PEP pep,PEPStop *stop)
957: {
960:   *stop = pep->stop;
961:   return 0;
962: }

964: /*@
965:    PEPSetScale - Specifies the scaling strategy to be used.

967:    Logically Collective on pep

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

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

984:    Notes:
985:    There are two non-exclusive scaling strategies, scalar and diagonal.

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

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

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: {
1011:   pep->scale = scale;
1012:   if (scale==PEP_SCALE_SCALAR || scale==PEP_SCALE_BOTH) {
1014:     if (alpha == PETSC_DEFAULT || alpha == PETSC_DECIDE) {
1015:       pep->sfactor = 0.0;
1016:       pep->sfactor_set = PETSC_FALSE;
1017:     } else {
1019:       pep->sfactor = alpha;
1020:       pep->sfactor_set = PETSC_TRUE;
1021:     }
1022:   }
1023:   if (scale==PEP_SCALE_DIAGONAL || scale==PEP_SCALE_BOTH) {
1024:     if (Dl) {
1027:       PetscObjectReference((PetscObject)Dl);
1028:       VecDestroy(&pep->Dl);
1029:       pep->Dl = Dl;
1030:     }
1031:     if (Dr) {
1034:       PetscObjectReference((PetscObject)Dr);
1035:       VecDestroy(&pep->Dr);
1036:       pep->Dr = Dr;
1037:     }
1040:     if (its==PETSC_DECIDE || its==PETSC_DEFAULT) pep->sits = 5;
1041:     else pep->sits = its;
1042:     if (lambda==PETSC_DECIDE || lambda==PETSC_DEFAULT) pep->slambda = 1.0;
1043:     else {
1045:       pep->slambda = lambda;
1046:     }
1047:   }
1048:   pep->state = PEP_STATE_INITIAL;
1049:   return 0;
1050: }

1052: /*@C
1053:    PEPGetScale - Gets the scaling strategy used by the PEP object, and the
1054:    associated parameters.

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

1058:    Input Parameter:
1059: .  pep - the eigensolver context

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

1069:    Level: intermediate

1071:    Note:
1072:    The user can specify NULL for any parameter that is not needed.

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

1077: .seealso: PEPSetScale(), PEPSetUp()
1078: @*/
1079: PetscErrorCode PEPGetScale(PEP pep,PEPScale *scale,PetscReal *alpha,Vec *Dl,Vec *Dr,PetscInt *its,PetscReal *lambda)
1080: {
1082:   if (scale)  *scale  = pep->scale;
1083:   if (alpha)  *alpha  = pep->sfactor;
1084:   if (Dl)     *Dl     = pep->Dl;
1085:   if (Dr)     *Dr     = pep->Dr;
1086:   if (its)    *its    = pep->sits;
1087:   if (lambda) *lambda = pep->slambda;
1088:   return 0;
1089: }

1091: /*@
1092:    PEPSetExtract - Specifies the extraction strategy to be used.

1094:    Logically Collective on pep

1096:    Input Parameters:
1097: +  pep     - the eigensolver context
1098: -  extract - extraction strategy

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

1103:    Level: intermediate

1105: .seealso: PEPGetExtract()
1106: @*/
1107: PetscErrorCode PEPSetExtract(PEP pep,PEPExtract extract)
1108: {
1111:   pep->extract = extract;
1112:   return 0;
1113: }

1115: /*@
1116:    PEPGetExtract - Gets the extraction strategy used by the PEP object.

1118:    Not Collective

1120:    Input Parameter:
1121: .  pep - the eigensolver context

1123:    Output Parameter:
1124: .  extract - extraction strategy

1126:    Level: intermediate

1128: .seealso: PEPSetExtract(), PEPExtract
1129: @*/
1130: PetscErrorCode PEPGetExtract(PEP pep,PEPExtract *extract)
1131: {
1134:   *extract = pep->extract;
1135:   return 0;
1136: }

1138: /*@
1139:    PEPSetRefine - Specifies the refinement type (and options) to be used
1140:    after the solve.

1142:    Logically Collective on pep

1144:    Input Parameters:
1145: +  pep    - the polynomial eigensolver context
1146: .  refine - refinement type
1147: .  npart  - number of partitions of the communicator
1148: .  tol    - the convergence tolerance
1149: .  its    - maximum number of refinement iterations
1150: -  scheme - which scheme to be used for solving the involved linear systems

1152:    Options Database Keys:
1153: +  -pep_refine <type> - refinement type, one of <none,simple,multiple>
1154: .  -pep_refine_partitions <n> - the number of partitions
1155: .  -pep_refine_tol <tol> - the tolerance
1156: .  -pep_refine_its <its> - number of iterations
1157: -  -pep_refine_scheme - to set the scheme for the linear solves

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

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

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

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

1182:    Level: intermediate

1184: .seealso: PEPGetRefine()
1185: @*/
1186: PetscErrorCode PEPSetRefine(PEP pep,PEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,PEPRefineScheme scheme)
1187: {
1188:   PetscMPIInt    size;

1196:   pep->refine = refine;
1197:   if (refine) {  /* process parameters only if not REFINE_NONE */
1198:     if (npart!=pep->npart) {
1199:       PetscSubcommDestroy(&pep->refinesubc);
1200:       KSPDestroy(&pep->refineksp);
1201:     }
1202:     if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
1203:       pep->npart = 1;
1204:     } else {
1205:       MPI_Comm_size(PetscObjectComm((PetscObject)pep),&size);
1207:       pep->npart = npart;
1208:     }
1209:     if (tol == PETSC_DEFAULT || tol == PETSC_DECIDE) {
1210:       pep->rtol = PETSC_DEFAULT;
1211:     } else {
1213:       pep->rtol = tol;
1214:     }
1215:     if (its==PETSC_DECIDE || its==PETSC_DEFAULT) {
1216:       pep->rits = PETSC_DEFAULT;
1217:     } else {
1219:       pep->rits = its;
1220:     }
1221:     pep->scheme = scheme;
1222:   }
1223:   pep->state = PEP_STATE_INITIAL;
1224:   return 0;
1225: }

1227: /*@C
1228:    PEPGetRefine - Gets the refinement strategy used by the PEP object, and the
1229:    associated parameters.

1231:    Not Collective

1233:    Input Parameter:
1234: .  pep - the polynomial eigensolver context

1236:    Output Parameters:
1237: +  refine - refinement type
1238: .  npart  - number of partitions of the communicator
1239: .  tol    - the convergence tolerance
1240: .  its    - maximum number of refinement iterations
1241: -  scheme - the scheme used for solving linear systems

1243:    Level: intermediate

1245:    Note:
1246:    The user can specify NULL for any parameter that is not needed.

1248: .seealso: PEPSetRefine()
1249: @*/
1250: PetscErrorCode PEPGetRefine(PEP pep,PEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,PEPRefineScheme *scheme)
1251: {
1253:   if (refine) *refine = pep->refine;
1254:   if (npart)  *npart  = pep->npart;
1255:   if (tol)    *tol    = pep->rtol;
1256:   if (its)    *its    = pep->rits;
1257:   if (scheme) *scheme = pep->scheme;
1258:   return 0;
1259: }

1261: /*@C
1262:    PEPSetOptionsPrefix - Sets the prefix used for searching for all
1263:    PEP options in the database.

1265:    Logically Collective on pep

1267:    Input Parameters:
1268: +  pep - the polynomial eigensolver context
1269: -  prefix - the prefix string to prepend to all PEP option requests

1271:    Notes:
1272:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1273:    The first character of all runtime options is AUTOMATICALLY the
1274:    hyphen.

1276:    For example, to distinguish between the runtime options for two
1277:    different PEP contexts, one could call
1278: .vb
1279:       PEPSetOptionsPrefix(pep1,"qeig1_")
1280:       PEPSetOptionsPrefix(pep2,"qeig2_")
1281: .ve

1283:    Level: advanced

1285: .seealso: PEPAppendOptionsPrefix(), PEPGetOptionsPrefix()
1286: @*/
1287: PetscErrorCode PEPSetOptionsPrefix(PEP pep,const char *prefix)
1288: {
1290:   if (!pep->st) PEPGetST(pep,&pep->st);
1291:   STSetOptionsPrefix(pep->st,prefix);
1292:   if (!pep->V) PEPGetBV(pep,&pep->V);
1293:   BVSetOptionsPrefix(pep->V,prefix);
1294:   if (!pep->ds) PEPGetDS(pep,&pep->ds);
1295:   DSSetOptionsPrefix(pep->ds,prefix);
1296:   if (!pep->rg) PEPGetRG(pep,&pep->rg);
1297:   RGSetOptionsPrefix(pep->rg,prefix);
1298:   PetscObjectSetOptionsPrefix((PetscObject)pep,prefix);
1299:   return 0;
1300: }

1302: /*@C
1303:    PEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1304:    PEP options in the database.

1306:    Logically Collective on pep

1308:    Input Parameters:
1309: +  pep - the polynomial eigensolver context
1310: -  prefix - the prefix string to prepend to all PEP option requests

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

1316:    Level: advanced

1318: .seealso: PEPSetOptionsPrefix(), PEPGetOptionsPrefix()
1319: @*/
1320: PetscErrorCode PEPAppendOptionsPrefix(PEP pep,const char *prefix)
1321: {
1323:   if (!pep->st) PEPGetST(pep,&pep->st);
1324:   STAppendOptionsPrefix(pep->st,prefix);
1325:   if (!pep->V) PEPGetBV(pep,&pep->V);
1326:   BVAppendOptionsPrefix(pep->V,prefix);
1327:   if (!pep->ds) PEPGetDS(pep,&pep->ds);
1328:   DSAppendOptionsPrefix(pep->ds,prefix);
1329:   if (!pep->rg) PEPGetRG(pep,&pep->rg);
1330:   RGAppendOptionsPrefix(pep->rg,prefix);
1331:   PetscObjectAppendOptionsPrefix((PetscObject)pep,prefix);
1332:   return 0;
1333: }

1335: /*@C
1336:    PEPGetOptionsPrefix - Gets the prefix used for searching for all
1337:    PEP options in the database.

1339:    Not Collective

1341:    Input Parameters:
1342: .  pep - the polynomial eigensolver context

1344:    Output Parameters:
1345: .  prefix - pointer to the prefix string used is returned

1347:    Note:
1348:    On the Fortran side, the user should pass in a string 'prefix' of
1349:    sufficient length to hold the prefix.

1351:    Level: advanced

1353: .seealso: PEPSetOptionsPrefix(), PEPAppendOptionsPrefix()
1354: @*/
1355: PetscErrorCode PEPGetOptionsPrefix(PEP pep,const char *prefix[])
1356: {
1359:   PetscObjectGetOptionsPrefix((PetscObject)pep,prefix);
1360:   return 0;
1361: }