Actual source code: nepopts.c

slepc-3.16.1 2021-11-17
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    NEP routines related to options that can be set via the command-line
 12:    or procedurally
 13: */

 15: #include <slepc/private/nepimpl.h>
 16: #include <petscdraw.h>

 18: /*@C
 19:    NEPMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
 20:    indicated by the user.

 22:    Collective on nep

 24:    Input Parameters:
 25: +  nep      - the nonlinear eigensolver context
 26: .  opt      - the command line option for this monitor
 27: .  name     - the monitor type one is seeking
 28: .  ctx      - an optional user context for the monitor, or NULL
 29: -  trackall - whether this monitor tracks all eigenvalues or not

 31:    Level: developer

 33: .seealso: NEPMonitorSet(), NEPSetTrackAll()
 34: @*/
 35: PetscErrorCode NEPMonitorSetFromOptions(NEP nep,const char opt[],const char name[],void *ctx,PetscBool trackall)
 36: {
 37:   PetscErrorCode       (*mfunc)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
 38:   PetscErrorCode       (*cfunc)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
 39:   PetscErrorCode       (*dfunc)(PetscViewerAndFormat**);
 40:   PetscViewerAndFormat *vf;
 41:   PetscViewer          viewer;
 42:   PetscViewerFormat    format;
 43:   PetscViewerType      vtype;
 44:   char                 key[PETSC_MAX_PATH_LEN];
 45:   PetscBool            flg;
 46:   PetscErrorCode       ierr;

 49:   PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->options,((PetscObject)nep)->prefix,opt,&viewer,&format,&flg);
 50:   if (!flg) return(0);

 52:   PetscViewerGetType(viewer,&vtype);
 53:   SlepcMonitorMakeKey_Internal(name,vtype,format,key);
 54:   PetscFunctionListFind(NEPMonitorList,key,&mfunc);
 55:   if (!mfunc) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Specified viewer and format not supported");
 56:   PetscFunctionListFind(NEPMonitorCreateList,key,&cfunc);
 57:   PetscFunctionListFind(NEPMonitorDestroyList,key,&dfunc);
 58:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
 59:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

 61:   (*cfunc)(viewer,format,ctx,&vf);
 62:   PetscObjectDereference((PetscObject)viewer);
 63:   NEPMonitorSet(nep,mfunc,vf,(PetscErrorCode(*)(void **))dfunc);
 64:   if (trackall) {
 65:     NEPSetTrackAll(nep,PETSC_TRUE);
 66:   }
 67:   return(0);
 68: }

 70: /*@
 71:    NEPSetFromOptions - Sets NEP options from the options database.
 72:    This routine must be called before NEPSetUp() if the user is to be
 73:    allowed to set the solver type.

 75:    Collective on nep

 77:    Input Parameters:
 78: .  nep - the nonlinear eigensolver context

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

 83:    Level: beginner
 84: @*/
 85: PetscErrorCode NEPSetFromOptions(NEP nep)
 86: {
 87:   PetscErrorCode  ierr;
 88:   char            type[256];
 89:   PetscBool       set,flg,flg1,flg2,flg3,flg4,flg5,bval;
 90:   PetscReal       r;
 91:   PetscScalar     s;
 92:   PetscInt        i,j,k;
 93:   NEPRefine       refine;
 94:   NEPRefineScheme scheme;

 98:   NEPRegisterAll();
 99:   PetscObjectOptionsBegin((PetscObject)nep);
100:     PetscOptionsFList("-nep_type","Nonlinear eigensolver method","NEPSetType",NEPList,(char*)(((PetscObject)nep)->type_name?((PetscObject)nep)->type_name:NEPRII),type,sizeof(type),&flg);
101:     if (flg) {
102:       NEPSetType(nep,type);
103:     } else if (!((PetscObject)nep)->type_name) {
104:       NEPSetType(nep,NEPRII);
105:     }

107:     PetscOptionsBoolGroupBegin("-nep_general","General nonlinear eigenvalue problem","NEPSetProblemType",&flg);
108:     if (flg) { NEPSetProblemType(nep,NEP_GENERAL); }
109:     PetscOptionsBoolGroupEnd("-nep_rational","Rational eigenvalue problem","NEPSetProblemType",&flg);
110:     if (flg) { NEPSetProblemType(nep,NEP_RATIONAL); }

112:     refine = nep->refine;
113:     PetscOptionsEnum("-nep_refine","Iterative refinement method","NEPSetRefine",NEPRefineTypes,(PetscEnum)refine,(PetscEnum*)&refine,&flg1);
114:     i = nep->npart;
115:     PetscOptionsInt("-nep_refine_partitions","Number of partitions of the communicator for iterative refinement","NEPSetRefine",nep->npart,&i,&flg2);
116:     r = nep->rtol;
117:     PetscOptionsReal("-nep_refine_tol","Tolerance for iterative refinement","NEPSetRefine",nep->rtol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/1000:nep->rtol,&r,&flg3);
118:     j = nep->rits;
119:     PetscOptionsInt("-nep_refine_its","Maximum number of iterations for iterative refinement","NEPSetRefine",nep->rits,&j,&flg4);
120:     scheme = nep->scheme;
121:     PetscOptionsEnum("-nep_refine_scheme","Scheme used for linear systems within iterative refinement","NEPSetRefine",NEPRefineSchemes,(PetscEnum)scheme,(PetscEnum*)&scheme,&flg5);
122:     if (flg1 || flg2 || flg3 || flg4 || flg5) { NEPSetRefine(nep,refine,i,r,j,scheme); }

124:     i = nep->max_it;
125:     PetscOptionsInt("-nep_max_it","Maximum number of iterations","NEPSetTolerances",nep->max_it,&i,&flg1);
126:     r = nep->tol;
127:     PetscOptionsReal("-nep_tol","Tolerance","NEPSetTolerances",SlepcDefaultTol(nep->tol),&r,&flg2);
128:     if (flg1 || flg2) { NEPSetTolerances(nep,r,i); }

130:     PetscOptionsBoolGroupBegin("-nep_conv_rel","Relative error convergence test","NEPSetConvergenceTest",&flg);
131:     if (flg) { NEPSetConvergenceTest(nep,NEP_CONV_REL); }
132:     PetscOptionsBoolGroup("-nep_conv_norm","Convergence test relative to the matrix norms","NEPSetConvergenceTest",&flg);
133:     if (flg) { NEPSetConvergenceTest(nep,NEP_CONV_NORM); }
134:     PetscOptionsBoolGroup("-nep_conv_abs","Absolute error convergence test","NEPSetConvergenceTest",&flg);
135:     if (flg) { NEPSetConvergenceTest(nep,NEP_CONV_ABS); }
136:     PetscOptionsBoolGroupEnd("-nep_conv_user","User-defined convergence test","NEPSetConvergenceTest",&flg);
137:     if (flg) { NEPSetConvergenceTest(nep,NEP_CONV_USER); }

139:     PetscOptionsBoolGroupBegin("-nep_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","NEPSetStoppingTest",&flg);
140:     if (flg) { NEPSetStoppingTest(nep,NEP_STOP_BASIC); }
141:     PetscOptionsBoolGroupEnd("-nep_stop_user","User-defined stopping test","NEPSetStoppingTest",&flg);
142:     if (flg) { NEPSetStoppingTest(nep,NEP_STOP_USER); }

144:     i = nep->nev;
145:     PetscOptionsInt("-nep_nev","Number of eigenvalues to compute","NEPSetDimensions",nep->nev,&i,&flg1);
146:     j = nep->ncv;
147:     PetscOptionsInt("-nep_ncv","Number of basis vectors","NEPSetDimensions",nep->ncv,&j,&flg2);
148:     k = nep->mpd;
149:     PetscOptionsInt("-nep_mpd","Maximum dimension of projected problem","NEPSetDimensions",nep->mpd,&k,&flg3);
150:     if (flg1 || flg2 || flg3) {
151:       NEPSetDimensions(nep,i,j,k);
152:     }

154:     PetscOptionsBoolGroupBegin("-nep_largest_magnitude","Compute largest eigenvalues in magnitude","NEPSetWhichEigenpairs",&flg);
155:     if (flg) { NEPSetWhichEigenpairs(nep,NEP_LARGEST_MAGNITUDE); }
156:     PetscOptionsBoolGroup("-nep_smallest_magnitude","Compute smallest eigenvalues in magnitude","NEPSetWhichEigenpairs",&flg);
157:     if (flg) { NEPSetWhichEigenpairs(nep,NEP_SMALLEST_MAGNITUDE); }
158:     PetscOptionsBoolGroup("-nep_largest_real","Compute eigenvalues with largest real parts","NEPSetWhichEigenpairs",&flg);
159:     if (flg) { NEPSetWhichEigenpairs(nep,NEP_LARGEST_REAL); }
160:     PetscOptionsBoolGroup("-nep_smallest_real","Compute eigenvalues with smallest real parts","NEPSetWhichEigenpairs",&flg);
161:     if (flg) { NEPSetWhichEigenpairs(nep,NEP_SMALLEST_REAL); }
162:     PetscOptionsBoolGroup("-nep_largest_imaginary","Compute eigenvalues with largest imaginary parts","NEPSetWhichEigenpairs",&flg);
163:     if (flg) { NEPSetWhichEigenpairs(nep,NEP_LARGEST_IMAGINARY); }
164:     PetscOptionsBoolGroup("-nep_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","NEPSetWhichEigenpairs",&flg);
165:     if (flg) { NEPSetWhichEigenpairs(nep,NEP_SMALLEST_IMAGINARY); }
166:     PetscOptionsBoolGroup("-nep_target_magnitude","Compute eigenvalues closest to target","NEPSetWhichEigenpairs",&flg);
167:     if (flg) { NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE); }
168:     PetscOptionsBoolGroup("-nep_target_real","Compute eigenvalues with real parts closest to target","NEPSetWhichEigenpairs",&flg);
169:     if (flg) { NEPSetWhichEigenpairs(nep,NEP_TARGET_REAL); }
170:     PetscOptionsBoolGroup("-nep_target_imaginary","Compute eigenvalues with imaginary parts closest to target","NEPSetWhichEigenpairs",&flg);
171:     if (flg) { NEPSetWhichEigenpairs(nep,NEP_TARGET_IMAGINARY); }
172:     PetscOptionsBoolGroupEnd("-nep_all","Compute all eigenvalues in a region","NEPSetWhichEigenpairs",&flg);
173:     if (flg) { NEPSetWhichEigenpairs(nep,NEP_ALL); }

