Actual source code: pepopts.c

slepc-3.8.2 2017-12-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, 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>       /*I "slepcpep.h" I*/
 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: .  name     - the monitor option name
 27: .  help     - message indicating what monitoring is done
 28: .  manual   - manual page for the monitor
 29: .  monitor  - the monitor function, whose context is a PetscViewerAndFormat
 30: -  trackall - whether this monitor tracks all eigenvalues or not

 32:    Level: developer

 34: .seealso: PEPMonitorSet(), PEPSetTrackAll(), PEPConvMonitorSetFromOptions()
 35: @*/
 36: PetscErrorCode PEPMonitorSetFromOptions(PEP pep,const char name[],const char help[],const char manual[],PetscErrorCode (*monitor)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscBool trackall)
 37: {
 38:   PetscErrorCode       ierr;
 39:   PetscBool            flg;
 40:   PetscViewer          viewer;
 41:   PetscViewerFormat    format;
 42:   PetscViewerAndFormat *vf;

 45:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,name,&viewer,&format,&flg);
 46:   if (flg) {
 47:     PetscViewerAndFormatCreate(viewer,format,&vf);
 48:     PetscObjectDereference((PetscObject)viewer);
 49:     PEPMonitorSet(pep,(PetscErrorCode (*)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
 50:     if (trackall) {
 51:       PEPSetTrackAll(pep,PETSC_TRUE);
 52:     }
 53:   }
 54:   return(0);
 55: }

 57: /*@C
 58:    PEPConvMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
 59:    indicated by the user (for monitors that only show iteration numbers of convergence).

 61:    Collective on PEP

 63:    Input Parameters:
 64: +  pep      - the polynomial eigensolver context
 65: .  name     - the monitor option name
 66: .  help     - message indicating what monitoring is done
 67: .  manual   - manual page for the monitor
 68: -  monitor  - the monitor function, whose context is a SlepcConvMonitor

 70:    Level: developer

 72: .seealso: PEPMonitorSet(), PEPMonitorSetFromOptions()
 73: @*/
 74: PetscErrorCode PEPConvMonitorSetFromOptions(PEP pep,const char name[],const char help[],const char manual[],PetscErrorCode (*monitor)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,SlepcConvMonitor))
 75: {
 76:   PetscErrorCode    ierr;
 77:   PetscBool         flg;
 78:   PetscViewer       viewer;
 79:   PetscViewerFormat format;
 80:   SlepcConvMonitor  ctx;

 83:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,name,&viewer,&format,&flg);
 84:   if (flg) {
 85:     SlepcConvMonitorCreate(viewer,format,&ctx);
 86:     PetscObjectDereference((PetscObject)viewer);
 87:     PEPMonitorSet(pep,(PetscErrorCode (*)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))monitor,ctx,(PetscErrorCode (*)(void**))SlepcConvMonitorDestroy);
 88:   }
 89:   return(0);
 90: }

 92: /*@
 93:    PEPSetFromOptions - Sets PEP options from the options database.
 94:    This routine must be called before PEPSetUp() if the user is to be
 95:    allowed to set the solver type.

 97:    Collective on PEP

 99:    Input Parameters:
100: .  pep - the polynomial eigensolver context

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

105:    Level: beginner
106: @*/
107: PetscErrorCode PEPSetFromOptions(PEP pep)
108: {
109:   PetscErrorCode  ierr;
110:   char            type[256];
111:   PetscBool       set,flg,flg1,flg2,flg3,flg4,flg5;
112:   PetscReal       r,t;
113:   PetscScalar     s;
114:   PetscInt        i,j,k;
115:   PetscDrawLG     lg;
116:   PEPScale        scale;
117:   PEPRefine       refine;
118:   PEPRefineScheme scheme;

122:   PEPRegisterAll();
123:   PetscObjectOptionsBegin((PetscObject)pep);
124:     PetscOptionsFList("-pep_type","Polynomial eigensolver method","PEPSetType",PEPList,(char*)(((PetscObject)pep)->type_name?((PetscObject)pep)->type_name:PEPTOAR),type,256,&flg);
125:     if (flg) {
126:       PEPSetType(pep,type);
127:     } else if (!((PetscObject)pep)->type_name) {
128:       PEPSetType(pep,PEPTOAR);
129:     }

131:     PetscOptionsBoolGroupBegin("-pep_general","General polynomial eigenvalue problem","PEPSetProblemType",&flg);
132:     if (flg) { PEPSetProblemType(pep,PEP_GENERAL); }
133:     PetscOptionsBoolGroup("-pep_hermitian","Hermitian polynomial eigenvalue problem","PEPSetProblemType",&flg);
134:     if (flg) { PEPSetProblemType(pep,PEP_HERMITIAN); }
135:     PetscOptionsBoolGroupEnd("-pep_gyroscopic","Gyroscopic polynomial eigenvalue problem","PEPSetProblemType",&flg);
136:     if (flg) { PEPSetProblemType(pep,PEP_GYROSCOPIC); }

138:     scale = pep->scale;
139:     PetscOptionsEnum("-pep_scale","Scaling strategy","PEPSetScale",PEPScaleTypes,(PetscEnum)scale,(PetscEnum*)&scale,&flg1);
140:     r = pep->sfactor;
141:     PetscOptionsReal("-pep_scale_factor","Scale factor","PEPSetScale",pep->sfactor,&r,&flg2);
142:     if (!flg2 && r==1.0) r = PETSC_DEFAULT;
143:     j = pep->sits;
144:     PetscOptionsInt("-pep_scale_its","Number of iterations in diagonal scaling","PEPSetScale",pep->sits,&j,&flg3);
145:     t = pep->slambda;
146:     PetscOptionsReal("-pep_scale_lambda","Estimate of eigenvalue (modulus) for diagonal scaling","PEPSetScale",pep->slambda,&t,&flg4);
147:     if (flg1 || flg2 || flg3 || flg4) { PEPSetScale(pep,scale,r,NULL,NULL,j,t); }

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

151:     refine = pep->refine;
152:     PetscOptionsEnum("-pep_refine","Iterative refinement method","PEPSetRefine",PEPRefineTypes,(PetscEnum)refine,(PetscEnum*)&refine,&flg1);
153:     i = pep->npart;
154:     PetscOptionsInt("-pep_refine_partitions","Number of partitions of the communicator for iterative refinement","PEPSetRefine",pep->npart,&i,&flg2);
155:     r = pep->rtol;
156:     PetscOptionsReal("-pep_refine_tol","Tolerance for iterative refinement","PEPSetRefine",pep->rtol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/1000:pep->rtol,&r,&flg3);
157:     j = pep->rits;
158:     PetscOptionsInt("-pep_refine_its","Maximum number of iterations for iterative refinement","PEPSetRefine",pep->rits,&j,&flg4);
159:     scheme = pep->scheme;
160:     PetscOptionsEnum("-pep_refine_scheme","Scheme used for linear systems within iterative refinement","PEPSetRefine",PEPRefineSchemes,(PetscEnum)scheme,(PetscEnum*)&scheme,&flg5);
161:     if (flg1 || flg2 || flg3 || flg4 || flg5) { PEPSetRefine(pep,refine,i,r,j,scheme); }

163:     i = pep->max_it? pep->max_it: PETSC_DEFAULT;
164:     PetscOptionsInt("-pep_max_it","Maximum number of iterations","PEPSetTolerances",pep->max_it,&i,&flg1);
165:     r = pep->tol;
166:     PetscOptionsReal("-pep_tol","Tolerance","PEPSetTolerances",pep->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:pep->tol,&r,&flg2);
167:     if (flg1 || flg2) { PEPSetTolerances(pep,r,i); }

169:     PetscOptionsBoolGroupBegin("-pep_conv_rel","Relative error convergence test","PEPSetConvergenceTest",&flg);
170:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_REL); }
171:     PetscOptionsBoolGroup("-pep_conv_norm","Convergence test relative to the matrix norms","PEPSetConvergenceTest",&flg);
172:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_NORM); }
173:     PetscOptionsBoolGroup("-pep_conv_abs","Absolute error convergence test","PEPSetConvergenceTest",&flg);
174:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_ABS); }
175:     PetscOptionsBoolGroupEnd("-pep_conv_user","User-defined convergence test","PEPSetConvergenceTest",&flg);
176:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_USER); }

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

