Actual source code: pepopts.c

slepc-3.6.1 2015-09-03
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-2015, 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*/

 29: /*@
 30:    PEPSetFromOptions - Sets PEP options from the options database.
 31:    This routine must be called before PEPSetUp() if the user is to be
 32:    allowed to set the solver type.

 34:    Collective on PEP

 36:    Input Parameters:
 37: .  pep - the polynomial eigensolver context

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

 42:    Level: beginner
 43: @*/
 44: PetscErrorCode PEPSetFromOptions(PEP pep)
 45: {
 46:   PetscErrorCode   ierr;
 47:   char             type[256],monfilename[PETSC_MAX_PATH_LEN];
 48:   PetscBool        flg,flg1,flg2,flg3,flg4;
 49:   PetscReal        r,t;
 50:   PetscScalar      s;
 51:   PetscInt         i,j,k;
 52:   PetscViewer      monviewer;
 53:   SlepcConvMonitor ctx;

 57:   PEPRegisterAll();
 58:   PetscObjectOptionsBegin((PetscObject)pep);
 59:     PetscOptionsFList("-pep_type","Polynomial Eigenvalue Problem method","PEPSetType",PEPList,(char*)(((PetscObject)pep)->type_name?((PetscObject)pep)->type_name:PEPTOAR),type,256,&flg);
 60:     if (flg) {
 61:       PEPSetType(pep,type);
 62:     } else if (!((PetscObject)pep)->type_name) {
 63:       PEPSetType(pep,PEPTOAR);
 64:     }

 66:     PetscOptionsBoolGroupBegin("-pep_general","general polynomial eigenvalue problem","PEPSetProblemType",&flg);
 67:     if (flg) { PEPSetProblemType(pep,PEP_GENERAL); }
 68:     PetscOptionsBoolGroup("-pep_hermitian","hermitian polynomial eigenvalue problem","PEPSetProblemType",&flg);
 69:     if (flg) { PEPSetProblemType(pep,PEP_HERMITIAN); }
 70:     PetscOptionsBoolGroupEnd("-pep_gyroscopic","gyroscopic polynomial eigenvalue problem","PEPSetProblemType",&flg);
 71:     if (flg) { PEPSetProblemType(pep,PEP_GYROSCOPIC); }

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

 75:     r = pep->sfactor;
 76:     PetscOptionsReal("-pep_scale_factor","Scale factor","PEPSetScale",pep->sfactor,&r,&flg1);
 77:     j = pep->sits;
 78:     PetscOptionsInt("-pep_scale_its","Number of iterations in diagonal scaling","PEPSetScale",pep->sits,&j,&flg2);
 79:     t = pep->slambda;
 80:     PetscOptionsReal("-pep_scale_lambda","Estimate of eigenvalue (modulus) for diagonal scaling","PEPSetScale",pep->slambda,&t,&flg3);
 81:     if (flg1 || flg2 || flg3) {
 82:       PEPSetScale(pep,pep->scale,r,NULL,NULL,j,t);
 83:     }

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

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

 89:     i = pep->npart;
 90:     PetscOptionsInt("-pep_refine_partitions","Number of partitions of the communicator for iterative refinement","PEPSetRefine",pep->npart,&i,&flg1);
 91:     r = pep->rtol;
 92:     PetscOptionsReal("-pep_refine_tol","Tolerance for iterative refinement","PEPSetRefine",pep->rtol,&r,&flg2);
 93:     j = pep->rits;
 94:     PetscOptionsInt("-pep_refine_its","Maximum number of iterations for iterative refinement","PEPSetRefine",pep->rits,&j,&flg3);
 95:     flg = pep->schur;
 96:     PetscOptionsBool("-pep_refine_schur","Use Schur complement for iterative refinement","PEPSetRefine",pep->schur,&flg,&flg4);
 97:     if (flg1 || flg2 || flg3 || flg4) {
 98:       PEPSetRefine(pep,pep->refine,i,r,j,flg);
 99:     }

101:     i = pep->max_it? pep->max_it: PETSC_DEFAULT;
102:     PetscOptionsInt("-pep_max_it","Maximum number of iterations","PEPSetTolerances",pep->max_it,&i,&flg1);
103:     r = pep->tol;
104:     PetscOptionsReal("-pep_tol","Tolerance","PEPSetTolerances",pep->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:pep->tol,&r,&flg2);
105:     if (flg1 || flg2) {
106:       PEPSetTolerances(pep,r,i);
107:     }

109:     PetscOptionsBoolGroupBegin("-pep_conv_eig","Relative error convergence test","PEPSetConvergenceTest",&flg);
110:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_EIG); }
111:     PetscOptionsBoolGroup("-pep_conv_linear","Convergence test related to the linearized eigenproblem","PEPSetConvergenceTest",&flg);
112:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_LINEAR); }
113:     PetscOptionsBoolGroupBegin("-pep_conv_norm","Convergence test related to the matrix norms","PEPSetConvergenceTest",&flg);
114:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_NORM); }
115:     PetscOptionsBoolGroup("-pep_conv_abs","Absolute error convergence test","PEPSetConvergenceTest",&flg);
116:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_ABS); }
117:     PetscOptionsBoolGroupEnd("-pep_conv_user","User-defined convergence test","PEPSetConvergenceTest",&flg);
118:     if (flg) { PEPSetConvergenceTest(pep,PEP_CONV_USER); }

