Actual source code: nepopts.c

  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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

 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;

 47:   PetscFunctionBegin;
 48:   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->options,((PetscObject)nep)->prefix,opt,&viewer,&format,&flg));
 49:   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);

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

 60:   PetscCall((*cfunc)(viewer,format,ctx,&vf));
 61:   PetscCall(PetscViewerDestroy(&viewer));
 62:   PetscCall(NEPMonitorSet(nep,mfunc,vf,(PetscCtxDestroyFn*)dfunc));
 63:   if (trackall) PetscCall(NEPSetTrackAll(nep,PETSC_TRUE));
 64:   PetscFunctionReturn(PETSC_SUCCESS);
 65: }

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

 72:    Collective

 74:    Input Parameters:
 75: .  nep - the nonlinear eigensolver context

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

 80:    Level: beginner

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

 94:   PetscFunctionBegin;
 96:   PetscCall(NEPRegisterAll());
 97:   PetscObjectOptionsBegin((PetscObject)nep);
 98:     PetscCall(PetscOptionsFList("-nep_type","Nonlinear eigensolver method","NEPSetType",NEPList,(char*)(((PetscObject)nep)->type_name?((PetscObject)nep)->type_name:NEPRII),type,sizeof(type),&flg));
 99:     if (flg) PetscCall(NEPSetType(nep,type));
100:     else if (!((PetscObject)nep)->type_name) PetscCall(NEPSetType(nep,NEPRII));

102:     PetscCall(PetscOptionsBoolGroupBegin("-nep_general","General nonlinear eigenvalue problem","NEPSetProblemType",&flg));
103:     if (flg) PetscCall(NEPSetProblemType(nep,NEP_GENERAL));
104:     PetscCall(PetscOptionsBoolGroupEnd("-nep_rational","Rational eigenvalue problem","NEPSetProblemType",&flg));
105:     if (flg) PetscCall(NEPSetProblemType(nep,NEP_RATIONAL));

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

119:     i = nep->max_it;
120:     PetscCall(PetscOptionsInt("-nep_max_it","Maximum number of iterations","NEPSetTolerances",nep->max_it,&i,&flg1));
121:     r = nep->tol;
122:     PetscCall(PetscOptionsReal("-nep_tol","Tolerance","NEPSetTolerances",SlepcDefaultTol(nep->tol),&r,&flg2));
123:     if (flg1 || flg2) PetscCall(NEPSetTolerances(nep,r,i));

125:     PetscCall(PetscOptionsBoolGroupBegin("-nep_conv_rel","Relative error convergence test","NEPSetConvergenceTest",&flg));
126:     if (flg) PetscCall(NEPSetConvergenceTest(nep,NEP_CONV_REL));
127:     PetscCall(PetscOptionsBoolGroup("-nep_conv_norm","Convergence test relative to the matrix norms","NEPSetConvergenceTest",&flg));
128:     if (flg) PetscCall(NEPSetConvergenceTest(nep,NEP_CONV_NORM));
129:     PetscCall(PetscOptionsBoolGroup("-nep_conv_abs","Absolute error convergence test","NEPSetConvergenceTest",&flg));
130:     if (flg) PetscCall(NEPSetConvergenceTest(nep,NEP_CONV_ABS));
131:     PetscCall(PetscOptionsBoolGroupEnd("-nep_conv_user","User-defined convergence test","NEPSetConvergenceTest",&flg));
132:     if (flg) PetscCall(NEPSetConvergenceTest(nep,NEP_CONV_USER));

134:     PetscCall(PetscOptionsBoolGroupBegin("-nep_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","NEPSetStoppingTest",&flg));
135:     if (flg) PetscCall(NEPSetStoppingTest(nep,NEP_STOP_BASIC));
136:     PetscCall(PetscOptionsBoolGroupEnd("-nep_stop_user","User-defined stopping test","NEPSetStoppingTest",&flg));
137:     if (flg) PetscCall(NEPSetStoppingTest(nep,NEP_STOP_USER));

139:     i = nep->nev;
140:     PetscCall(PetscOptionsInt("-nep_nev","Number of eigenvalues to compute","NEPSetDimensions",nep->nev,&i,&flg1));
141:     j = nep->ncv;
142:     PetscCall(PetscOptionsInt("-nep_ncv","Number of basis vectors","NEPSetDimensions",nep->ncv,&j,&flg2));
143:     k = nep->mpd;
144:     PetscCall(PetscOptionsInt("-nep_mpd","Maximum dimension of projected problem","NEPSetDimensions",nep->mpd,&k,&flg3));
145:     if (flg1 || flg2 || flg3) PetscCall(NEPSetDimensions(nep,i,j,k));

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

