Actual source code: pepopts.c

slepc-3.7.3 2016-09-29
Report Typos and Errors
  1: /*
  2:       PEP routines related to options that can be set via the command-line
  3:       or procedurally.

  5:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  6:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  7:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  9:    This file is part of SLEPc.

 11:    SLEPc is free software: you can redistribute it and/or modify it under  the
 12:    terms of version 3 of the GNU Lesser General Public License as published by
 13:    the Free Software Foundation.

 15:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 16:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 17:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 18:    more details.

 20:    You  should have received a copy of the GNU Lesser General  Public  License
 21:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 22:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 23: */

 25: #include <slepc/private/pepimpl.h>       /*I "slepcpep.h" I*/
 26: #include <petscdraw.h>

 30: /*@C
 31:    PEPMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
 32:    indicated by the user.

 34:    Collective on PEP

 36:    Input Parameters:
 37: +  pep      - the polynomial eigensolver context
 38: .  name     - the monitor option name
 39: .  help     - message indicating what monitoring is done
 40: .  manual   - manual page for the monitor
 41: .  monitor  - the monitor function, whose context is a PetscViewerAndFormat
 42: -  trackall - whether this monitor tracks all eigenvalues or not

 44:    Level: developer

 46: .seealso: PEPMonitorSet(), PEPSetTrackAll(), PEPConvMonitorSetFromOptions()
 47: @*/
 48: 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)
 49: {
 50:   PetscErrorCode       ierr;
 51:   PetscBool            flg;
 52:   PetscViewer          viewer;
 53:   PetscViewerFormat    format;
 54:   PetscViewerAndFormat *vf;

 57:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,name,&viewer,&format,&flg);
 58:   if (flg) {
 59:     PetscViewerAndFormatCreate(viewer,format,&vf);
 60:     PetscObjectDereference((PetscObject)viewer);
 61:     PEPMonitorSet(pep,(PetscErrorCode (*)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
 62:     if (trackall) {
 63:       PEPSetTrackAll(pep,PETSC_TRUE);
 64:     }
 65:   }
 66:   return(0);
 67: }

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

 75:    Collective on PEP

 77:    Input Parameters:
 78: +  pep      - the polynomial eigensolver context
 79: .  name     - the monitor option name
 80: .  help     - message indicating what monitoring is done
 81: .  manual   - manual page for the monitor
 82: -  monitor  - the monitor function, whose context is a SlepcConvMonitor

 84:    Level: developer

 86: .seealso: PEPMonitorSet(), PEPMonitorSetFromOptions()
 87: @*/
 88: PetscErrorCode PEPConvMonitorSetFromOptions(PEP pep,const char name[],const char help[],const char manual[],PetscErrorCode (*monitor)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,SlepcConvMonitor))
 89: {
 90:   PetscErrorCode    ierr;
 91:   PetscBool         flg;
 92:   PetscViewer       viewer;
 93:   PetscViewerFormat format;
 94:   SlepcConvMonitor  ctx;

 97:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)pep),((PetscObject)pep)->prefix,name,&viewer,&format,&flg);
 98:   if (flg) {
 99:     SlepcConvMonitorCreate(viewer,format,&ctx);
100:     PetscObjectDereference((PetscObject)viewer);
101:     PEPMonitorSet(pep,(PetscErrorCode (*)(PEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))monitor,ctx,(PetscErrorCode (*)(void**))SlepcConvMonitorDestroy);
102:   }
103:   return(0);
104: }

