Actual source code: nepopts.c
slepc-3.18.3 2023-03-24
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 on nep
24: Input Parameters:
25: + nep - the nonlinear eigensolver context
26: . opt - the command line option for this monitor
27: . name - the monitor type one is seeking
28: . ctx - an optional user context for the monitor, or NULL
29: - trackall - whether this monitor tracks all eigenvalues or not
31: Level: developer
33: .seealso: NEPMonitorSet(), NEPSetTrackAll()
34: @*/
35: PetscErrorCode NEPMonitorSetFromOptions(NEP nep,const char opt[],const char name[],void *ctx,PetscBool trackall)
36: {
37: PetscErrorCode (*mfunc)(NEP,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*);
38: PetscErrorCode (*cfunc)(PetscViewer,PetscViewerFormat,void*,PetscViewerAndFormat**);
39: PetscErrorCode (*dfunc)(PetscViewerAndFormat**);
40: PetscViewerAndFormat *vf;
41: PetscViewer viewer;
42: PetscViewerFormat format;
43: PetscViewerType vtype;
44: char key[PETSC_MAX_PATH_LEN];
45: PetscBool flg;
47: PetscOptionsGetViewer(PetscObjectComm((PetscObject)nep),((PetscObject)nep)->options,((PetscObject)nep)->prefix,opt,&viewer,&format,&flg);
48: if (!flg) return 0;
50: PetscViewerGetType(viewer,&vtype);
51: SlepcMonitorMakeKey_Internal(name,vtype,format,key);
52: PetscFunctionListFind(NEPMonitorList,key,&mfunc);
54: PetscFunctionListFind(NEPMonitorCreateList,key,&cfunc);
55: PetscFunctionListFind(NEPMonitorDestroyList,key,&dfunc);
56: if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
57: if (!dfunc) dfunc = PetscViewerAndFormatDestroy;
59: (*cfunc)(viewer,format,ctx,&vf);
60: PetscObjectDereference((PetscObject)viewer);
61: NEPMonitorSet(nep,mfunc,vf,(PetscErrorCode(*)(void **))dfunc);
62: if (trackall) NEPSetTrackAll(nep,PETSC_TRUE);
63: return 0;
64: }
66: /*@
67: NEPSetFromOptions - Sets NEP options from the options database.
68: This routine must be called before NEPSetUp() if the user is to be
69: allowed to set the solver type.
71: Collective on nep
73: Input Parameters:
74: . nep - the nonlinear eigensolver context
76: Notes:
77: To see all options, run your program with the -help option.
79: Level: beginner
81: .seealso: NEPSetOptionsPrefix()
82: @*/
83: PetscErrorCode NEPSetFromOptions(NEP nep)
84: {
85: char type[256];
86: PetscBool set,flg,flg1,flg2,flg3,flg4,flg5,bval;
87: PetscReal r;
88: PetscScalar s;
89: PetscInt i,j,k;
90: NEPRefine refine;
91: NEPRefineScheme scheme;
94: NEPRegisterAll();
95: PetscObjectOptionsBegin((PetscObject)nep);
96: PetscOptionsFList("-nep_type","Nonlinear eigensolver method","NEPSetType",NEPList,(char*)(((PetscObject)nep)->type_name?((PetscObject)nep)->type_name:NEPRII),type,sizeof(type),&flg);
97: if (flg) NEPSetType(nep,type);
98: else if (!((PetscObject)nep)->type_name) NEPSetType(nep,NEPRII);
100: PetscOptionsBoolGroupBegin("-nep_general","General nonlinear eigenvalue problem","NEPSetProblemType",&flg);
101: if (flg) NEPSetProblemType(nep,NEP_GENERAL);
102: PetscOptionsBoolGroupEnd("-nep_rational","Rational eigenvalue problem","NEPSetProblemType",&flg);
103: if (flg) NEPSetProblemType(nep,NEP_RATIONAL);
105: refine = nep->refine;
106: PetscOptionsEnum("-nep_refine","Iterative refinement method","NEPSetRefine",NEPRefineTypes,(PetscEnum)refine,(PetscEnum*)&refine,&flg1);
107: i = nep->npart;
108: PetscOptionsInt("-nep_refine_partitions","Number of partitions of the communicator for iterative refinement","NEPSetRefine",nep->npart,&i,&flg2);
109: r = nep->rtol;
110: PetscOptionsReal("-nep_refine_tol","Tolerance for iterative refinement","NEPSetRefine",nep->rtol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL/1000:nep->rtol,&r,&flg3);
111: j = nep->rits;
112: PetscOptionsInt("-nep_refine_its","Maximum number of iterations for iterative refinement","NEPSetRefine",nep->rits,&j,&flg4);
113: scheme = nep->scheme;
114: PetscOptionsEnum("-nep_refine_scheme","Scheme used for linear systems within iterative refinement","NEPSetRefine",NEPRefineSchemes,(PetscEnum)scheme,(PetscEnum*)&scheme,&flg5);
115: if (flg1 || flg2 || flg3 || flg4 || flg5) NEPSetRefine(nep,refine,i,r,j,scheme);
117: i = nep->max_it;
118: PetscOptionsInt("-nep_max_it","Maximum number of iterations","NEPSetTolerances",nep->max_it,&i,&flg1);
119: r = nep->tol;
120: PetscOptionsReal("-nep_tol","Tolerance","NEPSetTolerances",SlepcDefaultTol(nep->tol),&r,&flg2);
121: if (flg1 || flg2) NEPSetTolerances(nep,r,i);
123: PetscOptionsBoolGroupBegin("-nep_conv_rel","Relative error convergence test","NEPSetConvergenceTest",&flg);
124: if (flg) NEPSetConvergenceTest(nep,NEP_CONV_REL);
125: PetscOptionsBoolGroup("-nep_conv_norm","Convergence test relative to the matrix norms","NEPSetConvergenceTest",&flg);
126: if (flg) NEPSetConvergenceTest(nep,NEP_CONV_NORM);
127: PetscOptionsBoolGroup("-nep_conv_abs","Absolute error convergence test","NEPSetConvergenceTest",&flg);
128: if (flg) NEPSetConvergenceTest(nep,NEP_CONV_ABS);
129: PetscOptionsBoolGroupEnd("-nep_conv_user","User-defined convergence test","NEPSetConvergenceTest",&flg);
130: if (flg) NEPSetConvergenceTest(nep,NEP_CONV_USER);
132: PetscOptionsBoolGroupBegin("-nep_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","NEPSetStoppingTest",&flg);
133: if (flg) NEPSetStoppingTest(nep,NEP_STOP_BASIC);
134: PetscOptionsBoolGroupEnd("-nep_stop_user","User-defined stopping test","NEPSetStoppingTest",&flg);
135: if (flg) NEPSetStoppingTest(nep,NEP_STOP_USER);
137: i = nep->nev;
138: PetscOptionsInt("-nep_nev","Number of eigenvalues to compute","NEPSetDimensions",nep->nev,&i,&flg1);
139: j = nep->ncv;
140: PetscOptionsInt("-nep_ncv","Number of basis vectors","NEPSetDimensions",nep->ncv,&j,&flg2);
141: k = nep->mpd;
142: PetscOptionsInt("-nep_mpd","Maximum dimension of projected problem","NEPSetDimensions",nep->mpd,&k,&flg3);
143: if (flg1 || flg2 || flg3) NEPSetDimensions(nep,i,j,k);
145: PetscOptionsBoolGroupBegin("-nep_largest_magnitude","Compute largest eigenvalues in magnitude","NEPSetWhichEigenpairs",&flg);
146: if (flg) NEPSetWhichEigenpairs(nep,NEP_LARGEST_MAGNITUDE);
147: PetscOptionsBoolGroup("-nep_smallest_magnitude","Compute smallest eigenvalues in magnitude","NEPSetWhichEigenpairs",&flg);
148: if (flg) NEPSetWhichEigenpairs(nep,NEP_SMALLEST_MAGNITUDE);
149: PetscOptionsBoolGroup("-nep_largest_real","Compute eigenvalues with largest real parts","NEPSetWhichEigenpairs",&flg);
150: if (flg) NEPSetWhichEigenpairs(nep,NEP_LARGEST_REAL);
151: PetscOptionsBoolGroup("-nep_smallest_real","Compute eigenvalues with smallest real parts","NEPSetWhichEigenpairs",&flg);
152: if (flg) NEPSetWhichEigenpairs(nep,NEP_SMALLEST_REAL);
153: PetscOptionsBoolGroup("-nep_largest_imaginary","Compute eigenvalues with largest imaginary parts","NEPSetWhichEigenpairs",&flg);
154: if (flg) NEPSetWhichEigenpairs(nep,NEP_LARGEST_IMAGINARY);
155: PetscOptionsBoolGroup("-nep_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","NEPSetWhichEigenpairs",&flg);
156: if (flg) NEPSetWhichEigenpairs(nep,NEP_SMALLEST_IMAGINARY);
157: PetscOptionsBoolGroup("-nep_target_magnitude","Compute eigenvalues closest to target","NEPSetWhichEigenpairs",&flg);
158: if (flg) NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE);
159: PetscOptionsBoolGroup("-nep_target_real","Compute eigenvalues with real parts closest to target","NEPSetWhichEigenpairs",&flg);
160: if (flg) NEPSetWhichEigenpairs(nep,NEP_TARGET_REAL);
161: PetscOptionsBoolGroup("-nep_target_imaginary","Compute eigenvalues with imaginary parts closest to target","NEPSetWhichEigenpairs",&flg);
162: if (flg) NEPSetWhichEigenpairs(nep,NEP_TARGET_IMAGINARY);
163: PetscOptionsBoolGroupEnd("-nep_all","Compute all eigenvalues in a region","NEPSetWhichEigenpairs",&flg);
164: if (flg) NEPSetWhichEigenpairs(nep,NEP_ALL);
166: PetscOptionsScalar("-nep_target","Value of the target","NEPSetTarget",nep->target,&s,&flg);
167: if (flg) {
168: if (nep->which!=NEP_TARGET_REAL && nep->which!=NEP_TARGET_IMAGINARY) NEPSetWhichEigenpairs(nep,NEP_TARGET_MAGNITUDE);
169: NEPSetTarget(nep,s);
170: }
172: PetscOptionsBool("-nep_two_sided","Use two-sided variant (to compute left eigenvectors)","NEPSetTwoSided",nep->twosided,&bval,&flg);
173: if (flg) NEPSetTwoSided(nep,bval);
175: /* -----------------------------------------------------------------------*/
176: /*
177: Cancels all monitors hardwired into code before call to NEPSetFromOptions()
178: */
179: PetscOptionsBool("-nep_monitor_cancel","Remove any hardwired monitor routines","NEPMonitorCancel",PETSC_FALSE,&flg,&set);
180: if (set && flg) NEPMonitorCancel(nep);
181: NEPMonitorSetFromOptions(nep,"-nep_monitor","first_approximation",NULL,PETSC_FALSE);
182: NEPMonitorSetFromOptions(nep,"-nep_monitor_all","all_approximations",NULL,PETSC_TRUE);
183: NEPMonitorSetFromOptions(nep,"-nep_monitor_conv","convergence_history",NULL,PETSC_FALSE);
185: /* -----------------------------------------------------------------------*/
186: PetscOptionsName("-nep_view","Print detailed information on solver used","NEPView",NULL);
187: PetscOptionsName("-nep_view_vectors","View computed eigenvectors","NEPVectorsView",NULL);
188: PetscOptionsName("-nep_view_values","View computed eigenvalues","NEPValuesView",NULL);
189: PetscOptionsName("-nep_converged_reason","Print reason for convergence, and number of iterations","NEPConvergedReasonView",NULL);
190: PetscOptionsName("-nep_error_absolute","Print absolute errors of each eigenpair","NEPErrorView",NULL);
191: PetscOptionsName("-nep_error_relative","Print relative errors of each eigenpair","NEPErrorView",NULL);
193: PetscTryTypeMethod(nep,setfromoptions,PetscOptionsObject);
194: PetscObjectProcessOptionsHandlers((PetscObject)nep,PetscOptionsObject);
195: PetscOptionsEnd();
197: if (!nep->V) NEPGetBV(nep,&nep->V);
198: BVSetFromOptions(nep->V);
199: if (!nep->rg) NEPGetRG(nep,&nep->rg);
200: RGSetFromOptions(nep->rg);
201: if (nep->useds) {
202: if (!nep->ds) NEPGetDS(nep,&nep->ds);
203: DSSetFromOptions(nep->ds);
204: }
205: if (!nep->refineksp) NEPRefineGetKSP(nep,&nep->refineksp);
206: KSPSetFromOptions(nep->refineksp);
207: if (nep->fui==NEP_USER_INTERFACE_SPLIT) for (i=0;i<nep->nt;i++) FNSetFromOptions(nep->f[i]);
208: return 0;
209: }
211: /*@C
212: NEPGetTolerances - Gets the tolerance and maximum iteration count used
213: by the NEP convergence tests.
215: Not Collective
217: Input Parameter:
218: . nep - the nonlinear eigensolver context
220: Output Parameters:
221: + tol - the convergence tolerance
222: - maxits - maximum number of iterations
224: Notes:
225: The user can specify NULL for any parameter that is not needed.
227: Level: intermediate
229: .seealso: NEPSetTolerances()
230: @*/
231: PetscErrorCode NEPGetTolerances(NEP nep,PetscReal *tol,PetscInt *maxits)
232: {
234: if (tol) *tol = nep->tol;
235: if (maxits) *maxits = nep->max_it;
236: return 0;
237: }
239: /*@
240: NEPSetTolerances - Sets the tolerance and maximum iteration count used
241: by the NEP convergence tests.
243: Logically Collective on nep
245: Input Parameters:
246: + nep - the nonlinear eigensolver context
247: . tol - the convergence tolerance
248: - maxits - maximum number of iterations to use
250: Options Database Keys:
251: + -nep_tol <tol> - Sets the convergence tolerance
252: - -nep_max_it <maxits> - Sets the maximum number of iterations allowed
254: Notes:
255: Use PETSC_DEFAULT for either argument to assign a reasonably good value.
257: Level: intermediate
259: .seealso: NEPGetTolerances()
260: @*/
261: PetscErrorCode NEPSetTolerances(NEP nep,PetscReal tol,PetscInt maxits)
262: {
266: if (tol == PETSC_DEFAULT) {
267: nep->tol = PETSC_DEFAULT;
268: nep->state = NEP_STATE_INITIAL;
269: } else {
271: nep->tol = tol;
272: }
273: if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
274: nep->max_it = PETSC_DEFAULT;
275: nep->state = NEP_STATE_INITIAL;
276: } else {
278: nep->max_it = maxits;
279: }
280: return 0;
281: }
283: /*@C
284: NEPGetDimensions - Gets the number of eigenvalues to compute
285: and the dimension of the subspace.
287: Not Collective
289: Input Parameter:
290: . nep - the nonlinear eigensolver context
292: Output Parameters:
293: + nev - number of eigenvalues to compute
294: . ncv - the maximum dimension of the subspace to be used by the solver
295: - mpd - the maximum dimension allowed for the projected problem
297: Notes:
298: The user can specify NULL for any parameter that is not needed.
300: Level: intermediate
302: .seealso: NEPSetDimensions()
303: @*/
304: PetscErrorCode NEPGetDimensions(NEP nep,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
305: {
307: if (nev) *nev = nep->nev;
308: if (ncv) *ncv = nep->ncv;
309: if (mpd) *mpd = nep->mpd;
310: return 0;
311: }
313: /*@
314: NEPSetDimensions - Sets the number of eigenvalues to compute
315: and the dimension of the subspace.
317: Logically Collective on nep
319: Input Parameters:
320: + nep - the nonlinear eigensolver context
321: . nev - number of eigenvalues to compute
322: . ncv - the maximum dimension of the subspace to be used by the solver
323: - mpd - the maximum dimension allowed for the projected problem
325: Options Database Keys:
326: + -nep_nev <nev> - Sets the number of eigenvalues
327: . -nep_ncv <ncv> - Sets the dimension of the subspace
328: - -nep_mpd <mpd> - Sets the maximum projected dimension
330: Notes:
331: Use PETSC_DEFAULT for ncv and mpd to assign a reasonably good value, which is
332: dependent on the solution method.
334: The parameters ncv and mpd are intimately related, so that the user is advised
335: to set one of them at most. Normal usage is that
336: (a) in cases where nev is small, the user sets ncv (a reasonable default is 2*nev); and
337: (b) in cases where nev is large, the user sets mpd.
339: The value of ncv should always be between nev and (nev+mpd), typically
340: ncv=nev+mpd. If nev is not too large, mpd=nev is a reasonable choice, otherwise
341: a smaller value should be used.
343: Level: intermediate
345: .seealso: NEPGetDimensions()
346: @*/
347: PetscErrorCode NEPSetDimensions(NEP nep,PetscInt nev,PetscInt ncv,PetscInt mpd)
348: {
354: nep->nev = nev;
355: if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
356: nep->ncv = PETSC_DEFAULT;
357: } else {
359: nep->ncv = ncv;
360: }
361: if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
362: nep->mpd = PETSC_DEFAULT;
363: } else {
365: nep->mpd = mpd;
366: }
367: nep->state = NEP_STATE_INITIAL;
368: return 0;
369: }
371: /*@
372: NEPSetWhichEigenpairs - Specifies which portion of the spectrum is
373: to be sought.
375: Logically Collective on nep
377: Input Parameters:
378: + nep - eigensolver context obtained from NEPCreate()
379: - which - the portion of the spectrum to be sought
381: Possible values:
382: The parameter 'which' can have one of these values
384: + NEP_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
385: . NEP_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
386: . NEP_LARGEST_REAL - largest real parts
387: . NEP_SMALLEST_REAL - smallest real parts
388: . NEP_LARGEST_IMAGINARY - largest imaginary parts
389: . NEP_SMALLEST_IMAGINARY - smallest imaginary parts
390: . NEP_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
391: . NEP_TARGET_REAL - eigenvalues with real part closest to target
392: . NEP_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
393: . NEP_ALL - all eigenvalues contained in a given region
394: - NEP_WHICH_USER - user defined ordering set with NEPSetEigenvalueComparison()
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: Not all eigensolvers implemented in NEP account for all the possible values
410: stated above. If SLEPc is compiled for real numbers NEP_LARGEST_IMAGINARY
411: and NEP_SMALLEST_IMAGINARY use the absolute value of the imaginary part
412: for eigenvalue selection.
414: The target is a scalar value provided with NEPSetTarget().
416: NEP_ALL is intended for use in the context of the CISS solver for
417: computing all eigenvalues in a region.
419: Level: intermediate
421: .seealso: NEPGetWhichEigenpairs(), NEPSetTarget(), NEPSetEigenvalueComparison(), NEPWhich
422: @*/
423: PetscErrorCode NEPSetWhichEigenpairs(NEP nep,NEPWhich which)
424: {
427: switch (which) {
428: case NEP_LARGEST_MAGNITUDE:
429: case NEP_SMALLEST_MAGNITUDE:
430: case NEP_LARGEST_REAL:
431: case NEP_SMALLEST_REAL:
432: case NEP_LARGEST_IMAGINARY:
433: case NEP_SMALLEST_IMAGINARY:
434: case NEP_TARGET_MAGNITUDE:
435: case NEP_TARGET_REAL:
436: #if defined(PETSC_USE_COMPLEX)
437: case NEP_TARGET_IMAGINARY:
438: #endif
439: case NEP_ALL:
440: case NEP_WHICH_USER:
441: if (nep->which != which) {
442: nep->state = NEP_STATE_INITIAL;
443: nep->which = which;
444: }
445: break;
446: #if !defined(PETSC_USE_COMPLEX)
447: case NEP_TARGET_IMAGINARY:
448: SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"NEP_TARGET_IMAGINARY can be used only with complex scalars");
449: #endif
450: default:
451: SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
452: }
453: return 0;
454: }
456: /*@
457: NEPGetWhichEigenpairs - Returns which portion of the spectrum is to be
458: sought.
460: Not Collective
462: Input Parameter:
463: . nep - eigensolver context obtained from NEPCreate()
465: Output Parameter:
466: . which - the portion of the spectrum to be sought
468: Notes:
469: See NEPSetWhichEigenpairs() for possible values of 'which'.
471: Level: intermediate
473: .seealso: NEPSetWhichEigenpairs(), NEPWhich
474: @*/
475: PetscErrorCode NEPGetWhichEigenpairs(NEP nep,NEPWhich *which)
476: {
479: *which = nep->which;
480: return 0;
481: }
483: /*@C
484: NEPSetEigenvalueComparison - Specifies the eigenvalue comparison function
485: when NEPSetWhichEigenpairs() is set to NEP_WHICH_USER.
487: Logically Collective on nep
489: Input Parameters:
490: + nep - eigensolver context obtained from NEPCreate()
491: . func - a pointer to the comparison function
492: - ctx - a context pointer (the last parameter to the comparison function)
494: Calling Sequence of func:
495: $ func(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)
497: + ar - real part of the 1st eigenvalue
498: . ai - imaginary part of the 1st eigenvalue
499: . br - real part of the 2nd eigenvalue
500: . bi - imaginary part of the 2nd eigenvalue
501: . res - result of comparison
502: - ctx - optional context, as set by NEPSetEigenvalueComparison()
504: Note:
505: The returning parameter 'res' can be
506: + negative - if the 1st eigenvalue is preferred to the 2st one
507: . zero - if both eigenvalues are equally preferred
508: - positive - if the 2st eigenvalue is preferred to the 1st one
510: Level: advanced
512: .seealso: NEPSetWhichEigenpairs(), NEPWhich
513: @*/
514: PetscErrorCode NEPSetEigenvalueComparison(NEP nep,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)
515: {
517: nep->sc->comparison = func;
518: nep->sc->comparisonctx = ctx;
519: nep->which = NEP_WHICH_USER;
520: return 0;
521: }
523: /*@
524: NEPSetProblemType - Specifies the type of the nonlinear eigenvalue problem.
526: Logically Collective on nep
528: Input Parameters:
529: + nep - the nonlinear eigensolver context
530: - type - a known type of nonlinear eigenvalue problem
532: Options Database Keys:
533: + -nep_general - general problem with no particular structure
534: - -nep_rational - a rational eigenvalue problem defined in split form with all f_i rational
536: Notes:
537: Allowed values for the problem type are general (NEP_GENERAL), and rational
538: (NEP_RATIONAL).
540: This function is used to provide a hint to the NEP solver to exploit certain
541: properties of the nonlinear eigenproblem. This hint may be used or not,
542: depending on the solver. By default, no particular structure is assumed.
544: Level: intermediate
546: .seealso: NEPSetType(), NEPGetProblemType(), NEPProblemType
547: @*/
548: PetscErrorCode NEPSetProblemType(NEP nep,NEPProblemType type)
549: {
553: if (type != nep->problem_type) {
554: nep->problem_type = type;
555: nep->state = NEP_STATE_INITIAL;
556: }
557: return 0;
558: }
560: /*@
561: NEPGetProblemType - Gets the problem type from the NEP object.
563: Not Collective
565: Input Parameter:
566: . nep - the nonlinear eigensolver context
568: Output Parameter:
569: . type - the problem type
571: Level: intermediate
573: .seealso: NEPSetProblemType(), NEPProblemType
574: @*/
575: PetscErrorCode NEPGetProblemType(NEP nep,NEPProblemType *type)
576: {
579: *type = nep->problem_type;
580: return 0;
581: }
583: /*@
584: NEPSetTwoSided - Sets the solver to use a two-sided variant so that left
585: eigenvectors are also computed.
587: Logically Collective on nep
589: Input Parameters:
590: + nep - the eigensolver context
591: - twosided - whether the two-sided variant is to be used or not
593: Options Database Keys:
594: . -nep_two_sided <boolean> - Sets/resets the twosided flag
596: Notes:
597: If the user sets twosided=PETSC_TRUE then the solver uses a variant of
598: the algorithm that computes both right and left eigenvectors. This is
599: usually much more costly. This option is not available in all solvers.
601: When using two-sided solvers, the problem matrices must have both the
602: MatMult and MatMultTranspose operations defined.
604: Level: advanced
606: .seealso: NEPGetTwoSided(), NEPGetLeftEigenvector()
607: @*/
608: PetscErrorCode NEPSetTwoSided(NEP nep,PetscBool twosided)
609: {
612: if (twosided!=nep->twosided) {
613: nep->twosided = twosided;
614: nep->state = NEP_STATE_INITIAL;
615: }
616: return 0;
617: }
619: /*@
620: NEPGetTwoSided - Returns the flag indicating whether a two-sided variant
621: of the algorithm is being used or not.
623: Not Collective
625: Input Parameter:
626: . nep - the eigensolver context
628: Output Parameter:
629: . twosided - the returned flag
631: Level: advanced
633: .seealso: NEPSetTwoSided()
634: @*/
635: PetscErrorCode NEPGetTwoSided(NEP nep,PetscBool *twosided)
636: {
639: *twosided = nep->twosided;
640: return 0;
641: }
643: /*@C
644: NEPSetConvergenceTestFunction - Sets a function to compute the error estimate
645: used in the convergence test.
647: Logically Collective on nep
649: Input Parameters:
650: + nep - nonlinear eigensolver context obtained from NEPCreate()
651: . func - a pointer to the convergence test function
652: . ctx - context for private data for the convergence routine (may be null)
653: - destroy - a routine for destroying the context (may be null)
655: Calling Sequence of func:
656: $ func(NEP nep,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
658: + nep - nonlinear eigensolver context obtained from NEPCreate()
659: . eigr - real part of the eigenvalue
660: . eigi - imaginary part of the eigenvalue
661: . res - residual norm associated to the eigenpair
662: . errest - (output) computed error estimate
663: - ctx - optional context, as set by NEPSetConvergenceTestFunction()
665: Note:
666: If the error estimate returned by the convergence test function is less than
667: the tolerance, then the eigenvalue is accepted as converged.
669: Level: advanced
671: .seealso: NEPSetConvergenceTest(), NEPSetTolerances()
672: @*/
673: PetscErrorCode NEPSetConvergenceTestFunction(NEP nep,PetscErrorCode (*func)(NEP,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
674: {
676: if (nep->convergeddestroy) (*nep->convergeddestroy)(nep->convergedctx);
677: nep->convergeduser = func;
678: nep->convergeddestroy = destroy;
679: nep->convergedctx = ctx;
680: if (func == NEPConvergedRelative) nep->conv = NEP_CONV_REL;
681: else if (func == NEPConvergedNorm) nep->conv = NEP_CONV_NORM;
682: else if (func == NEPConvergedAbsolute) nep->conv = NEP_CONV_ABS;
683: else {
684: nep->conv = NEP_CONV_USER;
685: nep->converged = nep->convergeduser;
686: }
687: return 0;
688: }
690: /*@
691: NEPSetConvergenceTest - Specifies how to compute the error estimate
692: used in the convergence test.
694: Logically Collective on nep
696: Input Parameters:
697: + nep - nonlinear eigensolver context obtained from NEPCreate()
698: - conv - the type of convergence test
700: Options Database Keys:
701: + -nep_conv_abs - Sets the absolute convergence test
702: . -nep_conv_rel - Sets the convergence test relative to the eigenvalue
703: - -nep_conv_user - Selects the user-defined convergence test
705: Note:
706: The parameter 'conv' can have one of these values
707: + NEP_CONV_ABS - absolute error ||r||
708: . NEP_CONV_REL - error relative to the eigenvalue l, ||r||/|l|
709: . NEP_CONV_NORM - error relative matrix norms, ||r||/sum_i(|f_i(l)|*||A_i||)
710: - NEP_CONV_USER - function set by NEPSetConvergenceTestFunction()
712: Level: intermediate
714: .seealso: NEPGetConvergenceTest(), NEPSetConvergenceTestFunction(), NEPSetStoppingTest(), NEPConv
715: @*/
716: PetscErrorCode NEPSetConvergenceTest(NEP nep,NEPConv conv)
717: {
720: switch (conv) {
721: case NEP_CONV_ABS: nep->converged = NEPConvergedAbsolute; break;
722: case NEP_CONV_REL: nep->converged = NEPConvergedRelative; break;
723: case NEP_CONV_NORM: nep->converged = NEPConvergedNorm; break;
724: case NEP_CONV_USER:
726: nep->converged = nep->convergeduser;
727: break;
728: default:
729: SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
730: }
731: nep->conv = conv;
732: return 0;
733: }
735: /*@
736: NEPGetConvergenceTest - Gets the method used to compute the error estimate
737: used in the convergence test.
739: Not Collective
741: Input Parameters:
742: . nep - nonlinear eigensolver context obtained from NEPCreate()
744: Output Parameters:
745: . conv - the type of convergence test
747: Level: intermediate
749: .seealso: NEPSetConvergenceTest(), NEPConv
750: @*/
751: PetscErrorCode NEPGetConvergenceTest(NEP nep,NEPConv *conv)
752: {
755: *conv = nep->conv;
756: return 0;
757: }
759: /*@C
760: NEPSetStoppingTestFunction - Sets a function to decide when to stop the outer
761: iteration of the eigensolver.
763: Logically Collective on nep
765: Input Parameters:
766: + nep - nonlinear eigensolver context obtained from NEPCreate()
767: . func - pointer to the stopping test function
768: . ctx - context for private data for the stopping routine (may be null)
769: - destroy - a routine for destroying the context (may be null)
771: Calling Sequence of func:
772: $ func(NEP nep,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,NEPConvergedReason *reason,void *ctx)
774: + nep - nonlinear eigensolver context obtained from NEPCreate()
775: . its - current number of iterations
776: . max_it - maximum number of iterations
777: . nconv - number of currently converged eigenpairs
778: . nev - number of requested eigenpairs
779: . reason - (output) result of the stopping test
780: - ctx - optional context, as set by NEPSetStoppingTestFunction()
782: Note:
783: Normal usage is to first call the default routine NEPStoppingBasic() and then
784: set reason to NEP_CONVERGED_USER if some user-defined conditions have been
785: met. To let the eigensolver continue iterating, the result must be left as
786: NEP_CONVERGED_ITERATING.
788: Level: advanced
790: .seealso: NEPSetStoppingTest(), NEPStoppingBasic()
791: @*/
792: PetscErrorCode NEPSetStoppingTestFunction(NEP nep,PetscErrorCode (*func)(NEP,PetscInt,PetscInt,PetscInt,PetscInt,NEPConvergedReason*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
793: {
795: if (nep->stoppingdestroy) (*nep->stoppingdestroy)(nep->stoppingctx);
796: nep->stoppinguser = func;
797: nep->stoppingdestroy = destroy;
798: nep->stoppingctx = ctx;
799: if (func == NEPStoppingBasic) nep->stop = NEP_STOP_BASIC;
800: else {
801: nep->stop = NEP_STOP_USER;
802: nep->stopping = nep->stoppinguser;
803: }
804: return 0;
805: }
807: /*@
808: NEPSetStoppingTest - Specifies how to decide the termination of the outer
809: loop of the eigensolver.
811: Logically Collective on nep
813: Input Parameters:
814: + nep - nonlinear eigensolver context obtained from NEPCreate()
815: - stop - the type of stopping test
817: Options Database Keys:
818: + -nep_stop_basic - Sets the default stopping test
819: - -nep_stop_user - Selects the user-defined stopping test
821: Note:
822: The parameter 'stop' can have one of these values
823: + NEP_STOP_BASIC - default stopping test
824: - NEP_STOP_USER - function set by NEPSetStoppingTestFunction()
826: Level: advanced
828: .seealso: NEPGetStoppingTest(), NEPSetStoppingTestFunction(), NEPSetConvergenceTest(), NEPStop
829: @*/
830: PetscErrorCode NEPSetStoppingTest(NEP nep,NEPStop stop)
831: {
834: switch (stop) {
835: case NEP_STOP_BASIC: nep->stopping = NEPStoppingBasic; break;
836: case NEP_STOP_USER:
838: nep->stopping = nep->stoppinguser;
839: break;
840: default:
841: SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
842: }
843: nep->stop = stop;
844: return 0;
845: }
847: /*@
848: NEPGetStoppingTest - Gets the method used to decide the termination of the outer
849: loop of the eigensolver.
851: Not Collective
853: Input Parameters:
854: . nep - nonlinear eigensolver context obtained from NEPCreate()
856: Output Parameters:
857: . stop - the type of stopping test
859: Level: advanced
861: .seealso: NEPSetStoppingTest(), NEPStop
862: @*/
863: PetscErrorCode NEPGetStoppingTest(NEP nep,NEPStop *stop)
864: {
867: *stop = nep->stop;
868: return 0;
869: }
871: /*@
872: NEPSetTrackAll - Specifies if the solver must compute the residual of all
873: approximate eigenpairs or not.
875: Logically Collective on nep
877: Input Parameters:
878: + nep - the eigensolver context
879: - trackall - whether compute all residuals or not
881: Notes:
882: If the user sets trackall=PETSC_TRUE then the solver explicitly computes
883: the residual for each eigenpair approximation. Computing the residual is
884: usually an expensive operation and solvers commonly compute the associated
885: residual to the first unconverged eigenpair.
887: The option '-nep_monitor_all' automatically activates this option.
889: Level: developer
891: .seealso: NEPGetTrackAll()
892: @*/
893: PetscErrorCode NEPSetTrackAll(NEP nep,PetscBool trackall)
894: {
897: nep->trackall = trackall;
898: return 0;
899: }
901: /*@
902: NEPGetTrackAll - Returns the flag indicating whether all residual norms must
903: be computed or not.
905: Not Collective
907: Input Parameter:
908: . nep - the eigensolver context
910: Output Parameter:
911: . trackall - the returned flag
913: Level: developer
915: .seealso: NEPSetTrackAll()
916: @*/
917: PetscErrorCode NEPGetTrackAll(NEP nep,PetscBool *trackall)
918: {
921: *trackall = nep->trackall;
922: return 0;
923: }
925: /*@
926: NEPSetRefine - Specifies the refinement type (and options) to be used
927: after the solve.
929: Logically Collective on nep
931: Input Parameters:
932: + nep - the nonlinear eigensolver context
933: . refine - refinement type
934: . npart - number of partitions of the communicator
935: . tol - the convergence tolerance
936: . its - maximum number of refinement iterations
937: - scheme - which scheme to be used for solving the involved linear systems
939: Options Database Keys:
940: + -nep_refine <type> - refinement type, one of <none,simple,multiple>
941: . -nep_refine_partitions <n> - the number of partitions
942: . -nep_refine_tol <tol> - the tolerance
943: . -nep_refine_its <its> - number of iterations
944: - -nep_refine_scheme - to set the scheme for the linear solves
946: Notes:
947: By default, iterative refinement is disabled, since it may be very
948: costly. There are two possible refinement strategies, simple and multiple.
949: The simple approach performs iterative refinement on each of the
950: converged eigenpairs individually, whereas the multiple strategy works
951: with the invariant pair as a whole, refining all eigenpairs simultaneously.
952: The latter may be required for the case of multiple eigenvalues.
954: In some cases, especially when using direct solvers within the
955: iterative refinement method, it may be helpful for improved scalability
956: to split the communicator in several partitions. The npart parameter
957: indicates how many partitions to use (defaults to 1).
959: The tol and its parameters specify the stopping criterion. In the simple
960: method, refinement continues until the residual of each eigenpair is
961: below the tolerance (tol defaults to the NEP tol, but may be set to a
962: different value). In contrast, the multiple method simply performs its
963: refinement iterations (just one by default).
965: The scheme argument is used to change the way in which linear systems are
966: solved. Possible choices are explicit, mixed block elimination (MBE),
967: and Schur complement.
969: Level: intermediate
971: .seealso: NEPGetRefine()
972: @*/
973: PetscErrorCode NEPSetRefine(NEP nep,NEPRefine refine,PetscInt npart,PetscReal tol,PetscInt its,NEPRefineScheme scheme)
974: {
975: PetscMPIInt size;
983: nep->refine = refine;
984: if (refine) { /* process parameters only if not REFINE_NONE */
985: if (npart!=nep->npart) {
986: PetscSubcommDestroy(&nep->refinesubc);
987: KSPDestroy(&nep->refineksp);
988: }
989: if (npart == PETSC_DEFAULT || npart == PETSC_DECIDE) {
990: nep->npart = 1;
991: } else {
992: MPI_Comm_size(PetscObjectComm((PetscObject)nep),&size);
994: nep->npart = npart;
995: }
996: if (tol == PETSC_DEFAULT || tol == PETSC_DECIDE) {
997: nep->rtol = PETSC_DEFAULT;
998: } else {
1000: nep->rtol = tol;
1001: }
1002: if (its==PETSC_DECIDE || its==PETSC_DEFAULT) {
1003: nep->rits = PETSC_DEFAULT;
1004: } else {
1006: nep->rits = its;
1007: }
1008: nep->scheme = scheme;
1009: }
1010: nep->state = NEP_STATE_INITIAL;
1011: return 0;
1012: }
1014: /*@C
1015: NEPGetRefine - Gets the refinement strategy used by the NEP object, and the
1016: associated parameters.
1018: Not Collective
1020: Input Parameter:
1021: . nep - the nonlinear eigensolver context
1023: Output Parameters:
1024: + refine - refinement type
1025: . npart - number of partitions of the communicator
1026: . tol - the convergence tolerance
1027: . its - maximum number of refinement iterations
1028: - scheme - the scheme used for solving linear systems
1030: Level: intermediate
1032: Note:
1033: The user can specify NULL for any parameter that is not needed.
1035: .seealso: NEPSetRefine()
1036: @*/
1037: PetscErrorCode NEPGetRefine(NEP nep,NEPRefine *refine,PetscInt *npart,PetscReal *tol,PetscInt *its,NEPRefineScheme *scheme)
1038: {
1040: if (refine) *refine = nep->refine;
1041: if (npart) *npart = nep->npart;
1042: if (tol) *tol = nep->rtol;
1043: if (its) *its = nep->rits;
1044: if (scheme) *scheme = nep->scheme;
1045: return 0;
1046: }
1048: /*@C
1049: NEPSetOptionsPrefix - Sets the prefix used for searching for all
1050: NEP options in the database.
1052: Logically Collective on nep
1054: Input Parameters:
1055: + nep - the nonlinear eigensolver context
1056: - prefix - the prefix string to prepend to all NEP option requests
1058: Notes:
1059: A hyphen (-) must NOT be given at the beginning of the prefix name.
1060: The first character of all runtime options is AUTOMATICALLY the
1061: hyphen.
1063: For example, to distinguish between the runtime options for two
1064: different NEP contexts, one could call
1065: .vb
1066: NEPSetOptionsPrefix(nep1,"neig1_")
1067: NEPSetOptionsPrefix(nep2,"neig2_")
1068: .ve
1070: Level: advanced
1072: .seealso: NEPAppendOptionsPrefix(), NEPGetOptionsPrefix()
1073: @*/
1074: PetscErrorCode NEPSetOptionsPrefix(NEP nep,const char *prefix)
1075: {
1077: if (!nep->V) NEPGetBV(nep,&nep->V);
1078: BVSetOptionsPrefix(nep->V,prefix);
1079: if (!nep->ds) NEPGetDS(nep,&nep->ds);
1080: DSSetOptionsPrefix(nep->ds,prefix);
1081: if (!nep->rg) NEPGetRG(nep,&nep->rg);
1082: RGSetOptionsPrefix(nep->rg,prefix);
1083: PetscObjectSetOptionsPrefix((PetscObject)nep,prefix);
1084: return 0;
1085: }
1087: /*@C
1088: NEPAppendOptionsPrefix - Appends to the prefix used for searching for all
1089: NEP options in the database.
1091: Logically Collective on nep
1093: Input Parameters:
1094: + nep - the nonlinear eigensolver context
1095: - prefix - the prefix string to prepend to all NEP option requests
1097: Notes:
1098: A hyphen (-) must NOT be given at the beginning of the prefix name.
1099: The first character of all runtime options is AUTOMATICALLY the hyphen.
1101: Level: advanced
1103: .seealso: NEPSetOptionsPrefix(), NEPGetOptionsPrefix()
1104: @*/
1105: PetscErrorCode NEPAppendOptionsPrefix(NEP nep,const char *prefix)
1106: {
1108: if (!nep->V) NEPGetBV(nep,&nep->V);
1109: BVAppendOptionsPrefix(nep->V,prefix);
1110: if (!nep->ds) NEPGetDS(nep,&nep->ds);
1111: DSAppendOptionsPrefix(nep->ds,prefix);
1112: if (!nep->rg) NEPGetRG(nep,&nep->rg);
1113: RGAppendOptionsPrefix(nep->rg,prefix);
1114: PetscObjectAppendOptionsPrefix((PetscObject)nep,prefix);
1115: return 0;
1116: }
1118: /*@C
1119: NEPGetOptionsPrefix - Gets the prefix used for searching for all
1120: NEP options in the database.
1122: Not Collective
1124: Input Parameters:
1125: . nep - the nonlinear eigensolver context
1127: Output Parameters:
1128: . prefix - pointer to the prefix string used is returned
1130: Note:
1131: On the Fortran side, the user should pass in a string 'prefix' of
1132: sufficient length to hold the prefix.
1134: Level: advanced
1136: .seealso: NEPSetOptionsPrefix(), NEPAppendOptionsPrefix()
1137: @*/
1138: PetscErrorCode NEPGetOptionsPrefix(NEP nep,const char *prefix[])
1139: {
1142: PetscObjectGetOptionsPrefix((PetscObject)nep,prefix);
1143: return 0;
1144: }