Actual source code: pepopts.c

slepc-3.11.2 2019-07-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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)->options,((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)->options,((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,array[2]={0,0};
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:     PetscOptionsBoolGroup("-pep_hyperbolic","Hyperbolic polynomial eigenvalue problem","PEPSetProblemType",&flg);
136:     if (flg) { PEPSetProblemType(pep,PEP_HYPERBOLIC); }
137:     PetscOptionsBoolGroupEnd("-pep_gyroscopic","Gyroscopic polynomial eigenvalue problem","PEPSetProblemType",&flg);
138:     if (flg) { PEPSetProblemType(pep,PEP_GYROSCOPIC); }

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

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

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

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

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

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

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

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

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

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

222:     k = 2;
223:     PetscOptionsRealArray("-pep_interval","Computational interval (two real values separated with a comma without spaces)","PEPSetInterval",array,&k,&flg);
224:     if (flg) {
225:       if (k<2) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_SIZ,"Must pass two values in -pep_interval (comma-separated without spaces)");
226:       PEPSetWhichEigenpairs(pep,PEP_ALL);
227:       PEPSetInterval(pep,array[0],array[1]);
228:     }

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

259:     /* -----------------------------------------------------------------------*/
260:     PetscOptionsName("-pep_view","Print detailed information on solver used","PEPView",NULL);
261:     PetscOptionsName("-pep_view_vectors","View computed eigenvectors","PEPVectorsView",NULL);
262:     PetscOptionsName("-pep_view_values","View computed eigenvalues","PEPValuesView",NULL);
263:     PetscOptionsName("-pep_converged_reason","Print reason for convergence, and number of iterations","PEPReasonView",NULL);
264:     PetscOptionsName("-pep_error_absolute","Print absolute errors of each eigenpair","PEPErrorView",NULL);
265:     PetscOptionsName("-pep_error_relative","Print relative errors of each eigenpair","PEPErrorView",NULL);
266:     PetscOptionsName("-pep_error_backward","Print backward errors of each eigenpair","PEPErrorView",NULL);

268:     if (pep->ops->setfromoptions) {
269:       (*pep->ops->setfromoptions)(PetscOptionsObject,pep);
270:     }
271:     PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)pep);
272:   PetscOptionsEnd();

274:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
275:   BVSetFromOptions(pep->V);
276:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
277:   RGSetFromOptions(pep->rg);
278:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
279:   DSSetFromOptions(pep->ds);
280:   if (!pep->st) { PEPGetST(pep,&pep->st); }
281:   PEPSetDefaultST(pep);
282:   STSetFromOptions(pep->st);
283:   if (!pep->refineksp) { PEPRefineGetKSP(pep,&pep->refineksp); }
284:   KSPSetFromOptions(pep->refineksp);
285:   return(0);
286: }

288: /*@C
289:    PEPGetTolerances - Gets the tolerance and maximum iteration count used
290:    by the PEP convergence tests.

292:    Not Collective

294:    Input Parameter:
295: .  pep - the polynomial eigensolver context

297:    Output Parameters:
298: +  tol - the convergence tolerance
299: -  maxits - maximum number of iterations

301:    Notes:
302:    The user can specify NULL for any parameter that is not needed.

304:    Level: intermediate

306: .seealso: PEPSetTolerances()
307: @*/
308: PetscErrorCode PEPGetTolerances(PEP pep,PetscReal *tol,PetscInt *maxits)
309: {
312:   if (tol)    *tol    = pep->tol;
313:   if (maxits) *maxits = pep->max_it;
314:   return(0);
315: }