108: /*@
109:    PEPSetFromOptions - Sets PEP options from the options database.
110:    This routine must be called before PEPSetUp() if the user is to be
111:    allowed to set the solver type.

113:    Collective on PEP

115:    Input Parameters:
116: .  pep - the polynomial eigensolver context

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

121:    Level: beginner
122: @*/
123: PetscErrorCode PEPSetFromOptions(PEP pep)
124: {
126:   char           type[256];
127:   PetscBool      set,flg,flg1,flg2,flg3;
128:   PetscReal      r,t;
129:   PetscScalar    s;
130:   PetscInt       i,j,k;
131:   PetscDrawLG    lg;

135:   PEPRegisterAll();
136:   PetscObjectOptionsBegin((PetscObject)pep);
137:     PetscOptionsFList("-pep_type","Polynomial Eigenvalue Problem method","PEPSetType",PEPList,(char*)(((PetscObject)pep)->type_name?((PetscObject)pep)->type_name:PEPTOAR),type,256,&flg);
138:     if (flg) {
139:       PEPSetType(pep,type);
140:     } else if (!((PetscObject)pep)->type_name) {
141:       PEPSetType(pep,PEPTOAR);
142:     }

144:     PetscOptionsBoolGroupBegin("-pep_general","general polynomial eigenvalue problem","PEPSetProblemType",&flg);
145:     if (flg) { PEPSetProblemType(pep,PEP_GENERAL); }
146:     PetscOptionsBoolGroup("-pep_hermitian","hermitian polynomial eigenvalue problem","PEPSetProblemType",&flg);
147:     if (flg) { PEPSetProblemType(pep,PEP_HERMITIAN); }
148:     PetscOptionsBoolGroupEnd("-pep_gyroscopic","gyroscopic polynomial eigenvalue problem","PEPSetProblemType",&flg);
149:     if (flg) { PEPSetProblemType(pep,PEP_GYROSCOPIC); }

151:     PetscOptionsEnum("-pep_scale","Scaling strategy","PEPSetScale",PEPScaleTypes,(PetscEnum)pep->scale,(PetscEnum*)&pep->scale,NULL);

153:     r = pep->sfactor;
154:     PetscOptionsReal("-pep_scale_factor","Scale factor","PEPSetScale",pep->sfactor,&r,&flg1);
155:     j = pep->sits;
156:     PetscOptionsInt("-pep_scale_its","Number of iterations in diagonal scaling","PEPSetScale",pep->sits,&j,&flg2);
157:     t = pep->slambda;
158:     PetscOptionsReal("-pep_scale_lambda","Estimate of eigenvalue (modulus) for diagonal scaling","PEPSetScale",pep->slambda,&t,&flg3);
159:     if (flg1 || flg2 || flg3) {
160:       PEPSetScale(pep,pep->scale,r,NULL,NULL,j,t);
161:     }

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

165:     PetscOptionsEnum("-pep_refine","Iterative refinement method","PEPSetRefine",PEPRefineTypes,(PetscEnum)pep->refine,(PetscEnum*)&pep->refine,NULL);

167:     i = pep->npart;
168:     PetscOptionsInt("-pep_refine_partitions","Number of partitions of the communicator for iterative refinement","PEPSetRefine",pep->npart,&i,&flg1);
169:     r = pep->rtol;
170:     PetscOptionsReal("-pep_refine_tol","Tolerance for iterative refinement","PEPSetRefine",pep->rtol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/1000:pep->rtol,&r,&flg2);
171:     j = pep->rits;
172:     PetscOptionsInt("-pep_refine_its","Maximum number of iterations for iterative refinement","PEPSetRefine",pep->rits,&j,&flg3);
173:     if (flg1 || flg2 || flg3) {
174:       PEPSetRefine(pep,pep->refine,i,r,j,pep->scheme);
175:     }

177:     PetscOptionsEnum("-pep_refine_scheme","Scheme used for linear systems within iterative refinement","PEPSetRefine",PEPRefineSchemes,(PetscEnum)pep->scheme,(PetscEnum*)&pep->scheme,NULL);

179:     i = pep->max_it? pep->max_it: PETSC_DEFAULT;
180:     PetscOptionsInt("-pep_max_it","Maximum number of iterations","PEPSetTolerances",pep->max_it,&i,&flg1);
181:     r = pep->tol;
182:     PetscOptionsReal("-pep_tol","Tolerance","PEPSetTolerances",pep->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:pep->tol,&r,&flg2);
183:     if (flg1 || flg2) {
184:       PEPSetTolerances(pep,r,i);
185:     }

187:     PetscOptionsBoolGroupBegin("-pep_conv_rel","Relative error convergence test","PEPSetConvergenceTest",&flg);
188:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_REL); }
189:     PetscOptionsBoolGroupBegin("-pep_conv_norm","Convergence test relative to the matrix norms","PEPSetConvergenceTest",&flg);
190:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_NORM); }
191:     PetscOptionsBoolGroup("-pep_conv_abs","Absolute error convergence test","PEPSetConvergenceTest",&flg);
192:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_ABS); }
193:     PetscOptionsBoolGroupEnd("-pep_conv_user","User-defined convergence test","PEPSetConvergenceTest",&flg);
194:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_USER); }

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

201:     i = pep->nev;
202:     PetscOptionsInt("-pep_nev","Number of eigenvalues to compute","PEPSetDimensions",pep->nev,&i,&flg1);
203:     j = pep->ncv? pep->ncv: PETSC_DEFAULT;
204:     PetscOptionsInt("-pep_ncv","Number of basis vectors","PEPSetDimensions",pep->ncv,&j,&flg2);
205:     k = pep->mpd? pep->mpd: PETSC_DEFAULT;
206:     PetscOptionsInt("-pep_mpd","Maximum dimension of projected problem","PEPSetDimensions",pep->mpd,&k,&flg3);
207:     if (flg1 || flg2 || flg3) {
208:       PEPSetDimensions(pep,i,j,k);
209:     }

211:     PetscOptionsScalar("-pep_target","Value of the target","PEPSetTarget",pep->target,&s,&flg);
212:     if (flg) {
213:       PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE);
214:       PEPSetTarget(pep,s);
215:     }

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

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

249:     PetscOptionsBoolGroupBegin("-pep_largest_magnitude","compute largest eigenvalues in magnitude","PEPSetWhichEigenpairs",&flg);
250:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_LARGEST_MAGNITUDE); }
251:     PetscOptionsBoolGroup("-pep_smallest_magnitude","compute smallest eigenvalues in magnitude","PEPSetWhichEigenpairs",&flg);
252:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_SMALLEST_MAGNITUDE); }
253:     PetscOptionsBoolGroup("-pep_largest_real","compute largest real parts","PEPSetWhichEigenpairs",&flg);
254:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_LARGEST_REAL); }
255:     PetscOptionsBoolGroup("-pep_smallest_real","compute smallest real parts","PEPSetWhichEigenpairs",&flg);
256:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_SMALLEST_REAL); }
257:     PetscOptionsBoolGroup("-pep_largest_imaginary","compute largest imaginary parts","PEPSetWhichEigenpairs",&flg);
258:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_LARGEST_IMAGINARY); }
259:     PetscOptionsBoolGroup("-pep_smallest_imaginary","compute smallest imaginary parts","PEPSetWhichEigenpairs",&flg);
260:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_SMALLEST_IMAGINARY); }
261:     PetscOptionsBoolGroup("-pep_target_magnitude","compute nearest eigenvalues to target","PEPSetWhichEigenpairs",&flg);
262:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE); }
263:     PetscOptionsBoolGroup("-pep_target_real","compute eigenvalues with real parts close to target","PEPSetWhichEigenpairs",&flg);
264:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_TARGET_REAL); }
265:     PetscOptionsBoolGroupEnd("-pep_target_imaginary","compute eigenvalues with imaginary parts close to target","PEPSetWhichEigenpairs",&flg);
266:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_TARGET_IMAGINARY); }