183:     i = pep->nev;
184:     PetscOptionsInt("-pep_nev","Number of eigenvalues to compute","PEPSetDimensions",pep->nev,&i,&flg1);
185:     j = pep->ncv? pep->ncv: PETSC_DEFAULT;
186:     PetscOptionsInt("-pep_ncv","Number of basis vectors","PEPSetDimensions",pep->ncv,&j,&flg2);
187:     k = pep->mpd? pep->mpd: PETSC_DEFAULT;
188:     PetscOptionsInt("-pep_mpd","Maximum dimension of projected problem","PEPSetDimensions",pep->mpd,&k,&flg3);
189:     if (flg1 || flg2 || flg3) { PEPSetDimensions(pep,i,j,k); }

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

193:     PetscOptionsBoolGroupBegin("-pep_largest_magnitude","Compute largest eigenvalues in magnitude","PEPSetWhichEigenpairs",&flg);
194:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_LARGEST_MAGNITUDE); }
195:     PetscOptionsBoolGroup("-pep_smallest_magnitude","Compute smallest eigenvalues in magnitude","PEPSetWhichEigenpairs",&flg);
196:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_SMALLEST_MAGNITUDE); }
197:     PetscOptionsBoolGroup("-pep_largest_real","Compute eigenvalues with largest real parts","PEPSetWhichEigenpairs",&flg);
198:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_LARGEST_REAL); }
199:     PetscOptionsBoolGroup("-pep_smallest_real","Compute eigenvalues with smallest real parts","PEPSetWhichEigenpairs",&flg);
200:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_SMALLEST_REAL); }
201:     PetscOptionsBoolGroup("-pep_largest_imaginary","Compute eigenvalues with largest imaginary parts","PEPSetWhichEigenpairs",&flg);
202:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_LARGEST_IMAGINARY); }
203:     PetscOptionsBoolGroup("-pep_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","PEPSetWhichEigenpairs",&flg);
204:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_SMALLEST_IMAGINARY); }
205:     PetscOptionsBoolGroup("-pep_target_magnitude","Compute eigenvalues closest to target","PEPSetWhichEigenpairs",&flg);
206:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE); }
207:     PetscOptionsBoolGroup("-pep_target_real","Compute eigenvalues with real parts closest to target","PEPSetWhichEigenpairs",&flg);
208:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_TARGET_REAL); }
209:     PetscOptionsBoolGroupEnd("-pep_target_imaginary","Compute eigenvalues with imaginary parts closest to target","PEPSetWhichEigenpairs",&flg);
210:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_TARGET_IMAGINARY); }

212:     PetscOptionsScalar("-pep_target","Value of the target","PEPSetTarget",pep->target,&s,&flg);
213:     if (flg) {
214:       if (pep->which!=PEP_TARGET_REAL && pep->which!=PEP_TARGET_IMAGINARY) {
215:         PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE);
216:       }
217:       PEPSetTarget(pep,s);
218:     }