317: /*@
318:    PEPSetTolerances - Sets the tolerance and maximum iteration count used
319:    by the PEP convergence tests.

321:    Logically Collective on PEP

323:    Input Parameters:
324: +  pep - the polynomial eigensolver context
325: .  tol - the convergence tolerance
326: -  maxits - maximum number of iterations to use

328:    Options Database Keys:
329: +  -pep_tol <tol> - Sets the convergence tolerance
330: -  -pep_max_it <maxits> - Sets the maximum number of iterations allowed

332:    Notes:
333:    Use PETSC_DEFAULT for either argument to assign a reasonably good value.

335:    Level: intermediate

337: .seealso: PEPGetTolerances()
338: @*/
339: PetscErrorCode PEPSetTolerances(PEP pep,PetscReal tol,PetscInt maxits)
340: {
345:   if (tol == PETSC_DEFAULT) {
346:     pep->tol   = PETSC_DEFAULT;
347:     pep->state = PEP_STATE_INITIAL;
348:   } else {
349:     if (tol <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
350:     pep->tol = tol;
351:   }
352:   if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
353:     pep->max_it = 0;
354:     pep->state  = PEP_STATE_INITIAL;
355:   } else {
356:     if (maxits <= 0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
357:     pep->max_it = maxits;
358:   }
359:   return(0);
360: }

362: /*@C
363:    PEPGetDimensions - Gets the number of eigenvalues to compute
364:    and the dimension of the subspace.

366:    Not Collective

368:    Input Parameter:
369: .  pep - the polynomial eigensolver context

371:    Output Parameters:
372: +  nev - number of eigenvalues to compute
373: .  ncv - the maximum dimension of the subspace to be used by the solver
374: -  mpd - the maximum dimension allowed for the projected problem

376:    Notes:
377:    The user can specify NULL for any parameter that is not needed.

379:    Level: intermediate

381: .seealso: PEPSetDimensions()
382: @*/
383: PetscErrorCode PEPGetDimensions(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
384: {
387:   if (nev) *nev = pep->nev;
388:   if (ncv) *ncv = pep->ncv;
389:   if (mpd) *mpd = pep->mpd;
390:   return(0);
391: }

393: /*@
394:    PEPSetDimensions - Sets the number of eigenvalues to compute
395:    and the dimension of the subspace.

397:    Logically Collective on PEP

399:    Input Parameters:
400: +  pep - the polynomial eigensolver context
401: .  nev - number of eigenvalues to compute
402: .  ncv - the maximum dimension of the subspace to be used by the solver
403: -  mpd - the maximum dimension allowed for the projected problem

405:    Options Database Keys:
406: +  -pep_nev <nev> - Sets the number of eigenvalues
407: .  -pep_ncv <ncv> - Sets the dimension of the subspace
408: -  -pep_mpd <mpd> - Sets the maximum projected dimension

410:    Notes:
411:    Use PETSC_DEFAULT for ncv and mpd to assign a reasonably good value, which is
412:    dependent on the solution method.

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

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

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

426:    Level: intermediate

428: .seealso: PEPGetDimensions(), PEPSetInterval(), PEPSTOARSetDimensions()
429: @*/
430: PetscErrorCode PEPSetDimensions(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)
431: {
437:   if (nev<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
438:   pep->nev = nev;
439:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
440:     pep->ncv = 0;
441:   } else {
442:     if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
443:     pep->ncv = ncv;
444:   }
445:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
446:     pep->mpd = 0;
447:   } else {
448:     if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
449:     pep->mpd = mpd;
450:   }
451:   pep->state = PEP_STATE_INITIAL;
452:   return(0);
453: }

455: /*@
456:    PEPSetWhichEigenpairs - Specifies which portion of the spectrum is
457:    to be sought.

459:    Logically Collective on PEP

461:    Input Parameters:
462: +  pep   - eigensolver context obtained from PEPCreate()
463: -  which - the portion of the spectrum to be sought

465:    Possible values:
466:    The parameter 'which' can have one of these values

468: +     PEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
469: .     PEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
470: .     PEP_LARGEST_REAL - largest real parts
471: .     PEP_SMALLEST_REAL - smallest real parts
472: .     PEP_LARGEST_IMAGINARY - largest imaginary parts
473: .     PEP_SMALLEST_IMAGINARY - smallest imaginary parts
474: .     PEP_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
475: .     PEP_TARGET_REAL - eigenvalues with real part closest to target
476: .     PEP_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
477: .     PEP_ALL - all eigenvalues contained in a given interval
478: -     PEP_WHICH_USER - user defined ordering set with PEPSetEigenvalueComparison()

480:    Options Database Keys:
481: +   -pep_largest_magnitude - Sets largest eigenvalues in magnitude
482: .   -pep_smallest_magnitude - Sets smallest eigenvalues in magnitude
483: .   -pep_largest_real - Sets largest real parts
484: .   -pep_smallest_real - Sets smallest real parts
485: .   -pep_largest_imaginary - Sets largest imaginary parts
486: .   -pep_smallest_imaginary - Sets smallest imaginary parts
487: .   -pep_target_magnitude - Sets eigenvalues closest to target
488: .   -pep_target_real - Sets real parts closest to target
489: .   -pep_target_imaginary - Sets imaginary parts closest to target
490: -   -pep_all - Sets all eigenvalues in an interval

492:    Notes:
493:    Not all eigensolvers implemented in PEP account for all the possible values
494:    stated above. If SLEPc is compiled for real numbers PEP_LARGEST_IMAGINARY
495:    and PEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
496:    for eigenvalue selection.

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

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

503:    PEP_ALL is intended for use in combination with an interval (see
504:    PEPSetInterval()), when all eigenvalues within the interval are requested.
505:    In that case, the number of eigenvalues is unknown, so the nev parameter
506:    has a different sense, see PEPSetDimensions().

508:    Level: intermediate

510: .seealso: PEPGetWhichEigenpairs(), PEPSetTarget(), PEPSetInterval(),
511:           PEPSetDimensions(), PEPSetEigenvalueComparison(), PEPWhich
512: @*/
513: PetscErrorCode PEPSetWhichEigenpairs(PEP pep,PEPWhich which)
514: {
518:   switch (which) {
519:     case PEP_LARGEST_MAGNITUDE:
520:     case PEP_SMALLEST_MAGNITUDE:
521:     case PEP_LARGEST_REAL:
522:     case PEP_SMALLEST_REAL:
523:     case PEP_LARGEST_IMAGINARY:
524:     case PEP_SMALLEST_IMAGINARY:
525:     case PEP_TARGET_MAGNITUDE:
526:     case PEP_TARGET_REAL:
527: #if defined(PETSC_USE_COMPLEX)
528:     case PEP_TARGET_IMAGINARY:
529: #endif
530:     case PEP_ALL:
531:     case PEP_WHICH_USER:
532:       if (pep->which != which) {
533:         pep->state = PEP_STATE_INITIAL;
534:         pep->which = which;
535:       }
536:       break;
537:     default:
538:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
539:   }
540:   return(0);
541: }

543: /*@
544:     PEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
545:     sought.

547:     Not Collective

549:     Input Parameter:
550: .   pep - eigensolver context obtained from PEPCreate()

552:     Output Parameter:
553: .   which - the portion of the spectrum to be sought

555:     Notes:
556:     See PEPSetWhichEigenpairs() for possible values of 'which'.

558:     Level: intermediate

560: .seealso: PEPSetWhichEigenpairs(), PEPWhich
561: @*/
562: PetscErrorCode PEPGetWhichEigenpairs(PEP pep,PEPWhich *which)
563: {
567:   *which = pep->which;
568:   return(0);
569: }

571: /*@C
572:    PEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
573:    when PEPSetWhichEigenpairs() is set to PEP_WHICH_USER.

575:    Logically Collective on PEP

577:    Input Parameters:
578: +  pep  - eigensolver context obtained from PEPCreate()
579: .  func - a pointer to the comparison function
580: -  ctx  - a context pointer (the last parameter to the comparison function)

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

585: +   ar     - real part of the 1st eigenvalue
586: .   ai     - imaginary part of the 1st eigenvalue
587: .   br     - real part of the 2nd eigenvalue
588: .   bi     - imaginary part of the 2nd eigenvalue
589: .   res    - result of comparison
590: -   ctx    - optional context, as set by PEPSetEigenvalueComparison()

592:    Note:
593:    The returning parameter 'res' can be
594: +  negative - if the 1st eigenvalue is preferred to the 2st one
595: .  zero     - if both eigenvalues are equally preferred
596: -  positive - if the 2st eigenvalue is preferred to the 1st one

598:    Level: advanced

600: .seealso: PEPSetWhichEigenpairs(), PEPWhich
601: @*/
602: PetscErrorCode PEPSetEigenvalueComparison(PEP pep,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)
603: {
606:   pep->sc->comparison    = func;
607:   pep->sc->comparisonctx = ctx;
608:   pep->which             = PEP_WHICH_USER;
609:   return(0);
610: }

612: /*@
613:    PEPSetProblemType - Specifies the type of the polynomial eigenvalue problem.

615:    Logically Collective on PEP

617:    Input Parameters:
618: +  pep  - the polynomial eigensolver context
619: -  type - a known type of polynomial eigenvalue problem

621:    Options Database Keys:
622: +  -pep_general - general problem with no particular structure
623: .  -pep_hermitian - problem whose coefficient matrices are Hermitian
624: .  -pep_hyperbolic - Hermitian problem that satisfies the definition of hyperbolic
625: -  -pep_gyroscopic - problem with Hamiltonian structure

627:    Notes:
628:    Allowed values for the problem type are: general (PEP_GENERAL), Hermitian
629:    (PEP_HERMITIAN), hyperbolic (PEP_HYPERBOLIC), and gyroscopic (PEP_GYROSCOPIC).

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

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

640:    Level: intermediate

642: .seealso: PEPSetOperators(), PEPSetType(), PEPGetProblemType(), PEPProblemType
643: @*/
644: PetscErrorCode PEPSetProblemType(PEP pep,PEPProblemType type)
645: {
649:   if (type!=PEP_GENERAL && type!=PEP_HERMITIAN && type!=PEP_HYPERBOLIC && type!=PEP_GYROSCOPIC) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
650:   if (type != pep->problem_type) {
651:     pep->problem_type = type;
652:     pep->state = PEP_STATE_INITIAL;
653:   }
654:   return(0);
655: }

657: /*@
658:    PEPGetProblemType - Gets the problem type from the PEP object.

660:    Not Collective

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

665:    Output Parameter:
666: .  type - the problem type

668:    Level: intermediate

670: .seealso: PEPSetProblemType(), PEPProblemType
671: @*/
672: PetscErrorCode PEPGetProblemType(PEP pep,PEPProblemType *type)
673: {
677:   *type = pep->problem_type;
678:   return(0);
679: }

681: /*@
682:    PEPSetBasis - Specifies the type of polynomial basis used to describe the
683:    polynomial eigenvalue problem.

685:    Logically Collective on PEP

687:    Input Parameters:
688: +  pep   - the polynomial eigensolver context
689: -  basis - the type of polynomial basis

691:    Options Database Key:
692: .  -pep_basis <basis> - Select the basis type

694:    Notes:
695:    By default, the coefficient matrices passed via PEPSetOperators() are
696:    expressed in the monomial basis, i.e.
697:    P(lambda) = A_0 + lambda*A_1 + lambda^2*A_2 + ... + lambda^d*A_d.
698:    Other polynomial bases may have better numerical behaviour, but the user
699:    must then pass the coefficient matrices accordingly.

701:    Level: intermediate

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

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

717:    Not Collective

719:    Input Parameter:
720: .  pep - the polynomial eigensolver context

722:    Output Parameter:
723: .  basis - the polynomial basis

725:    Level: intermediate

727: .seealso: PEPSetBasis(), PEPBasis
728: @*/
729: PetscErrorCode PEPGetBasis(PEP pep,PEPBasis *basis)
730: {
734:   *basis = pep->basis;
735:   return(0);
736: }

738: /*@
739:    PEPSetTrackAll - Specifies if the solver must compute the residual of all
740:    approximate eigenpairs or not.

742:    Logically Collective on PEP

744:    Input Parameters:
745: +  pep      - the eigensolver context
746: -  trackall - whether compute all residuals or not

748:    Notes:
749:    If the user sets trackall=PETSC_TRUE then the solver explicitly computes
750:    the residual for each eigenpair approximation. Computing the residual is
751:    usually an expensive operation and solvers commonly compute the associated
752:    residual to the first unconverged eigenpair.
753:    The options '-pep_monitor_all' and '-pep_monitor_lg_all' automatically
754:    activate this option.

756:    Level: developer

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

769: /*@
770:    PEPGetTrackAll - Returns the flag indicating whether all residual norms must
771:    be computed or not.

773:    Not Collective

775:    Input Parameter:
776: .  pep - the eigensolver context

778:    Output Parameter:
779: .  trackall - the returned flag

781:    Level: developer

783: .seealso: PEPSetTrackAll()
784: @*/
785: PetscErrorCode PEPGetTrackAll(PEP pep,PetscBool *trackall)
786: {
790:   *trackall = pep->trackall;
791:   return(0);
792: }

794: /*@C
795:    PEPSetConvergenceTestFunction - Sets a function to compute the error estimate
796:    used in the convergence test.

798:    Logically Collective on PEP

800:    Input Parameters:
801: +  pep     - eigensolver context obtained from PEPCreate()
802: .  func    - a pointer to the convergence test function
803: .  ctx     - context for private data for the convergence routine (may be null)
804: -  destroy - a routine for destroying the context (may be null)

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

809: +   pep    - eigensolver context obtained from PEPCreate()
810: .   eigr   - real part of the eigenvalue
811: .   eigi   - imaginary part of the eigenvalue
812: .   res    - residual norm associated to the eigenpair
813: .   errest - (output) computed error estimate
814: -   ctx    - optional context, as set by PEPSetConvergenceTestFunction()

816:    Note:
817:    If the error estimate returned by the convergence test function is less than
818:    the tolerance, then the eigenvalue is accepted as converged.

820:    Level: advanced

822: .seealso: PEPSetConvergenceTest(), PEPSetTolerances()
823: @*/
824: PetscErrorCode PEPSetConvergenceTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
825: {

830:   if (pep->convergeddestroy) {
831:     (*pep->convergeddestroy)(pep->convergedctx);
832:   }
833:   pep->convergeduser    = func;
834:   pep->convergeddestroy = destroy;
835:   pep->convergedctx     = ctx;
836:   if (func == PEPConvergedRelative) pep->conv = PEP_CONV_REL;
837:   else if (func == PEPConvergedNorm) pep->conv = PEP_CONV_NORM;
838:   else if (func == PEPConvergedAbsolute) pep->conv = PEP_CONV_ABS;
839:   else {
840:     pep->conv      = PEP_CONV_USER;
841:     pep->converged = pep->convergeduser;
842:   }
843:   return(0);
844: }

846: /*@
847:    PEPSetConvergenceTest - Specifies how to compute the error estimate
848:    used in the convergence test.

850:    Logically Collective on PEP

852:    Input Parameters:
853: +  pep  - eigensolver context obtained from PEPCreate()
854: -  conv - the type of convergence test

856:    Options Database Keys:
857: +  -pep_conv_abs    - Sets the absolute convergence test
858: .  -pep_conv_rel    - Sets the convergence test relative to the eigenvalue
859: .  -pep_conv_norm   - Sets the convergence test relative to the matrix norms
860: -  -pep_conv_user   - Selects the user-defined convergence test

862:    Note:
863:    The parameter 'conv' can have one of these values
864: +     PEP_CONV_ABS    - absolute error ||r||
865: .     PEP_CONV_REL    - error relative to the eigenvalue l, ||r||/|l|
866: .     PEP_CONV_NORM   - error relative matrix norms, ||r||/sum_i(l^i*||A_i||)
867: -     PEP_CONV_USER   - function set by PEPSetConvergenceTestFunction()

869:    Level: intermediate

871: .seealso: PEPGetConvergenceTest(), PEPSetConvergenceTestFunction(), PEPSetStoppingTest(), PEPConv
872: @*/
873: PetscErrorCode PEPSetConvergenceTest(PEP pep,PEPConv conv)
874: {
878:   switch (conv) {
879:     case PEP_CONV_ABS:  pep->converged = PEPConvergedAbsolute; break;
880:     case PEP_CONV_REL:  pep->converged = PEPConvergedRelative; break;
881:     case PEP_CONV_NORM: pep->converged = PEPConvergedNorm; break;
882:     case PEP_CONV_USER:
883:       if (!pep->convergeduser) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetConvergenceTestFunction() first");
884:       pep->converged = pep->convergeduser;
885:       break;
886:     default:
887:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
888:   }
889:   pep->conv = conv;
890:   return(0);
891: }

893: /*@
894:    PEPGetConvergenceTest - Gets the method used to compute the error estimate
895:    used in the convergence test.

897:    Not Collective

899:    Input Parameters:
900: .  pep   - eigensolver context obtained from PEPCreate()

902:    Output Parameters:
903: .  conv  - the type of convergence test

905:    Level: intermediate

907: .seealso: PEPSetConvergenceTest(), PEPConv
908: @*/
909: PetscErrorCode PEPGetConvergenceTest(PEP pep,PEPConv *conv)
910: {
914:   *conv = pep->conv;
915:   return(0);
916: }

918: /*@C
919:    PEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
920:    iteration of the eigensolver.

922:    Logically Collective on PEP

924:    Input Parameters:
925: +  pep     - eigensolver context obtained from PEPCreate()
926: .  func    - pointer to the stopping test function
927: .  ctx     - context for private data for the stopping routine (may be null)
928: -  destroy - a routine for destroying the context (may be null)

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

933: +   pep    - eigensolver context obtained from PEPCreate()
934: .   its    - current number of iterations
935: .   max_it - maximum number of iterations
936: .   nconv  - number of currently converged eigenpairs
937: .   nev    - number of requested eigenpairs
938: .   reason - (output) result of the stopping test
939: -   ctx    - optional context, as set by PEPSetStoppingTestFunction()

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

947:    Level: advanced

949: .seealso: PEPSetStoppingTest(), PEPStoppingBasic()
950: @*/
951: PetscErrorCode PEPSetStoppingTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
952: {

957:   if (pep->stoppingdestroy) {
958:     (*pep->stoppingdestroy)(pep->stoppingctx);
959:   }
960:   pep->stoppinguser    = func;
961:   pep->stoppingdestroy = destroy;
962:   pep->stoppingctx     = ctx;
963:   if (func == PEPStoppingBasic) pep->stop = PEP_STOP_BASIC;
964:   else {
965:     pep->stop     = PEP_STOP_USER;
966:     pep->stopping = pep->stoppinguser;
967:   }
968:   return(0);
969: }

971: /*@
972:    PEPSetStoppingTest - Specifies how to decide the termination of the outer
973:    loop of the eigensolver.

975:    Logically Collective on PEP

977:    Input Parameters:
978: +  pep  - eigensolver context obtained from PEPCreate()
979: -  stop - the type of stopping test

981:    Options Database Keys:
982: +  -pep_stop_basic - Sets the default stopping test
983: -  -pep_stop_user  - Selects the user-defined stopping test

985:    Note:
986:    The parameter 'stop' can have one of these values
987: +     PEP_STOP_BASIC - default stopping test
988: -     PEP_STOP_USER  - function set by PEPSetStoppingTestFunction()

990:    Level: advanced

992: .seealso: PEPGetStoppingTest(), PEPSetStoppingTestFunction(), PEPSetConvergenceTest(), PEPStop
993: @*/
994: PetscErrorCode PEPSetStoppingTest(PEP pep,PEPStop stop)
995: {
999:   switch (stop) {
1000:     case PEP_STOP_BASIC: pep->stopping = PEPStoppingBasic; break;
1001:     case PEP_STOP_USER:
1002:       if (!pep->stoppinguser) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ORDER,"Must call PEPSetStoppingTestFunction() first");
1003:       pep->stopping = pep->stoppinguser;
1004:       break;
1005:     default:
1006:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
1007:   }
1008:   pep->stop = stop;
1009:   return(0);
1010: }