175:     PetscOptionsScalar("-nep_target","Value of the target","NEPSetTarget",nep->target,&s,&flg);
176:     if (flg) {
177:       if (nep->which!=NEP_TARGET_REAL && nep->which!=NEP_TARGET_IMAGINARY) {
178:         NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE);
179:       }
180:       NEPSetTarget(nep,s);
181:     }

183:     PetscOptionsBool("-nep_two_sided","Use two-sided variant (to compute left eigenvectors)","NEPSetTwoSided",nep->twosided,&bval,&flg);
184:     if (flg) { NEPSetTwoSided(nep,bval); }

186:     /* -----------------------------------------------------------------------*/
187:     /*
188:       Cancels all monitors hardwired into code before call to NEPSetFromOptions()
189:     */
190:     PetscOptionsBool("-nep_monitor_cancel","Remove any hardwired monitor routines","NEPMonitorCancel",PETSC_FALSE,&flg,&set);
191:     if (set && flg) { NEPMonitorCancel(nep); }
192:     NEPMonitorSetFromOptions(nep,"-nep_monitor","first_approximation",NULL,PETSC_FALSE);
193:     NEPMonitorSetFromOptions(nep,"-nep_monitor_all","all_approximations",NULL,PETSC_TRUE);
194:     NEPMonitorSetFromOptions(nep,"-nep_monitor_conv","convergence_history",NULL,PETSC_FALSE);