220:     /* -----------------------------------------------------------------------*/
221:     /*
222:       Cancels all monitors hardwired into code before call to PEPSetFromOptions()
223:     */
224:     PetscOptionsBool("-pep_monitor_cancel","Remove any hardwired monitor routines","PEPMonitorCancel",PETSC_FALSE,&flg,&set);
225:     if (set && flg) {
226:       PEPMonitorCancel(pep);
227:     }
228:     /*
229:       Text monitors
230:     */
231:     PEPMonitorSetFromOptions(pep,"-pep_monitor","Monitor first unconverged approximate eigenvalue and error estimate","PEPMonitorFirst",PEPMonitorFirst,PETSC_FALSE);
232:     PEPConvMonitorSetFromOptions(pep,"-pep_monitor_conv","Monitor approximate eigenvalues and error estimates as they converge","PEPMonitorConverged",PEPMonitorConverged);
233:     PEPMonitorSetFromOptions(pep,"-pep_monitor_all","Monitor approximate eigenvalues and error estimates","PEPMonitorAll",PEPMonitorAll,PETSC_TRUE);
234:     /*
235:       Line graph monitors
236:     */
237:     PetscOptionsBool("-pep_monitor_lg","Monitor first unconverged approximate error estimate graphically","PEPMonitorSet",PETSC_FALSE,&flg,&set);
238:     if (set && flg) {
239:       PEPMonitorLGCreate(PetscObjectComm((PetscObject)pep),NULL,"Error estimates",PETSC_DECIDE,PETSC_DECIDE,300,300,&lg);
240:       PEPMonitorSet(pep,PEPMonitorLG,lg,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
241:     }
242:     PetscOptionsBool("-pep_monitor_lg_all","Monitor error estimates graphically","PEPMonitorSet",PETSC_FALSE,&flg,&set);
243:     if (set && flg) {
244:       PEPMonitorLGCreate(PetscObjectComm((PetscObject)pep),NULL,"Error estimates",PETSC_DECIDE,PETSC_DECIDE,300,300,&lg);
245:       PEPMonitorSet(pep,PEPMonitorLGAll,lg,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
246:       PEPSetTrackAll(pep,PETSC_TRUE);
247:     }

249:     /* -----------------------------------------------------------------------*/
250:     PetscOptionsName("-pep_view","Print detailed information on solver used","PEPView",NULL);
251:     PetscOptionsName("-pep_view_vectors","View computed eigenvectors","PEPVectorsView",NULL);
252:     PetscOptionsName("-pep_view_values","View computed eigenvalues","PEPValuesView",NULL);
253:     PetscOptionsName("-pep_converged_reason","Print reason for convergence, and number of iterations","PEPReasonView",NULL);
254:     PetscOptionsName("-pep_error_absolute","Print absolute errors of each eigenpair","PEPErrorView",NULL);
255:     PetscOptionsName("-pep_error_relative","Print relative errors of each eigenpair","PEPErrorView",NULL);
256:     PetscOptionsName("-pep_error_backward","Print backward errors of each eigenpair","PEPErrorView",NULL);

258:     if (pep->ops->setfromoptions) {
259:       (*pep->ops->setfromoptions)(PetscOptionsObject,pep);
260:     }
261:     PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)pep);
262:   PetscOptionsEnd();

264:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
265:   BVSetFromOptions(pep->V);
266:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
267:   RGSetFromOptions(pep->rg);
268:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
269:   DSSetFromOptions(pep->ds);
270:   if (!pep->st) { PEPGetST(pep,&pep->st); }
271:   PEPSetDefaultST(pep);
272:   STSetFromOptions(pep->st);
273:   if (!pep->refineksp) { PEPRefineGetKSP(pep,&pep->refineksp); }
274:   KSPSetFromOptions(pep->refineksp);
275:   return(0);
276: }

278: /*@C
279:    PEPGetTolerances - Gets the tolerance and maximum iteration count used
280:    by the PEP convergence tests.

282:    Not Collective

284:    Input Parameter:
285: .  pep - the polynomial eigensolver context

287:    Output Parameters:
288: +  tol - the convergence tolerance
289: -  maxits - maximum number of iterations

291:    Notes:
292:    The user can specify NULL for any parameter that is not needed.

294:    Level: intermediate

296: .seealso: PEPSetTolerances()
297: @*/
298: PetscErrorCode PEPGetTolerances(PEP pep,PetscReal *tol,PetscInt *maxits)
299: {
302:   if (tol)    *tol    = pep->tol;
303:   if (maxits) *maxits = pep->max_it;
304:   return(0);
305: }