168:     PetscCall(PetscOptionsScalar("-nep_target","Value of the target","NEPSetTarget",nep->target,&s,&flg));
169:     if (flg) {
170:       if (nep->which!=NEP_TARGET_REAL && nep->which!=NEP_TARGET_IMAGINARY) PetscCall(NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE));
171:       PetscCall(NEPSetTarget(nep,s));
172:     }

174:     PetscCall(PetscOptionsBool("-nep_two_sided","Use two-sided variant (to compute left eigenvectors)","NEPSetTwoSided",nep->twosided,&bval,&flg));
175:     if (flg) PetscCall(NEPSetTwoSided(nep,bval));

177:     /* -----------------------------------------------------------------------*/
178:     /*
179:       Cancels all monitors hardwired into code before call to NEPSetFromOptions()
180:     */
181:     PetscCall(PetscOptionsBool("-nep_monitor_cancel","Remove any hardwired monitor routines","NEPMonitorCancel",PETSC_FALSE,&flg,&set));
182:     if (set && flg) PetscCall(NEPMonitorCancel(nep));
183:     PetscCall(NEPMonitorSetFromOptions(nep,"-nep_monitor","first_approximation",NULL,PETSC_FALSE));
184:     PetscCall(NEPMonitorSetFromOptions(nep,"-nep_monitor_all","all_approximations",NULL,PETSC_TRUE));
185:     PetscCall(NEPMonitorSetFromOptions(nep,"-nep_monitor_conv","convergence_history",NULL,PETSC_FALSE));

187:     /* -----------------------------------------------------------------------*/
188:     PetscCall(PetscOptionsName("-nep_view","Print detailed information on solver used","NEPView",&set));
189:     PetscCall(PetscOptionsName("-nep_view_vectors","View computed eigenvectors","NEPVectorsView",&set));
190:     PetscCall(PetscOptionsName("-nep_view_values","View computed eigenvalues","NEPValuesView",&set));
191:     PetscCall(PetscOptionsName("-nep_converged_reason","Print reason for convergence, and number of iterations","NEPConvergedReasonView",&set));
192:     PetscCall(PetscOptionsName("-nep_error_absolute","Print absolute errors of each eigenpair","NEPErrorView",&set));
193:     PetscCall(PetscOptionsName("-nep_error_relative","Print relative errors of each eigenpair","NEPErrorView",&set));

195:     PetscTryTypeMethod(nep,setfromoptions,PetscOptionsObject);
196:     PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)nep,PetscOptionsObject));
197:   PetscOptionsEnd();

199:   if (!nep->V) PetscCall(NEPGetBV(nep,&nep->V));
200:   PetscCall(BVSetFromOptions(nep->V));
201:   if (!nep->rg) PetscCall(NEPGetRG(nep,&nep->rg));
202:   PetscCall(RGSetFromOptions(nep->rg));
203:   if (nep->useds) {
204:     if (!nep->ds) PetscCall(NEPGetDS(nep,&nep->ds));
205:     PetscCall(NEPSetDSType(nep));
206:     PetscCall(DSSetFromOptions(nep->ds));
207:   }
208:   if (!nep->refineksp) PetscCall(NEPRefineGetKSP(nep,&nep->refineksp));
209:   PetscCall(KSPSetFromOptions(nep->refineksp));
210:   if (nep->fui==NEP_USER_INTERFACE_SPLIT) for (i=0;i<nep->nt;i++) PetscCall(FNSetFromOptions(nep->f[i]));
211:   PetscFunctionReturn(PETSC_SUCCESS);
212: }

214: /*@
215:    NEPGetTolerances - Gets the tolerance and maximum iteration count used
216:    by the NEP convergence tests.

218:    Not Collective

220:    Input Parameter:
221: .  nep - the nonlinear eigensolver context

223:    Output Parameters:
224: +  tol - the convergence tolerance
225: -  maxits - maximum number of iterations

227:    Notes:
228:    The user can specify NULL for any parameter that is not needed.

230:    Level: intermediate

232: .seealso: `NEPSetTolerances()`
233: @*/
234: PetscErrorCode NEPGetTolerances(NEP nep,PetscReal *tol,PetscInt *maxits)
235: {
236:   PetscFunctionBegin;
238:   if (tol)    *tol    = nep->tol;
239:   if (maxits) *maxits = nep->max_it;
240:   PetscFunctionReturn(PETSC_SUCCESS);
241: }

