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: [](ch:nep), `NEPMonitorSet()`, `NEPSetTrackAll()`
 34: @*/
 35: PetscErrorCode NEPMonitorSetFromOptions(NEP nep,const char opt[],const char name[],PetscCtx 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 configure the solver.

 72:    Collective

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

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

 80:    Level: beginner

 82: .seealso: [](ch:nep), `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(PetscOptionsBoolGroup("-nep_all","Compute all eigenvalues in a region","NEPSetWhichEigenpairs",&flg));
166:     if (flg) PetscCall(NEPSetWhichEigenpairs(nep,NEP_ALL));
167:     PetscCall(PetscOptionsBoolGroupEnd("-nep_which_user","Select the user-defined selection criterion","NEPSetWhichEigenpairs",&flg));
168:     if (flg) PetscCall(NEPSetWhichEigenpairs(nep,NEP_WHICH_USER));

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

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

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

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

197:     PetscTryTypeMethod(nep,setfromoptions,PetscOptionsObject);
198:     PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)nep,PetscOptionsObject));
199:   PetscOptionsEnd();

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

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

220:    Not Collective

222:    Input Parameter:
223: .  nep - the nonlinear eigensolver context

225:    Output Parameters:
226: +  tol - the convergence tolerance
227: -  maxits - maximum number of iterations

229:    Notes:
230:    The user can specify `NULL` for any parameter that is not needed.

232:    Level: intermediate

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

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

249:    Logically Collective

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

256:    Options Database Keys:
257: +  -nep_tol tol       - sets the convergence tolerance
258: -  -nep_max_it maxits - sets the maximum number of iterations allowed

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

266:    Level: intermediate

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

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

299:    Not Collective

301:    Input Parameter:
302: .  nep - the nonlinear eigensolver context

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

309:    Note:
310:    The user can specify `NULL` for any parameter that is not needed.

312:    Level: intermediate

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

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

330:    Logically Collective

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

338:    Options Database Keys:
339: +  -nep_nev nev - sets the number of eigenvalues
340: .  -nep_ncv ncv - sets the dimension of the subspace
341: -  -nep_mpd mpd - sets the maximum projected dimension

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

348:    The parameters `ncv` and `mpd` are intimately related, so that the user is advised
349:    to set one of them at most. Normal usage is\:

351:     1. in cases where `nev` is small, the user sets `ncv` (a reasonable default is `2*nev`).
352:     2. 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: [](ch:nep), `NEPGetDimensions()`
361: @*/
362: PetscErrorCode NEPSetDimensions(NEP nep,PetscInt nev,PetscInt ncv,PetscInt mpd)
363: {
364:   PetscFunctionBegin;
369:   if (nev != PETSC_CURRENT) {
370:     PetscCheck(nev>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
371:     nep->nev = nev;
372:   }
373:   if (ncv == PETSC_DETERMINE) {
374:     nep->ncv = PETSC_DETERMINE;
375:   } else if (ncv != PETSC_CURRENT) {
376:     PetscCheck(ncv>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
377:     nep->ncv = ncv;
378:   }
379:   if (mpd == PETSC_DETERMINE) {
380:     nep->mpd = PETSC_DETERMINE;
381:   } else if (mpd != PETSC_CURRENT) {
382:     PetscCheck(mpd>0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
383:     nep->mpd = mpd;
384:   }
385:   nep->state = NEP_STATE_INITIAL;
386:   PetscFunctionReturn(PETSC_SUCCESS);
387: }

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

393:    Logically Collective

395:    Input Parameters:
396: +  nep   - the nonlinear eigensolver context
397: -  which - the portion of the spectrum to be sought, see `NEPWhich` for possible values

399:    Options Database Keys:
400: +  -nep_largest_magnitude  - sets largest eigenvalues in magnitude
401: .  -nep_smallest_magnitude - sets smallest eigenvalues in magnitude
402: .  -nep_largest_real       - sets largest real parts
403: .  -nep_smallest_real      - sets smallest real parts
404: .  -nep_largest_imaginary  - sets largest imaginary parts
405: .  -nep_smallest_imaginary - sets smallest imaginary parts
406: .  -nep_target_magnitude   - sets eigenvalues closest to target
407: .  -nep_target_real        - sets real parts closest to target
408: .  -nep_target_imaginary   - sets imaginary parts closest to target
409: .  -nep_all                - sets all eigenvalues in a region
410: -  -nep_which_user         - select the user-defined selection criterion

412:    Notes:
413:    Not all eigensolvers implemented in `NEP` account for all the possible values
414:    of `which`. Also, some values make sense only for certain types of
415:    problems. If SLEPc is compiled for real numbers `NEP_LARGEST_IMAGINARY`
416:    and `NEP_SMALLEST_IMAGINARY` use the absolute value of the imaginary part
417:    for eigenvalue selection.

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

421:    The criterion `NEP_TARGET_IMAGINARY` is available only in case PETSc and
422:    SLEPc have been built with complex scalars.

424:    `NEP_ALL` is intended for use in the context of the `NEPCISS` solver for
425:    computing all eigenvalues in a region.

427:    Level: intermediate

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

465: /*@
466:     NEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
467:     sought.

469:     Not Collective

471:     Input Parameter:
472: .   nep - the nonlinear eigensolver context

474:     Output Parameter:
475: .   which - the portion of the spectrum to be sought

477:     Level: intermediate

479: .seealso: [](ch:nep), `NEPSetWhichEigenpairs()`, `NEPWhich`
480: @*/
481: PetscErrorCode NEPGetWhichEigenpairs(NEP nep,NEPWhich *which)
482: {
483:   PetscFunctionBegin;
485:   PetscAssertPointer(which,2);
486:   *which = nep->which;
487:   PetscFunctionReturn(PETSC_SUCCESS);
488: }

490: /*@C
491:    NEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
492:    when `NEPSetWhichEigenpairs()` is set to `NEP_WHICH_USER`.

494:    Logically Collective

496:    Input Parameters:
497: +  nep  - the nonlinear eigensolver context
498: .  comp - a pointer to the comparison function, see `SlepcEigenvalueComparisonFn` for the calling sequence
499: -  ctx  - a context pointer (the last parameter to the comparison function)

501:    Level: advanced

503: .seealso: [](ch:nep), `NEPSetWhichEigenpairs()`, `NEPWhich`
504: @*/
505: PetscErrorCode NEPSetEigenvalueComparison(NEP nep,SlepcEigenvalueComparisonFn *comp,PetscCtx ctx)
506: {
507:   PetscFunctionBegin;
509:   nep->sc->comparison    = comp;
510:   nep->sc->comparisonctx = ctx;
511:   nep->which             = NEP_WHICH_USER;
512:   PetscFunctionReturn(PETSC_SUCCESS);
513: }

515: /*@
516:    NEPSetProblemType - Specifies the type of the nonlinear eigenvalue problem.

518:    Logically Collective

520:    Input Parameters:
521: +  nep  - the nonlinear eigensolver context
522: -  type - a known type of nonlinear eigenvalue problem

524:    Options Database Keys:
525: +  -nep_general  - general problem with no particular structure
526: -  -nep_rational - a rational eigenvalue problem defined in split form with all $f_i$ rational

528:    Notes:
529:    See `NEPProblemType` for possible problem types.

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

535:    Level: intermediate

537: .seealso: [](ch:nep), `NEPSetType()`, `NEPGetProblemType()`, `NEPProblemType`
538: @*/
539: PetscErrorCode NEPSetProblemType(NEP nep,NEPProblemType type)
540: {
541:   PetscFunctionBegin;
544:   PetscCheck(type==NEP_GENERAL || type==NEP_RATIONAL,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
545:   if (type != nep->problem_type) {
546:     nep->problem_type = type;
547:     nep->state = NEP_STATE_INITIAL;
548:   }
549:   PetscFunctionReturn(PETSC_SUCCESS);
550: }

552: /*@
553:    NEPGetProblemType - Gets the problem type from the `NEP` object.

555:    Not Collective

557:    Input Parameter:
558: .  nep - the nonlinear eigensolver context

560:    Output Parameter:
561: .  type - the problem type

563:    Level: intermediate

565: .seealso: [](ch:nep), `NEPSetProblemType()`, `NEPProblemType`
566: @*/
567: PetscErrorCode NEPGetProblemType(NEP nep,NEPProblemType *type)
568: {
569:   PetscFunctionBegin;
571:   PetscAssertPointer(type,2);
572:   *type = nep->problem_type;
573:   PetscFunctionReturn(PETSC_SUCCESS);
574: }

576: /*@
577:    NEPSetTwoSided - Sets the solver to use a two-sided variant so that left
578:    eigenvectors are also computed.

580:    Logically Collective

582:    Input Parameters:
583: +  nep      - the nonlinear eigensolver context
584: -  twosided - whether the two-sided variant is to be used or not

586:    Options Database Key:
587: .  -nep_two_sided (true|false) - toggles the twosided flag

589:    Notes:
590:    If the user sets `twosided`=`PETSC_TRUE` then the solver uses a variant of
591:    the algorithm that computes both right and left eigenvectors. This is
592:    usually much more costly. This option is not available in all solvers,
593:    see table [](#tab:solversn).

595:    When using two-sided solvers, the problem matrices must have both the
596:    `MATOP_MULT` and `MATOP_MULT_TRANSPOSE` operations defined.

598:    Level: advanced

600: .seealso: [](ch:nep), `NEPGetTwoSided()`, `NEPGetLeftEigenvector()`
601: @*/
602: PetscErrorCode NEPSetTwoSided(NEP nep,PetscBool twosided)
603: {
604:   PetscFunctionBegin;
607:   if (twosided!=nep->twosided) {
608:     nep->twosided = twosided;
609:     nep->state    = NEP_STATE_INITIAL;
610:   }
611:   PetscFunctionReturn(PETSC_SUCCESS);
612: }

614: /*@
615:    NEPGetTwoSided - Returns the flag indicating whether a two-sided variant
616:    of the algorithm is being used or not.

618:    Not Collective

620:    Input Parameter:
621: .  nep - the nonlinear eigensolver context

623:    Output Parameter:
624: .  twosided - the returned flag

626:    Level: advanced

628: .seealso: [](ch:nep), `NEPSetTwoSided()`
629: @*/
630: PetscErrorCode NEPGetTwoSided(NEP nep,PetscBool *twosided)
631: {
632:   PetscFunctionBegin;
634:   PetscAssertPointer(twosided,2);
635:   *twosided = nep->twosided;
636:   PetscFunctionReturn(PETSC_SUCCESS);
637: }

639: /*@C
640:    NEPSetConvergenceTestFunction - Sets a function to compute the error estimate
641:    used in the convergence test.

643:    Logically Collective

645:    Input Parameters:
646: +  nep     - the nonlinear eigensolver context
647: .  conv    - convergence test function, see `NEPConvergenceTestFn` for the calling sequence
648: .  ctx     - context for private data for the convergence routine (may be `NULL`)
649: -  destroy - a routine for destroying the context (may be `NULL`), see `PetscCtxDestroyFn`
650:              for the calling sequence

652:    Notes:
653:    When this is called with a user-defined function, then the convergence
654:    criterion is set to `NEP_CONV_USER`, see `NEPSetConvergenceTest()`.

656:    If the error estimate returned by the convergence test function is less than
657:    the tolerance, then the eigenvalue is accepted as converged.

659:    Level: advanced

661: .seealso: [](ch:nep), `NEPSetConvergenceTest()`, `NEPSetTolerances()`
662: @*/
663: PetscErrorCode NEPSetConvergenceTestFunction(NEP nep,NEPConvergenceTestFn *conv,PetscCtx ctx,PetscCtxDestroyFn *destroy)
664: {
665:   PetscFunctionBegin;
667:   if (nep->convergeddestroy) PetscCall((*nep->convergeddestroy)(&nep->convergedctx));
668:   nep->convergeduser    = conv;
669:   nep->convergeddestroy = destroy;
670:   nep->convergedctx     = ctx;
671:   if (conv == NEPConvergedRelative) nep->conv = NEP_CONV_REL;
672:   else if (conv == NEPConvergedNorm) nep->conv = NEP_CONV_NORM;
673:   else if (conv == NEPConvergedAbsolute) nep->conv = NEP_CONV_ABS;
674:   else {
675:     nep->conv      = NEP_CONV_USER;
676:     nep->converged = nep->convergeduser;
677:   }
678:   PetscFunctionReturn(PETSC_SUCCESS);
679: }

681: /*@
682:    NEPSetConvergenceTest - Specifies how to compute the error estimate
683:    used in the convergence test.

685:    Logically Collective

687:    Input Parameters:
688: +  nep  - the nonlinear eigensolver context
689: -  conv - the type of convergence test, see `NEPConv` for possible values

691:    Options Database Keys:
692: +  -nep_conv_abs  - sets the absolute convergence test
693: .  -nep_conv_rel  - sets the convergence test relative to the eigenvalue
694: .  -nep_conv_norm - sets the convergence test relative to the matrix norms
695: -  -nep_conv_user - selects the user-defined convergence test

697:    Level: intermediate

699: .seealso: [](ch:nep), `NEPGetConvergenceTest()`, `NEPSetConvergenceTestFunction()`, `NEPSetStoppingTest()`, `NEPConv`
700: @*/
701: PetscErrorCode NEPSetConvergenceTest(NEP nep,NEPConv conv)
702: {
703:   PetscFunctionBegin;
706:   switch (conv) {
707:     case NEP_CONV_ABS:  nep->converged = NEPConvergedAbsolute; break;
708:     case NEP_CONV_REL:  nep->converged = NEPConvergedRelative; break;
709:     case NEP_CONV_NORM: nep->converged = NEPConvergedNorm; break;
710:     case NEP_CONV_USER:
711:       PetscCheck(nep->convergeduser,PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"Must call NEPSetConvergenceTestFunction() first");
712:       nep->converged = nep->convergeduser;
713:       break;
714:     default:
715:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
716:   }
717:   nep->conv = conv;
718:   PetscFunctionReturn(PETSC_SUCCESS);
719: }

721: /*@
722:    NEPGetConvergenceTest - Gets the method used to compute the error estimate
723:    used in the convergence test.

725:    Not Collective

727:    Input Parameter:
728: .  nep   - the nonlinear eigensolver context

730:    Output Parameter:
731: .  conv  - the type of convergence test

733:    Level: intermediate

735: .seealso: [](ch:nep), `NEPSetConvergenceTest()`, `NEPConv`
736: @*/
737: PetscErrorCode NEPGetConvergenceTest(NEP nep,NEPConv *conv)
738: {
739:   PetscFunctionBegin;
741:   PetscAssertPointer(conv,2);
742:   *conv = nep->conv;
743:   PetscFunctionReturn(PETSC_SUCCESS);
744: }

746: /*@C
747:    NEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
748:    iteration of the eigensolver.

750:    Logically Collective

752:    Input Parameters:
753: +  nep     - the nonlinear eigensolver context
754: .  stop    - the stopping test function, see `NEPStoppingTestFn` for the calling sequence
755: .  ctx     - context for private data for the stopping routine (may be `NULL`)
756: -  destroy - a routine for destroying the context (may be `NULL`), see `PetscCtxDestroyFn`
757:              for the calling sequence

759:    Note:
760:    When implementing a function for this, normal usage is to first call the
761:    default routine `NEPStoppingBasic()` and then set `reason` to `NEP_CONVERGED_USER`
762:    if some user-defined conditions have been met. To let the eigensolver continue
763:    iterating, the result must be left as `NEP_CONVERGED_ITERATING`.

765:    Level: advanced

767: .seealso: [](ch:nep), `NEPSetStoppingTest()`, `NEPStoppingBasic()`
768: @*/
769: PetscErrorCode NEPSetStoppingTestFunction(NEP nep,NEPStoppingTestFn *stop,PetscCtx ctx,PetscCtxDestroyFn *destroy)
770: {
771:   PetscFunctionBegin;
773:   if (nep->stoppingdestroy) PetscCall((*nep->stoppingdestroy)(&nep->stoppingctx));
774:   nep->stoppinguser    = stop;
775:   nep->stoppingdestroy = destroy;
776:   nep->stoppingctx     = ctx;
777:   if (stop == NEPStoppingBasic) nep->stop = NEP_STOP_BASIC;
778:   else {
779:     nep->stop     = NEP_STOP_USER;
780:     nep->stopping = nep->stoppinguser;
781:   }
782:   PetscFunctionReturn(PETSC_SUCCESS);
783: }

785: /*@
786:    NEPSetStoppingTest - Specifies how to decide the termination of the outer
787:    loop of the eigensolver.

789:    Logically Collective

791:    Input Parameters:
792: +  nep  - the nonlinear eigensolver context
793: -  stop - the type of stopping test, see `NEPStop`

795:    Options Database Keys:
796: +  -nep_stop_basic - sets the default stopping test
797: -  -nep_stop_user  - selects the user-defined stopping test

799:    Level: advanced

801: .seealso: [](ch:nep), `NEPGetStoppingTest()`, `NEPSetStoppingTestFunction()`, `NEPSetConvergenceTest()`, `NEPStop`
802: @*/
803: PetscErrorCode NEPSetStoppingTest(NEP nep,NEPStop stop)
804: {
805:   PetscFunctionBegin;
808:   switch (stop) {
809:     case NEP_STOP_BASIC: nep->stopping = NEPStoppingBasic; break;
810:     case NEP_STOP_USER:
811:       PetscCheck(nep->stoppinguser,PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"Must call NEPSetStoppingTestFunction() first");
812:       nep->stopping = nep->stoppinguser;
813:       break;
814:     default:
815:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
816:   }
817:   nep->stop = stop;
818:   PetscFunctionReturn(PETSC_SUCCESS);
819: }

821: /*@
822:    NEPGetStoppingTest - Gets the method used to decide the termination of the outer
823:    loop of the eigensolver.

825:    Not Collective

827:    Input Parameter:
828: .  nep   - the nonlinear eigensolver context

830:    Output Parameter:
831: .  stop  - the type of stopping test

833:    Level: advanced

835: .seealso: [](ch:nep), `NEPSetStoppingTest()`, `NEPStop`
836: @*/
837: PetscErrorCode NEPGetStoppingTest(NEP nep,NEPStop *stop)
838: {
839:   PetscFunctionBegin;
841:   PetscAssertPointer(stop,2);
842:   *stop = nep->stop;
843:   PetscFunctionReturn(PETSC_SUCCESS);
844: }

846: /*@
847:    NEPSetTrackAll - Specifies if the solver must compute the residual of all
848:    approximate eigenpairs or not.

850:    Logically Collective

852:    Input Parameters:
853: +  nep      - the nonlinear eigensolver context
854: -  trackall - whether compute all residuals or not

856:    Notes:
857:    If the user sets `trackall`=`PETSC_TRUE` then the solver explicitly computes
858:    the residual for each eigenpair approximation. Computing the residual is
859:    usually an expensive operation and solvers commonly compute the associated
860:    residual to the first unconverged eigenpair.

862:    The option `-nep_monitor_all` automatically activates this option.

864:    Level: developer

866: .seealso: [](ch:nep), `NEPGetTrackAll()`
867: @*/
868: PetscErrorCode NEPSetTrackAll(NEP nep,PetscBool trackall)
869: {
870:   PetscFunctionBegin;
873:   nep->trackall = trackall;
874:   PetscFunctionReturn(PETSC_SUCCESS);
875: }

877: /*@
878:    NEPGetTrackAll - Returns the flag indicating whether all residual norms must
879:    be computed or not.

881:    Not Collective

883:    Input Parameter:
884: .  nep - the nonlinear eigensolver context

886:    Output Parameter:
887: .  trackall - the returned flag

889:    Level: developer

891: .seealso: [](ch:nep), `NEPSetTrackAll()`
892: @*/
893: PetscErrorCode NEPGetTrackAll(NEP nep,PetscBool *trackall)
894: {
895:   PetscFunctionBegin;
897:   PetscAssertPointer(trackall,2);
898:   *trackall = nep->trackall;
899:   PetscFunctionReturn(PETSC_SUCCESS);
900: }

902: /*@
903:    NEPSetRefine - Specifies the refinement type (and options) to be used
904:    after the solve.

906:    Logically Collective

908:    Input Parameters:
909: +  nep    - the nonlinear eigensolver context
910: .  refine - refinement type, see `NEPRefine` for possible values
911: .  npart  - number of partitions of the communicator
912: .  tol    - the convergence tolerance
913: .  its    - maximum number of refinement iterations
914: -  scheme - which scheme to be used for solving the involved linear systems, see `NEPRefineScheme`
915:             for possible values

917:    Options Database Keys:
918: +  -nep_refine (none|simple|multiple)      - set the refinement type
919: .  -nep_refine_partitions npart            - set the number of partitions
920: .  -nep_refine_tol tol                     - set the tolerance
921: .  -nep_refine_its its                     - set the number of iterations
922: -  -nep_refine_scheme (schur|mbe|explicit) - set the scheme for the linear solves

924:    Notes:
925:    This function configures the parameters of Newton iterative refinement,
926:    see section [](#sec:refine) for a discussion of the different strategies
927:    in the context of polynomial eigenproblems.

929:    By default, iterative refinement is disabled, since it may be very
930:    costly. There are two possible refinement strategies, simple and multiple.
931:    The simple approach performs iterative refinement on each of the
932:    converged eigenpairs individually, whereas the multiple strategy works
933:    with the invariant pair as a whole, refining all eigenpairs simultaneously.
934:    The latter may be required for the case of multiple eigenvalues.

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

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

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

951:    Use `PETSC_CURRENT` to retain the current value of `npart`, `tol` or `its`. Use
952:    `PETSC_DETERMINE` to assign a default value.

954:    Level: intermediate

956: .seealso: [](ch:nep), [](#sec:refine), `NEPGetRefine()`
957: @*/
958: PetscErrorCode NEPSetRefine(NEP nep,NEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,NEPRefineScheme scheme)
959: {
960:   PetscMPIInt    size;

962:   PetscFunctionBegin;
969:   nep->refine = refine;
970:   if (refine) {  /* process parameters only if not REFINE_NONE */
971:     if (npart!=nep->npart) {
972:       PetscCall(PetscSubcommDestroy(&nep->refinesubc));
973:       PetscCall(KSPDestroy(&nep->refineksp));
974:     }
975:     if (npart == PETSC_DETERMINE) {
976:       nep->npart = 1;
977:     } else if (npart != PETSC_CURRENT) {
978:       PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)nep),&size));
979:       PetscCheck(npart>0 && npart<=size,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
980:       nep->npart = npart;
981:     }
982:     if (tol == (PetscReal)PETSC_DETERMINE) {
983:       nep->rtol = PETSC_DETERMINE;
984:     } else if (tol != (PetscReal)PETSC_CURRENT) {
985:       PetscCheck(tol>0.0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
986:       nep->rtol = tol;
987:     }
988:     if (its==PETSC_DETERMINE) {
989:       nep->rits = PETSC_DETERMINE;
990:     } else if (its != PETSC_CURRENT) {
991:       PetscCheck(its>=0,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be >= 0");
992:       nep->rits = its;
993:     }
994:     nep->scheme = scheme;
995:   }
996:   nep->state = NEP_STATE_INITIAL;
997:   PetscFunctionReturn(PETSC_SUCCESS);
998: }

1000: /*@
1001:    NEPGetRefine - Gets the refinement strategy used by the `NEP` object, and the
1002:    associated parameters.

1004:    Not Collective

1006:    Input Parameter:
1007: .  nep - the nonlinear eigensolver context

1009:    Output Parameters:
1010: +  refine - refinement type
1011: .  npart  - number of partitions of the communicator
1012: .  tol    - the convergence tolerance
1013: .  its    - maximum number of refinement iterations
1014: -  scheme - the scheme used for solving linear systems

1016:    Level: intermediate

1018:    Note:
1019:    The user can specify `NULL` for any parameter that is not needed.

1021: .seealso: [](ch:nep), `NEPSetRefine()`
1022: @*/
1023: PetscErrorCode NEPGetRefine(NEP nep,NEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,NEPRefineScheme *scheme)
1024: {
1025:   PetscFunctionBegin;
1027:   if (refine) *refine = nep->refine;
1028:   if (npart)  *npart  = nep->npart;
1029:   if (tol)    *tol    = nep->rtol;
1030:   if (its)    *its    = nep->rits;
1031:   if (scheme) *scheme = nep->scheme;
1032:   PetscFunctionReturn(PETSC_SUCCESS);
1033: }

1035: /*@
1036:    NEPSetOptionsPrefix - Sets the prefix used for searching for all
1037:    `NEP` options in the database.

1039:    Logically Collective

1041:    Input Parameters:
1042: +  nep    - the nonlinear eigensolver context
1043: -  prefix - the prefix string to prepend to all `NEP` option requests

1045:    Notes:
1046:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1047:    The first character of all runtime options is AUTOMATICALLY the
1048:    hyphen.

1050:    For example, to distinguish between the runtime options for two
1051:    different `NEP` contexts, one could call
1052: .vb
1053:    NEPSetOptionsPrefix(nep1,"neig1_")
1054:    NEPSetOptionsPrefix(nep2,"neig2_")
1055: .ve

1057:    Level: advanced

1059: .seealso: [](ch:nep), `NEPAppendOptionsPrefix()`, `NEPGetOptionsPrefix()`
1060: @*/
1061: PetscErrorCode NEPSetOptionsPrefix(NEP nep,const char prefix[])
1062: {
1063:   PetscFunctionBegin;
1065:   if (!nep->V) PetscCall(NEPGetBV(nep,&nep->V));
1066:   PetscCall(BVSetOptionsPrefix(nep->V,prefix));
1067:   if (!nep->ds) PetscCall(NEPGetDS(nep,&nep->ds));
1068:   PetscCall(DSSetOptionsPrefix(nep->ds,prefix));
1069:   if (!nep->rg) PetscCall(NEPGetRG(nep,&nep->rg));
1070:   PetscCall(RGSetOptionsPrefix(nep->rg,prefix));
1071:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)nep,prefix));
1072:   PetscFunctionReturn(PETSC_SUCCESS);
1073: }

1075: /*@
1076:    NEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1077:    `NEP` options in the database.

1079:    Logically Collective

1081:    Input Parameters:
1082: +  nep    - the nonlinear eigensolver context
1083: -  prefix - the prefix string to prepend to all `NEP` option requests

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

1089:    Level: advanced

1091: .seealso: [](ch:nep), `NEPSetOptionsPrefix()`, `NEPGetOptionsPrefix()`
1092: @*/
1093: PetscErrorCode NEPAppendOptionsPrefix(NEP nep,const char prefix[])
1094: {
1095:   PetscFunctionBegin;
1097:   if (!nep->V) PetscCall(NEPGetBV(nep,&nep->V));
1098:   PetscCall(BVAppendOptionsPrefix(nep->V,prefix));
1099:   if (!nep->ds) PetscCall(NEPGetDS(nep,&nep->ds));
1100:   PetscCall(DSAppendOptionsPrefix(nep->ds,prefix));
1101:   if (!nep->rg) PetscCall(NEPGetRG(nep,&nep->rg));
1102:   PetscCall(RGAppendOptionsPrefix(nep->rg,prefix));
1103:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)nep,prefix));
1104:   PetscFunctionReturn(PETSC_SUCCESS);
1105: }

1107: /*@
1108:    NEPGetOptionsPrefix - Gets the prefix used for searching for all
1109:    `NEP` options in the database.

1111:    Not Collective

1113:    Input Parameter:
1114: .  nep - the nonlinear eigensolver context

1116:    Output Parameter:
1117: .  prefix - pointer to the prefix string used is returned

1119:    Level: advanced

1121: .seealso: [](ch:nep), `NEPSetOptionsPrefix()`, `NEPAppendOptionsPrefix()`
1122: @*/
1123: PetscErrorCode NEPGetOptionsPrefix(NEP nep,const char *prefix[])
1124: {
1125:   PetscFunctionBegin;
1127:   PetscAssertPointer(prefix,2);
1128:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)nep,prefix));
1129:   PetscFunctionReturn(PETSC_SUCCESS);
1130: }