1012: /*@
1013:    PEPGetStoppingTest - Gets the method used to decide the termination of the outer
1014:    loop of the eigensolver.

1016:    Not Collective

1018:    Input Parameters:
1019: .  pep   - eigensolver context obtained from PEPCreate()

1021:    Output Parameters:
1022: .  stop  - the type of stopping test

1024:    Level: advanced

1026: .seealso: PEPSetStoppingTest(), PEPStop
1027: @*/
1028: PetscErrorCode PEPGetStoppingTest(PEP pep,PEPStop *stop)
1029: {
1033:   *stop = pep->stop;
1034:   return(0);
1035: }

1037: /*@
1038:    PEPSetScale - Specifies the scaling strategy to be used.

1040:    Logically Collective on PEP

1042:    Input Parameters:
1043: +  pep    - the eigensolver context
1044: .  scale  - scaling strategy
1045: .  alpha  - the scaling factor used in the scalar strategy
1046: .  Dl     - the left diagonal matrix of the diagonal scaling algorithm
1047: .  Dr     - the right diagonal matrix of the diagonal scaling algorithm
1048: .  its    - number of iterations of the diagonal scaling algorithm
1049: -  lambda - approximation to wanted eigenvalues (modulus)

1051:    Options Database Keys:
1052: +  -pep_scale <type> - scaling type, one of <none,scalar,diagonal,both>
1053: .  -pep_scale_factor <alpha> - the scaling factor
1054: .  -pep_scale_its <its> - number of iterations
1055: -  -pep_scale_lambda <lambda> - approximation to eigenvalues

1057:    Notes:
1058:    There are two non-exclusive scaling strategies: scalar and diagonal.

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

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

1076:    Level: intermediate

1078: .seealso: PEPGetScale()
1079: @*/
1080: PetscErrorCode PEPSetScale(PEP pep,PEPScale scale,PetscReal alpha,Vec Dl,Vec Dr,PetscInt its,PetscReal lambda)
1081: {

1087:   pep->scale = scale;
1088:   if (scale==PEP_SCALE_SCALAR || scale==PEP_SCALE_BOTH) {
1090:     if (alpha == PETSC_DEFAULT || alpha == PETSC_DECIDE) {
1091:       pep->sfactor = 0.0;
1092:       pep->sfactor_set = PETSC_FALSE;
1093:     } else {
1094:       if (alpha<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of alpha. Must be > 0");
1095:       pep->sfactor = alpha;
1096:       pep->sfactor_set = PETSC_TRUE;
1097:     }
1098:   }
1099:   if (scale==PEP_SCALE_DIAGONAL || scale==PEP_SCALE_BOTH) {
1100:     if (Dl) {
1103:       PetscObjectReference((PetscObject)Dl);
1104:       VecDestroy(&pep->Dl);
1105:       pep->Dl = Dl;
1106:     }
1107:     if (Dr) {
1110:       PetscObjectReference((PetscObject)Dr);
1111:       VecDestroy(&pep->Dr);
1112:       pep->Dr = Dr;
1113:     }
1116:     if (its==PETSC_DECIDE || its==PETSC_DEFAULT) pep->sits = 5;
1117:     else pep->sits = its;
1118:     if (lambda==PETSC_DECIDE || lambda==PETSC_DEFAULT) pep->slambda = 1.0;
1119:     else if (lambda<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of lambda. Must be > 0");
1120:     else pep->slambda = lambda;
1121:   }
1122:   pep->state = PEP_STATE_INITIAL;
1123:   return(0);
1124: }

1126: /*@C
1127:    PEPGetScale - Gets the scaling strategy used by the PEP object, and the
1128:    associated parameters.

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

1132:    Input Parameter:
1133: .  pep - the eigensolver context

1135:    Output Parameters:
1136: +  scale  - scaling strategy
1137: .  alpha  - the scaling factor used in the scalar strategy
1138: .  Dl     - the left diagonal matrix of the diagonal scaling algorithm
1139: .  Dr     - the right diagonal matrix of the diagonal scaling algorithm
1140: .  its    - number of iterations of the diagonal scaling algorithm
1141: -  lambda - approximation to wanted eigenvalues (modulus)

1143:    Level: intermediate

1145:    Note:
1146:    The user can specify NULL for any parameter that is not needed.

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

1151: .seealso: PEPSetScale(), PEPSetUp()
1152: @*/
1153: PetscErrorCode PEPGetScale(PEP pep,PEPScale *scale,PetscReal *alpha,Vec *Dl,Vec *Dr,PetscInt *its,PetscReal *lambda)
1154: {
1157:   if (scale)  *scale  = pep->scale;
1158:   if (alpha)  *alpha  = pep->sfactor;
1159:   if (Dl)     *Dl     = pep->Dl;
1160:   if (Dr)     *Dr     = pep->Dr;
1161:   if (its)    *its    = pep->sits;
1162:   if (lambda) *lambda = pep->slambda;
1163:   return(0);
1164: }

1166: /*@
1167:    PEPSetExtract - Specifies the extraction strategy to be used.

1169:    Logically Collective on PEP

1171:    Input Parameters:
1172: +  pep     - the eigensolver context
1173: -  extract - extraction strategy

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

1178:    Level: intermediate

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

1191: /*@
1192:    PEPGetExtract - Gets the extraction strategy used by the PEP object.

1194:    Not Collective

1196:    Input Parameter:
1197: .  pep - the eigensolver context

1199:    Output Parameter:
1200: .  extract - extraction strategy

1202:    Level: intermediate

1204: .seealso: PEPSetExtract()
1205: @*/
1206: PetscErrorCode PEPGetExtract(PEP pep,PEPExtract *extract)
1207: {
1210:   if (extract) *extract = pep->extract;
1211:   return(0);
1212: }

1214: /*@
1215:    PEPSetRefine - Specifies the refinement type (and options) to be used
1216:    after the solve.

1218:    Logically Collective on PEP

1220:    Input Parameters:
1221: +  pep    - the polynomial eigensolver context
1222: .  refine - refinement type
1223: .  npart  - number of partitions of the communicator
1224: .  tol    - the convergence tolerance
1225: .  its    - maximum number of refinement iterations
1226: -  scheme - which scheme to be used for solving the involved linear systems

1228:    Options Database Keys:
1229: +  -pep_refine <type> - refinement type, one of <none,simple,multiple>
1230: .  -pep_refine_partitions <n> - the number of partitions
1231: .  -pep_refine_tol <tol> - the tolerance
1232: .  -pep_refine_its <its> - number of iterations
1233: -  -pep_refine_scheme - to set the scheme for the linear solves

1235:    Notes:
1236:    By default, iterative refinement is disabled, since it may be very
1237:    costly. There are two possible refinement strategies: simple and multiple.
1238:    The simple approach performs iterative refinement on each of the
1239:    converged eigenpairs individually, whereas the multiple strategy works
1240:    with the invariant pair as a whole, refining all eigenpairs simultaneously.
1241:    The latter may be required for the case of multiple eigenvalues.

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

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

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

1258:    Level: intermediate

1260: .seealso: PEPGetRefine()
1261: @*/
1262: PetscErrorCode PEPSetRefine(PEP pep,PEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,PEPRefineScheme scheme)
1263: {
1265:   PetscMPIInt    size;

1274:   pep->refine = refine;
1275:   if (refine) {  /* process parameters only if not REFINE_NONE */
1276:     if (npart!=pep->npart) {
1277:       PetscSubcommDestroy(&pep->refinesubc);
1278:       KSPDestroy(&pep->refineksp);
1279:     }
1280:     if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
1281:       pep->npart = 1;
1282:     } else {
1283:       MPI_Comm_size(PetscObjectComm((PetscObject)pep),&size);
1284:       if (npart<1 || npart>size) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
1285:       pep->npart = npart;
1286:     }
1287:     if (tol == PETSC_DEFAULT || tol == PETSC_DECIDE) {
1288:       pep->rtol = PETSC_DEFAULT;
1289:     } else {
1290:       if (tol<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1291:       pep->rtol = tol;
1292:     }
1293:     if (its==PETSC_DECIDE || its==PETSC_DEFAULT) {
1294:       pep->rits = PETSC_DEFAULT;
1295:     } else {
1296:       if (its<0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be >= 0");
1297:       pep->rits = its;
1298:     }
1299:     pep->scheme = scheme;
1300:   }
1301:   pep->state = PEP_STATE_INITIAL;
1302:   return(0);
1303: }

1305: /*@C
1306:    PEPGetRefine - Gets the refinement strategy used by the PEP object, and the
1307:    associated parameters.

1309:    Not Collective

1311:    Input Parameter:
1312: .  pep - the polynomial eigensolver context

1314:    Output Parameters:
1315: +  refine - refinement type
1316: .  npart  - number of partitions of the communicator
1317: .  tol    - the convergence tolerance
1318: .  its    - maximum number of refinement iterations
1319: -  scheme - the scheme used for solving linear systems

1321:    Level: intermediate

1323:    Note:
1324:    The user can specify NULL for any parameter that is not needed.

1326: .seealso: PEPSetRefine()
1327: @*/
1328: PetscErrorCode PEPGetRefine(PEP pep,PEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,PEPRefineScheme *scheme)
1329: {
1332:   if (refine) *refine = pep->refine;
1333:   if (npart)  *npart  = pep->npart;
1334:   if (tol)    *tol    = pep->rtol;
1335:   if (its)    *its    = pep->rits;
1336:   if (scheme) *scheme = pep->scheme;
1337:   return(0);
1338: }

1340: /*@C
1341:    PEPSetOptionsPrefix - Sets the prefix used for searching for all
1342:    PEP options in the database.

1344:    Logically Collective on PEP

1346:    Input Parameters:
1347: +  pep - the polynomial eigensolver context
1348: -  prefix - the prefix string to prepend to all PEP option requests

1350:    Notes:
1351:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1352:    The first character of all runtime options is AUTOMATICALLY the
1353:    hyphen.

1355:    For example, to distinguish between the runtime options for two
1356:    different PEP contexts, one could call
1357: .vb
1358:       PEPSetOptionsPrefix(pep1,"qeig1_")
1359:       PEPSetOptionsPrefix(pep2,"qeig2_")
1360: .ve

1362:    Level: advanced

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

1372:   if (!pep->st) { PEPGetST(pep,&pep->st); }
1373:   STSetOptionsPrefix(pep->st,prefix);
1374:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
1375:   BVSetOptionsPrefix(pep->V,prefix);
1376:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
1377:   DSSetOptionsPrefix(pep->ds,prefix);
1378:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
1379:   RGSetOptionsPrefix(pep->rg,prefix);
1380:   PetscObjectSetOptionsPrefix((PetscObject)pep,prefix);
1381:   return(0);
1382: }

1384: /*@C
1385:    PEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1386:    PEP options in the database.

1388:    Logically Collective on PEP

1390:    Input Parameters:
1391: +  pep - the polynomial eigensolver context
1392: -  prefix - the prefix string to prepend to all PEP option requests

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

1398:    Level: advanced

1400: .seealso: PEPSetOptionsPrefix(), PEPGetOptionsPrefix()
1401: @*/
1402: PetscErrorCode PEPAppendOptionsPrefix(PEP pep,const char *prefix)
1403: {

1408:   if (!pep->st) { PEPGetST(pep,&pep->st); }
1409:   STAppendOptionsPrefix(pep->st,prefix);
1410:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
1411:   BVAppendOptionsPrefix(pep->V,prefix);
1412:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
1413:   DSAppendOptionsPrefix(pep->ds,prefix);
1414:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
1415:   RGAppendOptionsPrefix(pep->rg,prefix);
1416:   PetscObjectAppendOptionsPrefix((PetscObject)pep,prefix);
1417:   return(0);
1418: }

1420: /*@C
1421:    PEPGetOptionsPrefix - Gets the prefix used for searching for all
1422:    PEP options in the database.

1424:    Not Collective

1426:    Input Parameters:
1427: .  pep - the polynomial eigensolver context

1429:    Output Parameters:
1430: .  prefix - pointer to the prefix string used is returned

1432:    Note:
1433:    On the Fortran side, the user should pass in a string 'prefix' of
1434:    sufficient length to hold the prefix.

1436:    Level: advanced

1438: .seealso: PEPSetOptionsPrefix(), PEPAppendOptionsPrefix()
1439: @*/
1440: PetscErrorCode PEPGetOptionsPrefix(PEP pep,const char *prefix[])
1441: {

1447:   PetscObjectGetOptionsPrefix((PetscObject)pep,prefix);
1448:   return(0);
1449: }