1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, 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> /*I "slepcnep.h" I*/
16: #include <petscdraw.h>
18: /*@C
19: NEPMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
20: indicated by the user.
22: Collective on NEP 24: Input Parameters:
25: + nep - the nonlinear eigensolver context
26: . name - the monitor option name
27: . help - message indicating what monitoring is done
28: . manual - manual page for the monitor
29: . monitor - the monitor function, whose context is a PetscViewerAndFormat
30: - trackall - whether this monitor tracks all eigenvalues or not
32: Level: developer
34: .seealso: NEPMonitorSet(), NEPSetTrackAll(), NEPConvMonitorSetFromOptions()
35: @*/
36: PetscErrorCode NEPMonitorSetFromOptions(NEP nep,const char name[],const char help[],const char manual[],PetscErrorCode (*monitor)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscBool trackall) 37: {
38: PetscErrorCode ierr;
39: PetscBool flg;
40: PetscViewer viewer;
41: PetscViewerFormat format;
42: PetscViewerAndFormat *vf;
45: PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->options,((PetscObject)nep)->prefix,name,&viewer,&format,&flg);
46: if (flg) {
47: PetscViewerAndFormatCreate(viewer,format,&vf);
48: PetscObjectDereference((PetscObject)viewer);
49: NEPMonitorSet(nep,(PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
50: if (trackall) {
51: NEPSetTrackAll(nep,PETSC_TRUE);
52: }
53: }
54: return(0);
55: }
57: /*@C
58: NEPConvMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
59: indicated by the user (for monitors that only show iteration numbers of convergence).
61: Collective on NEP 63: Input Parameters:
64: + nep - the nonlinear eigensolver context
65: . name - the monitor option name
66: . help - message indicating what monitoring is done
67: . manual - manual page for the monitor
68: - monitor - the monitor function, whose context is a SlepcConvMonitor
70: Level: developer
72: .seealso: NEPMonitorSet(), NEPMonitorSetFromOptions()
73: @*/
74: PetscErrorCode NEPConvMonitorSetFromOptions(NEP nep,const char name[],const char help[],const char manual[],PetscErrorCode (*monitor)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,SlepcConvMonitor)) 75: {
76: PetscErrorCode ierr;
77: PetscBool flg;
78: PetscViewer viewer;
79: PetscViewerFormat format;
80: SlepcConvMonitor ctx;
83: PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->options,((PetscObject)nep)->prefix,name,&viewer,&format,&flg);
84: if (flg) {
85: SlepcConvMonitorCreate(viewer,format,&ctx);
86: PetscObjectDereference((PetscObject)viewer);
87: NEPMonitorSet(nep,(PetscErrorCode (*)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))monitor,ctx,(PetscErrorCode (*)(void**))SlepcConvMonitorDestroy);
88: }
89: return(0);
90: }
92: /*@
93: NEPSetFromOptions - Sets NEP options from the options database.
94: This routine must be called before NEPSetUp() if the user is to be
95: allowed to set the solver type.
97: Collective on NEP 99: Input Parameters:
100: . nep - the nonlinear eigensolver context
102: Notes:
103: To see all options, run your program with the -help option.
105: Level: beginner
106: @*/
107: PetscErrorCode NEPSetFromOptions(NEP nep)108: {
109: PetscErrorCode ierr;
110: char type[256];
111: PetscBool set,flg,flg1,flg2,flg3,flg4,flg5,bval;
112: PetscReal r;
113: PetscScalar s;
114: PetscInt i,j,k;
115: PetscDrawLG lg;
116: NEPRefine refine;
117: NEPRefineScheme scheme;
121: NEPRegisterAll();
122: PetscObjectOptionsBegin((PetscObject)nep);
123: PetscOptionsFList("-nep_type","Nonlinear eigensolver method","NEPSetType",NEPList,(char*)(((PetscObject)nep)->type_name?((PetscObject)nep)->type_name:NEPRII),type,256,&flg);
124: if (flg) {
125: NEPSetType(nep,type);
126: } else if (!((PetscObject)nep)->type_name) {
127: NEPSetType(nep,NEPRII);
128: }
130: PetscOptionsBoolGroupBegin("-nep_general","General nonlinear eigenvalue problem","NEPSetProblemType",&flg);
131: if (flg) { NEPSetProblemType(nep,NEP_GENERAL); }
132: PetscOptionsBoolGroupEnd("-nep_rational","Rational eigenvalue problem","NEPSetProblemType",&flg);
133: if (flg) { NEPSetProblemType(nep,NEP_RATIONAL); }
135: refine = nep->refine;
136: PetscOptionsEnum("-nep_refine","Iterative refinement method","NEPSetRefine",NEPRefineTypes,(PetscEnum)refine,(PetscEnum*)&refine,&flg1);
137: i = nep->npart;
138: PetscOptionsInt("-nep_refine_partitions","Number of partitions of the communicator for iterative refinement","NEPSetRefine",nep->npart,&i,&flg2);
139: r = nep->rtol;
140: PetscOptionsReal("-nep_refine_tol","Tolerance for iterative refinement","NEPSetRefine",nep->rtol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/1000:nep->rtol,&r,&flg3);
141: j = nep->rits;
142: PetscOptionsInt("-nep_refine_its","Maximum number of iterations for iterative refinement","NEPSetRefine",nep->rits,&j,&flg4);
143: scheme = nep->scheme;
144: PetscOptionsEnum("-nep_refine_scheme","Scheme used for linear systems within iterative refinement","NEPSetRefine",NEPRefineSchemes,(PetscEnum)scheme,(PetscEnum*)&scheme,&flg5);
145: if (flg1 || flg2 || flg3 || flg4 || flg5) { NEPSetRefine(nep,refine,i,r,j,scheme); }
147: i = nep->max_it? nep->max_it: PETSC_DEFAULT;
148: PetscOptionsInt("-nep_max_it","Maximum number of iterations","NEPSetTolerances",nep->max_it,&i,&flg1);
149: r = nep->tol;
150: PetscOptionsReal("-nep_tol","Tolerance","NEPSetTolerances",nep->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:nep->tol,&r,&flg2);
151: if (flg1 || flg2) { NEPSetTolerances(nep,r,i); }
153: PetscOptionsBoolGroupBegin("-nep_conv_rel","Relative error convergence test","NEPSetConvergenceTest",&flg);
154: if (flg) { NEPSetConvergenceTest(nep,NEP_CONV_REL); }
155: PetscOptionsBoolGroup("-nep_conv_norm","Convergence test relative to the matrix norms","NEPSetConvergenceTest",&flg);
156: if (flg) { NEPSetConvergenceTest(nep,NEP_CONV_NORM); }
157: PetscOptionsBoolGroup("-nep_conv_abs","Absolute error convergence test","NEPSetConvergenceTest",&flg);
158: if (flg) { NEPSetConvergenceTest(nep,NEP_CONV_ABS); }
159: PetscOptionsBoolGroupEnd("-nep_conv_user","User-defined convergence test","NEPSetConvergenceTest",&flg);
160: if (flg) { NEPSetConvergenceTest(nep,NEP_CONV_USER); }
162: PetscOptionsBoolGroupBegin("-nep_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","NEPSetStoppingTest",&flg);
163: if (flg) { NEPSetStoppingTest(nep,NEP_STOP_BASIC); }
164: PetscOptionsBoolGroupEnd("-nep_stop_user","User-defined stopping test","NEPSetStoppingTest",&flg);
165: if (flg) { NEPSetStoppingTest(nep,NEP_STOP_USER); }
167: i = nep->nev;
168: PetscOptionsInt("-nep_nev","Number of eigenvalues to compute","NEPSetDimensions",nep->nev,&i,&flg1);
169: j = nep->ncv? nep->ncv: PETSC_DEFAULT;
170: PetscOptionsInt("-nep_ncv","Number of basis vectors","NEPSetDimensions",nep->ncv,&j,&flg2);
171: k = nep->mpd? nep->mpd: PETSC_DEFAULT;
172: PetscOptionsInt("-nep_mpd","Maximum dimension of projected problem","NEPSetDimensions",nep->mpd,&k,&flg3);
173: if (flg1 || flg2 || flg3) {
174: NEPSetDimensions(nep,i,j,k);
175: }
177: PetscOptionsBoolGroupBegin("-nep_largest_magnitude","Compute largest eigenvalues in magnitude","NEPSetWhichEigenpairs",&flg);
178: if (flg) { NEPSetWhichEigenpairs(nep,NEP_LARGEST_MAGNITUDE); }
179: PetscOptionsBoolGroup("-nep_smallest_magnitude","Compute smallest eigenvalues in magnitude","NEPSetWhichEigenpairs",&flg);
180: if (flg) { NEPSetWhichEigenpairs(nep,NEP_SMALLEST_MAGNITUDE); }
181: PetscOptionsBoolGroup("-nep_largest_real","Compute eigenvalues with largest real parts","NEPSetWhichEigenpairs",&flg);
182: if (flg) { NEPSetWhichEigenpairs(nep,NEP_LARGEST_REAL); }
183: PetscOptionsBoolGroup("-nep_smallest_real","Compute eigenvalues with smallest real parts","NEPSetWhichEigenpairs",&flg);
184: if (flg) { NEPSetWhichEigenpairs(nep,NEP_SMALLEST_REAL); }
185: PetscOptionsBoolGroup("-nep_largest_imaginary","Compute eigenvalues with largest imaginary parts","NEPSetWhichEigenpairs",&flg);
186: if (flg) { NEPSetWhichEigenpairs(nep,NEP_LARGEST_IMAGINARY); }
187: PetscOptionsBoolGroup("-nep_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","NEPSetWhichEigenpairs",&flg);
188: if (flg) { NEPSetWhichEigenpairs(nep,NEP_SMALLEST_IMAGINARY); }
189: PetscOptionsBoolGroup("-nep_target_magnitude","Compute eigenvalues closest to target","NEPSetWhichEigenpairs",&flg);
190: if (flg) { NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE); }
191: PetscOptionsBoolGroup("-nep_target_real","Compute eigenvalues with real parts closest to target","NEPSetWhichEigenpairs",&flg);
192: if (flg) { NEPSetWhichEigenpairs(nep,NEP_TARGET_REAL); }
193: PetscOptionsBoolGroup("-nep_target_imaginary","Compute eigenvalues with imaginary parts closest to target","NEPSetWhichEigenpairs",&flg);
194: if (flg) { NEPSetWhichEigenpairs(nep,NEP_TARGET_IMAGINARY); }
195: PetscOptionsBoolGroupEnd("-nep_all","Compute all eigenvalues in a region","NEPSetWhichEigenpairs",&flg);
196: if (flg) { NEPSetWhichEigenpairs(nep,NEP_ALL); }
198: PetscOptionsScalar("-nep_target","Value of the target","NEPSetTarget",nep->target,&s,&flg);
199: if (flg) {
200: if (nep->which!=NEP_TARGET_REAL && nep->which!=NEP_TARGET_IMAGINARY) {
201: NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE);
202: }
203: NEPSetTarget(nep,s);
204: }
206: PetscOptionsBool("-nep_two_sided","Use two-sided variant (to compute left eigenvectors)","NEPSetTwoSided",nep->twosided,&bval,&flg);
207: if (flg) { NEPSetTwoSided(nep,bval); }
209: /* -----------------------------------------------------------------------*/
210: /*
211: Cancels all monitors hardwired into code before call to NEPSetFromOptions()
212: */
213: PetscOptionsBool("-nep_monitor_cancel","Remove any hardwired monitor routines","NEPMonitorCancel",PETSC_FALSE,&flg,&set);
214: if (set && flg) {
215: NEPMonitorCancel(nep);
216: }
217: /*
218: Text monitors
219: */
220: NEPMonitorSetFromOptions(nep,"-nep_monitor","Monitor first unconverged approximate eigenvalue and error estimate","NEPMonitorFirst",NEPMonitorFirst,PETSC_FALSE);
221: NEPConvMonitorSetFromOptions(nep,"-nep_monitor_conv","Monitor approximate eigenvalues and error estimates as they converge","NEPMonitorConverged",NEPMonitorConverged);
222: NEPMonitorSetFromOptions(nep,"-nep_monitor_all","Monitor approximate eigenvalues and error estimates","NEPMonitorAll",NEPMonitorAll,PETSC_TRUE);
223: /*
224: Line graph monitors
225: */
226: PetscOptionsBool("-nep_monitor_lg","Monitor first unconverged approximate error estimate graphically","NEPMonitorSet",PETSC_FALSE,&flg,&set);
227: if (set && flg) {
228: NEPMonitorLGCreate(PetscObjectComm((PetscObject)nep),NULL,"Error estimates",PETSC_DECIDE,PETSC_DECIDE,300,300,&lg);
229: NEPMonitorSet(nep,NEPMonitorLG,lg,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
230: }
231: PetscOptionsBool("-nep_monitor_lg_all","Monitor error estimates graphically","NEPMonitorSet",PETSC_FALSE,&flg,&set);
232: if (set && flg) {
233: NEPMonitorLGCreate(PetscObjectComm((PetscObject)nep),NULL,"Error estimates",PETSC_DECIDE,PETSC_DECIDE,300,300,&lg);
234: NEPMonitorSet(nep,NEPMonitorLGAll,lg,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
235: NEPSetTrackAll(nep,PETSC_TRUE);
236: }
238: /* -----------------------------------------------------------------------*/
239: PetscOptionsName("-nep_view","Print detailed information on solver used","NEPView",NULL);
240: PetscOptionsName("-nep_view_vectors","View computed eigenvectors","NEPVectorsView",NULL);
241: PetscOptionsName("-nep_view_values","View computed eigenvalues","NEPValuesView",NULL);
242: PetscOptionsName("-nep_converged_reason","Print reason for convergence, and number of iterations","NEPReasonView",NULL);
243: PetscOptionsName("-nep_error_absolute","Print absolute errors of each eigenpair","NEPErrorView",NULL);
244: PetscOptionsName("-nep_error_relative","Print relative errors of each eigenpair","NEPErrorView",NULL);
246: if (nep->ops->setfromoptions) {
247: (*nep->ops->setfromoptions)(PetscOptionsObject,nep);
248: }
249: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)nep);
250: PetscOptionsEnd();
252: if (!nep->V) { NEPGetBV(nep,&nep->V); }
253: BVSetFromOptions(nep->V);
254: if (!nep->rg) { NEPGetRG(nep,&nep->rg); }
255: RGSetFromOptions(nep->rg);
256: if (nep->useds) {
257: if (!nep->ds) { NEPGetDS(nep,&nep->ds); }
258: DSSetFromOptions(nep->ds);
259: }
260: if (!nep->refineksp) { NEPRefineGetKSP(nep,&nep->refineksp); }
261: KSPSetFromOptions(nep->refineksp);
262: if (nep->fui==NEP_USER_INTERFACE_SPLIT) for (i=0;i<nep->nt;i++) {FNSetFromOptions(nep->f[i]);}
263: return(0);
264: }
266: /*@C
267: NEPGetTolerances - Gets the tolerance and maximum iteration count used
268: by the NEP convergence tests.
270: Not Collective
272: Input Parameter:
273: . nep - the nonlinear eigensolver context
275: Output Parameters:
276: + tol - the convergence tolerance
277: - maxits - maximum number of iterations
279: Notes:
280: The user can specify NULL for any parameter that is not needed.
282: Level: intermediate
284: .seealso: NEPSetTolerances()
285: @*/
286: PetscErrorCode NEPGetTolerances(NEP nep,PetscReal *tol,PetscInt *maxits)287: {
290: if (tol) *tol = nep->tol;
291: if (maxits) *maxits = nep->max_it;
292: return(0);
293: }
295: /*@
296: NEPSetTolerances - Sets the tolerance and maximum iteration count used
297: by the NEP convergence tests.
299: Logically Collective on NEP301: Input Parameters:
302: + nep - the nonlinear eigensolver context
303: . tol - the convergence tolerance
304: - maxits - maximum number of iterations to use
306: Options Database Keys:
307: + -nep_tol <tol> - Sets the convergence tolerance
308: - -nep_max_it <maxits> - Sets the maximum number of iterations allowed
310: Notes:
311: Use PETSC_DEFAULT for either argument to assign a reasonably good value.
313: Level: intermediate
315: .seealso: NEPGetTolerances()
316: @*/
317: PetscErrorCode NEPSetTolerances(NEP nep,PetscReal tol,PetscInt maxits)318: {
323: if (tol == PETSC_DEFAULT) {
324: nep->tol = PETSC_DEFAULT;
325: nep->state = NEP_STATE_INITIAL;
326: } else {
327: if (tol <= 0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
328: nep->tol = tol;
329: }
330: if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
331: nep->max_it = 0;
332: nep->state = NEP_STATE_INITIAL;
333: } else {
334: if (maxits <= 0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
335: nep->max_it = maxits;
336: }
337: return(0);
338: }
340: /*@C
341: NEPGetDimensions - Gets the number of eigenvalues to compute
342: and the dimension of the subspace.
344: Not Collective
346: Input Parameter:
347: . nep - the nonlinear eigensolver context
349: Output Parameters:
350: + nev - number of eigenvalues to compute
351: . ncv - the maximum dimension of the subspace to be used by the solver
352: - mpd - the maximum dimension allowed for the projected problem
354: Notes:
355: The user can specify NULL for any parameter that is not needed.
357: Level: intermediate
359: .seealso: NEPSetDimensions()
360: @*/
361: PetscErrorCode NEPGetDimensions(NEP nep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)362: {
365: if (nev) *nev = nep->nev;
366: if (ncv) *ncv = nep->ncv;
367: if (mpd) *mpd = nep->mpd;
368: return(0);
369: }
371: /*@
372: NEPSetDimensions - Sets the number of eigenvalues to compute
373: and the dimension of the subspace.
375: Logically Collective on NEP377: Input Parameters:
378: + nep - the nonlinear eigensolver context
379: . nev - number of eigenvalues to compute
380: . ncv - the maximum dimension of the subspace to be used by the solver
381: - mpd - the maximum dimension allowed for the projected problem
383: Options Database Keys:
384: + -nep_nev <nev> - Sets the number of eigenvalues
385: . -nep_ncv <ncv> - Sets the dimension of the subspace
386: - -nep_mpd <mpd> - Sets the maximum projected dimension
388: Notes:
389: Use PETSC_DEFAULT for ncv and mpd to assign a reasonably good value, which is
390: dependent on the solution method.
392: The parameters ncv and mpd are intimately related, so that the user is advised
393: to set one of them at most. Normal usage is that
394: (a) in cases where nev is small, the user sets ncv (a reasonable default is 2*nev); and
395: (b) in cases where nev is large, the user sets mpd.
397: The value of ncv should always be between nev and (nev+mpd), typically
398: ncv=nev+mpd. If nev is not too large, mpd=nev is a reasonable choice, otherwise
399: a smaller value should be used.
401: Level: intermediate
403: .seealso: NEPGetDimensions()
404: @*/
405: PetscErrorCode NEPSetDimensions(NEP nep,PetscInt nev,PetscInt ncv,PetscInt mpd)406: {
412: if (nev<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
413: nep->nev = nev;
414: if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
415: nep->ncv = 0;
416: } else {
417: if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
418: nep->ncv = ncv;
419: }
420: if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
421: nep->mpd = 0;
422: } else {
423: if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
424: nep->mpd = mpd;
425: }
426: nep->state = NEP_STATE_INITIAL;
427: return(0);
428: }
430: /*@
431: NEPSetWhichEigenpairs - Specifies which portion of the spectrum is
432: to be sought.
434: Logically Collective on NEP436: Input Parameters:
437: + nep - eigensolver context obtained from NEPCreate()
438: - which - the portion of the spectrum to be sought
440: Possible values:
441: The parameter 'which' can have one of these values
443: + NEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
444: . NEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
445: . NEP_LARGEST_REAL - largest real parts
446: . NEP_SMALLEST_REAL - smallest real parts
447: . NEP_LARGEST_IMAGINARY - largest imaginary parts
448: . NEP_SMALLEST_IMAGINARY - smallest imaginary parts
449: . NEP_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
450: . NEP_TARGET_REAL - eigenvalues with real part closest to target
451: . NEP_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
452: . NEP_ALL - all eigenvalues contained in a given region
453: - NEP_WHICH_USER - user defined ordering set with NEPSetEigenvalueComparison()
455: Options Database Keys:
456: + -nep_largest_magnitude - Sets largest eigenvalues in magnitude
457: . -nep_smallest_magnitude - Sets smallest eigenvalues in magnitude
458: . -nep_largest_real - Sets largest real parts
459: . -nep_smallest_real - Sets smallest real parts
460: . -nep_largest_imaginary - Sets largest imaginary parts
461: . -nep_smallest_imaginary - Sets smallest imaginary parts
462: . -nep_target_magnitude - Sets eigenvalues closest to target
463: . -nep_target_real - Sets real parts closest to target
464: . -nep_target_imaginary - Sets imaginary parts closest to target
465: - -nep_all - Sets all eigenvalues in a region
467: Notes:
468: Not all eigensolvers implemented in NEP account for all the possible values
469: stated above. If SLEPc is compiled for real numbers NEP_LARGEST_IMAGINARY470: and NEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
471: for eigenvalue selection.
473: The target is a scalar value provided with NEPSetTarget().
475: NEP_ALL is intended for use in the context of the CISS solver for
476: computing all eigenvalues in a region.
478: Level: intermediate
480: .seealso: NEPGetWhichEigenpairs(), NEPSetTarget(), NEPSetEigenvalueComparison(), NEPWhich481: @*/
482: PetscErrorCode NEPSetWhichEigenpairs(NEP nep,NEPWhich which)483: {
487: switch (which) {
488: case NEP_LARGEST_MAGNITUDE:
489: case NEP_SMALLEST_MAGNITUDE:
490: case NEP_LARGEST_REAL:
491: case NEP_SMALLEST_REAL:
492: case NEP_LARGEST_IMAGINARY:
493: case NEP_SMALLEST_IMAGINARY:
494: case NEP_TARGET_MAGNITUDE:
495: case NEP_TARGET_REAL:
496: #if defined(PETSC_USE_COMPLEX)
497: case NEP_TARGET_IMAGINARY:
498: #endif
499: case NEP_ALL:
500: case NEP_WHICH_USER:
501: if (nep->which != which) {
502: nep->state = NEP_STATE_INITIAL;
503: nep->which = which;
504: }
505: break;
506: default:507: SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
508: }
509: return(0);
510: }
512: /*@
513: NEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
514: sought.
516: Not Collective
518: Input Parameter:
519: . nep - eigensolver context obtained from NEPCreate()
521: Output Parameter:
522: . which - the portion of the spectrum to be sought
524: Notes:
525: See NEPSetWhichEigenpairs() for possible values of 'which'.
527: Level: intermediate
529: .seealso: NEPSetWhichEigenpairs(), NEPWhich530: @*/
531: PetscErrorCode NEPGetWhichEigenpairs(NEP nep,NEPWhich *which)532: {
536: *which = nep->which;
537: return(0);
538: }
540: /*@C
541: NEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
542: when NEPSetWhichEigenpairs() is set to NEP_WHICH_USER.
544: Logically Collective on NEP546: Input Parameters:
547: + nep - eigensolver context obtained from NEPCreate()
548: . func - a pointer to the comparison function
549: - ctx - a context pointer (the last parameter to the comparison function)
551: Calling Sequence of func:
552: $ func(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)
554: + ar - real part of the 1st eigenvalue
555: . ai - imaginary part of the 1st eigenvalue
556: . br - real part of the 2nd eigenvalue
557: . bi - imaginary part of the 2nd eigenvalue
558: . res - result of comparison
559: - ctx - optional context, as set by NEPSetEigenvalueComparison()
561: Note:
562: The returning parameter 'res' can be
563: + negative - if the 1st eigenvalue is preferred to the 2st one
564: . zero - if both eigenvalues are equally preferred
565: - positive - if the 2st eigenvalue is preferred to the 1st one
567: Level: advanced
569: .seealso: NEPSetWhichEigenpairs(), NEPWhich570: @*/
571: PetscErrorCode NEPSetEigenvalueComparison(NEP nep,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)572: {
575: nep->sc->comparison = func;
576: nep->sc->comparisonctx = ctx;
577: nep->which = NEP_WHICH_USER;
578: return(0);
579: }
581: /*@
582: NEPSetProblemType - Specifies the type of the nonlinear eigenvalue problem.
584: Logically Collective on NEP586: Input Parameters:
587: + nep - the nonlinear eigensolver context
588: - type - a known type of nonlinear eigenvalue problem
590: Options Database Keys:
591: + -nep_general - general problem with no particular structure
592: - -nep_rational - a rational eigenvalue problem defined in split form with all f_i rational
594: Notes:
595: Allowed values for the problem type are: general (NEP_GENERAL), and rational
596: (NEP_RATIONAL).
598: This function is used to provide a hint to the NEP solver to exploit certain
599: properties of the nonlinear eigenproblem. This hint may be used or not,
600: depending on the solver. By default, no particular structure is assumed.
602: Level: intermediate
604: .seealso: NEPSetType(), NEPGetProblemType(), NEPProblemType605: @*/
606: PetscErrorCode NEPSetProblemType(NEP nep,NEPProblemType type)607: {
611: if (type!=NEP_GENERAL && type!=NEP_RATIONAL) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
612: if (type != nep->problem_type) {
613: nep->problem_type = type;
614: nep->state = NEP_STATE_INITIAL;
615: }
616: return(0);
617: }
619: /*@
620: NEPGetProblemType - Gets the problem type from the NEP object.
622: Not Collective
624: Input Parameter:
625: . nep - the nonlinear eigensolver context
627: Output Parameter:
628: . type - the problem type
630: Level: intermediate
632: .seealso: NEPSetProblemType(), NEPProblemType633: @*/
634: PetscErrorCode NEPGetProblemType(NEP nep,NEPProblemType *type)635: {
639: *type = nep->problem_type;
640: return(0);
641: }
643: /*@
644: NEPSetTwoSided - Sets the solver to use a two-sided variant so that left
645: eigenvectors are also computed.
647: Logically Collective on NEP649: Input Parameters:
650: + nep - the eigensolver context
651: - twosided - whether the two-sided variant is to be used or not
653: Options Database Keys:
654: . -nep_two_sided <boolean> - Sets/resets the twosided flag
656: Notes:
657: If the user sets twosided=PETSC_TRUE then the solver uses a variant of
658: the algorithm that computes both right and left eigenvectors. This is
659: usually much more costly. This option is not available in all solvers.
661: When using two-sided solvers, the problem matrices must have both the
662: MatMult and MatMultTranspose operations defined.
664: Level: advanced
666: .seealso: NEPGetTwoSided(), NEPGetLeftEigenvector()
667: @*/
668: PetscErrorCode NEPSetTwoSided(NEP nep,PetscBool twosided)669: {
673: if (twosided!=nep->twosided) {
674: nep->twosided = twosided;
675: nep->state = NEP_STATE_INITIAL;
676: }
677: return(0);
678: }
680: /*@
681: NEPGetTwoSided - Returns the flag indicating whether a two-sided variant
682: of the algorithm is being used or not.
684: Not Collective
686: Input Parameter:
687: . nep - the eigensolver context
689: Output Parameter:
690: . twosided - the returned flag
692: Level: advanced
694: .seealso: NEPSetTwoSided()
695: @*/
696: PetscErrorCode NEPGetTwoSided(NEP nep,PetscBool *twosided)697: {
701: *twosided = nep->twosided;
702: return(0);
703: }
705: /*@C
706: NEPSetConvergenceTestFunction - Sets a function to compute the error estimate
707: used in the convergence test.
709: Logically Collective on NEP711: Input Parameters:
712: + nep - nonlinear eigensolver context obtained from NEPCreate()
713: . func - a pointer to the convergence test function
714: . ctx - context for private data for the convergence routine (may be null)
715: - destroy - a routine for destroying the context (may be null)
717: Calling Sequence of func:
718: $ func(NEP nep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
720: + nep - nonlinear eigensolver context obtained from NEPCreate()
721: . eigr - real part of the eigenvalue
722: . eigi - imaginary part of the eigenvalue
723: . res - residual norm associated to the eigenpair
724: . errest - (output) computed error estimate
725: - ctx - optional context, as set by NEPSetConvergenceTestFunction()
727: Note:
728: If the error estimate returned by the convergence test function is less than
729: the tolerance, then the eigenvalue is accepted as converged.
731: Level: advanced
733: .seealso: NEPSetConvergenceTest(), NEPSetTolerances()
734: @*/
735: PetscErrorCode NEPSetConvergenceTestFunction(NEP nep,PetscErrorCode (*func)(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))736: {
741: if (nep->convergeddestroy) {
742: (*nep->convergeddestroy)(nep->convergedctx);
743: }
744: nep->convergeduser = func;
745: nep->convergeddestroy = destroy;
746: nep->convergedctx = ctx;
747: if (func == NEPConvergedRelative) nep->conv = NEP_CONV_REL;
748: else if (func == NEPConvergedNorm) nep->conv = NEP_CONV_NORM;
749: else if (func == NEPConvergedAbsolute) nep->conv = NEP_CONV_ABS;
750: else {
751: nep->conv = NEP_CONV_USER;
752: nep->converged = nep->convergeduser;
753: }
754: return(0);
755: }
757: /*@
758: NEPSetConvergenceTest - Specifies how to compute the error estimate
759: used in the convergence test.
761: Logically Collective on NEP763: Input Parameters:
764: + nep - nonlinear eigensolver context obtained from NEPCreate()
765: - conv - the type of convergence test
767: Options Database Keys:
768: + -nep_conv_abs - Sets the absolute convergence test
769: . -nep_conv_rel - Sets the convergence test relative to the eigenvalue
770: - -nep_conv_user - Selects the user-defined convergence test
772: Note:
773: The parameter 'conv' can have one of these values
774: + NEP_CONV_ABS - absolute error ||r||
775: . NEP_CONV_REL - error relative to the eigenvalue l, ||r||/|l|
776: . NEP_CONV_NORM - error relative matrix norms, ||r||/sum_i(|f_i(l)|*||A_i||)
777: - NEP_CONV_USER - function set by NEPSetConvergenceTestFunction()
779: Level: intermediate
781: .seealso: NEPGetConvergenceTest(), NEPSetConvergenceTestFunction(), NEPSetStoppingTest(), NEPConv782: @*/
783: PetscErrorCode NEPSetConvergenceTest(NEP nep,NEPConv conv)784: {
788: switch (conv) {
789: case NEP_CONV_ABS: nep->converged = NEPConvergedAbsolute; break;
790: case NEP_CONV_REL: nep->converged = NEPConvergedRelative; break;
791: case NEP_CONV_NORM: nep->converged = NEPConvergedNorm; break;
792: case NEP_CONV_USER:
793: if (!nep->convergeduser) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"Must call NEPSetConvergenceTestFunction() first");
794: nep->converged = nep->convergeduser;
795: break;
796: default:797: SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
798: }
799: nep->conv = conv;
800: return(0);
801: }
803: /*@
804: NEPGetConvergenceTest - Gets the method used to compute the error estimate
805: used in the convergence test.
807: Not Collective
809: Input Parameters:
810: . nep - nonlinear eigensolver context obtained from NEPCreate()
812: Output Parameters:
813: . conv - the type of convergence test
815: Level: intermediate
817: .seealso: NEPSetConvergenceTest(), NEPConv818: @*/
819: PetscErrorCode NEPGetConvergenceTest(NEP nep,NEPConv *conv)820: {
824: *conv = nep->conv;
825: return(0);
826: }
828: /*@C
829: NEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
830: iteration of the eigensolver.
832: Logically Collective on NEP834: Input Parameters:
835: + nep - nonlinear eigensolver context obtained from NEPCreate()
836: . func - pointer to the stopping test function
837: . ctx - context for private data for the stopping routine (may be null)
838: - destroy - a routine for destroying the context (may be null)
840: Calling Sequence of func:
841: $ func(NEP nep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,NEPConvergedReason *reason,void *ctx)
843: + nep - nonlinear eigensolver context obtained from NEPCreate()
844: . its - current number of iterations
845: . max_it - maximum number of iterations
846: . nconv - number of currently converged eigenpairs
847: . nev - number of requested eigenpairs
848: . reason - (output) result of the stopping test
849: - ctx - optional context, as set by NEPSetStoppingTestFunction()
851: Note:
852: Normal usage is to first call the default routine NEPStoppingBasic() and then
853: set reason to NEP_CONVERGED_USER if some user-defined conditions have been
854: met. To let the eigensolver continue iterating, the result must be left as
855: NEP_CONVERGED_ITERATING.
857: Level: advanced
859: .seealso: NEPSetStoppingTest(), NEPStoppingBasic()
860: @*/
861: PetscErrorCode NEPSetStoppingTestFunction(NEP nep,PetscErrorCode (*func)(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*),void* ctx,PetscErrorCode (*destroy)(void*))862: {
867: if (nep->stoppingdestroy) {
868: (*nep->stoppingdestroy)(nep->stoppingctx);
869: }
870: nep->stoppinguser = func;
871: nep->stoppingdestroy = destroy;
872: nep->stoppingctx = ctx;
873: if (func == NEPStoppingBasic) nep->stop = NEP_STOP_BASIC;
874: else {
875: nep->stop = NEP_STOP_USER;
876: nep->stopping = nep->stoppinguser;
877: }
878: return(0);
879: }
881: /*@
882: NEPSetStoppingTest - Specifies how to decide the termination of the outer
883: loop of the eigensolver.
885: Logically Collective on NEP887: Input Parameters:
888: + nep - nonlinear eigensolver context obtained from NEPCreate()
889: - stop - the type of stopping test
891: Options Database Keys:
892: + -nep_stop_basic - Sets the default stopping test
893: - -nep_stop_user - Selects the user-defined stopping test
895: Note:
896: The parameter 'stop' can have one of these values
897: + NEP_STOP_BASIC - default stopping test
898: - NEP_STOP_USER - function set by NEPSetStoppingTestFunction()
900: Level: advanced
902: .seealso: NEPGetStoppingTest(), NEPSetStoppingTestFunction(), NEPSetConvergenceTest(), NEPStop903: @*/
904: PetscErrorCode NEPSetStoppingTest(NEP nep,NEPStop stop)905: {
909: switch (stop) {
910: case NEP_STOP_BASIC: nep->stopping = NEPStoppingBasic; break;
911: case NEP_STOP_USER:
912: if (!nep->stoppinguser) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ORDER,"Must call NEPSetStoppingTestFunction() first");
913: nep->stopping = nep->stoppinguser;
914: break;
915: default:916: SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
917: }
918: nep->stop = stop;
919: return(0);
920: }
922: /*@
923: NEPGetStoppingTest - Gets the method used to decide the termination of the outer
924: loop of the eigensolver.
926: Not Collective
928: Input Parameters:
929: . nep - nonlinear eigensolver context obtained from NEPCreate()
931: Output Parameters:
932: . stop - the type of stopping test
934: Level: advanced
936: .seealso: NEPSetStoppingTest(), NEPStop937: @*/
938: PetscErrorCode NEPGetStoppingTest(NEP nep,NEPStop *stop)939: {
943: *stop = nep->stop;
944: return(0);
945: }
947: /*@
948: NEPSetTrackAll - Specifies if the solver must compute the residual of all
949: approximate eigenpairs or not.
951: Logically Collective on NEP953: Input Parameters:
954: + nep - the eigensolver context
955: - trackall - whether compute all residuals or not
957: Notes:
958: If the user sets trackall=PETSC_TRUE then the solver explicitly computes
959: the residual for each eigenpair approximation. Computing the residual is
960: usually an expensive operation and solvers commonly compute the associated
961: residual to the first unconverged eigenpair.
962: The options '-nep_monitor_all' and '-nep_monitor_lg_all' automatically
963: activate this option.
965: Level: developer
967: .seealso: NEPGetTrackAll()
968: @*/
969: PetscErrorCode NEPSetTrackAll(NEP nep,PetscBool trackall)970: {
974: nep->trackall = trackall;
975: return(0);
976: }
978: /*@
979: NEPGetTrackAll - Returns the flag indicating whether all residual norms must
980: be computed or not.
982: Not Collective
984: Input Parameter:
985: . nep - the eigensolver context
987: Output Parameter:
988: . trackall - the returned flag
990: Level: developer
992: .seealso: NEPSetTrackAll()
993: @*/
994: PetscErrorCode NEPGetTrackAll(NEP nep,PetscBool *trackall)995: {
999: *trackall = nep->trackall;
1000: return(0);
1001: }
1003: /*@
1004: NEPSetRefine - Specifies the refinement type (and options) to be used
1005: after the solve.
1007: Logically Collective on NEP1009: Input Parameters:
1010: + nep - the nonlinear eigensolver context
1011: . refine - refinement type
1012: . npart - number of partitions of the communicator
1013: . tol - the convergence tolerance
1014: . its - maximum number of refinement iterations
1015: - scheme - which scheme to be used for solving the involved linear systems
1017: Options Database Keys:
1018: + -nep_refine <type> - refinement type, one of <none,simple,multiple>
1019: . -nep_refine_partitions <n> - the number of partitions
1020: . -nep_refine_tol <tol> - the tolerance
1021: . -nep_refine_its <its> - number of iterations
1022: - -nep_refine_scheme - to set the scheme for the linear solves
1024: Notes:
1025: By default, iterative refinement is disabled, since it may be very
1026: costly. There are two possible refinement strategies: simple and multiple.
1027: The simple approach performs iterative refinement on each of the
1028: converged eigenpairs individually, whereas the multiple strategy works
1029: with the invariant pair as a whole, refining all eigenpairs simultaneously.
1030: The latter may be required for the case of multiple eigenvalues.
1032: In some cases, especially when using direct solvers within the
1033: iterative refinement method, it may be helpful for improved scalability
1034: to split the communicator in several partitions. The npart parameter
1035: indicates how many partitions to use (defaults to 1).
1037: The tol and its parameters specify the stopping criterion. In the simple
1038: method, refinement continues until the residual of each eigenpair is
1039: below the tolerance (tol defaults to the NEP tol, but may be set to a
1040: different value). In contrast, the multiple method simply performs its
1041: refinement iterations (just one by default).
1043: The scheme argument is used to change the way in which linear systems are
1044: solved. Possible choices are: explicit, mixed block elimination (MBE),
1045: and Schur complement.
1047: Level: intermediate
1049: .seealso: NEPGetRefine()
1050: @*/
1051: PetscErrorCode NEPSetRefine(NEP nep,NEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,NEPRefineScheme scheme)1052: {
1054: PetscMPIInt size;
1063: nep->refine = refine;
1064: if (refine) { /* process parameters only if not REFINE_NONE */
1065: if (npart!=nep->npart) {
1066: PetscSubcommDestroy(&nep->refinesubc);
1067: KSPDestroy(&nep->refineksp);
1068: }
1069: if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
1070: nep->npart = 1;
1071: } else {
1072: MPI_Comm_size(PetscObjectComm((PetscObject)nep),&size);
1073: if (npart<1 || npart>size) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of npart");
1074: nep->npart = npart;
1075: }
1076: if (tol == PETSC_DEFAULT || tol == PETSC_DECIDE) {
1077: nep->rtol = PETSC_DEFAULT;
1078: } else {
1079: if (tol<=0.0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
1080: nep->rtol = tol;
1081: }
1082: if (its==PETSC_DECIDE || its==PETSC_DEFAULT) {
1083: nep->rits = PETSC_DEFAULT;
1084: } else {
1085: if (its<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be >= 0");
1086: nep->rits = its;
1087: }
1088: nep->scheme = scheme;
1089: }
1090: nep->state = NEP_STATE_INITIAL;
1091: return(0);
1092: }
1094: /*@C
1095: NEPGetRefine - Gets the refinement strategy used by the NEP object, and the
1096: associated parameters.
1098: Not Collective
1100: Input Parameter:
1101: . nep - the nonlinear eigensolver context
1103: Output Parameters:
1104: + refine - refinement type
1105: . npart - number of partitions of the communicator
1106: . tol - the convergence tolerance
1107: - its - maximum number of refinement iterations
1108: - scheme - the scheme used for solving linear systems
1110: Level: intermediate
1112: Note:
1113: The user can specify NULL for any parameter that is not needed.
1115: .seealso: NEPSetRefine()
1116: @*/
1117: PetscErrorCode NEPGetRefine(NEP nep,NEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,NEPRefineScheme *scheme)1118: {
1121: if (refine) *refine = nep->refine;
1122: if (npart) *npart = nep->npart;
1123: if (tol) *tol = nep->rtol;
1124: if (its) *its = nep->rits;
1125: if (scheme) *scheme = nep->scheme;
1126: return(0);
1127: }
1129: /*@C
1130: NEPSetOptionsPrefix - Sets the prefix used for searching for all
1131: NEP options in the database.
1133: Logically Collective on NEP1135: Input Parameters:
1136: + nep - the nonlinear eigensolver context
1137: - prefix - the prefix string to prepend to all NEP option requests
1139: Notes:
1140: A hyphen (-) must NOT be given at the beginning of the prefix name.
1141: The first character of all runtime options is AUTOMATICALLY the
1142: hyphen.
1144: For example, to distinguish between the runtime options for two
1145: different NEP contexts, one could call
1146: .vb
1147: NEPSetOptionsPrefix(nep1,"neig1_")
1148: NEPSetOptionsPrefix(nep2,"neig2_")
1149: .ve
1151: Level: advanced
1153: .seealso: NEPAppendOptionsPrefix(), NEPGetOptionsPrefix()
1154: @*/
1155: PetscErrorCode NEPSetOptionsPrefix(NEP nep,const char *prefix)1156: {
1161: if (!nep->V) { NEPGetBV(nep,&nep->V); }
1162: BVSetOptionsPrefix(nep->V,prefix);
1163: if (!nep->ds) { NEPGetDS(nep,&nep->ds); }
1164: DSSetOptionsPrefix(nep->ds,prefix);
1165: if (!nep->rg) { NEPGetRG(nep,&nep->rg); }
1166: RGSetOptionsPrefix(nep->rg,prefix);
1167: PetscObjectSetOptionsPrefix((PetscObject)nep,prefix);
1168: return(0);
1169: }
1171: /*@C
1172: NEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1173: NEP options in the database.
1175: Logically Collective on NEP1177: Input Parameters:
1178: + nep - the nonlinear eigensolver context
1179: - prefix - the prefix string to prepend to all NEP option requests
1181: Notes:
1182: A hyphen (-) must NOT be given at the beginning of the prefix name.
1183: The first character of all runtime options is AUTOMATICALLY the hyphen.
1185: Level: advanced
1187: .seealso: NEPSetOptionsPrefix(), NEPGetOptionsPrefix()
1188: @*/
1189: PetscErrorCode NEPAppendOptionsPrefix(NEP nep,const char *prefix)1190: {
1195: if (!nep->V) { NEPGetBV(nep,&nep->V); }
1196: BVAppendOptionsPrefix(nep->V,prefix);
1197: if (!nep->ds) { NEPGetDS(nep,&nep->ds); }
1198: DSAppendOptionsPrefix(nep->ds,prefix);
1199: if (!nep->rg) { NEPGetRG(nep,&nep->rg); }
1200: RGAppendOptionsPrefix(nep->rg,prefix);
1201: PetscObjectAppendOptionsPrefix((PetscObject)nep,prefix);
1202: return(0);
1203: }
1205: /*@C
1206: NEPGetOptionsPrefix - Gets the prefix used for searching for all
1207: NEP options in the database.
1209: Not Collective
1211: Input Parameters:
1212: . nep - the nonlinear eigensolver context
1214: Output Parameters:
1215: . prefix - pointer to the prefix string used is returned
1217: Note:
1218: On the Fortran side, the user should pass in a string 'prefix' of
1219: sufficient length to hold the prefix.
1221: Level: advanced
1223: .seealso: NEPSetOptionsPrefix(), NEPAppendOptionsPrefix()
1224: @*/
1225: PetscErrorCode NEPGetOptionsPrefix(NEP nep,const char *prefix[])1226: {
1232: PetscObjectGetOptionsPrefix((PetscObject)nep,prefix);
1233: return(0);
1234: }