196:     /* -----------------------------------------------------------------------*/
197:     PetscOptionsName("-nep_view","Print detailed information on solver used","NEPView",NULL);
198:     PetscOptionsName("-nep_view_vectors","View computed eigenvectors","NEPVectorsView",NULL);
199:     PetscOptionsName("-nep_view_values","View computed eigenvalues","NEPValuesView",NULL);
200:     PetscOptionsName("-nep_converged_reason","Print reason for convergence, and number of iterations","NEPConvergedReasonView",NULL);
201:     PetscOptionsName("-nep_error_absolute","Print absolute errors of each eigenpair","NEPErrorView",NULL);
202:     PetscOptionsName("-nep_error_relative","Print relative errors of each eigenpair","NEPErrorView",NULL);

204:     if (nep->ops->setfromoptions) {
205:       (*nep->ops->setfromoptions)(PetscOptionsObject,nep);
206:     }
207:     PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)nep);
208:   PetscOptionsEnd();

210:   if (!nep->V) { NEPGetBV(nep,&nep->V); }
211:   BVSetFromOptions(nep->V);
212:   if (!nep->rg) { NEPGetRG(nep,&nep->rg); }
213:   RGSetFromOptions(nep->rg);
214:   if (nep->useds) {
215:     if (!nep->ds) { NEPGetDS(nep,&nep->ds); }
216:     DSSetFromOptions(nep->ds);
217:   }
218:   if (!nep->refineksp) { NEPRefineGetKSP(nep,&nep->refineksp); }
219:   KSPSetFromOptions(nep->refineksp);
220:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) for (i=0;i<nep->nt;i++) {FNSetFromOptions(nep->f[i]);}
221:   return(0);
222: }

224: /*@C
225:    NEPGetTolerances - Gets the tolerance and maximum iteration count used
226:    by the NEP convergence tests.

228:    Not Collective

230:    Input Parameter:
231: .  nep - the nonlinear eigensolver context

233:    Output Parameters:
234: +  tol - the convergence tolerance
235: -  maxits - maximum number of iterations

237:    Notes:
238:    The user can specify NULL for any parameter that is not needed.

240:    Level: intermediate

242: .seealso: NEPSetTolerances()
243: @*/
244: PetscErrorCode NEPGetTolerances(NEP nep,PetscReal *tol,PetscInt *maxits)
245: {
248:   if (tol)    *tol    = nep->tol;
249:   if (maxits) *maxits = nep->max_it;
250:   return(0);
251: }