243: /*@
244:    NEPSetTolerances - Sets the tolerance and maximum iteration count used
245:    by the NEP convergence tests.

247:    Logically Collective

249:    Input Parameters:
250: +  nep    - the nonlinear eigensolver context
251: .  tol    - the convergence tolerance
252: -  maxits - maximum number of iterations to use

254:    Options Database Keys:
255: +  -nep_tol <tol> - Sets the convergence tolerance
256: -  -nep_max_it <maxits> - Sets the maximum number of iterations allowed

258:    Notes:
259:    Use PETSC_CURRENT to retain the current value of any of the parameters.
260:    Use PETSC_DETERMINE for either argument to assign a default value computed
261:    internally (may be different in each solver).
262:    For maxits use PETSC_UMLIMITED to indicate there is no upper bound on this value.

264:    Level: intermediate

266: .seealso: `NEPGetTolerances()`
267: @*/
268: PetscErrorCode NEPSetTolerances(NEP nep,PetscReal tol,PetscInt maxits)
269: {
270:   PetscFunctionBegin;
274:   if (tol == (PetscReal)PETSC_DETERMINE) {
275:     nep->tol   = PETSC_DETERMINE;
276:     nep->state = NEP_STATE_INITIAL;
277:   } else if (tol != (PetscReal)PETSC_CURRENT) {
278:     PetscCheck(tol>0.0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
279:     nep->tol = tol;
280:   }
281:   if (maxits == PETSC_DETERMINE) {
282:     nep->max_it = PETSC_DETERMINE;
283:     nep->state  = NEP_STATE_INITIAL;
284:   } else if (maxits == PETSC_UNLIMITED) {
285:     nep->max_it = PETSC_INT_MAX;
286:   } else if (maxits != PETSC_CURRENT) {
287:     PetscCheck(maxits>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
288:     nep->max_it = maxits;
289:   }
290:   PetscFunctionReturn(PETSC_SUCCESS);
291: }

293: /*@
294:    NEPGetDimensions - Gets the number of eigenvalues to compute
295:    and the dimension of the subspace.

297:    Not Collective

299:    Input Parameter:
300: .  nep - the nonlinear eigensolver context

302:    Output Parameters:
303: +  nev - number of eigenvalues to compute
304: .  ncv - the maximum dimension of the subspace to be used by the solver
305: -  mpd - the maximum dimension allowed for the projected problem

307:    Notes:
308:    The user can specify NULL for any parameter that is not needed.

310:    Level: intermediate

312: .seealso: `NEPSetDimensions()`
313: @*/
314: PetscErrorCode NEPGetDimensions(NEP nep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
315: {
316:   PetscFunctionBegin;
318:   if (nev) *nev = nep->nev;
319:   if (ncv) *ncv = nep->ncv;
320:   if (mpd) *mpd = nep->mpd;
321:   PetscFunctionReturn(PETSC_SUCCESS);
322: }

324: /*@
325:    NEPSetDimensions - Sets the number of eigenvalues to compute
326:    and the dimension of the subspace.

328:    Logically Collective

330:    Input Parameters:
331: +  nep - the nonlinear eigensolver context
332: .  nev - number of eigenvalues to compute
333: .  ncv - the maximum dimension of the subspace to be used by the solver
334: -  mpd - the maximum dimension allowed for the projected problem

336:    Options Database Keys:
337: +  -nep_nev <nev> - Sets the number of eigenvalues
338: .  -nep_ncv <ncv> - Sets the dimension of the subspace
339: -  -nep_mpd <mpd> - Sets the maximum projected dimension

341:    Notes:
342:    Use PETSC_DETERMINE for ncv and mpd to assign a reasonably good value, which is
343:    dependent on the solution method. For any of the arguments, use PETSC_CURRENT
344:    to preserve the current value.

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

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

355:    Level: intermediate

357: .seealso: `NEPGetDimensions()`
358: @*/
359: PetscErrorCode NEPSetDimensions(NEP nep,PetscInt nev,PetscInt ncv,PetscInt mpd)
360: {
361:   PetscFunctionBegin;
366:   if (nev != PETSC_CURRENT) {
367:     PetscCheck(nev>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
368:     nep->nev = nev;
369:   }
370:   if (ncv == PETSC_DETERMINE) {
371:     nep->ncv = PETSC_DETERMINE;
372:   } else if (ncv != PETSC_CURRENT) {
373:     PetscCheck(ncv>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
374:     nep->ncv = ncv;
375:   }
376:   if (mpd == PETSC_DETERMINE) {
377:     nep->mpd = PETSC_DETERMINE;
378:   } else if (mpd != PETSC_CURRENT) {
379:     PetscCheck(mpd>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
380:     nep->mpd = mpd;
381:   }
382:   nep->state = NEP_STATE_INITIAL;
383:   PetscFunctionReturn(PETSC_SUCCESS);
384: }

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

390:     Logically Collective

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

396:     Options Database Keys:
397: +   -nep_largest_magnitude - Sets largest eigenvalues in magnitude
398: .   -nep_smallest_magnitude - Sets smallest eigenvalues in magnitude
399: .   -nep_largest_real - Sets largest real parts
400: .   -nep_smallest_real - Sets smallest real parts
401: .   -nep_largest_imaginary - Sets largest imaginary parts
402: .   -nep_smallest_imaginary - Sets smallest imaginary parts
403: .   -nep_target_magnitude - Sets eigenvalues closest to target
404: .   -nep_target_real - Sets real parts closest to target
405: .   -nep_target_imaginary - Sets imaginary parts closest to target
406: -   -nep_all - Sets all eigenvalues in a region

408:     Notes:
409:     The parameter 'which' can have one of these values

411: +     NEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
412: .     NEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
413: .     NEP_LARGEST_REAL - largest real parts
414: .     NEP_SMALLEST_REAL - smallest real parts
415: .     NEP_LARGEST_IMAGINARY - largest imaginary parts
416: .     NEP_SMALLEST_IMAGINARY - smallest imaginary parts
417: .     NEP_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
418: .     NEP_TARGET_REAL - eigenvalues with real part closest to target
419: .     NEP_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
420: .     NEP_ALL - all eigenvalues contained in a given region
421: -     NEP_WHICH_USER - user defined ordering set with NEPSetEigenvalueComparison()

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

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

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

433:     Level: intermediate

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

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

475:     Not Collective

477:     Input Parameter:
478: .   nep - eigensolver context obtained from NEPCreate()

480:     Output Parameter:
481: .   which - the portion of the spectrum to be sought

483:     Notes:
484:     See NEPSetWhichEigenpairs() for possible values of 'which'.

486:     Level: intermediate

488: .seealso: `NEPSetWhichEigenpairs()`, `NEPWhich`
489: @*/
490: PetscErrorCode NEPGetWhichEigenpairs(NEP nep,NEPWhich *which)
491: {
492:   PetscFunctionBegin;
494:   PetscAssertPointer(which,2);
495:   *which = nep->which;
496:   PetscFunctionReturn(PETSC_SUCCESS);
497: }

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

503:    Logically Collective

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

510:    Note:
511:    The returning parameter 'res' can be
512: +  negative - if the 1st eigenvalue is preferred to the 2st one
513: .  zero     - if both eigenvalues are equally preferred
514: -  positive - if the 2st eigenvalue is preferred to the 1st one

516:    Level: advanced

518: .seealso: `NEPSetWhichEigenpairs()`, `NEPWhich`
519: @*/
520: PetscErrorCode NEPSetEigenvalueComparison(NEP nep,SlepcEigenvalueComparisonFn *comp,void *ctx)
521: {
522:   PetscFunctionBegin;
524:   nep->sc->comparison    = comp;
525:   nep->sc->comparisonctx = ctx;
526:   nep->which             = NEP_WHICH_USER;
527:   PetscFunctionReturn(PETSC_SUCCESS);
528: }

530: /*@
531:    NEPSetProblemType - Specifies the type of the nonlinear eigenvalue problem.

533:    Logically Collective

535:    Input Parameters:
536: +  nep  - the nonlinear eigensolver context
537: -  type - a known type of nonlinear eigenvalue problem

539:    Options Database Keys:
540: +  -nep_general - general problem with no particular structure
541: -  -nep_rational - a rational eigenvalue problem defined in split form with all f_i rational

543:    Notes:
544:    Allowed values for the problem type are general (NEP_GENERAL), and rational
545:    (NEP_RATIONAL).

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

551:    Level: intermediate

553: .seealso: `NEPSetType()`, `NEPGetProblemType()`, `NEPProblemType`
554: @*/
555: PetscErrorCode NEPSetProblemType(NEP nep,NEPProblemType type)
556: {
557:   PetscFunctionBegin;
560:   PetscCheck(type==NEP_GENERAL || type==NEP_RATIONAL,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
561:   if (type != nep->problem_type) {
562:     nep->problem_type = type;
563:     nep->state = NEP_STATE_INITIAL;
564:   }
565:   PetscFunctionReturn(PETSC_SUCCESS);
566: }

568: /*@
569:    NEPGetProblemType - Gets the problem type from the NEP object.

571:    Not Collective

573:    Input Parameter:
574: .  nep - the nonlinear eigensolver context

576:    Output Parameter:
577: .  type - the problem type

579:    Level: intermediate

581: .seealso: `NEPSetProblemType()`, `NEPProblemType`
582: @*/
583: PetscErrorCode NEPGetProblemType(NEP nep,NEPProblemType *type)
584: {
585:   PetscFunctionBegin;
587:   PetscAssertPointer(type,2);
588:   *type = nep->problem_type;
589:   PetscFunctionReturn(PETSC_SUCCESS);
590: }

592: /*@
593:    NEPSetTwoSided - Sets the solver to use a two-sided variant so that left
594:    eigenvectors are also computed.

596:    Logically Collective

598:    Input Parameters:
599: +  nep      - the eigensolver context
600: -  twosided - whether the two-sided variant is to be used or not

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

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

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

613:    Level: advanced

615: .seealso: `NEPGetTwoSided()`, `NEPGetLeftEigenvector()`
616: @*/
617: PetscErrorCode NEPSetTwoSided(NEP nep,PetscBool twosided)
618: {
619:   PetscFunctionBegin;
622:   if (twosided!=nep->twosided) {
623:     nep->twosided = twosided;
624:     nep->state    = NEP_STATE_INITIAL;
625:   }
626:   PetscFunctionReturn(PETSC_SUCCESS);
627: }

629: /*@
630:    NEPGetTwoSided - Returns the flag indicating whether a two-sided variant
631:    of the algorithm is being used or not.

633:    Not Collective

635:    Input Parameter:
636: .  nep - the eigensolver context

638:    Output Parameter:
639: .  twosided - the returned flag

641:    Level: advanced

643: .seealso: `NEPSetTwoSided()`
644: @*/
645: PetscErrorCode NEPGetTwoSided(NEP nep,PetscBool *twosided)
646: {
647:   PetscFunctionBegin;
649:   PetscAssertPointer(twosided,2);
650:   *twosided = nep->twosided;
651:   PetscFunctionReturn(PETSC_SUCCESS);
652: }

654: /*@C
655:    NEPSetConvergenceTestFunction - Sets a function to compute the error estimate
656:    used in the convergence test.

658:    Logically Collective

660:    Input Parameters:
661: +  nep     - nonlinear eigensolver context obtained from NEPCreate()
662: .  conv    - convergence test function, see NEPConvergenceTestFn for the calling sequence
663: .  ctx     - context for private data for the convergence routine (may be NULL)
664: -  destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence

666:    Note:
667:    If the error estimate returned by the convergence test function is less than
668:    the tolerance, then the eigenvalue is accepted as converged.

670:    Level: advanced

672: .seealso: `NEPSetConvergenceTest()`, `NEPSetTolerances()`
673: @*/
674: PetscErrorCode NEPSetConvergenceTestFunction(NEP nep,NEPConvergenceTestFn *conv,void *ctx,PetscCtxDestroyFn *destroy)
675: {
676:   PetscFunctionBegin;
678:   if (nep->convergeddestroy) PetscCall((*nep->convergeddestroy)(&nep->convergedctx));
679:   nep->convergeduser    = conv;
680:   nep->convergeddestroy = destroy;
681:   nep->convergedctx     = ctx;
682:   if (conv == NEPConvergedRelative) nep->conv = NEP_CONV_REL;
683:   else if (conv == NEPConvergedNorm) nep->conv = NEP_CONV_NORM;
684:   else if (conv == NEPConvergedAbsolute) nep->conv = NEP_CONV_ABS;
685:   else {
686:     nep->conv      = NEP_CONV_USER;
687:     nep->converged = nep->convergeduser;
688:   }
689:   PetscFunctionReturn(PETSC_SUCCESS);
690: }

692: /*@
693:    NEPSetConvergenceTest - Specifies how to compute the error estimate
694:    used in the convergence test.

696:    Logically Collective

698:    Input Parameters:
699: +  nep  - nonlinear eigensolver context obtained from NEPCreate()
700: -  conv - the type of convergence test

702:    Options Database Keys:
703: +  -nep_conv_abs  - Sets the absolute convergence test
704: .  -nep_conv_rel  - Sets the convergence test relative to the eigenvalue
705: -  -nep_conv_user - Selects the user-defined convergence test

707:    Note:
708:    The parameter 'conv' can have one of these values
709: +     NEP_CONV_ABS  - absolute error ||r||
710: .     NEP_CONV_REL  - error relative to the eigenvalue l, ||r||/|l|
711: .     NEP_CONV_NORM - error relative matrix norms, ||r||/sum_i(|f_i(l)|*||A_i||)
712: -     NEP_CONV_USER - function set by NEPSetConvergenceTestFunction()

714:    Level: intermediate

716: .seealso: `NEPGetConvergenceTest()`, `NEPSetConvergenceTestFunction()`, `NEPSetStoppingTest()`, `NEPConv`
717: @*/
718: PetscErrorCode NEPSetConvergenceTest(NEP nep,NEPConv conv)
719: {
720:   PetscFunctionBegin;
723:   switch (conv) {
724:     case NEP_CONV_ABS:  nep->converged = NEPConvergedAbsolute; break;
725:     case NEP_CONV_REL:  nep->converged = NEPConvergedRelative; break;
726:     case NEP_CONV_NORM: nep->converged = NEPConvergedNorm; break;
727:     case NEP_CONV_USER:
728:       PetscCheck(nep->convergeduser,PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"Must call NEPSetConvergenceTestFunction() first");
729:       nep->converged = nep->convergeduser;
730:       break;
731:     default:
732:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
733:   }
734:   nep->conv = conv;
735:   PetscFunctionReturn(PETSC_SUCCESS);
736: }

738: /*@
739:    NEPGetConvergenceTest - Gets the method used to compute the error estimate
740:    used in the convergence test.

742:    Not Collective

744:    Input Parameters:
745: .  nep   - nonlinear eigensolver context obtained from NEPCreate()

747:    Output Parameters:
748: .  conv  - the type of convergence test

750:    Level: intermediate

752: .seealso: `NEPSetConvergenceTest()`, `NEPConv`
753: @*/
754: PetscErrorCode NEPGetConvergenceTest(NEP nep,NEPConv *conv)
755: {
756:   PetscFunctionBegin;
758:   PetscAssertPointer(conv,2);
759:   *conv = nep->conv;
760:   PetscFunctionReturn(PETSC_SUCCESS);
761: }

763: /*@C
764:    NEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
765:    iteration of the eigensolver.

767:    Logically Collective

769:    Input Parameters:
770: +  nep     - nonlinear eigensolver context obtained from NEPCreate()
771: .  stop    - the stopping test function, see NEPStoppingTestFn for the calling sequence
772: .  ctx     - context for private data for the stopping routine (may be NULL)
773: -  destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence

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

781:    Level: advanced

783: .seealso: `NEPSetStoppingTest()`, `NEPStoppingBasic()`
784: @*/
785: PetscErrorCode NEPSetStoppingTestFunction(NEP nep,NEPStoppingTestFn *stop,void *ctx,PetscCtxDestroyFn *destroy)
786: {
787:   PetscFunctionBegin;
789:   if (nep->stoppingdestroy) PetscCall((*nep->stoppingdestroy)(&nep->stoppingctx));
790:   nep->stoppinguser    = stop;
791:   nep->stoppingdestroy = destroy;
792:   nep->stoppingctx     = ctx;
793:   if (stop == NEPStoppingBasic) nep->stop = NEP_STOP_BASIC;
794:   else {
795:     nep->stop     = NEP_STOP_USER;
796:     nep->stopping = nep->stoppinguser;
797:   }
798:   PetscFunctionReturn(PETSC_SUCCESS);
799: }

801: /*@
802:    NEPSetStoppingTest - Specifies how to decide the termination of the outer
803:    loop of the eigensolver.

805:    Logically Collective

807:    Input Parameters:
808: +  nep  - nonlinear eigensolver context obtained from NEPCreate()
809: -  stop - the type of stopping test

811:    Options Database Keys:
812: +  -nep_stop_basic - Sets the default stopping test
813: -  -nep_stop_user  - Selects the user-defined stopping test

815:    Note:
816:    The parameter 'stop' can have one of these values
817: +     NEP_STOP_BASIC - default stopping test
818: -     NEP_STOP_USER  - function set by NEPSetStoppingTestFunction()

820:    Level: advanced

822: .seealso: `NEPGetStoppingTest()`, `NEPSetStoppingTestFunction()`, `NEPSetConvergenceTest()`, `NEPStop`
823: @*/
824: PetscErrorCode NEPSetStoppingTest(NEP nep,NEPStop stop)
825: {
826:   PetscFunctionBegin;
829:   switch (stop) {
830:     case NEP_STOP_BASIC: nep->stopping = NEPStoppingBasic; break;
831:     case NEP_STOP_USER:
832:       PetscCheck(nep->stoppinguser,PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"Must call NEPSetStoppingTestFunction() first");
833:       nep->stopping = nep->stoppinguser;
834:       break;
835:     default:
836:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
837:   }
838:   nep->stop = stop;
839:   PetscFunctionReturn(PETSC_SUCCESS);
840: }

842: /*@
843:    NEPGetStoppingTest - Gets the method used to decide the termination of the outer
844:    loop of the eigensolver.

846:    Not Collective

848:    Input Parameters:
849: .  nep   - nonlinear eigensolver context obtained from NEPCreate()

851:    Output Parameters:
852: .  stop  - the type of stopping test

854:    Level: advanced

856: .seealso: `NEPSetStoppingTest()`, `NEPStop`
857: @*/
858: PetscErrorCode NEPGetStoppingTest(NEP nep,NEPStop *stop)
859: {
860:   PetscFunctionBegin;
862:   PetscAssertPointer(stop,2);
863:   *stop = nep->stop;
864:   PetscFunctionReturn(PETSC_SUCCESS);
865: }

867: /*@
868:    NEPSetTrackAll - Specifies if the solver must compute the residual of all
869:    approximate eigenpairs or not.

871:    Logically Collective

873:    Input Parameters:
874: +  nep      - the eigensolver context
875: -  trackall - whether compute all residuals or not

877:    Notes:
878:    If the user sets trackall=PETSC_TRUE then the solver explicitly computes
879:    the residual for each eigenpair approximation. Computing the residual is
880:    usually an expensive operation and solvers commonly compute the associated
881:    residual to the first unconverged eigenpair.

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

885:    Level: developer

887: .seealso: `NEPGetTrackAll()`
888: @*/
889: PetscErrorCode NEPSetTrackAll(NEP nep,PetscBool trackall)
890: {
891:   PetscFunctionBegin;
894:   nep->trackall = trackall;
895:   PetscFunctionReturn(PETSC_SUCCESS);
896: }

898: /*@
899:    NEPGetTrackAll - Returns the flag indicating whether all residual norms must
900:    be computed or not.

902:    Not Collective

904:    Input Parameter:
905: .  nep - the eigensolver context

907:    Output Parameter:
908: .  trackall - the returned flag

910:    Level: developer

912: .seealso: `NEPSetTrackAll()`
913: @*/
914: PetscErrorCode NEPGetTrackAll(NEP nep,PetscBool *trackall)
915: {
916:   PetscFunctionBegin;
918:   PetscAssertPointer(trackall,2);
919:   *trackall = nep->trackall;
920:   PetscFunctionReturn(PETSC_SUCCESS);
921: }

923: /*@
924:    NEPSetRefine - Specifies the refinement type (and options) to be used
925:    after the solve.

927:    Logically Collective

929:    Input Parameters:
930: +  nep    - the nonlinear eigensolver context
931: .  refine - refinement type
932: .  npart  - number of partitions of the communicator
933: .  tol    - the convergence tolerance
934: .  its    - maximum number of refinement iterations
935: -  scheme - which scheme to be used for solving the involved linear systems

937:    Options Database Keys:
938: +  -nep_refine <type> - refinement type, one of <none,simple,multiple>
939: .  -nep_refine_partitions <n> - the number of partitions
940: .  -nep_refine_tol <tol> - the tolerance
941: .  -nep_refine_its <its> - number of iterations
942: -  -nep_refine_scheme - to set the scheme for the linear solves

944:    Notes:
945:    By default, iterative refinement is disabled, since it may be very
946:    costly. There are two possible refinement strategies, simple and multiple.
947:    The simple approach performs iterative refinement on each of the
948:    converged eigenpairs individually, whereas the multiple strategy works
949:    with the invariant pair as a whole, refining all eigenpairs simultaneously.
950:    The latter may be required for the case of multiple eigenvalues.

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

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

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

967:    Use PETSC_CURRENT to retain the current value of npart, tol or its. Use
968:    PETSC_DETERMINE to assign a default value.

970:    Level: intermediate

972: .seealso: `NEPGetRefine()`
973: @*/
974: PetscErrorCode NEPSetRefine(NEP nep,NEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,NEPRefineScheme scheme)
975: {
976:   PetscMPIInt    size;

978:   PetscFunctionBegin;
985:   nep->refine = refine;
986:   if (refine) {  /* process parameters only if not REFINE_NONE */
987:     if (npart!=nep->npart) {
988:       PetscCall(PetscSubcommDestroy(&nep->refinesubc));
989:       PetscCall(KSPDestroy(&nep->refineksp));
990:     }
991:     if (npart == PETSC_DETERMINE) {
992:       nep->npart = 1;
993:     } else if (npart != PETSC_CURRENT) {
994:       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)nep),&size));
995:       PetscCheck(npart>0 && npart<=size,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
996:       nep->npart = npart;
997:     }
998:     if (tol == (PetscReal)PETSC_DETERMINE) {
999:       nep->rtol = PETSC_DETERMINE;
1000:     } else if (tol != (PetscReal)PETSC_CURRENT) {
1001:       PetscCheck(tol>0.0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1002:       nep->rtol = tol;
1003:     }
1004:     if (its==PETSC_DETERMINE) {
1005:       nep->rits = PETSC_DETERMINE;
1006:     } else if (its != PETSC_CURRENT) {
1007:       PetscCheck(its>=0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be >= 0");
1008:       nep->rits = its;
1009:     }
1010:     nep->scheme = scheme;
1011:   }
1012:   nep->state = NEP_STATE_INITIAL;
1013:   PetscFunctionReturn(PETSC_SUCCESS);
1014: }

1016: /*@
1017:    NEPGetRefine - Gets the refinement strategy used by the NEP object, and the
1018:    associated parameters.

1020:    Not Collective

1022:    Input Parameter:
1023: .  nep - the nonlinear eigensolver context

1025:    Output Parameters:
1026: +  refine - refinement type
1027: .  npart  - number of partitions of the communicator
1028: .  tol    - the convergence tolerance
1029: .  its    - maximum number of refinement iterations
1030: -  scheme - the scheme used for solving linear systems

1032:    Level: intermediate

1034:    Note:
1035:    The user can specify NULL for any parameter that is not needed.

1037: .seealso: `NEPSetRefine()`
1038: @*/
1039: PetscErrorCode NEPGetRefine(NEP nep,NEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,NEPRefineScheme *scheme)
1040: {
1041:   PetscFunctionBegin;
1043:   if (refine) *refine = nep->refine;
1044:   if (npart)  *npart  = nep->npart;
1045:   if (tol)    *tol    = nep->rtol;
1046:   if (its)    *its    = nep->rits;
1047:   if (scheme) *scheme = nep->scheme;
1048:   PetscFunctionReturn(PETSC_SUCCESS);
1049: }

1051: /*@
1052:    NEPSetOptionsPrefix - Sets the prefix used for searching for all
1053:    NEP options in the database.

1055:    Logically Collective

1057:    Input Parameters:
1058: +  nep - the nonlinear eigensolver context
1059: -  prefix - the prefix string to prepend to all NEP option requests

1061:    Notes:
1062:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1063:    The first character of all runtime options is AUTOMATICALLY the
1064:    hyphen.

1066:    For example, to distinguish between the runtime options for two
1067:    different NEP contexts, one could call
1068: .vb
1069:       NEPSetOptionsPrefix(nep1,"neig1_")
1070:       NEPSetOptionsPrefix(nep2,"neig2_")
1071: .ve

1073:    Level: advanced

1075: .seealso: `NEPAppendOptionsPrefix()`, `NEPGetOptionsPrefix()`
1076: @*/
1077: PetscErrorCode NEPSetOptionsPrefix(NEP nep,const char *prefix)
1078: {
1079:   PetscFunctionBegin;
1081:   if (!nep->V) PetscCall(NEPGetBV(nep,&nep->V));
1082:   PetscCall(BVSetOptionsPrefix(nep->V,prefix));
1083:   if (!nep->ds) PetscCall(NEPGetDS(nep,&nep->ds));
1084:   PetscCall(DSSetOptionsPrefix(nep->ds,prefix));
1085:   if (!nep->rg) PetscCall(NEPGetRG(nep,&nep->rg));
1086:   PetscCall(RGSetOptionsPrefix(nep->rg,prefix));
1087:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)nep,prefix));
1088:   PetscFunctionReturn(PETSC_SUCCESS);
1089: }

1091: /*@
1092:    NEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1093:    NEP options in the database.

1095:    Logically Collective

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 hyphen.

1105:    Level: advanced

1107: .seealso: `NEPSetOptionsPrefix()`, `NEPGetOptionsPrefix()`
1108: @*/
1109: PetscErrorCode NEPAppendOptionsPrefix(NEP nep,const char *prefix)
1110: {
1111:   PetscFunctionBegin;
1113:   if (!nep->V) PetscCall(NEPGetBV(nep,&nep->V));
1114:   PetscCall(BVAppendOptionsPrefix(nep->V,prefix));
1115:   if (!nep->ds) PetscCall(NEPGetDS(nep,&nep->ds));
1116:   PetscCall(DSAppendOptionsPrefix(nep->ds,prefix));
1117:   if (!nep->rg) PetscCall(NEPGetRG(nep,&nep->rg));
1118:   PetscCall(RGAppendOptionsPrefix(nep->rg,prefix));
1119:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)nep,prefix));
1120:   PetscFunctionReturn(PETSC_SUCCESS);
1121: }

1123: /*@
1124:    NEPGetOptionsPrefix - Gets the prefix used for searching for all
1125:    NEP options in the database.

1127:    Not Collective

1129:    Input Parameters:
1130: .  nep - the nonlinear eigensolver context

1132:    Output Parameters:
1133: .  prefix - pointer to the prefix string used is returned

1135:    Note:
1136:    On the Fortran side, the user should pass in a string 'prefix' of
1137:    sufficient length to hold the prefix.

1139:    Level: advanced

1141: .seealso: `NEPSetOptionsPrefix()`, `NEPAppendOptionsPrefix()`
1142: @*/
1143: PetscErrorCode NEPGetOptionsPrefix(NEP nep,const char *prefix[])
1144: {
1145:   PetscFunctionBegin;
1147:   PetscAssertPointer(prefix,2);
1148:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)nep,prefix));
1149:   PetscFunctionReturn(PETSC_SUCCESS);
1150: }