Actual source code: nepopts.c

slepc-main 2024-12-17
Report Typos and Errors
  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:    Calling sequence of comp:
511: $  PetscErrorCode comp(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)
512: +   ar     - real part of the 1st eigenvalue
513: .   ai     - imaginary part of the 1st eigenvalue
514: .   br     - real part of the 2nd eigenvalue
515: .   bi     - imaginary part of the 2nd eigenvalue
516: .   res    - result of comparison
517: -   ctx    - optional context, as set by NEPSetEigenvalueComparison()

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

525:    Level: advanced

527: .seealso: NEPSetWhichEigenpairs(), NEPWhich
528: @*/
529: PetscErrorCode NEPSetEigenvalueComparison(NEP nep,PetscErrorCode (*comp)(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx),void* ctx)
530: {
531:   PetscFunctionBegin;
533:   nep->sc->comparison    = comp;
534:   nep->sc->comparisonctx = ctx;
535:   nep->which             = NEP_WHICH_USER;
536:   PetscFunctionReturn(PETSC_SUCCESS);
537: }

539: /*@
540:    NEPSetProblemType - Specifies the type of the nonlinear eigenvalue problem.

542:    Logically Collective

544:    Input Parameters:
545: +  nep  - the nonlinear eigensolver context
546: -  type - a known type of nonlinear eigenvalue problem

548:    Options Database Keys:
549: +  -nep_general - general problem with no particular structure
550: -  -nep_rational - a rational eigenvalue problem defined in split form with all f_i rational

552:    Notes:
553:    Allowed values for the problem type are general (NEP_GENERAL), and rational
554:    (NEP_RATIONAL).

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

560:    Level: intermediate

562: .seealso: NEPSetType(), NEPGetProblemType(), NEPProblemType
563: @*/
564: PetscErrorCode NEPSetProblemType(NEP nep,NEPProblemType type)
565: {
566:   PetscFunctionBegin;
569:   PetscCheck(type==NEP_GENERAL || type==NEP_RATIONAL,PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
570:   if (type != nep->problem_type) {
571:     nep->problem_type = type;
572:     nep->state = NEP_STATE_INITIAL;
573:   }
574:   PetscFunctionReturn(PETSC_SUCCESS);
575: }

577: /*@
578:    NEPGetProblemType - Gets the problem type from the NEP object.

580:    Not Collective

582:    Input Parameter:
583: .  nep - the nonlinear eigensolver context

585:    Output Parameter:
586: .  type - the problem type

588:    Level: intermediate

590: .seealso: NEPSetProblemType(), NEPProblemType
591: @*/
592: PetscErrorCode NEPGetProblemType(NEP nep,NEPProblemType *type)
593: {
594:   PetscFunctionBegin;
596:   PetscAssertPointer(type,2);
597:   *type = nep->problem_type;
598:   PetscFunctionReturn(PETSC_SUCCESS);
599: }

601: /*@
602:    NEPSetTwoSided - Sets the solver to use a two-sided variant so that left
603:    eigenvectors are also computed.

605:    Logically Collective

607:    Input Parameters:
608: +  nep      - the eigensolver context
609: -  twosided - whether the two-sided variant is to be used or not

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

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

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

622:    Level: advanced

624: .seealso: NEPGetTwoSided(), NEPGetLeftEigenvector()
625: @*/
626: PetscErrorCode NEPSetTwoSided(NEP nep,PetscBool twosided)
627: {
628:   PetscFunctionBegin;
631:   if (twosided!=nep->twosided) {
632:     nep->twosided = twosided;
633:     nep->state    = NEP_STATE_INITIAL;
634:   }
635:   PetscFunctionReturn(PETSC_SUCCESS);
636: }

638: /*@
639:    NEPGetTwoSided - Returns the flag indicating whether a two-sided variant
640:    of the algorithm is being used or not.

642:    Not Collective

644:    Input Parameter:
645: .  nep - the eigensolver context

647:    Output Parameter:
648: .  twosided - the returned flag

650:    Level: advanced

652: .seealso: NEPSetTwoSided()
653: @*/
654: PetscErrorCode NEPGetTwoSided(NEP nep,PetscBool *twosided)
655: {
656:   PetscFunctionBegin;
658:   PetscAssertPointer(twosided,2);
659:   *twosided = nep->twosided;
660:   PetscFunctionReturn(PETSC_SUCCESS);
661: }

663: /*@C
664:    NEPSetConvergenceTestFunction - Sets a function to compute the error estimate
665:    used in the convergence test.

667:    Logically Collective

669:    Input Parameters:
670: +  nep     - nonlinear eigensolver context obtained from NEPCreate()
671: .  conv    - convergence test function, see NEPConvergenceTestFn for the calling sequence
672: .  ctx     - context for private data for the convergence routine (may be NULL)
673: -  destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence

675:    Note:
676:    If the error estimate returned by the convergence test function is less than
677:    the tolerance, then the eigenvalue is accepted as converged.

679:    Level: advanced

681: .seealso: NEPSetConvergenceTest(), NEPSetTolerances()
682: @*/
683: PetscErrorCode NEPSetConvergenceTestFunction(NEP nep,NEPConvergenceTestFn *conv,void* ctx,PetscCtxDestroyFn *destroy)
684: {
685:   PetscFunctionBegin;
687:   if (nep->convergeddestroy) PetscCall((*nep->convergeddestroy)(&nep->convergedctx));
688:   nep->convergeduser    = conv;
689:   nep->convergeddestroy = destroy;
690:   nep->convergedctx     = ctx;
691:   if (conv == NEPConvergedRelative) nep->conv = NEP_CONV_REL;
692:   else if (conv == NEPConvergedNorm) nep->conv = NEP_CONV_NORM;
693:   else if (conv == NEPConvergedAbsolute) nep->conv = NEP_CONV_ABS;
694:   else {
695:     nep->conv      = NEP_CONV_USER;
696:     nep->converged = nep->convergeduser;
697:   }
698:   PetscFunctionReturn(PETSC_SUCCESS);
699: }

701: /*@
702:    NEPSetConvergenceTest - Specifies how to compute the error estimate
703:    used in the convergence test.

705:    Logically Collective

707:    Input Parameters:
708: +  nep  - nonlinear eigensolver context obtained from NEPCreate()
709: -  conv - the type of convergence test

711:    Options Database Keys:
712: +  -nep_conv_abs  - Sets the absolute convergence test
713: .  -nep_conv_rel  - Sets the convergence test relative to the eigenvalue
714: -  -nep_conv_user - Selects the user-defined convergence test

716:    Note:
717:    The parameter 'conv' can have one of these values
718: +     NEP_CONV_ABS  - absolute error ||r||
719: .     NEP_CONV_REL  - error relative to the eigenvalue l, ||r||/|l|
720: .     NEP_CONV_NORM - error relative matrix norms, ||r||/sum_i(|f_i(l)|*||A_i||)
721: -     NEP_CONV_USER - function set by NEPSetConvergenceTestFunction()

723:    Level: intermediate

725: .seealso: NEPGetConvergenceTest(), NEPSetConvergenceTestFunction(), NEPSetStoppingTest(), NEPConv
726: @*/
727: PetscErrorCode NEPSetConvergenceTest(NEP nep,NEPConv conv)
728: {
729:   PetscFunctionBegin;
732:   switch (conv) {
733:     case NEP_CONV_ABS:  nep->converged = NEPConvergedAbsolute; break;
734:     case NEP_CONV_REL:  nep->converged = NEPConvergedRelative; break;
735:     case NEP_CONV_NORM: nep->converged = NEPConvergedNorm; break;
736:     case NEP_CONV_USER:
737:       PetscCheck(nep->convergeduser,PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"Must call NEPSetConvergenceTestFunction() first");
738:       nep->converged = nep->convergeduser;
739:       break;
740:     default:
741:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
742:   }
743:   nep->conv = conv;
744:   PetscFunctionReturn(PETSC_SUCCESS);
745: }

747: /*@
748:    NEPGetConvergenceTest - Gets the method used to compute the error estimate
749:    used in the convergence test.

751:    Not Collective

753:    Input Parameters:
754: .  nep   - nonlinear eigensolver context obtained from NEPCreate()

756:    Output Parameters:
757: .  conv  - the type of convergence test

759:    Level: intermediate

761: .seealso: NEPSetConvergenceTest(), NEPConv
762: @*/
763: PetscErrorCode NEPGetConvergenceTest(NEP nep,NEPConv *conv)
764: {
765:   PetscFunctionBegin;
767:   PetscAssertPointer(conv,2);
768:   *conv = nep->conv;
769:   PetscFunctionReturn(PETSC_SUCCESS);
770: }

772: /*@C
773:    NEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
774:    iteration of the eigensolver.

776:    Logically Collective

778:    Input Parameters:
779: +  nep     - nonlinear eigensolver context obtained from NEPCreate()
780: .  stop    - the stopping test function, see NEPStoppingTestFn for the calling sequence
781: .  ctx     - context for private data for the stopping routine (may be NULL)
782: -  destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence

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

790:    Level: advanced

792: .seealso: NEPSetStoppingTest(), NEPStoppingBasic()
793: @*/
794: PetscErrorCode NEPSetStoppingTestFunction(NEP nep,NEPStoppingTestFn *stop,void* ctx,PetscCtxDestroyFn *destroy)
795: {
796:   PetscFunctionBegin;
798:   if (nep->stoppingdestroy) PetscCall((*nep->stoppingdestroy)(&nep->stoppingctx));
799:   nep->stoppinguser    = stop;
800:   nep->stoppingdestroy = destroy;
801:   nep->stoppingctx     = ctx;
802:   if (stop == NEPStoppingBasic) nep->stop = NEP_STOP_BASIC;
803:   else {
804:     nep->stop     = NEP_STOP_USER;
805:     nep->stopping = nep->stoppinguser;
806:   }
807:   PetscFunctionReturn(PETSC_SUCCESS);
808: }

810: /*@
811:    NEPSetStoppingTest - Specifies how to decide the termination of the outer
812:    loop of the eigensolver.

814:    Logically Collective

816:    Input Parameters:
817: +  nep  - nonlinear eigensolver context obtained from NEPCreate()
818: -  stop - the type of stopping test

820:    Options Database Keys:
821: +  -nep_stop_basic - Sets the default stopping test
822: -  -nep_stop_user  - Selects the user-defined stopping test

824:    Note:
825:    The parameter 'stop' can have one of these values
826: +     NEP_STOP_BASIC - default stopping test
827: -     NEP_STOP_USER  - function set by NEPSetStoppingTestFunction()

829:    Level: advanced

831: .seealso: NEPGetStoppingTest(), NEPSetStoppingTestFunction(), NEPSetConvergenceTest(), NEPStop
832: @*/
833: PetscErrorCode NEPSetStoppingTest(NEP nep,NEPStop stop)
834: {
835:   PetscFunctionBegin;
838:   switch (stop) {
839:     case NEP_STOP_BASIC: nep->stopping = NEPStoppingBasic; break;
840:     case NEP_STOP_USER:
841:       PetscCheck(nep->stoppinguser,PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"Must call NEPSetStoppingTestFunction() first");
842:       nep->stopping = nep->stoppinguser;
843:       break;
844:     default:
845:       SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
846:   }
847:   nep->stop = stop;
848:   PetscFunctionReturn(PETSC_SUCCESS);
849: }

851: /*@
852:    NEPGetStoppingTest - Gets the method used to decide the termination of the outer
853:    loop of the eigensolver.

855:    Not Collective

857:    Input Parameters:
858: .  nep   - nonlinear eigensolver context obtained from NEPCreate()

860:    Output Parameters:
861: .  stop  - the type of stopping test

863:    Level: advanced

865: .seealso: NEPSetStoppingTest(), NEPStop
866: @*/
867: PetscErrorCode NEPGetStoppingTest(NEP nep,NEPStop *stop)
868: {
869:   PetscFunctionBegin;
871:   PetscAssertPointer(stop,2);
872:   *stop = nep->stop;
873:   PetscFunctionReturn(PETSC_SUCCESS);
874: }

876: /*@
877:    NEPSetTrackAll - Specifies if the solver must compute the residual of all
878:    approximate eigenpairs or not.

880:    Logically Collective

882:    Input Parameters:
883: +  nep      - the eigensolver context
884: -  trackall - whether compute all residuals or not

886:    Notes:
887:    If the user sets trackall=PETSC_TRUE then the solver explicitly computes
888:    the residual for each eigenpair approximation. Computing the residual is
889:    usually an expensive operation and solvers commonly compute the associated
890:    residual to the first unconverged eigenpair.

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

894:    Level: developer

896: .seealso: NEPGetTrackAll()
897: @*/
898: PetscErrorCode NEPSetTrackAll(NEP nep,PetscBool trackall)
899: {
900:   PetscFunctionBegin;
903:   nep->trackall = trackall;
904:   PetscFunctionReturn(PETSC_SUCCESS);
905: }

907: /*@
908:    NEPGetTrackAll - Returns the flag indicating whether all residual norms must
909:    be computed or not.

911:    Not Collective

913:    Input Parameter:
914: .  nep - the eigensolver context

916:    Output Parameter:
917: .  trackall - the returned flag

919:    Level: developer

921: .seealso: NEPSetTrackAll()
922: @*/
923: PetscErrorCode NEPGetTrackAll(NEP nep,PetscBool *trackall)
924: {
925:   PetscFunctionBegin;
927:   PetscAssertPointer(trackall,2);
928:   *trackall = nep->trackall;
929:   PetscFunctionReturn(PETSC_SUCCESS);
930: }

932: /*@
933:    NEPSetRefine - Specifies the refinement type (and options) to be used
934:    after the solve.

936:    Logically Collective

938:    Input Parameters:
939: +  nep    - the nonlinear eigensolver context
940: .  refine - refinement type
941: .  npart  - number of partitions of the communicator
942: .  tol    - the convergence tolerance
943: .  its    - maximum number of refinement iterations
944: -  scheme - which scheme to be used for solving the involved linear systems

946:    Options Database Keys:
947: +  -nep_refine <type> - refinement type, one of <none,simple,multiple>
948: .  -nep_refine_partitions <n> - the number of partitions
949: .  -nep_refine_tol <tol> - the tolerance
950: .  -nep_refine_its <its> - number of iterations
951: -  -nep_refine_scheme - to set the scheme for the linear solves

953:    Notes:
954:    By default, iterative refinement is disabled, since it may be very
955:    costly. There are two possible refinement strategies, simple and multiple.
956:    The simple approach performs iterative refinement on each of the
957:    converged eigenpairs individually, whereas the multiple strategy works
958:    with the invariant pair as a whole, refining all eigenpairs simultaneously.
959:    The latter may be required for the case of multiple eigenvalues.

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

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

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

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

979:    Level: intermediate

981: .seealso: NEPGetRefine()
982: @*/
983: PetscErrorCode NEPSetRefine(NEP nep,NEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,NEPRefineScheme scheme)
984: {
985:   PetscMPIInt    size;

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

1025: /*@
1026:    NEPGetRefine - Gets the refinement strategy used by the NEP object, and the
1027:    associated parameters.

1029:    Not Collective

1031:    Input Parameter:
1032: .  nep - the nonlinear eigensolver context

1034:    Output Parameters:
1035: +  refine - refinement type
1036: .  npart  - number of partitions of the communicator
1037: .  tol    - the convergence tolerance
1038: .  its    - maximum number of refinement iterations
1039: -  scheme - the scheme used for solving linear systems

1041:    Level: intermediate

1043:    Note:
1044:    The user can specify NULL for any parameter that is not needed.

1046: .seealso: NEPSetRefine()
1047: @*/
1048: PetscErrorCode NEPGetRefine(NEP nep,NEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,NEPRefineScheme *scheme)
1049: {
1050:   PetscFunctionBegin;
1052:   if (refine) *refine = nep->refine;
1053:   if (npart)  *npart  = nep->npart;
1054:   if (tol)    *tol    = nep->rtol;
1055:   if (its)    *its    = nep->rits;
1056:   if (scheme) *scheme = nep->scheme;
1057:   PetscFunctionReturn(PETSC_SUCCESS);
1058: }

1060: /*@
1061:    NEPSetOptionsPrefix - Sets the prefix used for searching for all
1062:    NEP options in the database.

1064:    Logically Collective

1066:    Input Parameters:
1067: +  nep - the nonlinear eigensolver context
1068: -  prefix - the prefix string to prepend to all NEP option requests

1070:    Notes:
1071:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1072:    The first character of all runtime options is AUTOMATICALLY the
1073:    hyphen.

1075:    For example, to distinguish between the runtime options for two
1076:    different NEP contexts, one could call
1077: .vb
1078:       NEPSetOptionsPrefix(nep1,"neig1_")
1079:       NEPSetOptionsPrefix(nep2,"neig2_")
1080: .ve

1082:    Level: advanced

1084: .seealso: NEPAppendOptionsPrefix(), NEPGetOptionsPrefix()
1085: @*/
1086: PetscErrorCode NEPSetOptionsPrefix(NEP nep,const char *prefix)
1087: {
1088:   PetscFunctionBegin;
1090:   if (!nep->V) PetscCall(NEPGetBV(nep,&nep->V));
1091:   PetscCall(BVSetOptionsPrefix(nep->V,prefix));
1092:   if (!nep->ds) PetscCall(NEPGetDS(nep,&nep->ds));
1093:   PetscCall(DSSetOptionsPrefix(nep->ds,prefix));
1094:   if (!nep->rg) PetscCall(NEPGetRG(nep,&nep->rg));
1095:   PetscCall(RGSetOptionsPrefix(nep->rg,prefix));
1096:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)nep,prefix));
1097:   PetscFunctionReturn(PETSC_SUCCESS);
1098: }