120:     i = pep->nev;
121:     PetscOptionsInt("-pep_nev","Number of eigenvalues to compute","PEPSetDimensions",pep->nev,&i,&flg1);
122:     j = pep->ncv? pep->ncv: PETSC_DEFAULT;
123:     PetscOptionsInt("-pep_ncv","Number of basis vectors","PEPSetDimensions",pep->ncv,&j,&flg2);
124:     k = pep->mpd? pep->mpd: PETSC_DEFAULT;
125:     PetscOptionsInt("-pep_mpd","Maximum dimension of projected problem","PEPSetDimensions",pep->mpd,&k,&flg3);
126:     if (flg1 || flg2 || flg3) {
127:       PEPSetDimensions(pep,i,j,k);
128:     }

130:     PetscOptionsScalar("-pep_target","Value of the target","PEPSetTarget",pep->target,&s,&flg);
131:     if (flg) {
132:       PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE);
133:       PEPSetTarget(pep,s);
134:     }

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

138:     /* -----------------------------------------------------------------------*/
139:     /*
140:       Cancels all monitors hardwired into code before call to PEPSetFromOptions()
141:     */
142:     flg  = PETSC_FALSE;
143:     PetscOptionsBool("-pep_monitor_cancel","Remove any hardwired monitor routines","PEPMonitorCancel",flg,&flg,NULL);
144:     if (flg) {
145:       PEPMonitorCancel(pep);
146:     }
147:     /*
148:       Prints approximate eigenvalues and error estimates at each iteration
149:     */
150:     PetscOptionsString("-pep_monitor","Monitor first unconverged approximate eigenvalue and error estimate","PEPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
151:     if (flg) {
152:       PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pep),monfilename,&monviewer);
153:       PEPMonitorSet(pep,PEPMonitorFirst,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
154:     }
155:     PetscOptionsString("-pep_monitor_conv","Monitor approximate eigenvalues and error estimates as they converge","PEPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
156:     if (flg) {
157:       PetscNew(&ctx);
158:       PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pep),monfilename,&ctx->viewer);
159:       PEPMonitorSet(pep,PEPMonitorConverged,ctx,(PetscErrorCode (*)(void**))SlepcConvMonitorDestroy);
160:     }
161:     PetscOptionsString("-pep_monitor_all","Monitor approximate eigenvalues and error estimates","PEPMonitorSet","stdout",monfilename,PETSC_MAX_PATH_LEN,&flg);
162:     if (flg) {
163:       PetscViewerASCIIOpen(PetscObjectComm((PetscObject)pep),monfilename,&monviewer);
164:       PEPMonitorSet(pep,PEPMonitorAll,monviewer,(PetscErrorCode (*)(void**))PetscViewerDestroy);
165:       PEPSetTrackAll(pep,PETSC_TRUE);
166:     }
167:     flg = PETSC_FALSE;
168:     PetscOptionsBool("-pep_monitor_lg","Monitor first unconverged approximate error estimate graphically","PEPMonitorSet",flg,&flg,NULL);
169:     if (flg) {
170:       PEPMonitorSet(pep,PEPMonitorLG,NULL,NULL);
171:     }
172:     flg = PETSC_FALSE;
173:     PetscOptionsBool("-pep_monitor_lg_all","Monitor error estimates graphically","PEPMonitorSet",flg,&flg,NULL);
174:     if (flg) {
175:       PEPMonitorSet(pep,PEPMonitorLGAll,NULL,NULL);
176:       PEPSetTrackAll(pep,PETSC_TRUE);
177:     }
178:   /* -----------------------------------------------------------------------*/

180:     PetscOptionsBoolGroupBegin("-pep_largest_magnitude","compute largest eigenvalues in magnitude","PEPSetWhichEigenpairs",&flg);
181:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_LARGEST_MAGNITUDE); }
182:     PetscOptionsBoolGroup("-pep_smallest_magnitude","compute smallest eigenvalues in magnitude","PEPSetWhichEigenpairs",&flg);
183:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_SMALLEST_MAGNITUDE); }
184:     PetscOptionsBoolGroup("-pep_largest_real","compute largest real parts","PEPSetWhichEigenpairs",&flg);
185:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_LARGEST_REAL); }
186:     PetscOptionsBoolGroup("-pep_smallest_real","compute smallest real parts","PEPSetWhichEigenpairs",&flg);
187:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_SMALLEST_REAL); }
188:     PetscOptionsBoolGroup("-pep_largest_imaginary","compute largest imaginary parts","PEPSetWhichEigenpairs",&flg);
189:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_LARGEST_IMAGINARY); }
190:     PetscOptionsBoolGroup("-pep_smallest_imaginary","compute smallest imaginary parts","PEPSetWhichEigenpairs",&flg);
191:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_SMALLEST_IMAGINARY); }
192:     PetscOptionsBoolGroup("-pep_target_magnitude","compute nearest eigenvalues to target","PEPSetWhichEigenpairs",&flg);
193:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_TARGET_MAGNITUDE); }
194:     PetscOptionsBoolGroup("-pep_target_real","compute eigenvalues with real parts close to target","PEPSetWhichEigenpairs",&flg);
195:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_TARGET_REAL); }
196:     PetscOptionsBoolGroupEnd("-pep_target_imaginary","compute eigenvalues with imaginary parts close to target","PEPSetWhichEigenpairs",&flg);
197:     if (flg) { PEPSetWhichEigenpairs(pep,PEP_TARGET_IMAGINARY); }