268:     PetscOptionsName("-pep_view","Print detailed information on solver used","PEPView",NULL);
269:     PetscOptionsName("-pep_view_vectors","View computed eigenvectors","PEPVectorsView",NULL);
270:     PetscOptionsName("-pep_view_values","View computed eigenvalues","PEPValuesView",NULL);
271:     PetscOptionsName("-pep_converged_reason","Print reason for convergence, and number of iterations","PEPReasonView",NULL);
272:     PetscOptionsName("-pep_error_absolute","Print absolute errors of each eigenpair","PEPErrorView",NULL);
273:     PetscOptionsName("-pep_error_relative","Print relative errors of each eigenpair","PEPErrorView",NULL);
274:     PetscOptionsName("-pep_error_backward","Print backward errors of each eigenpair","PEPErrorView",NULL);

276:     if (pep->ops->setfromoptions) {
277:       (*pep->ops->setfromoptions)(PetscOptionsObject,pep);
278:     }
279:     PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)pep);
280:   PetscOptionsEnd();

282:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
283:   BVSetFromOptions(pep->V);
284:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
285:   RGSetFromOptions(pep->rg);
286:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
287:   DSSetFromOptions(pep->ds);
288:   if (!pep->st) { PEPGetST(pep,&pep->st); }
289:   STSetFromOptions(pep->st);
290:   if (!pep->refineksp) { PEPRefineGetKSP(pep,&pep->refineksp); }
291:   KSPSetFromOptions(pep->refineksp);
292:   return(0);
293: }

297: /*@
298:    PEPGetTolerances - Gets the tolerance and maximum iteration count used
299:    by the PEP convergence tests.

301:    Not Collective

303:    Input Parameter:
304: .  pep - the polynomial eigensolver context

306:    Output Parameters:
307: +  tol - the convergence tolerance
308: -  maxits - maximum number of iterations

310:    Notes:
311:    The user can specify NULL for any parameter that is not needed.

313:    Level: intermediate

315: .seealso: PEPSetTolerances()
316: @*/
317: PetscErrorCode PEPGetTolerances(PEP pep,PetscReal *tol,PetscInt *maxits)
318: {
321:   if (tol)    *tol    = pep->tol;
322:   if (maxits) *maxits = pep->max_it;
323:   return(0);
324: }