1100: /*@
1101:    NEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1102:    NEP options in the database.

1104:    Logically Collective

1106:    Input Parameters:
1107: +  nep - the nonlinear eigensolver context
1108: -  prefix - the prefix string to prepend to all NEP option requests

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

1114:    Level: advanced

1116: .seealso: NEPSetOptionsPrefix(), NEPGetOptionsPrefix()
1117: @*/
1118: PetscErrorCode NEPAppendOptionsPrefix(NEP nep,const char *prefix)
1119: {
1120:   PetscFunctionBegin;
1122:   if (!nep->V) PetscCall(NEPGetBV(nep,&nep->V));
1123:   PetscCall(BVAppendOptionsPrefix(nep->V,prefix));
1124:   if (!nep->ds) PetscCall(NEPGetDS(nep,&nep->ds));
1125:   PetscCall(DSAppendOptionsPrefix(nep->ds,prefix));
1126:   if (!nep->rg) PetscCall(NEPGetRG(nep,&nep->rg));
1127:   PetscCall(RGAppendOptionsPrefix(nep->rg,prefix));
1128:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)nep,prefix));
1129:   PetscFunctionReturn(PETSC_SUCCESS);
1130: }

1132: /*@
1133:    NEPGetOptionsPrefix - Gets the prefix used for searching for all
1134:    NEP options in the database.

1136:    Not Collective

1138:    Input Parameters:
1139: .  nep - the nonlinear eigensolver context

1141:    Output Parameters:
1142: .  prefix - pointer to the prefix string used is returned

1144:    Note:
1145:    On the Fortran side, the user should pass in a string 'prefix' of
1146:    sufficient length to hold the prefix.

1148:    Level: advanced

1150: .seealso: NEPSetOptionsPrefix(), NEPAppendOptionsPrefix()
1151: @*/
1152: PetscErrorCode NEPGetOptionsPrefix(NEP nep,const char *prefix[])
1153: {
1154:   PetscFunctionBegin;
1156:   PetscAssertPointer(prefix,2);
1157:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)nep,prefix));
1158:   PetscFunctionReturn(PETSC_SUCCESS);
1159: }