199:     PetscOptionsName("-pep_view","Print detailed information on solver used","PEPView",0);
200:     PetscOptionsName("-pep_view_vectors","View computed eigenvectors","PEPVectorsView",0);
201:     PetscOptionsName("-pep_view_values","View computed eigenvalues","PEPValuesView",0);
202:     PetscOptionsName("-pep_converged_reason","Print reason for convergence, and number of iterations","PEPReasonView",0);
203:     PetscOptionsName("-pep_error_absolute","Print absolute errors of each eigenpair","PEPErrorView",0);
204:     PetscOptionsName("-pep_error_relative","Print relative errors of each eigenpair","PEPErrorView",0);
205:     PetscOptionsName("-pep_error_backward","Print backward errors of each eigenpair","PEPErrorView",0);

207:     if (pep->ops->setfromoptions) {
208:       (*pep->ops->setfromoptions)(PetscOptionsObject,pep);
209:     }
210:     PetscObjectProcessOptionsHandlers((PetscObject)pep);
211:   PetscOptionsEnd();

213:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
214:   BVSetFromOptions(pep->V);
215:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
216:   RGSetFromOptions(pep->rg);
217:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
218:   DSSetFromOptions(pep->ds);
219:   if (!pep->st) { PEPGetST(pep,&pep->st); }
220:   STSetFromOptions(pep->st);
221:   if (!pep->refineksp) { PEPRefineGetKSP(pep,&pep->refineksp); }
222:   KSPSetFromOptions(pep->refineksp);
223:   PetscRandomSetFromOptions(pep->rand);
224:   return(0);
225: }

229: /*@
230:    PEPGetTolerances - Gets the tolerance and maximum iteration count used
231:    by the PEP convergence tests.

233:    Not Collective

235:    Input Parameter:
236: .  pep - the polynomial eigensolver context

238:    Output Parameters:
239: +  tol - the convergence tolerance
240: -  maxits - maximum number of iterations

242:    Notes:
243:    The user can specify NULL for any parameter that is not needed.

245:    Level: intermediate

247: .seealso: PEPSetTolerances()
248: @*/
249: PetscErrorCode PEPGetTolerances(PEP pep,PetscReal *tol,PetscInt *maxits)
250: {
253:   if (tol)    *tol    = pep->tol;
254:   if (maxits) *maxits = pep->max_it;
255:   return(0);
256: }

260: /*@
261:    PEPSetTolerances - Sets the tolerance and maximum iteration count used
262:    by the PEP convergence tests.

264:    Logically Collective on PEP

266:    Input Parameters:
267: +  pep - the polynomial eigensolver context
268: .  tol - the convergence tolerance
269: -  maxits - maximum number of iterations to use

271:    Options Database Keys:
272: +  -pep_tol <tol> - Sets the convergence tolerance
273: -  -pep_max_it <maxits> - Sets the maximum number of iterations allowed

275:    Notes:
276:    Use PETSC_DEFAULT for either argument to assign a reasonably good value.

278:    Level: intermediate

280: .seealso: PEPGetTolerances()
281: @*/
282: PetscErrorCode PEPSetTolerances(PEP pep,PetscReal tol,PetscInt maxits)
283: {
288:   if (tol == PETSC_DEFAULT) {
289:     pep->tol   = PETSC_DEFAULT;
290:     pep->state = PEP_STATE_INITIAL;
291:   } else {
292:     if (tol <= 0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
293:     pep->tol = tol;
294:   }
295:   if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
296:     pep->max_it = 0;
297:     pep->state  = PEP_STATE_INITIAL;
298:   } else {
299:     if (maxits <= 0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
300:     pep->max_it = maxits;
301:   }
302:   return(0);
303: }

307: /*@
308:    PEPGetDimensions - Gets the number of eigenvalues to compute
309:    and the dimension of the subspace.

311:    Not Collective

313:    Input Parameter:
314: .  pep - the polynomial eigensolver context

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

321:    Notes:
322:    The user can specify NULL for any parameter that is not needed.

324:    Level: intermediate

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

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

344:    Logically Collective on PEP

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

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

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

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

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

370:    Level: intermediate

372: .seealso: PEPGetDimensions()
373: @*/
374: PetscErrorCode PEPSetDimensions(PEP pep,PetscInt nev,PetscInt ncv,PetscInt mpd)
375: {
381:   if (nev<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
382:   pep->nev = nev;
383:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
384:     pep->ncv = 0;
385:   } else {
386:     if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
387:     pep->ncv = ncv;
388:   }
389:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
390:     pep->mpd = 0;
391:   } else {
392:     if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
393:     pep->mpd = mpd;
394:   }
395:   pep->state = PEP_STATE_INITIAL;
396:   return(0);
397: }

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