253: /*@
254:    NEPSetTolerances - Sets the tolerance and maximum iteration count used
255:    by the NEP convergence tests.

257:    Logically Collective on nep

259:    Input Parameters:
260: +  nep    - the nonlinear eigensolver context
261: .  tol    - the convergence tolerance
262: -  maxits - maximum number of iterations to use

264:    Options Database Keys:
265: +  -nep_tol <tol> - Sets the convergence tolerance
266: -  -nep_max_it <maxits> - Sets the maximum number of iterations allowed

268:    Notes:
269:    Use PETSC_DEFAULT for either argument to assign a reasonably good value.

271:    Level: intermediate

273: .seealso: NEPGetTolerances()
274: @*/
275: PetscErrorCode NEPSetTolerances(NEP nep,PetscReal tol,PetscInt maxits)
276: {
281:   if (tol == PETSC_DEFAULT) {
282:     nep->tol   = PETSC_DEFAULT;
283:     nep->state = NEP_STATE_INITIAL;
284:   } else {
285:     if (tol <= 0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
286:     nep->tol = tol;
287:   }
288:   if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
289:     nep->max_it = PETSC_DEFAULT;
290:     nep->state  = NEP_STATE_INITIAL;
291:   } else {
292:     if (maxits <= 0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
293:     nep->max_it = maxits;
294:   }
295:   return(0);
296: }

298: /*@C
299:    NEPGetDimensions - Gets the number of eigenvalues to compute
300:    and the dimension of the subspace.

302:    Not Collective

304:    Input Parameter:
305: .  nep - the nonlinear eigensolver context

307:    Output Parameters:
308: +  nev - number of eigenvalues to compute
309: .  ncv - the maximum dimension of the subspace to be used by the solver
310: -  mpd - the maximum dimension allowed for the projected problem

312:    Notes:
313:    The user can specify NULL for any parameter that is not needed.

315:    Level: intermediate

317: .seealso: NEPSetDimensions()
318: @*/
319: PetscErrorCode NEPGetDimensions(NEP nep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
320: {
323:   if (nev) *nev = nep->nev;
324:   if (ncv) *ncv = nep->ncv;
325:   if (mpd) *mpd = nep->mpd;
326:   return(0);
327: }

329: /*@
330:    NEPSetDimensions - Sets the number of eigenvalues to compute
331:    and the dimension of the subspace.

333:    Logically Collective on nep

335:    Input Parameters:
336: +  nep - the nonlinear eigensolver context
337: .  nev - number of eigenvalues to compute
338: .  ncv - the maximum dimension of the subspace to be used by the solver
339: -  mpd - the maximum dimension allowed for the projected problem

341:    Options Database Keys:
342: +  -nep_nev <nev> - Sets the number of eigenvalues
343: .  -nep_ncv <ncv> - Sets the dimension of the subspace
344: -  -nep_mpd <mpd> - Sets the maximum projected dimension

346:    Notes:
347:    Use PETSC_DEFAULT for ncv and mpd to assign a reasonably good value, which is
348:    dependent on the solution method.

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

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

359:    Level: intermediate

361: .seealso: NEPGetDimensions()
362: @*/
363: PetscErrorCode NEPSetDimensions(NEP nep,PetscInt nev,PetscInt ncv,PetscInt mpd)
364: {
370:   if (nev<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
371:   nep->nev = nev;
372:   if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
373:     nep->ncv = PETSC_DEFAULT;
374:   } else {
375:     if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
376:     nep->ncv = ncv;
377:   }
378:   if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
379:     nep->mpd = PETSC_DEFAULT;
380:   } else {
381:     if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
382:     nep->mpd = mpd;
383:   }
384:   nep->state = NEP_STATE_INITIAL;
385:   return(0);
386: }

388: /*@
389:     NEPSetWhichEigenpairs - Specifies which portion of the spectrum is
390:     to be sought.

392:     Logically Collective on nep

394:     Input Parameters:
395: +   nep   - eigensolver context obtained from NEPCreate()
396: -   which - the portion of the spectrum to be sought

398:     Possible values:
399:     The parameter 'which' can have one of these values

401: +     NEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
402: .     NEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
403: .     NEP_LARGEST_REAL - largest real parts
404: .     NEP_SMALLEST_REAL - smallest real parts
405: .     NEP_LARGEST_IMAGINARY - largest imaginary parts
406: .     NEP_SMALLEST_IMAGINARY - smallest imaginary parts
407: .     NEP_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
408: .     NEP_TARGET_REAL - eigenvalues with real part closest to target
409: .     NEP_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
410: .     NEP_ALL - all eigenvalues contained in a given region
411: -     NEP_WHICH_USER - user defined ordering set with NEPSetEigenvalueComparison()

413:     Options Database Keys:
414: +   -nep_largest_magnitude - Sets largest eigenvalues in magnitude
415: .   -nep_smallest_magnitude - Sets smallest eigenvalues in magnitude
416: .   -nep_largest_real - Sets largest real parts
417: .   -nep_smallest_real - Sets smallest real parts
418: .   -nep_largest_imaginary - Sets largest imaginary parts
419: .   -nep_smallest_imaginary - Sets smallest imaginary parts
420: .   -nep_target_magnitude - Sets eigenvalues closest to target
421: .   -nep_target_real - Sets real parts closest to target
422: .   -nep_target_imaginary - Sets imaginary parts closest to target
423: -   -nep_all - Sets all eigenvalues in a region

425:     Notes:
426:     Not all eigensolvers implemented in NEP account for all the possible values
427:     stated above. If SLEPc is compiled for real numbers NEP_LARGEST_IMAGINARY
428:     and NEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
429:     for eigenvalue selection.

431:     The target is a scalar value provided with NEPSetTarget().

433:     NEP_ALL is intended for use in the context of the CISS solver for
434:     computing all eigenvalues in a region.

436:     Level: intermediate

438: .seealso: NEPGetWhichEigenpairs(), NEPSetTarget(), NEPSetEigenvalueComparison(), NEPWhich
439: @*/
440: PetscErrorCode NEPSetWhichEigenpairs(NEP nep,NEPWhich which)
441: {
445:   switch (which) {
446:     case NEP_LARGEST_MAGNITUDE:
447:     case NEP_SMALLEST_MAGNITUDE:
448:     case NEP_LARGEST_REAL:
449:     case NEP_SMALLEST_REAL:
450:     case NEP_LARGEST_IMAGINARY:
451:     case NEP_SMALLEST_IMAGINARY:
452:     case NEP_TARGET_MAGNITUDE:
453:     case NEP_TARGET_REAL:
454: #if defined(PETSC_USE_COMPLEX)
455:     case NEP_TARGET_IMAGINARY:
456: #endif
457:     case NEP_ALL:
458:     case NEP_WHICH_USER:
459:       if (nep->which != which) {
460:         nep->state = NEP_STATE_INITIAL;
461:         nep->which = which;
462:       }
463:       break;
464: #if !defined(PETSC_USE_COMPLEX)
465:     case NEP_TARGET_IMAGINARY:
466:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NEP_TARGET_IMAGINARY can be used only with complex scalars");
467: #endif
468:     default:
469:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
470:   }
471:   return(0);
472: }

474: /*@
475:     NEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
476:     sought.

478:     Not Collective

480:     Input Parameter:
481: .   nep - eigensolver context obtained from NEPCreate()

483:     Output Parameter:
484: .   which - the portion of the spectrum to be sought

486:     Notes:
487:     See NEPSetWhichEigenpairs() for possible values of 'which'.

489:     Level: intermediate

491: .seealso: NEPSetWhichEigenpairs(), NEPWhich
492: @*/
493: PetscErrorCode NEPGetWhichEigenpairs(NEP nep,NEPWhich *which)
494: {
498:   *which = nep->which;
499:   return(0);
500: }

502: /*@C
503:    NEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
504:    when NEPSetWhichEigenpairs() is set to NEP_WHICH_USER.

506:    Logically Collective on nep

508:    Input Parameters:
509: +  nep  - eigensolver context obtained from NEPCreate()
510: .  func - a pointer to the comparison function
511: -  ctx  - a context pointer (the last parameter to the comparison function)

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

516: +   ar     - real part of the 1st eigenvalue
517: .   ai     - imaginary part of the 1st eigenvalue
518: .   br     - real part of the 2nd eigenvalue
519: .   bi     - imaginary part of the 2nd eigenvalue
520: .   res    - result of comparison
521: -   ctx    - optional context, as set by NEPSetEigenvalueComparison()

523:    Note:
524:    The returning parameter 'res' can be
525: +  negative - if the 1st eigenvalue is preferred to the 2st one
526: .  zero     - if both eigenvalues are equally preferred
527: -  positive - if the 2st eigenvalue is preferred to the 1st one

529:    Level: advanced

531: .seealso: NEPSetWhichEigenpairs(), NEPWhich
532: @*/
533: PetscErrorCode NEPSetEigenvalueComparison(NEP nep,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)
534: {
537:   nep->sc->comparison    = func;
538:   nep->sc->comparisonctx = ctx;
539:   nep->which             = NEP_WHICH_USER;
540:   return(0);
541: }

543: /*@
544:    NEPSetProblemType - Specifies the type of the nonlinear eigenvalue problem.

546:    Logically Collective on nep

548:    Input Parameters:
549: +  nep  - the nonlinear eigensolver context
550: -  type - a known type of nonlinear eigenvalue problem

552:    Options Database Keys:
553: +  -nep_general - general problem with no particular structure
554: -  -nep_rational - a rational eigenvalue problem defined in split form with all f_i rational

556:    Notes:
557:    Allowed values for the problem type are: general (NEP_GENERAL), and rational
558:    (NEP_RATIONAL).

560:    This function is used to provide a hint to the NEP solver to exploit certain
561:    properties of the nonlinear eigenproblem. This hint may be used or not,
562:    depending on the solver. By default, no particular structure is assumed.

564:    Level: intermediate

566: .seealso: NEPSetType(), NEPGetProblemType(), NEPProblemType
567: @*/
568: PetscErrorCode NEPSetProblemType(NEP nep,NEPProblemType type)
569: {
573:   if (type!=NEP_GENERAL && type!=NEP_RATIONAL) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
574:   if (type != nep->problem_type) {
575:     nep->problem_type = type;
576:     nep->state = NEP_STATE_INITIAL;
577:   }
578:   return(0);
579: }

581: /*@
582:    NEPGetProblemType - Gets the problem type from the NEP object.

584:    Not Collective

586:    Input Parameter:
587: .  nep - the nonlinear eigensolver context

589:    Output Parameter:
590: .  type - the problem type

592:    Level: intermediate

594: .seealso: NEPSetProblemType(), NEPProblemType
595: @*/
596: PetscErrorCode NEPGetProblemType(NEP nep,NEPProblemType *type)
597: {
601:   *type = nep->problem_type;
602:   return(0);
603: }

605: /*@
606:    NEPSetTwoSided - Sets the solver to use a two-sided variant so that left
607:    eigenvectors are also computed.

609:    Logically Collective on nep

611:    Input Parameters:
612: +  nep      - the eigensolver context
613: -  twosided - whether the two-sided variant is to be used or not

615:    Options Database Keys:
616: .  -nep_two_sided <boolean> - Sets/resets the twosided flag

618:    Notes:
619:    If the user sets twosided=PETSC_TRUE then the solver uses a variant of
620:    the algorithm that computes both right and left eigenvectors. This is
621:    usually much more costly. This option is not available in all solvers.

623:    When using two-sided solvers, the problem matrices must have both the
624:    MatMult and MatMultTranspose operations defined.

626:    Level: advanced

628: .seealso: NEPGetTwoSided(), NEPGetLeftEigenvector()
629: @*/
630: PetscErrorCode NEPSetTwoSided(NEP nep,PetscBool twosided)
631: {
635:   if (twosided!=nep->twosided) {
636:     nep->twosided = twosided;
637:     nep->state    = NEP_STATE_INITIAL;
638:   }
639:   return(0);
640: }

642: /*@
643:    NEPGetTwoSided - Returns the flag indicating whether a two-sided variant
644:    of the algorithm is being used or not.

646:    Not Collective

648:    Input Parameter:
649: .  nep - the eigensolver context

651:    Output Parameter:
652: .  twosided - the returned flag

654:    Level: advanced

656: .seealso: NEPSetTwoSided()
657: @*/
658: PetscErrorCode NEPGetTwoSided(NEP nep,PetscBool *twosided)
659: {
663:   *twosided = nep->twosided;
664:   return(0);
665: }

667: /*@C
668:    NEPSetConvergenceTestFunction - Sets a function to compute the error estimate
669:    used in the convergence test.

671:    Logically Collective on nep

673:    Input Parameters:
674: +  nep     - nonlinear eigensolver context obtained from NEPCreate()
675: .  func    - a pointer to the convergence test function
676: .  ctx     - context for private data for the convergence routine (may be null)
677: -  destroy - a routine for destroying the context (may be null)

679:    Calling Sequence of func:
680: $   func(NEP nep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)

682: +   nep    - nonlinear eigensolver context obtained from NEPCreate()
683: .   eigr   - real part of the eigenvalue
684: .   eigi   - imaginary part of the eigenvalue
685: .   res    - residual norm associated to the eigenpair
686: .   errest - (output) computed error estimate
687: -   ctx    - optional context, as set by NEPSetConvergenceTestFunction()

689:    Note:
690:    If the error estimate returned by the convergence test function is less than
691:    the tolerance, then the eigenvalue is accepted as converged.

693:    Level: advanced

695: .seealso: NEPSetConvergenceTest(), NEPSetTolerances()
696: @*/
697: PetscErrorCode NEPSetConvergenceTestFunction(NEP nep,PetscErrorCode (*func)(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
698: {

703:   if (nep->convergeddestroy) {
704:     (*nep->convergeddestroy)(nep->convergedctx);
705:   }
706:   nep->convergeduser    = func;
707:   nep->convergeddestroy = destroy;
708:   nep->convergedctx     = ctx;
709:   if (func == NEPConvergedRelative) nep->conv = NEP_CONV_REL;
710:   else if (func == NEPConvergedNorm) nep->conv = NEP_CONV_NORM;
711:   else if (func == NEPConvergedAbsolute) nep->conv = NEP_CONV_ABS;
712:   else {
713:     nep->conv      = NEP_CONV_USER;
714:     nep->converged = nep->convergeduser;
715:   }
716:   return(0);
717: }

719: /*@
720:    NEPSetConvergenceTest - Specifies how to compute the error estimate
721:    used in the convergence test.

723:    Logically Collective on nep

725:    Input Parameters:
726: +  nep  - nonlinear eigensolver context obtained from NEPCreate()
727: -  conv - the type of convergence test

729:    Options Database Keys:
730: +  -nep_conv_abs  - Sets the absolute convergence test
731: .  -nep_conv_rel  - Sets the convergence test relative to the eigenvalue
732: -  -nep_conv_user - Selects the user-defined convergence test

734:    Note:
735:    The parameter 'conv' can have one of these values
736: +     NEP_CONV_ABS  - absolute error ||r||
737: .     NEP_CONV_REL  - error relative to the eigenvalue l, ||r||/|l|
738: .     NEP_CONV_NORM - error relative matrix norms, ||r||/sum_i(|f_i(l)|*||A_i||)
739: -     NEP_CONV_USER - function set by NEPSetConvergenceTestFunction()

741:    Level: intermediate

743: .seealso: NEPGetConvergenceTest(), NEPSetConvergenceTestFunction(), NEPSetStoppingTest(), NEPConv
744: @*/
745: PetscErrorCode NEPSetConvergenceTest(NEP nep,NEPConv conv)
746: {
750:   switch (conv) {
751:     case NEP_CONV_ABS:  nep->converged = NEPConvergedAbsolute; break;
752:     case NEP_CONV_REL:  nep->converged = NEPConvergedRelative; break;
753:     case NEP_CONV_NORM: nep->converged = NEPConvergedNorm; break;
754:     case NEP_CONV_USER:
755:       if (!nep->convergeduser) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"Must call NEPSetConvergenceTestFunction() first");
756:       nep->converged = nep->convergeduser;
757:       break;
758:     default:
759:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
760:   }
761:   nep->conv = conv;
762:   return(0);
763: }

765: /*@
766:    NEPGetConvergenceTest - Gets the method used to compute the error estimate
767:    used in the convergence test.

769:    Not Collective

771:    Input Parameters:
772: .  nep   - nonlinear eigensolver context obtained from NEPCreate()

774:    Output Parameters:
775: .  conv  - the type of convergence test

777:    Level: intermediate

779: .seealso: NEPSetConvergenceTest(), NEPConv
780: @*/
781: PetscErrorCode NEPGetConvergenceTest(NEP nep,NEPConv *conv)
782: {
786:   *conv = nep->conv;
787:   return(0);
788: }

790: /*@C
791:    NEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
792:    iteration of the eigensolver.

794:    Logically Collective on nep

796:    Input Parameters:
797: +  nep     - nonlinear eigensolver context obtained from NEPCreate()
798: .  func    - pointer to the stopping test function
799: .  ctx     - context for private data for the stopping routine (may be null)
800: -  destroy - a routine for destroying the context (may be null)

802:    Calling Sequence of func:
803: $   func(NEP nep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,NEPConvergedReason *reason,void *ctx)

805: +   nep    - nonlinear eigensolver context obtained from NEPCreate()
806: .   its    - current number of iterations
807: .   max_it - maximum number of iterations
808: .   nconv  - number of currently converged eigenpairs
809: .   nev    - number of requested eigenpairs
810: .   reason - (output) result of the stopping test
811: -   ctx    - optional context, as set by NEPSetStoppingTestFunction()

813:    Note:
814:    Normal usage is to first call the default routine NEPStoppingBasic() and then
815:    set reason to NEP_CONVERGED_USER if some user-defined conditions have been
816:    met. To let the eigensolver continue iterating, the result must be left as
817:    NEP_CONVERGED_ITERATING.

819:    Level: advanced

821: .seealso: NEPSetStoppingTest(), NEPStoppingBasic()
822: @*/
823: PetscErrorCode NEPSetStoppingTestFunction(NEP nep,PetscErrorCode (*func)(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
824: {

829:   if (nep->stoppingdestroy) {
830:     (*nep->stoppingdestroy)(nep->stoppingctx);
831:   }
832:   nep->stoppinguser    = func;
833:   nep->stoppingdestroy = destroy;
834:   nep->stoppingctx     = ctx;
835:   if (func == NEPStoppingBasic) nep->stop = NEP_STOP_BASIC;
836:   else {
837:     nep->stop     = NEP_STOP_USER;
838:     nep->stopping = nep->stoppinguser;
839:   }
840:   return(0);
841: }

843: /*@
844:    NEPSetStoppingTest - Specifies how to decide the termination of the outer
845:    loop of the eigensolver.

847:    Logically Collective on nep

849:    Input Parameters:
850: +  nep  - nonlinear eigensolver context obtained from NEPCreate()
851: -  stop - the type of stopping test

853:    Options Database Keys:
854: +  -nep_stop_basic - Sets the default stopping test
855: -  -nep_stop_user  - Selects the user-defined stopping test

857:    Note:
858:    The parameter 'stop' can have one of these values
859: +     NEP_STOP_BASIC - default stopping test
860: -     NEP_STOP_USER  - function set by NEPSetStoppingTestFunction()

862:    Level: advanced

864: .seealso: NEPGetStoppingTest(), NEPSetStoppingTestFunction(), NEPSetConvergenceTest(), NEPStop
865: @*/
866: PetscErrorCode NEPSetStoppingTest(NEP nep,NEPStop stop)
867: {
871:   switch (stop) {
872:     case NEP_STOP_BASIC: nep->stopping = NEPStoppingBasic; break;
873:     case NEP_STOP_USER:
874:       if (!nep->stoppinguser) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"Must call NEPSetStoppingTestFunction() first");
875:       nep->stopping = nep->stoppinguser;
876:       break;
877:     default:
878:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
879:   }
880:   nep->stop = stop;
881:   return(0);
882: }

884: /*@
885:    NEPGetStoppingTest - Gets the method used to decide the termination of the outer
886:    loop of the eigensolver.

888:    Not Collective

890:    Input Parameters:
891: .  nep   - nonlinear eigensolver context obtained from NEPCreate()

893:    Output Parameters:
894: .  stop  - the type of stopping test

896:    Level: advanced

898: .seealso: NEPSetStoppingTest(), NEPStop
899: @*/
900: PetscErrorCode NEPGetStoppingTest(NEP nep,NEPStop *stop)
901: {
905:   *stop = nep->stop;
906:   return(0);
907: }

909: /*@
910:    NEPSetTrackAll - Specifies if the solver must compute the residual of all
911:    approximate eigenpairs or not.

913:    Logically Collective on nep

915:    Input Parameters:
916: +  nep      - the eigensolver context
917: -  trackall - whether compute all residuals or not

919:    Notes:
920:    If the user sets trackall=PETSC_TRUE then the solver explicitly computes
921:    the residual for each eigenpair approximation. Computing the residual is
922:    usually an expensive operation and solvers commonly compute the associated
923:    residual to the first unconverged eigenpair.

925:    The option '-nep_monitor_all' automatically activates this option.

927:    Level: developer

929: .seealso: NEPGetTrackAll()
930: @*/
931: PetscErrorCode NEPSetTrackAll(NEP nep,PetscBool trackall)
932: {
936:   nep->trackall = trackall;
937:   return(0);
938: }

940: /*@
941:    NEPGetTrackAll - Returns the flag indicating whether all residual norms must
942:    be computed or not.

944:    Not Collective

946:    Input Parameter:
947: .  nep - the eigensolver context

949:    Output Parameter:
950: .  trackall - the returned flag

952:    Level: developer

954: .seealso: NEPSetTrackAll()
955: @*/
956: PetscErrorCode NEPGetTrackAll(NEP nep,PetscBool *trackall)
957: {
961:   *trackall = nep->trackall;
962:   return(0);
963: }

965: /*@
966:    NEPSetRefine - Specifies the refinement type (and options) to be used
967:    after the solve.

969:    Logically Collective on nep

971:    Input Parameters:
972: +  nep    - the nonlinear eigensolver context
973: .  refine - refinement type
974: .  npart  - number of partitions of the communicator
975: .  tol    - the convergence tolerance
976: .  its    - maximum number of refinement iterations
977: -  scheme - which scheme to be used for solving the involved linear systems

979:    Options Database Keys:
980: +  -nep_refine <type> - refinement type, one of <none,simple,multiple>
981: .  -nep_refine_partitions <n> - the number of partitions
982: .  -nep_refine_tol <tol> - the tolerance
983: .  -nep_refine_its <its> - number of iterations
984: -  -nep_refine_scheme - to set the scheme for the linear solves

986:    Notes:
987:    By default, iterative refinement is disabled, since it may be very
988:    costly. There are two possible refinement strategies: simple and multiple.
989:    The simple approach performs iterative refinement on each of the
990:    converged eigenpairs individually, whereas the multiple strategy works
991:    with the invariant pair as a whole, refining all eigenpairs simultaneously.
992:    The latter may be required for the case of multiple eigenvalues.

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

999:    The tol and its parameters specify the stopping criterion. In the simple
1000:    method, refinement continues until the residual of each eigenpair is
1001:    below the tolerance (tol defaults to the NEP tol, but may be set to a
1002:    different value). In contrast, the multiple method simply performs its
1003:    refinement iterations (just one by default).

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

1009:    Level: intermediate

1011: .seealso: NEPGetRefine()
1012: @*/
1013: PetscErrorCode NEPSetRefine(NEP nep,NEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,NEPRefineScheme scheme)
1014: {
1016:   PetscMPIInt    size;

1025:   nep->refine = refine;
1026:   if (refine) {  /* process parameters only if not REFINE_NONE */
1027:     if (npart!=nep->npart) {
1028:       PetscSubcommDestroy(&nep->refinesubc);
1029:       KSPDestroy(&nep->refineksp);
1030:     }
1031:     if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
1032:       nep->npart = 1;
1033:     } else {
1034:       MPI_Comm_size(PetscObjectComm((PetscObject)nep),&size);
1035:       if (npart<1 || npart>size) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
1036:       nep->npart = npart;
1037:     }
1038:     if (tol == PETSC_DEFAULT || tol == PETSC_DECIDE) {
1039:       nep->rtol = PETSC_DEFAULT;
1040:     } else {
1041:       if (tol<=0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1042:       nep->rtol = tol;
1043:     }
1044:     if (its==PETSC_DECIDE || its==PETSC_DEFAULT) {
1045:       nep->rits = PETSC_DEFAULT;
1046:     } else {
1047:       if (its<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be >= 0");
1048:       nep->rits = its;
1049:     }
1050:     nep->scheme = scheme;
1051:   }
1052:   nep->state = NEP_STATE_INITIAL;
1053:   return(0);
1054: }

1056: /*@C
1057:    NEPGetRefine - Gets the refinement strategy used by the NEP object, and the
1058:    associated parameters.

1060:    Not Collective

1062:    Input Parameter:
1063: .  nep - the nonlinear eigensolver context

1065:    Output Parameters:
1066: +  refine - refinement type
1067: .  npart  - number of partitions of the communicator
1068: .  tol    - the convergence tolerance
1069: .  its    - maximum number of refinement iterations
1070: -  scheme - the scheme used for solving linear systems

1072:    Level: intermediate

1074:    Note:
1075:    The user can specify NULL for any parameter that is not needed.

1077: .seealso: NEPSetRefine()
1078: @*/
1079: PetscErrorCode NEPGetRefine(NEP nep,NEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,NEPRefineScheme *scheme)
1080: {
1083:   if (refine) *refine = nep->refine;
1084:   if (npart)  *npart  = nep->npart;
1085:   if (tol)    *tol    = nep->rtol;
1086:   if (its)    *its    = nep->rits;
1087:   if (scheme) *scheme = nep->scheme;
1088:   return(0);
1089: }

1091: /*@C
1092:    NEPSetOptionsPrefix - Sets the prefix used for searching for all
1093:    NEP options in the database.

1095:    Logically Collective on nep

1097:    Input Parameters:
1098: +  nep - the nonlinear eigensolver context
1099: -  prefix - the prefix string to prepend to all NEP option requests

1101:    Notes:
1102:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1103:    The first character of all runtime options is AUTOMATICALLY the
1104:    hyphen.

1106:    For example, to distinguish between the runtime options for two
1107:    different NEP contexts, one could call
1108: .vb
1109:       NEPSetOptionsPrefix(nep1,"neig1_")
1110:       NEPSetOptionsPrefix(nep2,"neig2_")
1111: .ve

1113:    Level: advanced

1115: .seealso: NEPAppendOptionsPrefix(), NEPGetOptionsPrefix()
1116: @*/
1117: PetscErrorCode NEPSetOptionsPrefix(NEP nep,const char *prefix)
1118: {

1123:   if (!nep->V) { NEPGetBV(nep,&nep->V); }
1124:   BVSetOptionsPrefix(nep->V,prefix);
1125:   if (!nep->ds) { NEPGetDS(nep,&nep->ds); }
1126:   DSSetOptionsPrefix(nep->ds,prefix);
1127:   if (!nep->rg) { NEPGetRG(nep,&nep->rg); }
1128:   RGSetOptionsPrefix(nep->rg,prefix);
1129:   PetscObjectSetOptionsPrefix((PetscObject)nep,prefix);
1130:   return(0);
1131: }

1133: /*@C
1134:    NEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1135:    NEP options in the database.

1137:    Logically Collective on nep

1139:    Input Parameters:
1140: +  nep - the nonlinear eigensolver context
1141: -  prefix - the prefix string to prepend to all NEP option requests

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

1147:    Level: advanced

1149: .seealso: NEPSetOptionsPrefix(), NEPGetOptionsPrefix()
1150: @*/
1151: PetscErrorCode NEPAppendOptionsPrefix(NEP nep,const char *prefix)
1152: {

1157:   if (!nep->V) { NEPGetBV(nep,&nep->V); }
1158:   BVAppendOptionsPrefix(nep->V,prefix);
1159:   if (!nep->ds) { NEPGetDS(nep,&nep->ds); }
1160:   DSAppendOptionsPrefix(nep->ds,prefix);
1161:   if (!nep->rg) { NEPGetRG(nep,&nep->rg); }
1162:   RGAppendOptionsPrefix(nep->rg,prefix);
1163:   PetscObjectAppendOptionsPrefix((PetscObject)nep,prefix);
1164:   return(0);
1165: }

1167: /*@C
1168:    NEPGetOptionsPrefix - Gets the prefix used for searching for all
1169:    NEP options in the database.

1171:    Not Collective

1173:    Input Parameters:
1174: .  nep - the nonlinear eigensolver context

1176:    Output Parameters:
1177: .  prefix - pointer to the prefix string used is returned

1179:    Note:
1180:    On the Fortran side, the user should pass in a string 'prefix' of
1181:    sufficient length to hold the prefix.

1183:    Level: advanced

1185: .seealso: NEPSetOptionsPrefix(), NEPAppendOptionsPrefix()
1186: @*/
1187: PetscErrorCode NEPGetOptionsPrefix(NEP nep,const char *prefix[])
1188: {

1194:   PetscObjectGetOptionsPrefix((PetscObject)nep,prefix);
1195:   return(0);
1196: }