328: /*@
329:    PEPSetTolerances - Sets the tolerance and maximum iteration count used
330:    by the PEP convergence tests.

332:    Logically Collective on PEP

334:    Input Parameters:
335: +  pep - the polynomial eigensolver context
336: .  tol - the convergence tolerance
337: -  maxits - maximum number of iterations to use

339:    Options Database Keys:
340: +  -pep_tol <tol> - Sets the convergence tolerance
341: -  -pep_max_it <maxits> - Sets the maximum number of iterations allowed

343:    Notes:
344:    Use PETSC_DEFAULT for either argument to assign a reasonably good value.

346:    Level: intermediate

348: .seealso: PEPGetTolerances()
349: @*/
350: PetscErrorCode PEPSetTolerances(PEP pep,PetscReal tol,PetscInt maxits)
351: {
356:   if (tol == PETSC_DEFAULT) {
357:     pep->tol   = PETSC_DEFAULT;
358:     pep->state = PEP_STATE_INITIAL;
359:   } else {
360:     if (tol <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
361:     pep->tol = tol;
362:   }
363:   if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
364:     pep->max_it = 0;
365:     pep->state  = PEP_STATE_INITIAL;
366:   } else {
367:     if (maxits <= 0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
368:     pep->max_it = maxits;
369:   }
370:   return(0);
371: }

375: /*@
376:    PEPGetDimensions - Gets the number of eigenvalues to compute
377:    and the dimension of the subspace.

379:    Not Collective

381:    Input Parameter:
382: .  pep - the polynomial eigensolver context

384:    Output Parameters:
385: +  nev - number of eigenvalues to compute
386: .  ncv - the maximum dimension of the subspace to be used by the solver
387: -  mpd - the maximum dimension allowed for the projected problem

389:    Notes:
390:    The user can specify NULL for any parameter that is not needed.

392:    Level: intermediate

394: .seealso: PEPSetDimensions()
395: @*/
396: PetscErrorCode PEPGetDimensions(PEP pep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
397: {
400:   if (nev) *nev = pep->nev;
401:   if (ncv) *ncv = pep->ncv;
402:   if (mpd) *mpd = pep->mpd;
403:   return(0);
404: }

408: /*@
409:    PEPSetDimensions - Sets the number of eigenvalues to compute
410:    and the dimension of the subspace.

412:    Logically Collective on PEP

414:    Input Parameters:
415: +  pep - the polynomial eigensolver context
416: .  nev - number of eigenvalues to compute
417: .  ncv - the maximum dimension of the subspace to be used by the solver
418: -  mpd - the maximum dimension allowed for the projected problem

420:    Options Database Keys:
421: +  -pep_nev <nev> - Sets the number of eigenvalues
422: .  -pep_ncv <ncv> - Sets the dimension of the subspace
423: -  -pep_mpd <mpd> - Sets the maximum projected dimension

425:    Notes:
426:    Use PETSC_DEFAULT for ncv and mpd to assign a reasonably good value, which is
427:    dependent on the solution method.

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

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

438:    Level: intermediate

440: .seealso: PEPGetDimensions()
441: @*/
442: PetscErrorCode PEPSetDimensions(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)
443: {
449:   if (nev<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
450:   pep->nev = nev;
451:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
452:     pep->ncv = 0;
453:   } else {
454:     if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
455:     pep->ncv = ncv;
456:   }
457:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
458:     pep->mpd = 0;
459:   } else {
460:     if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
461:     pep->mpd = mpd;
462:   }
463:   pep->state = PEP_STATE_INITIAL;
464:   return(0);
465: }

469: /*@
470:    PEPSetWhichEigenpairs - Specifies which portion of the spectrum is
471:    to be sought.

473:    Logically Collective on PEP

475:    Input Parameters:
476: +  pep   - eigensolver context obtained from PEPCreate()
477: -  which - the portion of the spectrum to be sought

479:    Possible values:
480:    The parameter 'which' can have one of these values

482: +     PEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
483: .     PEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
484: .     PEP_LARGEST_REAL - largest real parts
485: .     PEP_SMALLEST_REAL - smallest real parts
486: .     PEP_LARGEST_IMAGINARY - largest imaginary parts
487: .     PEP_SMALLEST_IMAGINARY - smallest imaginary parts
488: .     PEP_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
489: .     PEP_TARGET_REAL - eigenvalues with real part closest to target
490: .     PEP_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
491: -     PEP_WHICH_USER - user defined ordering set with PEPSetEigenvalueComparison()

493:    Options Database Keys:
494: +   -pep_largest_magnitude - Sets largest eigenvalues in magnitude
495: .   -pep_smallest_magnitude - Sets smallest eigenvalues in magnitude
496: .   -pep_largest_real - Sets largest real parts
497: .   -pep_smallest_real - Sets smallest real parts
498: .   -pep_largest_imaginary - Sets largest imaginary parts
499: .   -pep_smallest_imaginary - Sets smallest imaginary parts
500: .   -pep_target_magnitude - Sets eigenvalues closest to target
501: .   -pep_target_real - Sets real parts closest to target
502: -   -pep_target_imaginary - Sets imaginary parts closest to target

504:    Notes:
505:    Not all eigensolvers implemented in PEP account for all the possible values
506:    stated above. If SLEPc is compiled for real numbers PEP_LARGEST_IMAGINARY
507:    and PEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
508:    for eigenvalue selection.

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

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

515:    Level: intermediate

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

550: /*@
551:     PEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
552:     sought.

554:     Not Collective

556:     Input Parameter:
557: .   pep - eigensolver context obtained from PEPCreate()

559:     Output Parameter:
560: .   which - the portion of the spectrum to be sought

562:     Notes:
563:     See PEPSetWhichEigenpairs() for possible values of 'which'.

565:     Level: intermediate

567: .seealso: PEPSetWhichEigenpairs(), PEPWhich
568: @*/
569: PetscErrorCode PEPGetWhichEigenpairs(PEP pep,PEPWhich *which)
570: {
574:   *which = pep->which;
575:   return(0);
576: }

580: /*@C
581:    PEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
582:    when PEPSetWhichEigenpairs() is set to PEP_WHICH_USER.

584:    Logically Collective on PEP

586:    Input Parameters:
587: +  pep  - eigensolver context obtained from PEPCreate()
588: .  func - a pointer to the comparison function
589: -  ctx  - a context pointer (the last parameter to the comparison function)

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

594: +   ar     - real part of the 1st eigenvalue
595: .   ai     - imaginary part of the 1st eigenvalue
596: .   br     - real part of the 2nd eigenvalue
597: .   bi     - imaginary part of the 2nd eigenvalue
598: .   res    - result of comparison
599: -   ctx    - optional context, as set by PEPSetEigenvalueComparison()

601:    Note:
602:    The returning parameter 'res' can be
603: +  negative - if the 1st eigenvalue is preferred to the 2st one
604: .  zero     - if both eigenvalues are equally preferred
605: -  positive - if the 2st eigenvalue is preferred to the 1st one

607:    Level: advanced

609: .seealso: PEPSetWhichEigenpairs(), PEPWhich
610: @*/
611: PetscErrorCode PEPSetEigenvalueComparison(PEP pep,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)
612: {
615:   pep->sc->comparison    = func;
616:   pep->sc->comparisonctx = ctx;
617:   pep->which             = PEP_WHICH_USER;
618:   return(0);
619: }

623: /*@
624:    PEPSetProblemType - Specifies the type of the polynomial eigenvalue problem.

626:    Logically Collective on PEP

628:    Input Parameters:
629: +  pep  - the polynomial eigensolver context
630: -  type - a known type of polynomial eigenvalue problem

632:    Options Database Keys:
633: +  -pep_general - general problem with no particular structure
634: .  -pep_hermitian - problem whose coefficient matrices are Hermitian
635: -  -pep_gyroscopic - problem with Hamiltonian structure

637:    Notes:
638:    Allowed values for the problem type are: general (PEP_GENERAL), Hermitian
639:    (PEP_HERMITIAN), and gyroscopic (PEP_GYROSCOPIC).

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

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

648:    Level: intermediate

650: .seealso: PEPSetOperators(), PEPSetType(), PEPGetProblemType(), PEPProblemType
651: @*/
652: PetscErrorCode PEPSetProblemType(PEP pep,PEPProblemType type)
653: {
657:   if (type!=PEP_GENERAL && type!=PEP_HERMITIAN && type!=PEP_GYROSCOPIC) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
658:   if (type != pep->problem_type) {
659:     pep->problem_type = type;
660:     pep->state = PEP_STATE_INITIAL;
661:   }
662:   return(0);
663: }

667: /*@
668:    PEPGetProblemType - Gets the problem type from the PEP object.

670:    Not Collective

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

675:    Output Parameter:
676: .  type - name of PEP problem type

678:    Level: intermediate

680: .seealso: PEPSetProblemType(), PEPProblemType
681: @*/
682: PetscErrorCode PEPGetProblemType(PEP pep,PEPProblemType *type)
683: {
687:   *type = pep->problem_type;
688:   return(0);
689: }

693: /*@
694:    PEPSetBasis - Specifies the type of polynomial basis used to describe the
695:    polynomial eigenvalue problem.

697:    Logically Collective on PEP

699:    Input Parameters:
700: +  pep   - the polynomial eigensolver context
701: -  basis - the type of polynomial basis

703:    Options Database Key:
704: .  -pep_basis <basis> - Select the basis type

706:    Notes:
707:    By default, the coefficient matrices passed via PEPSetOperators() are
708:    expressed in the monomial basis, i.e. 
709:    P(lambda) = A_0 + lambda*A_1 + lambda^2*A_2 + ... + lambda^d*A_d.
710:    Other polynomial bases may have better numerical behaviour, but the user
711:    must then pass the coefficient matrices accordingly.

713:    Level: intermediate

715: .seealso: PEPSetOperators(), PEPGetBasis(), PEPBasis
716: @*/
717: PetscErrorCode PEPSetBasis(PEP pep,PEPBasis basis)
718: {
722:   pep->basis = basis;
723:   return(0);
724: }

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

731:    Not Collective

733:    Input Parameter:
734: .  pep - the polynomial eigensolver context

736:    Output Parameter:
737: .  basis - the polynomial basis

739:    Level: intermediate

741: .seealso: PEPSetBasis(), PEPBasis
742: @*/
743: PetscErrorCode PEPGetBasis(PEP pep,PEPBasis *basis)
744: {
748:   *basis = pep->basis;
749:   return(0);
750: }

754: /*@
755:    PEPSetTrackAll - Specifies if the solver must compute the residual of all
756:    approximate eigenpairs or not.

758:    Logically Collective on PEP

760:    Input Parameters:
761: +  pep      - the eigensolver context
762: -  trackall - whether compute all residuals or not

764:    Notes:
765:    If the user sets trackall=PETSC_TRUE then the solver explicitly computes
766:    the residual for each eigenpair approximation. Computing the residual is
767:    usually an expensive operation and solvers commonly compute the associated
768:    residual to the first unconverged eigenpair.
769:    The options '-pep_monitor_all' and '-pep_monitor_lg_all' automatically
770:    activate this option.

772:    Level: developer

774: .seealso: PEPGetTrackAll()
775: @*/
776: PetscErrorCode PEPSetTrackAll(PEP pep,PetscBool trackall)
777: {
781:   pep->trackall = trackall;
782:   return(0);
783: }

787: /*@
788:    PEPGetTrackAll - Returns the flag indicating whether all residual norms must
789:    be computed or not.

791:    Not Collective

793:    Input Parameter:
794: .  pep - the eigensolver context

796:    Output Parameter:
797: .  trackall - the returned flag

799:    Level: developer

801: .seealso: PEPSetTrackAll()
802: @*/
803: PetscErrorCode PEPGetTrackAll(PEP pep,PetscBool *trackall)
804: {
808:   *trackall = pep->trackall;
809:   return(0);
810: }

814: /*@C
815:    PEPSetConvergenceTestFunction - Sets a function to compute the error estimate
816:    used in the convergence test.

818:    Logically Collective on PEP

820:    Input Parameters:
821: +  pep     - eigensolver context obtained from PEPCreate()
822: .  func    - a pointer to the convergence test function
823: .  ctx     - context for private data for the convergence routine (may be null)
824: -  destroy - a routine for destroying the context (may be null)

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

829: +   pep    - eigensolver context obtained from PEPCreate()
830: .   eigr   - real part of the eigenvalue
831: .   eigi   - imaginary part of the eigenvalue
832: .   res    - residual norm associated to the eigenpair
833: .   errest - (output) computed error estimate
834: -   ctx    - optional context, as set by PEPSetConvergenceTestFunction()

836:    Note:
837:    If the error estimate returned by the convergence test function is less than
838:    the tolerance, then the eigenvalue is accepted as converged.

840:    Level: advanced

842: .seealso: PEPSetConvergenceTest(), PEPSetTolerances()
843: @*/
844: PetscErrorCode PEPSetConvergenceTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
845: {

850:   if (pep->convergeddestroy) {
851:     (*pep->convergeddestroy)(pep->convergedctx);
852:   }
853:   pep->converged        = func;
854:   pep->convergeddestroy = destroy;
855:   pep->convergedctx     = ctx;
856:   if (func == PEPConvergedRelative) pep->conv = PEP_CONV_REL;
857:   else if (func == PEPConvergedNorm) pep->conv = PEP_CONV_NORM;
858:   else if (func == PEPConvergedAbsolute) pep->conv = PEP_CONV_ABS;
859:   else pep->conv = PEP_CONV_USER;
860:   return(0);
861: }

865: /*@
866:    PEPSetConvergenceTest - Specifies how to compute the error estimate
867:    used in the convergence test.

869:    Logically Collective on PEP

871:    Input Parameters:
872: +  pep  - eigensolver context obtained from PEPCreate()
873: -  conv - the type of convergence test

875:    Options Database Keys:
876: +  -pep_conv_abs    - Sets the absolute convergence test
877: .  -pep_conv_rel    - Sets the convergence test relative to the eigenvalue
878: .  -pep_conv_norm   - Sets the convergence test relative to the matrix norms
879: -  -pep_conv_user   - Selects the user-defined convergence test

881:    Note:
882:    The parameter 'conv' can have one of these values
883: +     PEP_CONV_ABS    - absolute error ||r||
884: .     PEP_CONV_REL    - error relative to the eigenvalue l, ||r||/|l|
885: .     PEP_CONV_NORM   - error relative matrix norms, ||r||/sum_i(l^i*||A_i||)
886: -     PEP_CONV_USER   - function set by PEPSetConvergenceTestFunction()

888:    Level: intermediate

890: .seealso: PEPGetConvergenceTest(), PEPSetConvergenceTestFunction(), PEPSetStoppingTest(), PEPConv
891: @*/
892: PetscErrorCode PEPSetConvergenceTest(PEP pep,PEPConv conv)
893: {
897:   switch (conv) {
898:     case PEP_CONV_ABS:    pep->converged = PEPConvergedAbsolute; break;
899:     case PEP_CONV_REL:    pep->converged = PEPConvergedRelative; break;
900:     case PEP_CONV_NORM:   pep->converged = PEPConvergedNorm; break;
901:     case PEP_CONV_USER: break;
902:     default:
903:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
904:   }
905:   pep->conv = conv;
906:   return(0);
907: }

911: /*@
912:    PEPGetConvergenceTest - Gets the method used to compute the error estimate
913:    used in the convergence test.

915:    Not Collective

917:    Input Parameters:
918: .  pep   - eigensolver context obtained from PEPCreate()

920:    Output Parameters:
921: .  conv  - the type of convergence test

923:    Level: intermediate

925: .seealso: PEPSetConvergenceTest(), PEPConv
926: @*/
927: PetscErrorCode PEPGetConvergenceTest(PEP pep,PEPConv *conv)
928: {
932:   *conv = pep->conv;
933:   return(0);
934: }

938: /*@C
939:    PEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
940:    iteration of the eigensolver.

942:    Logically Collective on PEP

944:    Input Parameters:
945: +  pep     - eigensolver context obtained from PEPCreate()
946: .  func    - pointer to the stopping test function
947: .  ctx     - context for private data for the stopping routine (may be null)
948: -  destroy - a routine for destroying the context (may be null)

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

953: +   pep    - eigensolver context obtained from PEPCreate()
954: .   its    - current number of iterations
955: .   max_it - maximum number of iterations
956: .   nconv  - number of currently converged eigenpairs
957: .   nev    - number of requested eigenpairs
958: .   reason - (output) result of the stopping test
959: -   ctx    - optional context, as set by PEPSetStoppingTestFunction()

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

967:    Level: advanced

969: .seealso: PEPSetStoppingTest(), PEPStoppingBasic()
970: @*/
971: PetscErrorCode PEPSetStoppingTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscInt,PetscInt,PetscInt,PetscInt,PEPConvergedReason*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
972: {

977:   if (pep->stoppingdestroy) {
978:     (*pep->stoppingdestroy)(pep->stoppingctx);
979:   }
980:   pep->stopping        = func;
981:   pep->stoppingdestroy = destroy;
982:   pep->stoppingctx     = ctx;
983:   if (func == PEPStoppingBasic) pep->stop = PEP_STOP_BASIC;
984:   else pep->stop = PEP_STOP_USER;
985:   return(0);
986: }

990: /*@
991:    PEPSetStoppingTest - Specifies how to decide the termination of the outer
992:    loop of the eigensolver.

994:    Logically Collective on PEP

996:    Input Parameters:
997: +  pep  - eigensolver context obtained from PEPCreate()
998: -  stop - the type of stopping test

1000:    Options Database Keys:
1001: +  -pep_stop_basic - Sets the default stopping test
1002: -  -pep_stop_user  - Selects the user-defined stopping test

1004:    Note:
1005:    The parameter 'stop' can have one of these values
1006: +     PEP_STOP_BASIC - default stopping test
1007: -     PEP_STOP_USER  - function set by PEPSetStoppingTestFunction()

1009:    Level: advanced

1011: .seealso: PEPGetStoppingTest(), PEPSetStoppingTestFunction(), PEPSetConvergenceTest(), PEPStop
1012: @*/
1013: PetscErrorCode PEPSetStoppingTest(PEP pep,PEPStop stop)
1014: {
1018:   switch (stop) {
1019:     case PEP_STOP_BASIC: pep->stopping = PEPStoppingBasic; break;
1020:     case PEP_STOP_USER:  break;
1021:     default:
1022:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
1023:   }
1024:   pep->stop = stop;
1025:   return(0);
1026: }

1030: /*@
1031:    PEPGetStoppingTest - Gets the method used to decide the termination of the outer
1032:    loop of the eigensolver.

1034:    Not Collective

1036:    Input Parameters:
1037: .  pep   - eigensolver context obtained from PEPCreate()

1039:    Output Parameters:
1040: .  stop  - the type of stopping test

1042:    Level: advanced

1044: .seealso: PEPSetStoppingTest(), PEPStop
1045: @*/
1046: PetscErrorCode PEPGetStoppingTest(PEP pep,PEPStop *stop)
1047: {
1051:   *stop = pep->stop;
1052:   return(0);
1053: }

1057: /*@
1058:    PEPSetScale - Specifies the scaling strategy to be used.

1060:    Logically Collective on PEP

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

1071:    Options Database Keys:
1072: +  -pep_scale <type> - scaling type, one of <none,scalar,diagonal,both>
1073: .  -pep_scale_factor <alpha> - the scaling factor
1074: .  -pep_scale_its <its> - number of iterations
1075: -  -pep_scale_lambda <lambda> - approximation to eigenvalues

1077:    Notes:
1078:    There are two non-exclusive scaling strategies: scalar and diagonal.

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

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

1096:    Level: intermediate

1098: .seealso: PEPGetScale()
1099: @*/
1100: PetscErrorCode PEPSetScale(PEP pep,PEPScale scale,PetscReal alpha,Vec Dl,Vec Dr,PetscInt its,PetscReal lambda)
1101: {

1107:   pep->scale = scale;
1108:   if (scale==PEP_SCALE_SCALAR || scale==PEP_SCALE_BOTH) {
1110:     if (alpha == PETSC_DEFAULT || alpha == PETSC_DECIDE) {
1111:       pep->sfactor = 0.0;
1112:       pep->sfactor_set = PETSC_FALSE;
1113:     } else {
1114:       if (alpha<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of alpha. Must be > 0");
1115:       pep->sfactor = alpha;
1116:       pep->sfactor_set = PETSC_TRUE;
1117:     }
1118:   }
1119:   if (scale==PEP_SCALE_DIAGONAL || scale==PEP_SCALE_BOTH) {
1120:     if (Dl) {
1123:       PetscObjectReference((PetscObject)Dl);
1124:       VecDestroy(&pep->Dl);
1125:       pep->Dl = Dl;
1126:     }
1127:     if (Dr) {
1130:       PetscObjectReference((PetscObject)Dr);
1131:       VecDestroy(&pep->Dr);
1132:       pep->Dr = Dr;
1133:     }
1136:     if (its==PETSC_DECIDE || its==PETSC_DEFAULT) pep->sits = 5;
1137:     else pep->sits = its;
1138:     if (lambda==PETSC_DECIDE || lambda==PETSC_DEFAULT) pep->slambda = 1.0;
1139:     else if (lambda<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of lambda. Must be > 0");
1140:     else pep->slambda = lambda;
1141:   }
1142:   pep->state = PEP_STATE_INITIAL;
1143:   return(0);
1144: }

1148: /*@
1149:    PEPGetScale - Gets the scaling strategy used by the PEP object, and the
1150:    associated parameters.

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

1154:    Input Parameter:
1155: .  pep - the eigensolver context

1157:    Output Parameters:
1158: +  scale  - scaling strategy
1159: .  alpha  - the scaling factor used in the scalar strategy
1160: .  Dl     - the left diagonal matrix of the diagonal scaling algorithm
1161: .  Dr     - the right diagonal matrix of the diagonal scaling algorithm
1162: .  its    - number of iterations of the diagonal scaling algorithm
1163: -  lambda - approximation to wanted eigenvalues (modulus)

1165:    Level: intermediate

1167:    Note:
1168:    The user can specify NULL for any parameter that is not needed.

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

1173: .seealso: PEPSetScale(), PEPSetUp()
1174: @*/
1175: PetscErrorCode PEPGetScale(PEP pep,PEPScale *scale,PetscReal *alpha,Vec *Dl,Vec *Dr,PetscInt *its,PetscReal *lambda)
1176: {
1179:   if (scale)  *scale  = pep->scale;
1180:   if (alpha)  *alpha  = pep->sfactor;
1181:   if (Dl)     *Dl     = pep->Dl;
1182:   if (Dr)     *Dr     = pep->Dr;
1183:   if (its)    *its    = pep->sits;
1184:   if (lambda) *lambda = pep->slambda;
1185:   return(0);
1186: }

1190: /*@
1191:    PEPSetExtract - Specifies the extraction strategy to be used.

1193:    Logically Collective on PEP

1195:    Input Parameters:
1196: +  pep     - the eigensolver context
1197: -  extract - extraction strategy

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

1202:    Level: intermediate

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

1217: /*@
1218:    PEPGetExtract - Gets the extraction strategy used by the PEP object.

1220:    Not Collective

1222:    Input Parameter:
1223: .  pep - the eigensolver context

1225:    Output Parameter:
1226: .  extract - extraction strategy

1228:    Level: intermediate

1230: .seealso: PEPSetExtract()
1231: @*/
1232: PetscErrorCode PEPGetExtract(PEP pep,PEPExtract *extract)
1233: {
1236:   if (extract) *extract = pep->extract;
1237:   return(0);
1238: }

1242: /*@
1243:    PEPSetRefine - Specifies the refinement type (and options) to be used
1244:    after the solve.

1246:    Logically Collective on PEP

1248:    Input Parameters:
1249: +  pep    - the polynomial eigensolver context
1250: .  refine - refinement type
1251: .  npart  - number of partitions of the communicator
1252: .  tol    - the convergence tolerance
1253: .  its    - maximum number of refinement iterations
1254: -  scheme - which scheme to be used for solving the involved linear systems

1256:    Options Database Keys:
1257: +  -pep_refine <type> - refinement type, one of <none,simple,multiple>
1258: .  -pep_refine_partitions <n> - the number of partitions
1259: .  -pep_refine_tol <tol> - the tolerance
1260: .  -pep_refine_its <its> - number of iterations
1261: -  -pep_refine_scheme - to set the scheme for the linear solves

1263:    Notes:
1264:    By default, iterative refinement is disabled, since it may be very
1265:    costly. There are two possible refinement strategies: simple and multiple.
1266:    The simple approach performs iterative refinement on each of the
1267:    converged eigenpairs individually, whereas the multiple strategy works
1268:    with the invariant pair as a whole, refining all eigenpairs simultaneously.
1269:    The latter may be required for the case of multiple eigenvalues.

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

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

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

1286:    Level: intermediate

1288: .seealso: PEPGetRefine()
1289: @*/
1290: PetscErrorCode PEPSetRefine(PEP pep,PEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,PEPRefineScheme scheme)
1291: {
1293:   PetscMPIInt    size;

1302:   pep->refine = refine;
1303:   if (refine) {  /* process parameters only if not REFINE_NONE */
1304:     if (npart!=pep->npart) {
1305:       PetscSubcommDestroy(&pep->refinesubc);
1306:       KSPDestroy(&pep->refineksp);
1307:     }
1308:     if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
1309:       pep->npart = 1;
1310:     } else {
1311:       MPI_Comm_size(PetscObjectComm((PetscObject)pep),&size);
1312:       if (npart<1 || npart>size) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
1313:       pep->npart = npart;
1314:     }
1315:     if (tol == PETSC_DEFAULT || tol == PETSC_DECIDE) {
1316:       pep->rtol = PetscMax(pep->tol/1000,PETSC_MACHINE_EPSILON);
1317:     } else {
1318:       if (tol<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1319:       pep->rtol = tol;
1320:     }
1321:     if (its==PETSC_DECIDE || its==PETSC_DEFAULT) {
1322:       pep->rits = PETSC_DEFAULT;
1323:     } else {
1324:       if (its<0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be >= 0");
1325:       pep->rits = its;
1326:     }
1327:     pep->scheme = scheme;
1328:   }
1329:   pep->state = PEP_STATE_INITIAL;
1330:   return(0);
1331: }

1335: /*@
1336:    PEPGetRefine - Gets the refinement strategy used by the PEP object, and the
1337:    associated parameters.

1339:    Not Collective

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

1344:    Output Parameters:
1345: +  refine - refinement type
1346: .  npart  - number of partitions of the communicator
1347: .  tol    - the convergence tolerance
1348: .  its    - maximum number of refinement iterations
1349: -  scheme - the scheme used for solving linear systems

1351:    Level: intermediate

1353:    Note:
1354:    The user can specify NULL for any parameter that is not needed.

1356: .seealso: PEPSetRefine()
1357: @*/
1358: PetscErrorCode PEPGetRefine(PEP pep,PEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,PEPRefineScheme *scheme)
1359: {
1362:   if (refine) *refine = pep->refine;
1363:   if (npart)  *npart  = pep->npart;
1364:   if (tol)    *tol    = pep->rtol;
1365:   if (its)    *its    = pep->rits;
1366:   if (scheme) *scheme = pep->scheme;
1367:   return(0);
1368: }

1372: /*@C
1373:    PEPSetOptionsPrefix - Sets the prefix used for searching for all
1374:    PEP options in the database.

1376:    Logically Collective on PEP

1378:    Input Parameters:
1379: +  pep - the polynomial eigensolver context
1380: -  prefix - the prefix string to prepend to all PEP option requests

1382:    Notes:
1383:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1384:    The first character of all runtime options is AUTOMATICALLY the
1385:    hyphen.

1387:    For example, to distinguish between the runtime options for two
1388:    different PEP contexts, one could call
1389: .vb
1390:       PEPSetOptionsPrefix(pep1,"qeig1_")
1391:       PEPSetOptionsPrefix(pep2,"qeig2_")
1392: .ve

1394:    Level: advanced

1396: .seealso: PEPAppendOptionsPrefix(), PEPGetOptionsPrefix()
1397: @*/
1398: PetscErrorCode PEPSetOptionsPrefix(PEP pep,const char *prefix)
1399: {

1404:   if (!pep->st) { PEPGetST(pep,&pep->st); }
1405:   STSetOptionsPrefix(pep->st,prefix);
1406:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
1407:   BVSetOptionsPrefix(pep->V,prefix);
1408:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
1409:   DSSetOptionsPrefix(pep->ds,prefix);
1410:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
1411:   RGSetOptionsPrefix(pep->rg,prefix);
1412:   PetscObjectSetOptionsPrefix((PetscObject)pep,prefix);
1413:   return(0);
1414: }

1418: /*@C
1419:    PEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1420:    PEP options in the database.

1422:    Logically Collective on PEP

1424:    Input Parameters:
1425: +  pep - the polynomial eigensolver context
1426: -  prefix - the prefix string to prepend to all PEP option requests

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

1432:    Level: advanced

1434: .seealso: PEPSetOptionsPrefix(), PEPGetOptionsPrefix()
1435: @*/
1436: PetscErrorCode PEPAppendOptionsPrefix(PEP pep,const char *prefix)
1437: {
1439:   PetscBool      flg;
1440:   EPS            eps;

1444:   if (!pep->st) { PEPGetST(pep,&pep->st); }
1445:   STAppendOptionsPrefix(pep->st,prefix);
1446:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
1447:   BVSetOptionsPrefix(pep->V,prefix);
1448:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
1449:   DSSetOptionsPrefix(pep->ds,prefix);
1450:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
1451:   RGSetOptionsPrefix(pep->rg,prefix);
1452:   PetscObjectAppendOptionsPrefix((PetscObject)pep,prefix);
1453:   PetscObjectTypeCompare((PetscObject)pep,PEPLINEAR,&flg);
1454:   if (flg) {
1455:     PEPLinearGetEPS(pep,&eps);
1456:     EPSSetOptionsPrefix(eps,((PetscObject)pep)->prefix);
1457:     EPSAppendOptionsPrefix(eps,"pep_");
1458:   }
1459:   return(0);
1460: }

1464: /*@C
1465:    PEPGetOptionsPrefix - Gets the prefix used for searching for all
1466:    PEP options in the database.

1468:    Not Collective

1470:    Input Parameters:
1471: .  pep - the polynomial eigensolver context

1473:    Output Parameters:
1474: .  prefix - pointer to the prefix string used is returned

1476:    Note:
1477:    On the Fortran side, the user should pass in a string 'prefix' of
1478:    sufficient length to hold the prefix.

1480:    Level: advanced

1482: .seealso: PEPSetOptionsPrefix(), PEPAppendOptionsPrefix()
1483: @*/
1484: PetscErrorCode PEPGetOptionsPrefix(PEP pep,const char *prefix[])
1485: {

1491:   PetscObjectGetOptionsPrefix((PetscObject)pep,prefix);
1492:   return(0);
1493: }