405:    Logically Collective on PEP

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

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

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

425:    Options Database Keys:
426: +   -pep_largest_magnitude - Sets largest eigenvalues in magnitude
427: .   -pep_smallest_magnitude - Sets smallest eigenvalues in magnitude
428: .   -pep_largest_real - Sets largest real parts
429: .   -pep_smallest_real - Sets smallest real parts
430: .   -pep_largest_imaginary - Sets largest imaginary parts
431: .   -pep_smallest_imaginary - Sets smallest imaginary parts
432: .   -pep_target_magnitude - Sets eigenvalues closest to target
433: .   -pep_target_real - Sets real parts closest to target
434: -   -pep_target_imaginary - Sets imaginary parts closest to target

436:    Notes:
437:    Not all eigensolvers implemented in PEP account for all the possible values
438:    stated above. If SLEPc is compiled for real numbers PEP_LARGEST_IMAGINARY
439:    and PEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
440:    for eigenvalue selection.

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

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

447:    Level: intermediate

449: .seealso: PEPGetWhichEigenpairs(), PEPSetTarget(), PEPSetEigenvalueComparison(), PEPWhich
450: @*/
451: PetscErrorCode PEPSetWhichEigenpairs(PEP pep,PEPWhich which)
452: {
456:   switch (which) {
457:     case PEP_LARGEST_MAGNITUDE:
458:     case PEP_SMALLEST_MAGNITUDE:
459:     case PEP_LARGEST_REAL:
460:     case PEP_SMALLEST_REAL:
461:     case PEP_LARGEST_IMAGINARY:
462:     case PEP_SMALLEST_IMAGINARY:
463:     case PEP_TARGET_MAGNITUDE:
464:     case PEP_TARGET_REAL:
465: #if defined(PETSC_USE_COMPLEX)
466:     case PEP_TARGET_IMAGINARY:
467: #endif
468:     case PEP_WHICH_USER:
469:       if (pep->which != which) {
470:         pep->state = PEP_STATE_INITIAL;
471:         pep->which = which;
472:       }
473:       break;
474:     default:
475:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
476:   }
477:   return(0);
478: }

482: /*@
483:     PEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
484:     sought.

486:     Not Collective

488:     Input Parameter:
489: .   pep - eigensolver context obtained from PEPCreate()

491:     Output Parameter:
492: .   which - the portion of the spectrum to be sought

494:     Notes:
495:     See PEPSetWhichEigenpairs() for possible values of 'which'.

497:     Level: intermediate

499: .seealso: PEPSetWhichEigenpairs(), PEPWhich
500: @*/
501: PetscErrorCode PEPGetWhichEigenpairs(PEP pep,PEPWhich *which)
502: {
506:   *which = pep->which;
507:   return(0);
508: }

512: /*@C
513:    PEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
514:    when PEPSetWhichEigenpairs() is set to PEP_WHICH_USER.

516:    Logically Collective on PEP

518:    Input Parameters:
519: +  pep  - eigensolver context obtained from PEPCreate()
520: .  func - a pointer to the comparison function
521: -  ctx  - a context pointer (the last parameter to the comparison function)

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

526: +   ar     - real part of the 1st eigenvalue
527: .   ai     - imaginary part of the 1st eigenvalue
528: .   br     - real part of the 2nd eigenvalue
529: .   bi     - imaginary part of the 2nd eigenvalue
530: .   res    - result of comparison
531: -   ctx    - optional context, as set by PEPSetEigenvalueComparison()

533:    Note:
534:    The returning parameter 'res' can be:
535: +  negative - if the 1st eigenvalue is preferred to the 2st one
536: .  zero     - if both eigenvalues are equally preferred
537: -  positive - if the 2st eigenvalue is preferred to the 1st one

539:    Level: advanced

541: .seealso: PEPSetWhichEigenpairs(), PEPWhich
542: @*/
543: PetscErrorCode PEPSetEigenvalueComparison(PEP pep,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)
544: {
547:   pep->sc->comparison    = func;
548:   pep->sc->comparisonctx = ctx;
549:   pep->which             = PEP_WHICH_USER;
550:   return(0);
551: }