307: /*@
308:    PEPSetTolerances - Sets the tolerance and maximum iteration count used
309:    by the PEP convergence tests.

311:    Logically Collective on PEP

313:    Input Parameters:
314: +  pep - the polynomial eigensolver context
315: .  tol - the convergence tolerance
316: -  maxits - maximum number of iterations to use

318:    Options Database Keys:
319: +  -pep_tol <tol> - Sets the convergence tolerance
320: -  -pep_max_it <maxits> - Sets the maximum number of iterations allowed

322:    Notes:
323:    Use PETSC_DEFAULT for either argument to assign a reasonably good value.

325:    Level: intermediate

327: .seealso: PEPGetTolerances()
328: @*/
329: PetscErrorCode PEPSetTolerances(PEP pep,PetscReal tol,PetscInt maxits)
330: {
335:   if (tol == PETSC_DEFAULT) {
336:     pep->tol   = PETSC_DEFAULT;
337:     pep->state = PEP_STATE_INITIAL;
338:   } else {
339:     if (tol <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
340:     pep->tol = tol;
341:   }
342:   if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
343:     pep->max_it = 0;
344:     pep->state  = PEP_STATE_INITIAL;
345:   } else {
346:     if (maxits <= 0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
347:     pep->max_it = maxits;
348:   }
349:   return(0);
350: }

352: /*@C
353:    PEPGetDimensions - Gets the number of eigenvalues to compute
354:    and the dimension of the subspace.

356:    Not Collective

358:    Input Parameter:
359: .  pep - the polynomial eigensolver context

361:    Output Parameters:
362: +  nev - number of eigenvalues to compute
363: .  ncv - the maximum dimension of the subspace to be used by the solver
364: -  mpd - the maximum dimension allowed for the projected problem

366:    Notes:
367:    The user can specify NULL for any parameter that is not needed.

369:    Level: intermediate

371: .seealso: PEPSetDimensions()
372: @*/
373: PetscErrorCode PEPGetDimensions(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
374: {
377:   if (nev) *nev = pep->nev;
378:   if (ncv) *ncv = pep->ncv;
379:   if (mpd) *mpd = pep->mpd;
380:   return(0);
381: }

383: /*@
384:    PEPSetDimensions - Sets the number of eigenvalues to compute
385:    and the dimension of the subspace.

387:    Logically Collective on PEP

389:    Input Parameters:
390: +  pep - the polynomial eigensolver context
391: .  nev - number of eigenvalues to compute
392: .  ncv - the maximum dimension of the subspace to be used by the solver
393: -  mpd - the maximum dimension allowed for the projected problem

395:    Options Database Keys:
396: +  -pep_nev <nev> - Sets the number of eigenvalues
397: .  -pep_ncv <ncv> - Sets the dimension of the subspace
398: -  -pep_mpd <mpd> - Sets the maximum projected dimension

400:    Notes:
401:    Use PETSC_DEFAULT for ncv and mpd to assign a reasonably good value, which is
402:    dependent on the solution method.

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

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

413:    Level: intermediate

415: .seealso: PEPGetDimensions()
416: @*/
417: PetscErrorCode PEPSetDimensions(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)
418: {
424:   if (nev<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
425:   pep->nev = nev;
426:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
427:     pep->ncv = 0;
428:   } else {
429:     if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
430:     pep->ncv = ncv;
431:   }
432:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
433:     pep->mpd = 0;
434:   } else {
435:     if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
436:     pep->mpd = mpd;
437:   }
438:   pep->state = PEP_STATE_INITIAL;
439:   return(0);
440: }

442: /*@
443:    PEPSetWhichEigenpairs - Specifies which portion of the spectrum is
444:    to be sought.

446:    Logically Collective on PEP

448:    Input Parameters:
449: +  pep   - eigensolver context obtained from PEPCreate()
450: -  which - the portion of the spectrum to be sought

452:    Possible values:
453:    The parameter 'which' can have one of these values

455: +     PEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
456: .     PEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
457: .     PEP_LARGEST_REAL - largest real parts
458: .     PEP_SMALLEST_REAL - smallest real parts
459: .     PEP_LARGEST_IMAGINARY - largest imaginary parts
460: .     PEP_SMALLEST_IMAGINARY - smallest imaginary parts
461: .     PEP_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
462: .     PEP_TARGET_REAL - eigenvalues with real part closest to target
463: .     PEP_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
464: -     PEP_WHICH_USER - user defined ordering set with PEPSetEigenvalueComparison()

466:    Options Database Keys:
467: +   -pep_largest_magnitude - Sets largest eigenvalues in magnitude
468: .   -pep_smallest_magnitude - Sets smallest eigenvalues in magnitude
469: .   -pep_largest_real - Sets largest real parts
470: .   -pep_smallest_real - Sets smallest real parts
471: .   -pep_largest_imaginary - Sets largest imaginary parts
472: .   -pep_smallest_imaginary - Sets smallest imaginary parts
473: .   -pep_target_magnitude - Sets eigenvalues closest to target
474: .   -pep_target_real - Sets real parts closest to target
475: -   -pep_target_imaginary - Sets imaginary parts closest to target

477:    Notes:
478:    Not all eigensolvers implemented in PEP account for all the possible values
479:    stated above. If SLEPc is compiled for real numbers PEP_LARGEST_IMAGINARY
480:    and PEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
481:    for eigenvalue selection.

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

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

488:    Level: intermediate

490: .seealso: PEPGetWhichEigenpairs(), PEPSetTarget(), PEPSetEigenvalueComparison(), PEPWhich
491: @*/
492: PetscErrorCode PEPSetWhichEigenpairs(PEP pep,PEPWhich which)
493: {
497:   switch (which) {
498:     case PEP_LARGEST_MAGNITUDE:
499:     case PEP_SMALLEST_MAGNITUDE:
500:     case PEP_LARGEST_REAL:
501:     case PEP_SMALLEST_REAL:
502:     case PEP_LARGEST_IMAGINARY:
503:     case PEP_SMALLEST_IMAGINARY:
504:     case PEP_TARGET_MAGNITUDE:
505:     case PEP_TARGET_REAL:
506: #if defined(PETSC_USE_COMPLEX)
507:     case PEP_TARGET_IMAGINARY:
508: #endif
509:     case PEP_WHICH_USER:
510:       if (pep->which != which) {
511:         pep->state = PEP_STATE_INITIAL;
512:         pep->which = which;
513:       }
514:       break;
515:     default:
516:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
517:   }
518:   return(0);
519: }

521: /*@
522:     PEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
523:     sought.

525:     Not Collective

527:     Input Parameter:
528: .   pep - eigensolver context obtained from PEPCreate()

530:     Output Parameter:
531: .   which - the portion of the spectrum to be sought

533:     Notes:
534:     See PEPSetWhichEigenpairs() for possible values of 'which'.

536:     Level: intermediate

538: .seealso: PEPSetWhichEigenpairs(), PEPWhich
539: @*/
540: PetscErrorCode PEPGetWhichEigenpairs(PEP pep,PEPWhich *which)
541: {
545:   *which = pep->which;
546:   return(0);
547: }

549: /*@C
550:    PEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
551:    when PEPSetWhichEigenpairs() is set to PEP_WHICH_USER.

553:    Logically Collective on PEP

555:    Input Parameters:
556: +  pep  - eigensolver context obtained from PEPCreate()
557: .  func - a pointer to the comparison function
558: -  ctx  - a context pointer (the last parameter to the comparison function)

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

563: +   ar     - real part of the 1st eigenvalue
564: .   ai     - imaginary part of the 1st eigenvalue
565: .   br     - real part of the 2nd eigenvalue
566: .   bi     - imaginary part of the 2nd eigenvalue
567: .   res    - result of comparison
568: -   ctx    - optional context, as set by PEPSetEigenvalueComparison()

570:    Note:
571:    The returning parameter 'res' can be
572: +  negative - if the 1st eigenvalue is preferred to the 2st one
573: .  zero     - if both eigenvalues are equally preferred
574: -  positive - if the 2st eigenvalue is preferred to the 1st one

576:    Level: advanced

578: .seealso: PEPSetWhichEigenpairs(), PEPWhich
579: @*/
580: PetscErrorCode PEPSetEigenvalueComparison(PEP pep,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)
581: {
584:   pep->sc->comparison    = func;
585:   pep->sc->comparisonctx = ctx;
586:   pep->which             = PEP_WHICH_USER;
587:   return(0);
588: }

590: /*@
591:    PEPSetProblemType - Specifies the type of the polynomial eigenvalue problem.

593:    Logically Collective on PEP

595:    Input Parameters:
596: +  pep  - the polynomial eigensolver context
597: -  type - a known type of polynomial eigenvalue problem

599:    Options Database Keys:
600: +  -pep_general - general problem with no particular structure
601: .  -pep_hermitian - problem whose coefficient matrices are Hermitian
602: -  -pep_gyroscopic - problem with Hamiltonian structure

604:    Notes:
605:    Allowed values for the problem type are: general (PEP_GENERAL), Hermitian
606:    (PEP_HERMITIAN), and gyroscopic (PEP_GYROSCOPIC).

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

611:    If the problem matrices are Hermitian (symmetric in the real case) or
612:    Hermitian/skew-Hermitian then the solver can exploit this fact to perform
613:    less operations or provide better stability.

615:    Level: intermediate

617: .seealso: PEPSetOperators(), PEPSetType(), PEPGetProblemType(), PEPProblemType
618: @*/
619: PetscErrorCode PEPSetProblemType(PEP pep,PEPProblemType type)
620: {
624:   if (type!=PEP_GENERAL && type!=PEP_HERMITIAN && type!=PEP_GYROSCOPIC) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
625:   if (type != pep->problem_type) {
626:     pep->problem_type = type;
627:     pep->state = PEP_STATE_INITIAL;
628:   }
629:   return(0);
630: }

632: /*@
633:    PEPGetProblemType - Gets the problem type from the PEP object.

635:    Not Collective

637:    Input Parameter:
638: .  pep - the polynomial eigensolver context

640:    Output Parameter:
641: .  type - the problem type

643:    Level: intermediate

645: .seealso: PEPSetProblemType(), PEPProblemType
646: @*/
647: PetscErrorCode PEPGetProblemType(PEP pep,PEPProblemType *type)
648: {
652:   *type = pep->problem_type;
653:   return(0);
654: }

656: /*@
657:    PEPSetBasis - Specifies the type of polynomial basis used to describe the
658:    polynomial eigenvalue problem.

660:    Logically Collective on PEP

662:    Input Parameters:
663: +  pep   - the polynomial eigensolver context
664: -  basis - the type of polynomial basis

666:    Options Database Key:
667: .  -pep_basis <basis> - Select the basis type

669:    Notes:
670:    By default, the coefficient matrices passed via PEPSetOperators() are
671:    expressed in the monomial basis, i.e.
672:    P(lambda) = A_0 + lambda*A_1 + lambda^2*A_2 + ... + lambda^d*A_d.
673:    Other polynomial bases may have better numerical behaviour, but the user
674:    must then pass the coefficient matrices accordingly.

676:    Level: intermediate

678: .seealso: PEPSetOperators(), PEPGetBasis(), PEPBasis
679: @*/
680: PetscErrorCode PEPSetBasis(PEP pep,PEPBasis basis)
681: {
685:   pep->basis = basis;
686:   return(0);
687: }

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

692:    Not Collective

694:    Input Parameter:
695: .  pep - the polynomial eigensolver context

697:    Output Parameter:
698: .  basis - the polynomial basis

700:    Level: intermediate

702: .seealso: PEPSetBasis(), PEPBasis
703: @*/
704: PetscErrorCode PEPGetBasis(PEP pep,PEPBasis *basis)
705: {
709:   *basis = pep->basis;
710:   return(0);
711: }

713: /*@
714:    PEPSetTrackAll - Specifies if the solver must compute the residual of all
715:    approximate eigenpairs or not.

717:    Logically Collective on PEP

719:    Input Parameters:
720: +  pep      - the eigensolver context
721: -  trackall - whether compute all residuals or not

723:    Notes:
724:    If the user sets trackall=PETSC_TRUE then the solver explicitly computes
725:    the residual for each eigenpair approximation. Computing the residual is
726:    usually an expensive operation and solvers commonly compute the associated
727:    residual to the first unconverged eigenpair.
728:    The options '-pep_monitor_all' and '-pep_monitor_lg_all' automatically
729:    activate this option.

731:    Level: developer

733: .seealso: PEPGetTrackAll()
734: @*/
735: PetscErrorCode PEPSetTrackAll(PEP pep,PetscBool trackall)
736: {
740:   pep->trackall = trackall;
741:   return(0);
742: }

744: /*@
745:    PEPGetTrackAll - Returns the flag indicating whether all residual norms must
746:    be computed or not.

748:    Not Collective

750:    Input Parameter:
751: .  pep - the eigensolver context

753:    Output Parameter:
754: .  trackall - the returned flag

756:    Level: developer

758: .seealso: PEPSetTrackAll()
759: @*/
760: PetscErrorCode PEPGetTrackAll(PEP pep,PetscBool *trackall)
761: {
765:   *trackall = pep->trackall;
766:   return(0);
767: }

769: /*@C
770:    PEPSetConvergenceTestFunction - Sets a function to compute the error estimate
771:    used in the convergence test.

773:    Logically Collective on PEP

775:    Input Parameters:
776: +  pep     - eigensolver context obtained from PEPCreate()
777: .  func    - a pointer to the convergence test function
778: .  ctx     - context for private data for the convergence routine (may be null)
779: -  destroy - a routine for destroying the context (may be null)

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

784: +   pep    - eigensolver context obtained from PEPCreate()
785: .   eigr   - real part of the eigenvalue
786: .   eigi   - imaginary part of the eigenvalue
787: .   res    - residual norm associated to the eigenpair
788: .   errest - (output) computed error estimate
789: -   ctx    - optional context, as set by PEPSetConvergenceTestFunction()

791:    Note:
792:    If the error estimate returned by the convergence test function is less than
793:    the tolerance, then the eigenvalue is accepted as converged.

795:    Level: advanced

797: .seealso: PEPSetConvergenceTest(), PEPSetTolerances()
798: @*/
799: PetscErrorCode PEPSetConvergenceTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
800: {

805:   if (pep->convergeddestroy) {
806:     (*pep->convergeddestroy)(pep->convergedctx);
807:   }
808:   pep->convergeduser    = func;
809:   pep->convergeddestroy = destroy;
810:   pep->convergedctx     = ctx;
811:   if (func == PEPConvergedRelative) pep->conv = PEP_CONV_REL;
812:   else if (func == PEPConvergedNorm) pep->conv = PEP_CONV_NORM;
813:   else if (func == PEPConvergedAbsolute) pep->conv = PEP_CONV_ABS;
814:   else {
815:     pep->conv      = PEP_CONV_USER;
816:     pep->converged = pep->convergeduser;
817:   }
818:   return(0);
819: }

821: /*@
822:    PEPSetConvergenceTest - Specifies how to compute the error estimate
823:    used in the convergence test.

825:    Logically Collective on PEP

827:    Input Parameters:
828: +  pep  - eigensolver context obtained from PEPCreate()
829: -  conv - the type of convergence test

831:    Options Database Keys:
832: +  -pep_conv_abs    - Sets the absolute convergence test
833: .  -pep_conv_rel    - Sets the convergence test relative to the eigenvalue
834: .  -pep_conv_norm   - Sets the convergence test relative to the matrix norms
835: -  -pep_conv_user   - Selects the user-defined convergence test

837:    Note:
838:    The parameter 'conv' can have one of these values
839: +     PEP_CONV_ABS    - absolute error ||r||
840: .     PEP_CONV_REL    - error relative to the eigenvalue l, ||r||/|l|
841: .     PEP_CONV_NORM   - error relative matrix norms, ||r||/sum_i(l^i*||A_i||)
842: -     PEP_CONV_USER   - function set by PEPSetConvergenceTestFunction()

844:    Level: intermediate

846: .seealso: PEPGetConvergenceTest(), PEPSetConvergenceTestFunction(), PEPSetStoppingTest(), PEPConv
847: @*/
848: PetscErrorCode PEPSetConvergenceTest(PEP pep,PEPConv conv)
849: {
853:   switch (conv) {
854:     case PEP_CONV_ABS:  pep->converged = PEPConvergedAbsolute; break;
855:     case PEP_CONV_REL:  pep->converged = PEPConvergedRelative; break;
856:     case PEP_CONV_NORM: pep->converged = PEPConvergedNorm; break;
857:     case PEP_CONV_USER:
858:       if (!pep->convergeduser) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetConvergenceTestFunction() first");
859:       pep->converged = pep->convergeduser;
860:       break;
861:     default:
862:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
863:   }
864:   pep->conv = conv;
865:   return(0);
866: }

868: /*@
869:    PEPGetConvergenceTest - Gets the method used to compute the error estimate
870:    used in the convergence test.

872:    Not Collective

874:    Input Parameters:
875: .  pep   - eigensolver context obtained from PEPCreate()

877:    Output Parameters:
878: .  conv  - the type of convergence test

880:    Level: intermediate

882: .seealso: PEPSetConvergenceTest(), PEPConv
883: @*/
884: PetscErrorCode PEPGetConvergenceTest(PEP pep,PEPConv *conv)
885: {
889:   *conv = pep->conv;
890:   return(0);
891: }

893: /*@C
894:    PEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
895:    iteration of the eigensolver.

897:    Logically Collective on PEP

899:    Input Parameters:
900: +  pep     - eigensolver context obtained from PEPCreate()
901: .  func    - pointer to the stopping test function
902: .  ctx     - context for private data for the stopping routine (may be null)
903: -  destroy - a routine for destroying the context (may be null)

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

908: +   pep    - eigensolver context obtained from PEPCreate()
909: .   its    - current number of iterations
910: .   max_it - maximum number of iterations
911: .   nconv  - number of currently converged eigenpairs
912: .   nev    - number of requested eigenpairs
913: .   reason - (output) result of the stopping test
914: -   ctx    - optional context, as set by PEPSetStoppingTestFunction()

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

922:    Level: advanced

924: .seealso: PEPSetStoppingTest(), PEPStoppingBasic()
925: @*/
926: PetscErrorCode PEPSetStoppingTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
927: {

932:   if (pep->stoppingdestroy) {
933:     (*pep->stoppingdestroy)(pep->stoppingctx);
934:   }
935:   pep->stoppinguser    = func;
936:   pep->stoppingdestroy = destroy;
937:   pep->stoppingctx     = ctx;
938:   if (func == PEPStoppingBasic) pep->stop = PEP_STOP_BASIC;
939:   else {
940:     pep->stop     = PEP_STOP_USER;
941:     pep->stopping = pep->stoppinguser;
942:   }
943:   return(0);
944: }

946: /*@
947:    PEPSetStoppingTest - Specifies how to decide the termination of the outer
948:    loop of the eigensolver.

950:    Logically Collective on PEP

952:    Input Parameters:
953: +  pep  - eigensolver context obtained from PEPCreate()
954: -  stop - the type of stopping test

956:    Options Database Keys:
957: +  -pep_stop_basic - Sets the default stopping test
958: -  -pep_stop_user  - Selects the user-defined stopping test

960:    Note:
961:    The parameter 'stop' can have one of these values
962: +     PEP_STOP_BASIC - default stopping test
963: -     PEP_STOP_USER  - function set by PEPSetStoppingTestFunction()

965:    Level: advanced

967: .seealso: PEPGetStoppingTest(), PEPSetStoppingTestFunction(), PEPSetConvergenceTest(), PEPStop
968: @*/
969: PetscErrorCode PEPSetStoppingTest(PEP pep,PEPStop stop)
970: {
974:   switch (stop) {
975:     case PEP_STOP_BASIC: pep->stopping = PEPStoppingBasic; break;
976:     case PEP_STOP_USER:
977:       if (!pep->stoppinguser) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetStoppingTestFunction() first");
978:       pep->stopping = pep->stoppinguser;
979:       break;
980:     default:
981:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
982:   }
983:   pep->stop = stop;
984:   return(0);
985: }

987: /*@
988:    PEPGetStoppingTest - Gets the method used to decide the termination of the outer
989:    loop of the eigensolver.

991:    Not Collective

993:    Input Parameters:
994: .  pep   - eigensolver context obtained from PEPCreate()

996:    Output Parameters:
997: .  stop  - the type of stopping test

999:    Level: advanced

1001: .seealso: PEPSetStoppingTest(), PEPStop
1002: @*/
1003: PetscErrorCode PEPGetStoppingTest(PEP pep,PEPStop *stop)
1004: {
1008:   *stop = pep->stop;
1009:   return(0);
1010: }

1012: /*@C
1013:    PEPSetScale - Specifies the scaling strategy to be used.

1015:    Logically Collective on PEP

1017:    Input Parameters:
1018: +  pep    - the eigensolver context
1019: .  scale  - scaling strategy
1020: .  alpha  - the scaling factor used in the scalar strategy
1021: .  Dl     - the left diagonal matrix of the diagonal scaling algorithm
1022: .  Dr     - the right diagonal matrix of the diagonal scaling algorithm
1023: .  its    - number of iterations of the diagonal scaling algorithm
1024: -  lambda - approximation to wanted eigenvalues (modulus)

1026:    Options Database Keys:
1027: +  -pep_scale <type> - scaling type, one of <none,scalar,diagonal,both>
1028: .  -pep_scale_factor <alpha> - the scaling factor
1029: .  -pep_scale_its <its> - number of iterations
1030: -  -pep_scale_lambda <lambda> - approximation to eigenvalues

1032:    Notes:
1033:    There are two non-exclusive scaling strategies: scalar and diagonal.

1035:    In the scalar strategy, scaling is applied to the eigenvalue, that is,
1036:    mu = lambda/alpha is the new eigenvalue and all matrices are scaled
1037:    accordingly. After solving the scaled problem, the original lambda is
1038:    recovered. Parameter 'alpha' must be positive. Use PETSC_DECIDE to let
1039:    the solver compute a reasonable scaling factor.

1041:    In the diagonal strategy, the solver works implicitly with matrix Dl*A*Dr,
1042:    where Dl and Dr are appropriate diagonal matrices. This improves the accuracy
1043:    of the computed results in some cases. The user may provide the Dr and Dl
1044:    matrices represented as Vec objects storing diagonal elements. If not
1045:    provided, these matrices are computed internally. This option requires
1046:    that the polynomial coefficient matrices are of MATAIJ type.
1047:    The parameter 'its' is the number of iterations performed by the method.
1048:    Parameter 'lambda' must be positive. Use PETSC_DECIDE or set lambda = 1.0 if
1049:    no information about eigenvalues is available.

1051:    Level: intermediate

1053: .seealso: PEPGetScale()
1054: @*/
1055: PetscErrorCode PEPSetScale(PEP pep,PEPScale scale,PetscReal alpha,Vec Dl,Vec Dr,PetscInt its,PetscReal lambda)
1056: {

1062:   pep->scale = scale;
1063:   if (scale==PEP_SCALE_SCALAR || scale==PEP_SCALE_BOTH) {
1065:     if (alpha == PETSC_DEFAULT || alpha == PETSC_DECIDE) {
1066:       pep->sfactor = 0.0;
1067:       pep->sfactor_set = PETSC_FALSE;
1068:     } else {
1069:       if (alpha<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of alpha. Must be > 0");
1070:       pep->sfactor = alpha;
1071:       pep->sfactor_set = PETSC_TRUE;
1072:     }
1073:   }
1074:   if (scale==PEP_SCALE_DIAGONAL || scale==PEP_SCALE_BOTH) {
1075:     if (Dl) {
1078:       PetscObjectReference((PetscObject)Dl);
1079:       VecDestroy(&pep->Dl);
1080:       pep->Dl = Dl;
1081:     }
1082:     if (Dr) {
1085:       PetscObjectReference((PetscObject)Dr);
1086:       VecDestroy(&pep->Dr);
1087:       pep->Dr = Dr;
1088:     }
1091:     if (its==PETSC_DECIDE || its==PETSC_DEFAULT) pep->sits = 5;
1092:     else pep->sits = its;
1093:     if (lambda==PETSC_DECIDE || lambda==PETSC_DEFAULT) pep->slambda = 1.0;
1094:     else if (lambda<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of lambda. Must be > 0");
1095:     else pep->slambda = lambda;
1096:   }
1097:   pep->state = PEP_STATE_INITIAL;
1098:   return(0);
1099: }

1101: /*@C
1102:    PEPGetScale - Gets the scaling strategy used by the PEP object, and the
1103:    associated parameters.

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

1107:    Input Parameter:
1108: .  pep - the eigensolver context

1110:    Output Parameters:
1111: +  scale  - scaling strategy
1112: .  alpha  - the scaling factor used in the scalar strategy
1113: .  Dl     - the left diagonal matrix of the diagonal scaling algorithm
1114: .  Dr     - the right diagonal matrix of the diagonal scaling algorithm
1115: .  its    - number of iterations of the diagonal scaling algorithm
1116: -  lambda - approximation to wanted eigenvalues (modulus)

1118:    Level: intermediate

1120:    Note:
1121:    The user can specify NULL for any parameter that is not needed.

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

1126: .seealso: PEPSetScale(), PEPSetUp()
1127: @*/
1128: PetscErrorCode PEPGetScale(PEP pep,PEPScale *scale,PetscReal *alpha,Vec *Dl,Vec *Dr,PetscInt *its,PetscReal *lambda)
1129: {
1132:   if (scale)  *scale  = pep->scale;
1133:   if (alpha)  *alpha  = pep->sfactor;
1134:   if (Dl)     *Dl     = pep->Dl;
1135:   if (Dr)     *Dr     = pep->Dr;
1136:   if (its)    *its    = pep->sits;
1137:   if (lambda) *lambda = pep->slambda;
1138:   return(0);
1139: }

1141: /*@
1142:    PEPSetExtract - Specifies the extraction strategy to be used.

1144:    Logically Collective on PEP

1146:    Input Parameters:
1147: +  pep     - the eigensolver context
1148: -  extract - extraction strategy

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

1153:    Level: intermediate

1155: .seealso: PEPGetExtract()
1156: @*/
1157: PetscErrorCode PEPSetExtract(PEP pep,PEPExtract extract)
1158: {
1162:   pep->extract = extract;
1163:   return(0);
1164: }

1166: /*@
1167:    PEPGetExtract - Gets the extraction strategy used by the PEP object.

1169:    Not Collective

1171:    Input Parameter:
1172: .  pep - the eigensolver context

1174:    Output Parameter:
1175: .  extract - extraction strategy

1177:    Level: intermediate

1179: .seealso: PEPSetExtract()
1180: @*/
1181: PetscErrorCode PEPGetExtract(PEP pep,PEPExtract *extract)
1182: {
1185:   if (extract) *extract = pep->extract;
1186:   return(0);
1187: }

1189: /*@
1190:    PEPSetRefine - Specifies the refinement type (and options) to be used
1191:    after the solve.

1193:    Logically Collective on PEP

1195:    Input Parameters:
1196: +  pep    - the polynomial eigensolver context
1197: .  refine - refinement type
1198: .  npart  - number of partitions of the communicator
1199: .  tol    - the convergence tolerance
1200: .  its    - maximum number of refinement iterations
1201: -  scheme - which scheme to be used for solving the involved linear systems

1203:    Options Database Keys:
1204: +  -pep_refine <type> - refinement type, one of <none,simple,multiple>
1205: .  -pep_refine_partitions <n> - the number of partitions
1206: .  -pep_refine_tol <tol> - the tolerance
1207: .  -pep_refine_its <its> - number of iterations
1208: -  -pep_refine_scheme - to set the scheme for the linear solves

1210:    Notes:
1211:    By default, iterative refinement is disabled, since it may be very
1212:    costly. There are two possible refinement strategies: simple and multiple.
1213:    The simple approach performs iterative refinement on each of the
1214:    converged eigenpairs individually, whereas the multiple strategy works
1215:    with the invariant pair as a whole, refining all eigenpairs simultaneously.
1216:    The latter may be required for the case of multiple eigenvalues.

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

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

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

1233:    Level: intermediate

1235: .seealso: PEPGetRefine()
1236: @*/
1237: PetscErrorCode PEPSetRefine(PEP pep,PEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,PEPRefineScheme scheme)
1238: {
1240:   PetscMPIInt    size;

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

1280: /*@C
1281:    PEPGetRefine - Gets the refinement strategy used by the PEP object, and the
1282:    associated parameters.

1284:    Not Collective

1286:    Input Parameter:
1287: .  pep - the polynomial eigensolver context

1289:    Output Parameters:
1290: +  refine - refinement type
1291: .  npart  - number of partitions of the communicator
1292: .  tol    - the convergence tolerance
1293: .  its    - maximum number of refinement iterations
1294: -  scheme - the scheme used for solving linear systems

1296:    Level: intermediate

1298:    Note:
1299:    The user can specify NULL for any parameter that is not needed.

1301: .seealso: PEPSetRefine()
1302: @*/
1303: PetscErrorCode PEPGetRefine(PEP pep,PEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,PEPRefineScheme *scheme)
1304: {
1307:   if (refine) *refine = pep->refine;
1308:   if (npart)  *npart  = pep->npart;
1309:   if (tol)    *tol    = pep->rtol;
1310:   if (its)    *its    = pep->rits;
1311:   if (scheme) *scheme = pep->scheme;
1312:   return(0);
1313: }

1315: /*@C
1316:    PEPSetOptionsPrefix - Sets the prefix used for searching for all
1317:    PEP options in the database.

1319:    Logically Collective on PEP

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

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

1330:    For example, to distinguish between the runtime options for two
1331:    different PEP contexts, one could call
1332: .vb
1333:       PEPSetOptionsPrefix(pep1,"qeig1_")
1334:       PEPSetOptionsPrefix(pep2,"qeig2_")
1335: .ve

1337:    Level: advanced

1339: .seealso: PEPAppendOptionsPrefix(), PEPGetOptionsPrefix()
1340: @*/
1341: PetscErrorCode PEPSetOptionsPrefix(PEP pep,const char *prefix)
1342: {

1347:   if (!pep->st) { PEPGetST(pep,&pep->st); }
1348:   STSetOptionsPrefix(pep->st,prefix);
1349:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
1350:   BVSetOptionsPrefix(pep->V,prefix);
1351:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
1352:   DSSetOptionsPrefix(pep->ds,prefix);
1353:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
1354:   RGSetOptionsPrefix(pep->rg,prefix);
1355:   PetscObjectSetOptionsPrefix((PetscObject)pep,prefix);
1356:   return(0);
1357: }

1359: /*@C
1360:    PEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1361:    PEP options in the database.

1363:    Logically Collective on PEP

1365:    Input Parameters:
1366: +  pep - the polynomial eigensolver context
1367: -  prefix - the prefix string to prepend to all PEP option requests

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

1373:    Level: advanced

1375: .seealso: PEPSetOptionsPrefix(), PEPGetOptionsPrefix()
1376: @*/
1377: PetscErrorCode PEPAppendOptionsPrefix(PEP pep,const char *prefix)
1378: {

1383:   if (!pep->st) { PEPGetST(pep,&pep->st); }
1384:   STAppendOptionsPrefix(pep->st,prefix);
1385:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
1386:   BVAppendOptionsPrefix(pep->V,prefix);
1387:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
1388:   DSAppendOptionsPrefix(pep->ds,prefix);
1389:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
1390:   RGAppendOptionsPrefix(pep->rg,prefix);
1391:   PetscObjectAppendOptionsPrefix((PetscObject)pep,prefix);
1392:   return(0);
1393: }

1395: /*@C
1396:    PEPGetOptionsPrefix - Gets the prefix used for searching for all
1397:    PEP options in the database.

1399:    Not Collective

1401:    Input Parameters:
1402: .  pep - the polynomial eigensolver context

1404:    Output Parameters:
1405: .  prefix - pointer to the prefix string used is returned

1407:    Note:
1408:    On the Fortran side, the user should pass in a string 'prefix' of
1409:    sufficient length to hold the prefix.

1411:    Level: advanced

1413: .seealso: PEPSetOptionsPrefix(), PEPAppendOptionsPrefix()
1414: @*/
1415: PetscErrorCode PEPGetOptionsPrefix(PEP pep,const char *prefix[])
1416: {

1422:   PetscObjectGetOptionsPrefix((PetscObject)pep,prefix);
1423:   return(0);
1424: }