Actual source code: nepopts.c

slepc-3.16.0 2021-09-30
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:   PetscFunctionListFind(NEPMonitorCreateList,key,&cfunc);
 56:   PetscFunctionListFind(NEPMonitorDestroyList,key,&dfunc);
 57:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
 58:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

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

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

 74:    Collective on nep

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

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

227:    Not Collective

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

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

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

239:    Level: intermediate

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

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

256:    Logically Collective on nep

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

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

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

270:    Level: intermediate

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

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

301:    Not Collective

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

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

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

314:    Level: intermediate

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

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

332:    Logically Collective on nep

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

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

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

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

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

358:    Level: intermediate

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

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

391:     Logically Collective on nep

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

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

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

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

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

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

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

435:     Level: intermediate

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

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

477:     Not Collective

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

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

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

488:     Level: intermediate

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

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

505:    Logically Collective on nep

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

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

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

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

528:    Level: advanced

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

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

545:    Logically Collective on nep

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

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

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

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

563:    Level: intermediate

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

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

583:    Not Collective

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

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

591:    Level: intermediate

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

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

608:    Logically Collective on nep

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

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

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

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

625:    Level: advanced

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

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

645:    Not Collective

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

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

653:    Level: advanced

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

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

670:    Logically Collective on nep

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

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

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

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

692:    Level: advanced

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

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

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

722:    Logically Collective on nep

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

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

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

740:    Level: intermediate

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

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

768:    Not Collective

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

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

776:    Level: intermediate

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

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

793:    Logically Collective on nep

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

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

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

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

818:    Level: advanced

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

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

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

846:    Logically Collective on nep

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

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

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

861:    Level: advanced

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

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

887:    Not Collective

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

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

895:    Level: advanced

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

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

912:    Logically Collective on nep

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

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

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

926:    Level: developer

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

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

943:    Not Collective

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

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

951:    Level: developer

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

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

968:    Logically Collective on nep

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

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

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

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

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

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

1008:    Level: intermediate

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

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

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

1059:    Not Collective

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

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

1071:    Level: intermediate

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

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

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

1094:    Logically Collective on nep

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

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

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

1112:    Level: advanced

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

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

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

1136:    Logically Collective on nep

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

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

1146:    Level: advanced

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

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

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

1170:    Not Collective

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

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

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

1182:    Level: advanced

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

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