555: /*@
556:    PEPSetProblemType - Specifies the type of the polynomial eigenvalue problem.

558:    Logically Collective on PEP

560:    Input Parameters:
561: +  pep  - the polynomial eigensolver context
562: -  type - a known type of polynomial eigenvalue problem

564:    Options Database Keys:
565: +  -pep_general - general problem with no particular structure
566: .  -pep_hermitian - problem whose coefficient matrices are Hermitian
567: -  -pep_gyroscopic - problem with Hamiltonian structure

569:    Notes:
570:    Allowed values for the problem type are: general (PEP_GENERAL), Hermitian
571:    (PEP_HERMITIAN), and gyroscopic (PEP_GYROSCOPIC).

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

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

580:    Level: intermediate

582: .seealso: PEPSetOperators(), PEPSetType(), PEPGetProblemType(), PEPProblemType
583: @*/
584: PetscErrorCode PEPSetProblemType(PEP pep,PEPProblemType type)
585: {
589:   if (type!=PEP_GENERAL && type!=PEP_HERMITIAN && type!=PEP_GYROSCOPIC) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
590:   pep->problem_type = type;
591:   return(0);
592: }

596: /*@
597:    PEPGetProblemType - Gets the problem type from the PEP object.

599:    Not Collective

601:    Input Parameter:
602: .  pep - the polynomial eigensolver context

604:    Output Parameter:
605: .  type - name of PEP problem type

607:    Level: intermediate

609: .seealso: PEPSetProblemType(), PEPProblemType
610: @*/
611: PetscErrorCode PEPGetProblemType(PEP pep,PEPProblemType *type)
612: {
616:   *type = pep->problem_type;
617:   return(0);
618: }

622: /*@
623:    PEPSetBasis - Specifies the type of polynomial basis used to describe the
624:    polynomial eigenvalue problem.

626:    Logically Collective on PEP

628:    Input Parameters:
629: +  pep   - the polynomial eigensolver context
630: -  basis - the type of polynomial basis

632:    Options Database Key:
633: .  -pep_basis <basis> - Select the basis type

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

642:    Level: intermediate

644: .seealso: PEPSetOperators(), PEPGetBasis(), PEPBasis
645: @*/
646: PetscErrorCode PEPSetBasis(PEP pep,PEPBasis basis)
647: {
651:   pep->basis = basis;
652:   return(0);
653: }

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

660:    Not Collective

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

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

668:    Level: intermediate

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

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

687:    Logically Collective on PEP

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

693:    Notes:
694:    If the user sets trackall=PETSC_TRUE then the solver explicitly computes
695:    the residual for each eigenpair approximation. Computing the residual is
696:    usually an expensive operation and solvers commonly compute the associated
697:    residual to the first unconverged eigenpair.
698:    The options '-pep_monitor_all' and '-pep_monitor_lg_all' automatically
699:    activate this option.

701:    Level: developer

703: .seealso: PEPGetTrackAll()
704: @*/
705: PetscErrorCode PEPSetTrackAll(PEP pep,PetscBool trackall)
706: {
710:   pep->trackall = trackall;
711:   return(0);
712: }

716: /*@
717:    PEPGetTrackAll - Returns the flag indicating whether all residual norms must
718:    be computed or not.

720:    Not Collective

722:    Input Parameter:
723: .  pep - the eigensolver context

725:    Output Parameter:
726: .  trackall - the returned flag

728:    Level: developer

730: .seealso: PEPSetTrackAll()
731: @*/
732: PetscErrorCode PEPGetTrackAll(PEP pep,PetscBool *trackall)
733: {
737:   *trackall = pep->trackall;
738:   return(0);
739: }

743: /*@C
744:    PEPSetConvergenceTestFunction - Sets a function to compute the error estimate
745:    used in the convergence test.

747:    Logically Collective on PEP

749:    Input Parameters:
750: +  pep     - eigensolver context obtained from PEPCreate()
751: .  func    - a pointer to the convergence test function
752: .  ctx     - [optional] context for private data for the convergence routine
753: -  destroy - [optional] destructor for the context (may be NULL;
754:              PETSC_NULL_FUNCTION in Fortran)

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

759: +   pep    - eigensolver context obtained from PEPCreate()
760: .   eigr   - real part of the eigenvalue
761: .   eigi   - imaginary part of the eigenvalue
762: .   res    - residual norm associated to the eigenpair
763: .   errest - (output) computed error estimate
764: -   ctx    - optional context, as set by PEPSetConvergenceTest()

766:    Note:
767:    If the error estimate returned by the convergence test function is less than
768:    the tolerance, then the eigenvalue is accepted as converged.

770:    Level: advanced

772: .seealso: PEPSetConvergenceTest(), PEPSetTolerances()
773: @*/
774: PetscErrorCode PEPSetConvergenceTestFunction(PEP pep,PetscErrorCode (*func)(PEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
775: {

780:   if (pep->convergeddestroy) {
781:     (*pep->convergeddestroy)(pep->convergedctx);
782:   }
783:   pep->converged        = func;
784:   pep->convergeddestroy = destroy;
785:   pep->convergedctx     = ctx;
786:   if (func == PEPConvergedEigRelative) pep->conv = PEP_CONV_EIG;
787:   else if (func == PEPConvergedLinear) pep->conv = PEP_CONV_LINEAR;
788:   else if (func == PEPConvergedNorm) pep->conv = PEP_CONV_NORM;
789:   else if (func == PEPConvergedAbsolute) pep->conv = PEP_CONV_ABS;
790:   else pep->conv = PEP_CONV_USER;
791:   return(0);
792: }

796: /*@
797:    PEPSetConvergenceTest - Specifies how to compute the error estimate
798:    used in the convergence test.

800:    Logically Collective on PEP

802:    Input Parameters:
803: +  pep  - eigensolver context obtained from PEPCreate()
804: -  conv - the type of convergence test

806:    Options Database Keys:
807: +  -pep_conv_abs    - Sets the absolute convergence test
808: .  -pep_conv_eig    - Sets the convergence test relative to the eigenvalue
809: .  -pep_conv_linear - Sets the convergence test related to the linearized eigenproblem
810: -  -pep_conv_user   - Selects the user-defined convergence test

812:    Note:
813:    The parameter 'conv' can have one of these values
814: +     PEP_CONV_ABS    - absolute error ||r||
815: .     PEP_CONV_EIG    - error relative to the eigenvalue l, ||r||/|l|
816: .     PEP_CONV_LINEAR - error related to the linearized eigenproblem
817: .     PEP_CONV_NORM   - error relative matrix norms, ||r||/sum_i(l^i*||A_i||)
818: -     PEP_CONV_USER   - function set by PEPSetConvergenceTestFunction()

820:    Level: intermediate

822: .seealso: PEPGetConvergenceTest(), PEPSetConvergenceTestFunction(), PEPConv
823: @*/
824: PetscErrorCode PEPSetConvergenceTest(PEP pep,PEPConv conv)
825: {
829:   switch (conv) {
830:     case PEP_CONV_ABS:    pep->converged = PEPConvergedAbsolute; break;
831:     case PEP_CONV_EIG:    pep->converged = PEPConvergedEigRelative; break;
832:     case PEP_CONV_LINEAR: pep->converged = PEPConvergedLinear; break;
833:     case PEP_CONV_NORM:   pep->converged = PEPConvergedNorm; break;
834:     case PEP_CONV_USER: break;
835:     default:
836:       SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
837:   }
838:   pep->conv = conv;
839:   return(0);
840: }

844: /*@
845:    PEPGetConvergenceTest - Gets the method used to compute the error estimate
846:    used in the convergence test.

848:    Not Collective

850:    Input Parameters:
851: .  pep   - eigensolver context obtained from PEPCreate()

853:    Output Parameters:
854: .  conv  - the type of convergence test

856:    Level: intermediate

858: .seealso: PEPSetConvergenceTest(), PEPConv
859: @*/
860: PetscErrorCode PEPGetConvergenceTest(PEP pep,PEPConv *conv)
861: {
865:   *conv = pep->conv;
866:   return(0);
867: }

871: /*@
872:    PEPSetScale - Specifies the scaling strategy to be used.

874:    Logically Collective on PEP

876:    Input Parameters:
877: +  pep    - the eigensolver context
878: .  scale  - scaling strategy
879: .  alpha  - the scaling factor used in the scalar strategy
880: .  Dl     - the left diagonal matrix of the diagonal scaling algorithm
881: .  Dr     - the right diagonal matrix of the diagonal scaling algorithm
882: .  its    - number of iterations of the diagonal scaling algorithm
883: -  lambda - approximation to wanted eigenvalues (modulus)

885:    Options Database Keys:
886: +  -pep_scale <type> - scaling type, one of <none,scalar,diagonal,both>
887: .  -pep_scale_factor <alpha> - the scaling factor
888: .  -pep_scale_its <its> - number of iterations
889: -  -pep_scale_lambda <lambda> - approximation to eigenvalues

891:    Notes:
892:    There are two non-exclusive scaling strategies: scalar and diagonal.

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

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

910:    Level: intermediate

912: .seealso: PEPGetScale()
913: @*/
914: PetscErrorCode PEPSetScale(PEP pep,PEPScale scale,PetscReal alpha,Vec Dl,Vec Dr,PetscInt its,PetscReal lambda)
915: {

921:   pep->scale = scale;
922:   if (scale==PEP_SCALE_SCALAR || scale==PEP_SCALE_BOTH) {
924:     if (alpha == PETSC_DEFAULT || alpha == PETSC_DECIDE) {
925:       pep->sfactor = 0.0;
926:       pep->sfactor_set = PETSC_FALSE;
927:     } else {
928:       if (alpha<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of alpha. Must be > 0");
929:       pep->sfactor = alpha;
930:       pep->sfactor_set = PETSC_TRUE;
931:     }
932:   }
933:   if (scale==PEP_SCALE_DIAGONAL || scale==PEP_SCALE_BOTH) {
934:     if (Dl) {
937:       PetscObjectReference((PetscObject)Dl);
938:       VecDestroy(&pep->Dl);
939:       pep->Dl = Dl;
940:     }
941:     if (Dr) {
944:       PetscObjectReference((PetscObject)Dr);
945:       VecDestroy(&pep->Dr);
946:       pep->Dr = Dr;
947:     }
950:     if (its==PETSC_DECIDE || its==PETSC_DEFAULT) pep->sits = 5;
951:     else pep->sits = its;
952:     if (lambda==PETSC_DECIDE || lambda==PETSC_DEFAULT) pep->slambda = 1.0;
953:     else if (lambda<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of lambda. Must be > 0");
954:     else pep->slambda = lambda;
955:   }
956:   pep->state = PEP_STATE_INITIAL;
957:   return(0);
958: }

962: /*@
963:    PEPGetScale - Gets the scaling strategy used by the PEP object, and the
964:    associated parameters.

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

968:    Input Parameter:
969: .  pep - the eigensolver context

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

979:    Level: intermediate

981:    Note:
982:    The user can specify NULL for any parameter that is not needed.

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

987: .seealso: PEPSetScale(), PEPSetUp()
988: @*/
989: PetscErrorCode PEPGetScale(PEP pep,PEPScale *scale,PetscReal *alpha,Vec *Dl,Vec *Dr,PetscInt *its,PetscReal *lambda)
990: {
993:   if (scale)  *scale  = pep->scale;
994:   if (alpha)  *alpha  = pep->sfactor;
995:   if (Dl)     *Dl     = pep->Dl;
996:   if (Dr)     *Dr     = pep->Dr;
997:   if (its)    *its    = pep->sits;
998:   if (lambda) *lambda = pep->slambda;
999:   return(0);
1000: }

1004: /*@
1005:    PEPSetExtract - Specifies the extraction strategy to be used.

1007:    Logically Collective on PEP

1009:    Input Parameters:
1010: +  pep     - the eigensolver context
1011: -  extract - extraction strategy

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

1016:    Level: intermediate

1018: .seealso: PEPGetExtract()
1019: @*/
1020: PetscErrorCode PEPSetExtract(PEP pep,PEPExtract extract)
1021: {
1025:   pep->extract = extract;
1026:   return(0);
1027: }

1031: /*@
1032:    PEPGetExtract - Gets the extraction strategy used by the PEP object.

1034:    Not Collective

1036:    Input Parameter:
1037: .  pep - the eigensolver context

1039:    Output Parameter:
1040: .  extract - extraction strategy

1042:    Level: intermediate

1044: .seealso: PEPSetExtract()
1045: @*/
1046: PetscErrorCode PEPGetExtract(PEP pep,PEPExtract *extract)
1047: {
1050:   if (extract) *extract = pep->extract;
1051:   return(0);
1052: }

1056: /*@
1057:    PEPSetRefine - Specifies the refinement type (and options) to be used
1058:    after the solve.

1060:    Logically Collective on PEP

1062:    Input Parameters:
1063: +  pep    - the polynomial eigensolver context
1064: .  refine - refinement type
1065: .  npart  - number of partitions of the communicator
1066: .  tol    - the convergence tolerance
1067: .  its    - maximum number of refinement iterations
1068: -  schur  - boolean flag to activate the Schur complement approach

1070:    Options Database Keys:
1071: +  -pep_refine <type> - refinement type, one of <none,simple,multiple>
1072: .  -pep_refine_partitions <n> - the number of partitions
1073: .  -pep_refine_tol <tol> - the tolerance
1074: .  -pep_refine_its <its> - number of iterations
1075: -  -pep_refine_schur - to set the Schur complement approach

1077:    Notes:
1078:    By default, iterative refinement is disabled, since it may be very
1079:    costly. There are two possible refinement strategies: simple and multiple.
1080:    The simple approach performs iterative refinement on each of the
1081:    converged eigenpairs individually, whereas the multiple strategy works
1082:    with the invariant pair as a whole, refining all eigenpairs simultaneously.
1083:    The latter may be required for the case of multiple eigenvalues.

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

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

1096:    The schur flag is used to change the way in which linear systems are
1097:    solved, so that a Schur complement approach is used instead of explicitly
1098:    building the coefficient matrix.

1100:    Level: intermediate

1102: .seealso: PEPGetRefine()
1103: @*/
1104: PetscErrorCode PEPSetRefine(PEP pep,PEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,PetscBool schur)
1105: {
1107:   PetscMPIInt    size;

1116:   pep->refine = refine;
1117:   if (refine) {  /* process parameters only if not REFINE_NONE */
1118:     if (npart!=pep->npart) {
1119:       PetscSubcommDestroy(&pep->refinesubc);
1120:       KSPDestroy(&pep->refineksp);
1121:     }
1122:     if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
1123:       pep->npart = 1;
1124:     } else {
1125:       MPI_Comm_size(PetscObjectComm((PetscObject)pep),&size);
1126:       if (npart<1 || npart>size) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
1127:       pep->npart = npart;
1128:     }
1129:     if (tol == PETSC_DEFAULT || tol == PETSC_DECIDE) {
1130:       pep->rtol = pep->tol;
1131:     } else {
1132:       if (tol<=0.0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1133:       pep->rtol = tol;
1134:     }
1135:     if (its==PETSC_DECIDE || its==PETSC_DEFAULT) {
1136:       pep->rits = PETSC_DEFAULT;
1137:     } else {
1138:       if (its<0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be >= 0");
1139:       pep->rits = its;
1140:     }
1141:     pep->schur = schur;
1142:   }
1143:   pep->state = PEP_STATE_INITIAL;
1144:   return(0);
1145: }

1149: /*@
1150:    PEPGetRefine - Gets the refinement strategy used by the PEP object, and the
1151:    associated parameters.

1153:    Not Collective

1155:    Input Parameter:
1156: .  pep - the polynomial eigensolver context

1158:    Output Parameters:
1159: +  refine - refinement type
1160: .  npart  - number of partitions of the communicator
1161: .  tol    - the convergence tolerance
1162: .  its    - maximum number of refinement iterations
1163: -  schur  - whether the Schur complement approach is being used

1165:    Level: intermediate

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

1170: .seealso: PEPSetRefine()
1171: @*/
1172: PetscErrorCode PEPGetRefine(PEP pep,PEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,PetscBool *schur)
1173: {
1176:   if (refine) *refine = pep->refine;
1177:   if (npart)  *npart  = pep->npart;
1178:   if (tol)    *tol    = pep->rtol;
1179:   if (its)    *its    = pep->rits;
1180:   if (schur)  *schur  = pep->schur;
1181:   return(0);
1182: }

1186: /*@C
1187:    PEPSetOptionsPrefix - Sets the prefix used for searching for all
1188:    PEP options in the database.

1190:    Logically Collective on PEP

1192:    Input Parameters:
1193: +  pep - the polynomial eigensolver context
1194: -  prefix - the prefix string to prepend to all PEP option requests

1196:    Notes:
1197:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1198:    The first character of all runtime options is AUTOMATICALLY the
1199:    hyphen.

1201:    For example, to distinguish between the runtime options for two
1202:    different PEP contexts, one could call
1203: .vb
1204:       PEPSetOptionsPrefix(pep1,"qeig1_")
1205:       PEPSetOptionsPrefix(pep2,"qeig2_")
1206: .ve

1208:    Level: advanced

1210: .seealso: PEPAppendOptionsPrefix(), PEPGetOptionsPrefix()
1211: @*/
1212: PetscErrorCode PEPSetOptionsPrefix(PEP pep,const char *prefix)
1213: {

1218:   if (!pep->st) { PEPGetST(pep,&pep->st); }
1219:   STSetOptionsPrefix(pep->st,prefix);
1220:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
1221:   BVSetOptionsPrefix(pep->V,prefix);
1222:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
1223:   DSSetOptionsPrefix(pep->ds,prefix);
1224:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
1225:   RGSetOptionsPrefix(pep->rg,prefix);
1226:   PetscObjectSetOptionsPrefix((PetscObject)pep,prefix);
1227:   return(0);
1228: }

1232: /*@C
1233:    PEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1234:    PEP options in the database.

1236:    Logically Collective on PEP

1238:    Input Parameters:
1239: +  pep - the polynomial eigensolver context
1240: -  prefix - the prefix string to prepend to all PEP option requests

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

1246:    Level: advanced

1248: .seealso: PEPSetOptionsPrefix(), PEPGetOptionsPrefix()
1249: @*/
1250: PetscErrorCode PEPAppendOptionsPrefix(PEP pep,const char *prefix)
1251: {
1253:   PetscBool      flg;
1254:   EPS            eps;

1258:   if (!pep->st) { PEPGetST(pep,&pep->st); }
1259:   STAppendOptionsPrefix(pep->st,prefix);
1260:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
1261:   BVSetOptionsPrefix(pep->V,prefix);
1262:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
1263:   DSSetOptionsPrefix(pep->ds,prefix);
1264:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
1265:   RGSetOptionsPrefix(pep->rg,prefix);
1266:   PetscObjectAppendOptionsPrefix((PetscObject)pep,prefix);
1267:   PetscObjectTypeCompare((PetscObject)pep,PEPLINEAR,&flg);
1268:   if (flg) {
1269:     PEPLinearGetEPS(pep,&eps);
1270:     EPSSetOptionsPrefix(eps,((PetscObject)pep)->prefix);
1271:     EPSAppendOptionsPrefix(eps,"pep_");
1272:   }
1273:   return(0);
1274: }

1278: /*@C
1279:    PEPGetOptionsPrefix - Gets the prefix used for searching for all
1280:    PEP options in the database.

1282:    Not Collective

1284:    Input Parameters:
1285: .  pep - the polynomial eigensolver context

1287:    Output Parameters:
1288: .  prefix - pointer to the prefix string used is returned

1290:    Notes: On the fortran side, the user should pass in a string 'prefix' of
1291:    sufficient length to hold the prefix.

1293:    Level: advanced

1295: .seealso: PEPSetOptionsPrefix(), PEPAppendOptionsPrefix()
1296: @*/
1297: PetscErrorCode PEPGetOptionsPrefix(PEP pep,const char *prefix[])
1298: {

1304:   PetscObjectGetOptionsPrefix((PetscObject)pep,prefix);
1305:   return(0);
1306: }