Actual source code: epsopts.c
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: EPS routines related to options that can be set via the command-line
12: or procedurally.
13: */
15: #include <slepc/private/epsimpl.h>
16: #include <petscdraw.h>
18: /*@C
19: EPSMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
20: indicated by the user.
22: Collective
24: Input Parameters:
25: + eps - the linear eigensolver context
26: . opt - the command line option for this monitor
27: . name - the monitor type one is seeking
28: . ctx - an optional user context for the monitor, or `NULL`
29: - trackall - whether this monitor tracks all eigenvalues or not
31: Level: developer
33: .seealso: [](ch:eps), `EPSMonitorSet()`, `EPSSetTrackAll()`
34: @*/
35: PetscErrorCode EPSMonitorSetFromOptions(EPS eps,const char opt[],const char name[],PetscCtx ctx,PetscBool trackall)
36: {
37: PetscErrorCode (*mfunc)(EPS,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)eps),((PetscObject)eps)->options,((PetscObject)eps)->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(EPSMonitorList,key,&mfunc));
54: PetscCheck(mfunc,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Specified viewer and format not supported");
55: PetscCall(PetscFunctionListFind(EPSMonitorCreateList,key,&cfunc));
56: PetscCall(PetscFunctionListFind(EPSMonitorDestroyList,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(EPSMonitorSet(eps,mfunc,vf,(PetscCtxDestroyFn*)dfunc));
63: if (trackall) PetscCall(EPSSetTrackAll(eps,PETSC_TRUE));
64: PetscFunctionReturn(PETSC_SUCCESS);
65: }
67: /*@
68: EPSSetFromOptions - Sets `EPS` options from the options database.
69: This routine must be called before `EPSSetUp()` if the user is to be
70: allowed to configure the solver.
72: Collective
74: Input Parameter:
75: . eps - the linear eigensolver context
77: Note:
78: To see all options, run your program with the `-help` option.
80: Level: beginner
82: .seealso: [](ch:eps), `EPSSetOptionsPrefix()`
83: @*/
84: PetscErrorCode EPSSetFromOptions(EPS eps)
85: {
86: char type[256];
87: PetscBool set,flg,flg1,flg2,flg3,bval;
88: PetscReal r,array[2]={0,0};
89: PetscScalar s;
90: PetscInt i,j,k;
91: EPSBalance bal;
93: PetscFunctionBegin;
95: PetscCall(EPSRegisterAll());
96: PetscObjectOptionsBegin((PetscObject)eps);
97: PetscCall(PetscOptionsFList("-eps_type","Eigensolver method","EPSSetType",EPSList,(char*)(((PetscObject)eps)->type_name?((PetscObject)eps)->type_name:EPSKRYLOVSCHUR),type,sizeof(type),&flg));
98: if (flg) PetscCall(EPSSetType(eps,type));
99: else if (!((PetscObject)eps)->type_name) PetscCall(EPSSetType(eps,EPSKRYLOVSCHUR));
101: PetscCall(PetscOptionsBoolGroupBegin("-eps_hermitian","Hermitian eigenvalue problem","EPSSetProblemType",&flg));
102: if (flg) PetscCall(EPSSetProblemType(eps,EPS_HEP));
103: PetscCall(PetscOptionsBoolGroup("-eps_gen_hermitian","Generalized Hermitian eigenvalue problem","EPSSetProblemType",&flg));
104: if (flg) PetscCall(EPSSetProblemType(eps,EPS_GHEP));
105: PetscCall(PetscOptionsBoolGroup("-eps_non_hermitian","Non-Hermitian eigenvalue problem","EPSSetProblemType",&flg));
106: if (flg) PetscCall(EPSSetProblemType(eps,EPS_NHEP));
107: PetscCall(PetscOptionsBoolGroup("-eps_gen_non_hermitian","Generalized non-Hermitian eigenvalue problem","EPSSetProblemType",&flg));
108: if (flg) PetscCall(EPSSetProblemType(eps,EPS_GNHEP));
109: PetscCall(PetscOptionsBoolGroup("-eps_pos_gen_non_hermitian","Generalized non-Hermitian eigenvalue problem with positive semi-definite B","EPSSetProblemType",&flg));
110: if (flg) PetscCall(EPSSetProblemType(eps,EPS_PGNHEP));
111: PetscCall(PetscOptionsBoolGroup("-eps_gen_indefinite","Generalized Hermitian-indefinite eigenvalue problem","EPSSetProblemType",&flg));
112: if (flg) PetscCall(EPSSetProblemType(eps,EPS_GHIEP));
113: PetscCall(PetscOptionsBoolGroup("-eps_bse","Structured Bethe-Salpeter eigenvalue problem","EPSSetProblemType",&flg));
114: if (flg) PetscCall(EPSSetProblemType(eps,EPS_BSE));
115: PetscCall(PetscOptionsBoolGroup("-eps_hamiltonian","Structured Hamiltonian eigenvalue problem","EPSSetProblemType",&flg));
116: if (flg) PetscCall(EPSSetProblemType(eps,EPS_HAMILT));
117: PetscCall(PetscOptionsBoolGroupEnd("-eps_lrep","Structured Linear Response eigenvalue problem","EPSSetProblemType",&flg));
118: if (flg) PetscCall(EPSSetProblemType(eps,EPS_LREP));
120: PetscCall(PetscOptionsBoolGroupBegin("-eps_ritz","Rayleigh-Ritz extraction","EPSSetExtraction",&flg));
121: if (flg) PetscCall(EPSSetExtraction(eps,EPS_RITZ));
122: PetscCall(PetscOptionsBoolGroup("-eps_harmonic","Harmonic Ritz extraction","EPSSetExtraction",&flg));
123: if (flg) PetscCall(EPSSetExtraction(eps,EPS_HARMONIC));
124: PetscCall(PetscOptionsBoolGroup("-eps_harmonic_relative","Relative harmonic Ritz extraction","EPSSetExtraction",&flg));
125: if (flg) PetscCall(EPSSetExtraction(eps,EPS_HARMONIC_RELATIVE));
126: PetscCall(PetscOptionsBoolGroup("-eps_harmonic_right","Right harmonic Ritz extraction","EPSSetExtraction",&flg));
127: if (flg) PetscCall(EPSSetExtraction(eps,EPS_HARMONIC_RIGHT));
128: PetscCall(PetscOptionsBoolGroup("-eps_harmonic_largest","Largest harmonic Ritz extraction","EPSSetExtraction",&flg));
129: if (flg) PetscCall(EPSSetExtraction(eps,EPS_HARMONIC_LARGEST));
130: PetscCall(PetscOptionsBoolGroup("-eps_refined","Refined Ritz extraction","EPSSetExtraction",&flg));
131: if (flg) PetscCall(EPSSetExtraction(eps,EPS_REFINED));
132: PetscCall(PetscOptionsBoolGroupEnd("-eps_refined_harmonic","Refined harmonic Ritz extraction","EPSSetExtraction",&flg));
133: if (flg) PetscCall(EPSSetExtraction(eps,EPS_REFINED_HARMONIC));
135: bal = eps->balance;
136: PetscCall(PetscOptionsEnum("-eps_balance","Balancing method","EPSSetBalance",EPSBalanceTypes,(PetscEnum)bal,(PetscEnum*)&bal,&flg1));
137: j = eps->balance_its;
138: PetscCall(PetscOptionsInt("-eps_balance_its","Number of iterations in balancing","EPSSetBalance",eps->balance_its,&j,&flg2));
139: r = eps->balance_cutoff;
140: PetscCall(PetscOptionsReal("-eps_balance_cutoff","Cutoff value in balancing","EPSSetBalance",eps->balance_cutoff,&r,&flg3));
141: if (flg1 || flg2 || flg3) PetscCall(EPSSetBalance(eps,bal,j,r));
143: i = eps->max_it;
144: PetscCall(PetscOptionsInt("-eps_max_it","Maximum number of iterations","EPSSetTolerances",eps->max_it,&i,&flg1));
145: r = eps->tol;
146: PetscCall(PetscOptionsReal("-eps_tol","Tolerance","EPSSetTolerances",SlepcDefaultTol(eps->tol),&r,&flg2));
147: if (flg1 || flg2) PetscCall(EPSSetTolerances(eps,r,i));
149: r = eps->thres;
150: PetscCall(PetscOptionsReal("-eps_threshold_absolute","Absolute threshold","EPSSetThreshold",r,&r,&flg));
151: if (flg) PetscCall(EPSSetThreshold(eps,r,PETSC_FALSE));
152: PetscCall(PetscOptionsReal("-eps_threshold_relative","Relative threshold","EPSSetThreshold",r,&r,&flg));
153: if (flg) PetscCall(EPSSetThreshold(eps,r,PETSC_TRUE));
155: PetscCall(PetscOptionsBoolGroupBegin("-eps_conv_rel","Relative error convergence test","EPSSetConvergenceTest",&flg));
156: if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_REL));
157: PetscCall(PetscOptionsBoolGroup("-eps_conv_norm","Convergence test relative to the eigenvalue and the matrix norms","EPSSetConvergenceTest",&flg));
158: if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_NORM));
159: PetscCall(PetscOptionsBoolGroup("-eps_conv_abs","Absolute error convergence test","EPSSetConvergenceTest",&flg));
160: if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_ABS));
161: PetscCall(PetscOptionsBoolGroupEnd("-eps_conv_user","User-defined convergence test","EPSSetConvergenceTest",&flg));
162: if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_USER));
164: PetscCall(PetscOptionsBoolGroupBegin("-eps_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","EPSSetStoppingTest",&flg));
165: if (flg) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_BASIC));
166: PetscCall(PetscOptionsBoolGroup("-eps_stop_threshold","Stop iteration if a converged eigenvalue is below/above the threshold","EPSSetStoppingTest",&flg));
167: if (flg) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_THRESHOLD));
168: PetscCall(PetscOptionsBoolGroupEnd("-eps_stop_user","User-defined stopping test","EPSSetStoppingTest",&flg));
169: if (flg) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_USER));
171: i = eps->nev;
172: PetscCall(PetscOptionsInt("-eps_nev","Number of eigenvalues to compute","EPSSetDimensions",eps->nev,&i,&flg1));
173: if (!flg1) i = PETSC_CURRENT;
174: j = eps->ncv;
175: PetscCall(PetscOptionsInt("-eps_ncv","Number of basis vectors","EPSSetDimensions",eps->ncv,&j,&flg2));
176: k = eps->mpd;
177: PetscCall(PetscOptionsInt("-eps_mpd","Maximum dimension of projected problem","EPSSetDimensions",eps->mpd,&k,&flg3));
178: if (flg1 || flg2 || flg3) PetscCall(EPSSetDimensions(eps,i,j,k));
180: PetscCall(PetscOptionsBoolGroupBegin("-eps_largest_magnitude","Compute largest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg));
181: if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE));
182: PetscCall(PetscOptionsBoolGroup("-eps_smallest_magnitude","Compute smallest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg));
183: if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_MAGNITUDE));
184: PetscCall(PetscOptionsBoolGroup("-eps_largest_real","Compute eigenvalues with largest real parts","EPSSetWhichEigenpairs",&flg));
185: if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL));
186: PetscCall(PetscOptionsBoolGroup("-eps_smallest_real","Compute eigenvalues with smallest real parts","EPSSetWhichEigenpairs",&flg));
187: if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
188: PetscCall(PetscOptionsBoolGroup("-eps_largest_imaginary","Compute eigenvalues with largest imaginary parts","EPSSetWhichEigenpairs",&flg));
189: if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_IMAGINARY));
190: PetscCall(PetscOptionsBoolGroup("-eps_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","EPSSetWhichEigenpairs",&flg));
191: if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_IMAGINARY));
192: PetscCall(PetscOptionsBoolGroup("-eps_target_magnitude","Compute eigenvalues closest to target","EPSSetWhichEigenpairs",&flg));
193: if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
194: PetscCall(PetscOptionsBoolGroup("-eps_target_real","Compute eigenvalues with real parts closest to target","EPSSetWhichEigenpairs",&flg));
195: if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_REAL));
196: PetscCall(PetscOptionsBoolGroup("-eps_target_imaginary","Compute eigenvalues with imaginary parts closest to target","EPSSetWhichEigenpairs",&flg));
197: if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_IMAGINARY));
198: PetscCall(PetscOptionsBoolGroup("-eps_all","Compute all eigenvalues in an interval or a region","EPSSetWhichEigenpairs",&flg));
199: if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_ALL));
200: PetscCall(PetscOptionsBoolGroupEnd("-eps_which_user","Select the user-defined selection criterion","EPSSetWhichEigenpairs",&flg));
201: if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_WHICH_USER));
203: PetscCall(PetscOptionsScalar("-eps_target","Value of the target","EPSSetTarget",eps->target,&s,&flg));
204: if (flg) {
205: if (eps->which!=EPS_TARGET_REAL && eps->which!=EPS_TARGET_IMAGINARY) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
206: PetscCall(EPSSetTarget(eps,s));
207: }
209: k = 2;
210: PetscCall(PetscOptionsRealArray("-eps_interval","Computational interval (two real values separated with a comma without spaces)","EPSSetInterval",array,&k,&flg));
211: if (flg) {
212: PetscCheck(k>1,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_SIZ,"Must pass two values in -eps_interval (comma-separated without spaces)");
213: PetscCall(EPSSetWhichEigenpairs(eps,EPS_ALL));
214: PetscCall(EPSSetInterval(eps,array[0],array[1]));
215: }
217: PetscCall(PetscOptionsBool("-eps_true_residual","Compute true residuals explicitly","EPSSetTrueResidual",eps->trueres,&eps->trueres,NULL));
218: PetscCall(PetscOptionsBool("-eps_purify","Postprocess eigenvectors for purification","EPSSetPurify",eps->purify,&bval,&flg));
219: if (flg) PetscCall(EPSSetPurify(eps,bval));
220: PetscCall(PetscOptionsBool("-eps_two_sided","Use two-sided variant (to compute left eigenvectors)","EPSSetTwoSided",eps->twosided,&bval,&flg));
221: if (flg) PetscCall(EPSSetTwoSided(eps,bval));
223: /* -----------------------------------------------------------------------*/
224: /*
225: Cancels all monitors hardwired into code before call to EPSSetFromOptions()
226: */
227: PetscCall(PetscOptionsBool("-eps_monitor_cancel","Remove any hardwired monitor routines","EPSMonitorCancel",PETSC_FALSE,&flg,&set));
228: if (set && flg) PetscCall(EPSMonitorCancel(eps));
229: PetscCall(EPSMonitorSetFromOptions(eps,"-eps_monitor","first_approximation",NULL,PETSC_FALSE));
230: PetscCall(EPSMonitorSetFromOptions(eps,"-eps_monitor_all","all_approximations",NULL,PETSC_TRUE));
231: PetscCall(EPSMonitorSetFromOptions(eps,"-eps_monitor_conv","convergence_history",NULL,PETSC_FALSE));
233: /* -----------------------------------------------------------------------*/
234: PetscCall(PetscOptionsName("-eps_view","Print detailed information on solver used","EPSView",&set));
235: PetscCall(PetscOptionsName("-eps_view_vectors","View computed eigenvectors","EPSVectorsView",&set));
236: PetscCall(PetscOptionsName("-eps_view_values","View computed eigenvalues","EPSValuesView",&set));
237: PetscCall(PetscOptionsName("-eps_converged_reason","Print reason for convergence, and number of iterations","EPSConvergedReasonView",&set));
238: PetscCall(PetscOptionsName("-eps_error_absolute","Print absolute errors of each eigenpair","EPSErrorView",&set));
239: PetscCall(PetscOptionsName("-eps_error_relative","Print relative errors of each eigenpair","EPSErrorView",&set));
240: PetscCall(PetscOptionsName("-eps_error_backward","Print backward errors of each eigenpair","EPSErrorView",&set));
242: PetscTryTypeMethod(eps,setfromoptions,PetscOptionsObject);
243: PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)eps,PetscOptionsObject));
244: PetscOptionsEnd();
246: if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
247: PetscCall(BVSetFromOptions(eps->V));
248: if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
249: PetscCall(RGSetFromOptions(eps->rg));
250: if (eps->useds) {
251: if (!eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
252: PetscCall(EPSSetDSType(eps));
253: PetscCall(DSSetFromOptions(eps->ds));
254: }
255: if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
256: PetscCall(EPSSetDefaultST(eps));
257: PetscCall(STSetFromOptions(eps->st));
258: PetscFunctionReturn(PETSC_SUCCESS);
259: }
261: /*@
262: EPSGetTolerances - Gets the tolerance and maximum iteration count used
263: by the `EPS` convergence tests.
265: Not Collective
267: Input Parameter:
268: . eps - the linear eigensolver context
270: Output Parameters:
271: + tol - the convergence tolerance
272: - maxits - maximum number of iterations
274: Notes:
275: The user can specify `NULL` for any parameter that is not needed.
277: Level: intermediate
279: .seealso: [](ch:eps), `EPSSetTolerances()`
280: @*/
281: PetscErrorCode EPSGetTolerances(EPS eps,PetscReal *tol,PetscInt *maxits)
282: {
283: PetscFunctionBegin;
285: if (tol) *tol = eps->tol;
286: if (maxits) *maxits = eps->max_it;
287: PetscFunctionReturn(PETSC_SUCCESS);
288: }
290: /*@
291: EPSSetTolerances - Sets the tolerance and maximum iteration count used
292: by the `EPS` convergence tests.
294: Logically Collective
296: Input Parameters:
297: + eps - the linear eigensolver context
298: . tol - the convergence tolerance
299: - maxits - maximum number of iterations to use
301: Options Database Keys:
302: + -eps_tol tol - sets the convergence tolerance
303: - -eps_max_it maxits - sets the maximum number of iterations allowed
305: Note:
306: Use `PETSC_CURRENT` to retain the current value of any of the parameters.
307: Use `PETSC_DETERMINE` for either argument to assign a default value computed
308: internally (may be different in each solver).
309: For `maxits` use `PETSC_UNLIMITED` to indicate there is no upper bound on this value.
311: Level: intermediate
313: .seealso: [](ch:eps), `EPSGetTolerances()`
314: @*/
315: PetscErrorCode EPSSetTolerances(EPS eps,PetscReal tol,PetscInt maxits)
316: {
317: PetscFunctionBegin;
321: if (tol == (PetscReal)PETSC_DETERMINE) {
322: eps->tol = PETSC_DETERMINE;
323: eps->state = EPS_STATE_INITIAL;
324: } else if (tol != (PetscReal)PETSC_CURRENT) {
325: PetscCheck(tol>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
326: eps->tol = tol;
327: }
328: if (maxits == PETSC_DETERMINE) {
329: eps->max_it = PETSC_DETERMINE;
330: eps->state = EPS_STATE_INITIAL;
331: } else if (maxits == PETSC_UNLIMITED) {
332: eps->max_it = PETSC_INT_MAX;
333: } else if (maxits != PETSC_CURRENT) {
334: PetscCheck(maxits>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
335: eps->max_it = maxits;
336: }
337: PetscFunctionReturn(PETSC_SUCCESS);
338: }
340: /*@
341: EPSGetDimensions - Gets the number of eigenvalues to compute
342: and the dimension of the subspace.
344: Not Collective
346: Input Parameter:
347: . eps - the linear 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: Level: intermediate
356: .seealso: [](ch:eps), `EPSSetDimensions()`
357: @*/
358: PetscErrorCode EPSGetDimensions(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
359: {
360: PetscFunctionBegin;
362: if (nev) *nev = eps->nev? eps->nev: 1;
363: if (ncv) *ncv = eps->ncv;
364: if (mpd) *mpd = eps->mpd;
365: PetscFunctionReturn(PETSC_SUCCESS);
366: }
368: /*@
369: EPSSetDimensions - Sets the number of eigenvalues to compute
370: and the dimension of the subspace.
372: Logically Collective
374: Input Parameters:
375: + eps - the linear eigensolver context
376: . nev - number of eigenvalues to compute
377: . ncv - the maximum dimension of the subspace to be used by the solver
378: - mpd - the maximum dimension allowed for the projected problem
380: Options Database Keys:
381: + -eps_nev nev - sets the number of eigenvalues
382: . -eps_ncv ncv - sets the dimension of the subspace
383: - -eps_mpd mpd - sets the maximum projected dimension
385: Notes:
386: Use `PETSC_DETERMINE` for `ncv` and `mpd` to assign a reasonably good value, which is
387: dependent on the solution method. For any of the arguments, use `PETSC_CURRENT`
388: to preserve the current value.
390: The parameters `ncv` and `mpd` are intimately related, so that the user is advised
391: to set one of them at most. Normal usage is\:
393: 1. In cases where `nev` is small, the user sets `ncv` (a reasonable default is `2*nev`).
394: 2. In cases where `nev` is large, the user sets `mpd`.
396: The value of `ncv` should always be between `nev` and `(nev+mpd)`, typically
397: `ncv=nev+mpd`. If `nev` is not too large, `mpd=nev` is a reasonable choice, otherwise
398: a smaller value should be used.
400: When computing all eigenvalues in an interval, see `EPSSetInterval()`, these
401: parameters lose relevance, and tuning must be done with
402: `EPSKrylovSchurSetDimensions()`.
404: Level: intermediate
406: .seealso: [](ch:eps), `EPSGetDimensions()`, `EPSSetInterval()`, `EPSKrylovSchurSetDimensions()`
407: @*/
408: PetscErrorCode EPSSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
409: {
410: PetscFunctionBegin;
415: if (nev != PETSC_CURRENT) {
416: PetscCheck(nev>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
417: eps->nev = nev;
418: }
419: if (ncv == PETSC_DETERMINE) {
420: eps->ncv = PETSC_DETERMINE;
421: } else if (ncv != PETSC_CURRENT) {
422: PetscCheck(ncv>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
423: eps->ncv = ncv;
424: }
425: if (mpd == PETSC_DETERMINE) {
426: eps->mpd = PETSC_DETERMINE;
427: } else if (mpd != PETSC_CURRENT) {
428: PetscCheck(mpd>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
429: eps->mpd = mpd;
430: }
431: eps->state = EPS_STATE_INITIAL;
432: PetscFunctionReturn(PETSC_SUCCESS);
433: }
435: /*@
436: EPSSetWhichEigenpairs - Specifies which portion of the spectrum is
437: to be sought.
439: Logically Collective
441: Input Parameters:
442: + eps - the linear eigensolver context
443: - which - the portion of the spectrum to be sought, see `EPSWhich` for possible values
445: Options Database Keys:
446: + -eps_largest_magnitude - sets largest eigenvalues in magnitude
447: . -eps_smallest_magnitude - sets smallest eigenvalues in magnitude
448: . -eps_largest_real - sets largest real parts
449: . -eps_smallest_real - sets smallest real parts
450: . -eps_largest_imaginary - sets largest imaginary parts
451: . -eps_smallest_imaginary - sets smallest imaginary parts
452: . -eps_target_magnitude - sets eigenvalues closest to target
453: . -eps_target_real - sets real parts closest to target
454: . -eps_target_imaginary - sets imaginary parts closest to target
455: . -eps_all - sets all eigenvalues in an interval or region
456: - -eps_which_user - select the user-defined selection criterion
458: Notes:
459: Not all eigensolvers implemented in `EPS` account for all the possible values
460: of `which`. Also, some values make sense only for certain types of
461: problems. If SLEPc is compiled for real numbers `EPS_LARGEST_IMAGINARY`
462: and `EPS_SMALLEST_IMAGINARY` use the absolute value of the imaginary part
463: for eigenvalue selection.
465: The target is a scalar value provided with `EPSSetTarget()`.
467: The criterion `EPS_TARGET_IMAGINARY` is available only in case PETSc and
468: SLEPc have been built with complex scalars.
470: `EPS_ALL` is intended for use in combination with an interval (see
471: `EPSSetInterval()`), when all eigenvalues within the interval are requested,
472: or in the context of the `EPSCISS` solver for computing all eigenvalues in a region.
474: Level: intermediate
476: .seealso: [](ch:eps), `EPSGetWhichEigenpairs()`, `EPSSetTarget()`, `EPSSetInterval()`, `EPSSetDimensions()`, `EPSSetEigenvalueComparison()`, `EPSWhich`
477: @*/
478: PetscErrorCode EPSSetWhichEigenpairs(EPS eps,EPSWhich which)
479: {
480: PetscFunctionBegin;
483: switch (which) {
484: case EPS_LARGEST_MAGNITUDE:
485: case EPS_SMALLEST_MAGNITUDE:
486: case EPS_LARGEST_REAL:
487: case EPS_SMALLEST_REAL:
488: case EPS_LARGEST_IMAGINARY:
489: case EPS_SMALLEST_IMAGINARY:
490: case EPS_TARGET_MAGNITUDE:
491: case EPS_TARGET_REAL:
492: #if defined(PETSC_USE_COMPLEX)
493: case EPS_TARGET_IMAGINARY:
494: #endif
495: case EPS_ALL:
496: case EPS_WHICH_USER:
497: if (eps->which != which) {
498: eps->state = EPS_STATE_INITIAL;
499: eps->which = which;
500: }
501: break;
502: #if !defined(PETSC_USE_COMPLEX)
503: case EPS_TARGET_IMAGINARY:
504: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"EPS_TARGET_IMAGINARY can be used only with complex scalars");
505: #endif
506: default:
507: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
508: }
509: PetscFunctionReturn(PETSC_SUCCESS);
510: }
512: /*@
513: EPSGetWhichEigenpairs - Returns which portion of the spectrum is to be
514: sought.
516: Not Collective
518: Input Parameter:
519: . eps - the linear eigensolver context
521: Output Parameter:
522: . which - the portion of the spectrum to be sought, see `EPSWhich` for possible values
524: Level: intermediate
526: .seealso: [](ch:eps), `EPSSetWhichEigenpairs()`, `EPSWhich`
527: @*/
528: PetscErrorCode EPSGetWhichEigenpairs(EPS eps,EPSWhich *which)
529: {
530: PetscFunctionBegin;
532: PetscAssertPointer(which,2);
533: *which = eps->which;
534: PetscFunctionReturn(PETSC_SUCCESS);
535: }
537: /*@
538: EPSSetThreshold - Sets the threshold used in the threshold stopping test.
540: Logically Collective
542: Input Parameters:
543: + eps - the linear eigensolver context
544: . thres - the threshold value
545: - rel - whether the threshold is relative or not
547: Options Database Keys:
548: + -eps_threshold_absolute thres - sets an absolute threshold
549: - -eps_threshold_relative thres - sets a relative threshold
551: Notes:
552: This function internally calls `EPSSetStoppingTest()` to set a special stopping
553: test based on the threshold, where eigenvalues are computed in sequence until
554: one of the computed eigenvalues is below the threshold `thres` (in magnitude).
555: This is the interpretation in case of searching for largest eigenvalues in magnitude,
556: see `EPSSetWhichEigenpairs()`.
558: If the solver is configured to compute smallest magnitude eigenvalues, then the
559: threshold must be interpreted in the opposite direction, i.e., the computation
560: will stop when one of the computed values is above the threshold (in magnitude).
562: The threshold can also be used when computing largest/smallest real eigenvalues
563: (i.e, rightmost or leftmost), in which case the threshold is allowed to be
564: negative. The solver will stop when one of the computed eigenvalues is above
565: or below the threshold (considering the real part of the eigenvalue). This mode
566: is allowed only in problem types whose eigenvalues are always real (e.g., `EPS_HEP`).
568: In the case of largest magnitude eigenvalues, the threshold can be made relative
569: with respect to the dominant eigenvalue. Otherwise, the argument `rel` should be
570: `PETSC_FALSE`.
572: An additional use case is with target magnitude selection of eigenvalues (e.g.,
573: with shift-and-invert), but this must be used with caution to avoid unexpected
574: behavior. With an absolute threshold, the solver will assume that leftmost
575: eigenvalues are being computed (e.g., with `target`=0 for a problem with real
576: positive eigenvalues). In case of a relative threshold, a value of `thres`<1
577: implies that the wanted eigenvalues are the largest ones, and otherwise the
578: solver assumes that smallest eigenvalues are being computed.
580: The test against the threshold is done for converged eigenvalues, which
581: implies that the final number of converged eigenvalues will be at least
582: one more than the actual number of values below/above the threshold.
584: Since the number of computed eigenvalues is not known a priori, the solver
585: will need to reallocate the basis of vectors internally, to have enough room
586: to accommodate all the eigenvectors. Hence, this option must be used with
587: caution to avoid out-of-memory problems. The recommendation is to set the value
588: of `ncv` to be larger than the estimated number of eigenvalues, to minimize the
589: number of reallocations.
591: If a number of wanted eigenvalues has been set with `EPSSetDimensions()`
592: it is also taken into account and the solver will stop when one of the two
593: conditions (threshold or number of converged values) is met.
595: Use `EPSSetStoppingTest()` to return to the usual computation of a fixed number
596: of eigenvalues.
598: Level: advanced
600: .seealso: [](ch:eps), `EPSGetThreshold()`, `EPSSetStoppingTest()`, `EPSSetDimensions()`, `EPSSetWhichEigenpairs()`, `EPSSetProblemType()`
601: @*/
602: PetscErrorCode EPSSetThreshold(EPS eps,PetscReal thres,PetscBool rel)
603: {
604: PetscFunctionBegin;
608: if (eps->thres != thres || eps->threlative != rel) {
609: eps->thres = thres;
610: eps->threlative = rel;
611: eps->state = EPS_STATE_INITIAL;
612: PetscCall(EPSSetStoppingTest(eps,EPS_STOP_THRESHOLD));
613: }
614: PetscFunctionReturn(PETSC_SUCCESS);
615: }
617: /*@
618: EPSGetThreshold - Gets the threshold used by the threshold stopping test.
620: Not Collective
622: Input Parameter:
623: . eps - the linear eigensolver context
625: Output Parameters:
626: + thres - the threshold
627: - rel - whether the threshold is relative or not
629: Level: advanced
631: .seealso: [](ch:eps), `EPSSetThreshold()`
632: @*/
633: PetscErrorCode EPSGetThreshold(EPS eps,PetscReal *thres,PetscBool *rel)
634: {
635: PetscFunctionBegin;
637: if (thres) *thres = eps->thres;
638: if (rel) *rel = eps->threlative;
639: PetscFunctionReturn(PETSC_SUCCESS);
640: }
642: /*@C
643: EPSSetEigenvalueComparison - Specifies the eigenvalue comparison function
644: when `EPSSetWhichEigenpairs()` is set to `EPS_WHICH_USER`.
646: Logically Collective
648: Input Parameters:
649: + eps - the linear eigensolver context
650: . func - the comparison function, see `SlepcEigenvalueComparisonFn` for the calling sequence
651: - ctx - a context pointer (the last parameter to the comparison function)
653: Level: advanced
655: .seealso: [](ch:eps), `EPSSetWhichEigenpairs()`, `EPSWhich`
656: @*/
657: PetscErrorCode EPSSetEigenvalueComparison(EPS eps,SlepcEigenvalueComparisonFn *func,PetscCtx ctx)
658: {
659: PetscFunctionBegin;
661: eps->sc->comparison = func;
662: eps->sc->comparisonctx = ctx;
663: eps->which = EPS_WHICH_USER;
664: PetscFunctionReturn(PETSC_SUCCESS);
665: }
667: /*@C
668: EPSSetArbitrarySelection - Specifies a function intended to look for
669: eigenvalues according to an arbitrary selection criterion. This criterion
670: can be based on a computation involving the current eigenvector approximation.
672: Logically Collective
674: Input Parameters:
675: + eps - the linear eigensolver context
676: . func - the arbitrary selection function, see `SlepcArbitrarySelectionFn` for a calling sequence
677: - ctx - a context pointer (the last parameter to the arbitrary selection function)
679: Notes:
680: This provides a mechanism to select eigenpairs by evaluating a user-defined
681: function. When a function has been provided, the default selection based on
682: sorting the eigenvalues is replaced by the sorting of the results of this
683: function (with the same sorting criterion given in `EPSSetWhichEigenpairs()`).
685: For instance, suppose you want to compute those eigenvectors that maximize
686: a certain computable expression. Then implement the computation using
687: the arguments `xr` and `xi`, and return the result in `rr`. Then set the standard
688: sorting by magnitude so that the eigenpair with largest value of `rr` is
689: selected.
691: This evaluation function is collective, that is, all processes call it and
692: it can use collective operations; furthermore, the computed result must
693: be the same in all processes.
695: The result of `func` is expressed as a complex number so that it is possible to
696: use the standard eigenvalue sorting functions, but normally only `rr` is used.
697: Set `ri` to zero unless it is meaningful in your application.
699: Level: advanced
701: .seealso: [](ch:eps), `EPSSetWhichEigenpairs()`, `EPSSetArbitrarySelectionContextDestroy()`
702: @*/
703: PetscErrorCode EPSSetArbitrarySelection(EPS eps,SlepcArbitrarySelectionFn *func,PetscCtx ctx)
704: {
705: PetscFunctionBegin;
707: if (eps->arbitrarydestroy) PetscCall((*eps->arbitrarydestroy)(&eps->arbitraryctx));
708: eps->arbitrary = func;
709: eps->arbitraryctx = ctx;
710: eps->state = EPS_STATE_INITIAL;
711: PetscFunctionReturn(PETSC_SUCCESS);
712: }
714: /*@C
715: EPSSetArbitrarySelectionContextDestroy - Set a context destroy function for the
716: context used in the arbitrary selection.
718: Logically Collective
720: Input Parameters:
721: + eps - the linear eigensolver context
722: - destroy - context destroy function, see `PetscCtxDestroyFn` for its calling sequence
724: Level: advanced
726: .seealso: [](ch:eps), `EPSSetArbitrarySelection()`
727: @*/
728: PetscErrorCode EPSSetArbitrarySelectionContextDestroy(EPS eps,PetscCtxDestroyFn *destroy)
729: {
730: PetscFunctionBegin;
732: eps->arbitrarydestroy = destroy;
733: PetscFunctionReturn(PETSC_SUCCESS);
734: }
736: /*@C
737: EPSSetConvergenceTestFunction - Sets a function to compute the error estimate
738: used in the convergence test.
740: Logically Collective
742: Input Parameters:
743: + eps - the linear eigensolver context
744: . func - convergence test function, see `EPSConvergenceTestFn` for the calling sequence
745: . ctx - context for private data for the convergence routine (may be `NULL`)
746: - destroy - a routine for destroying the context (may be `NULL`), see `PetscCtxDestroyFn`
747: for the calling sequence
749: Notes:
750: When this is called with a user-defined function, then the convergence
751: criterion is set to `EPS_CONV_USER`, see `EPSSetConvergenceTest()`.
753: If the error estimate returned by the convergence test function is less than
754: the tolerance, then the eigenvalue is accepted as converged.
756: Level: advanced
758: .seealso: [](ch:eps), `EPSSetConvergenceTest()`, `EPSSetTolerances()`
759: @*/
760: PetscErrorCode EPSSetConvergenceTestFunction(EPS eps,EPSConvergenceTestFn *func,PetscCtx ctx,PetscCtxDestroyFn *destroy)
761: {
762: PetscFunctionBegin;
764: if (eps->convergeddestroy) PetscCall((*eps->convergeddestroy)(&eps->convergedctx));
765: eps->convergeduser = func;
766: eps->convergeddestroy = destroy;
767: eps->convergedctx = ctx;
768: if (func == EPSConvergedRelative) eps->conv = EPS_CONV_REL;
769: else if (func == EPSConvergedNorm) eps->conv = EPS_CONV_NORM;
770: else if (func == EPSConvergedAbsolute) eps->conv = EPS_CONV_ABS;
771: else {
772: eps->conv = EPS_CONV_USER;
773: eps->converged = eps->convergeduser;
774: }
775: PetscFunctionReturn(PETSC_SUCCESS);
776: }
778: /*@
779: EPSSetConvergenceTest - Specifies how to compute the error estimate
780: used in the convergence test.
782: Logically Collective
784: Input Parameters:
785: + eps - the linear eigensolver context
786: - conv - the type of convergence test, see `EPSConv` for possible values
788: Options Database Keys:
789: + -eps_conv_abs - sets the absolute convergence test
790: . -eps_conv_rel - sets the convergence test relative to the eigenvalue
791: . -eps_conv_norm - sets the convergence test relative to the matrix norms
792: - -eps_conv_user - selects the user-defined convergence test
794: Level: intermediate
796: .seealso: [](ch:eps), `EPSGetConvergenceTest()`, `EPSSetConvergenceTestFunction()`, `EPSSetStoppingTest()`, `EPSConv`
797: @*/
798: PetscErrorCode EPSSetConvergenceTest(EPS eps,EPSConv conv)
799: {
800: PetscFunctionBegin;
803: switch (conv) {
804: case EPS_CONV_ABS: eps->converged = EPSConvergedAbsolute; break;
805: case EPS_CONV_REL: eps->converged = EPSConvergedRelative; break;
806: case EPS_CONV_NORM: eps->converged = EPSConvergedNorm; break;
807: case EPS_CONV_USER:
808: PetscCheck(eps->convergeduser,PetscObjectComm((PetscObject)eps),PETSC_ERR_ORDER,"Must call EPSSetConvergenceTestFunction() first");
809: eps->converged = eps->convergeduser;
810: break;
811: default:
812: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
813: }
814: eps->conv = conv;
815: PetscFunctionReturn(PETSC_SUCCESS);
816: }
818: /*@
819: EPSGetConvergenceTest - Gets the method used to compute the error estimate
820: used in the convergence test.
822: Not Collective
824: Input Parameter:
825: . eps - the linear eigensolver context
827: Output Parameter:
828: . conv - the type of convergence test
830: Level: intermediate
832: .seealso: [](ch:eps), `EPSSetConvergenceTest()`, `EPSConv`
833: @*/
834: PetscErrorCode EPSGetConvergenceTest(EPS eps,EPSConv *conv)
835: {
836: PetscFunctionBegin;
838: PetscAssertPointer(conv,2);
839: *conv = eps->conv;
840: PetscFunctionReturn(PETSC_SUCCESS);
841: }
843: /*@C
844: EPSSetStoppingTestFunction - Sets a function to decide when to stop the outer
845: iteration of the eigensolver.
847: Logically Collective
849: Input Parameters:
850: + eps - the linear eigensolver context
851: . func - stopping test function, see `EPSStoppingTestFn` for the calling sequence
852: . ctx - context for private data for the stopping routine (may be `NULL`)
853: - destroy - a routine for destroying the context (may be `NULL`), see `PetscCtxDestroyFn`
854: for the calling sequence
856: Note:
857: When implementing a function for this, normal usage is to first call the
858: default routine `EPSStoppingBasic()` and then set `reason` to `EPS_CONVERGED_USER`
859: if some user-defined conditions have been met. To let the eigensolver continue
860: iterating, the result must be left as `EPS_CONVERGED_ITERATING`.
862: Level: advanced
864: .seealso: [](ch:eps), `EPSSetStoppingTest()`, `EPSStoppingBasic()`
865: @*/
866: PetscErrorCode EPSSetStoppingTestFunction(EPS eps,EPSStoppingTestFn *func,PetscCtx ctx,PetscCtxDestroyFn *destroy)
867: {
868: PetscFunctionBegin;
870: if (eps->stoppingdestroy) PetscCall((*eps->stoppingdestroy)(&eps->stoppingctx));
871: eps->stoppinguser = func;
872: eps->stoppingdestroy = destroy;
873: eps->stoppingctx = ctx;
874: if (func == EPSStoppingBasic) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_BASIC));
875: else if (func == EPSStoppingThreshold) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_THRESHOLD));
876: else {
877: eps->stop = EPS_STOP_USER;
878: eps->stopping = eps->stoppinguser;
879: }
880: PetscFunctionReturn(PETSC_SUCCESS);
881: }
883: /*@
884: EPSSetStoppingTest - Specifies how to decide the termination of the outer
885: loop of the eigensolver.
887: Logically Collective
889: Input Parameters:
890: + eps - the linear eigensolver context
891: - stop - the type of stopping test, see `EPSStop`
893: Options Database Keys:
894: + -eps_stop_basic - sets the default stopping test
895: . -eps_stop_threshold - sets the threshold stopping test
896: - -eps_stop_user - selects the user-defined stopping test
898: Level: advanced
900: .seealso: [](ch:eps), `EPSGetStoppingTest()`, `EPSSetStoppingTestFunction()`, `EPSSetConvergenceTest()`, `EPSStop`
901: @*/
902: PetscErrorCode EPSSetStoppingTest(EPS eps,EPSStop stop)
903: {
904: PetscFunctionBegin;
907: switch (stop) {
908: case EPS_STOP_BASIC: eps->stopping = EPSStoppingBasic; break;
909: case EPS_STOP_THRESHOLD: eps->stopping = EPSStoppingThreshold; break;
910: case EPS_STOP_USER:
911: PetscCheck(eps->stoppinguser,PetscObjectComm((PetscObject)eps),PETSC_ERR_ORDER,"Must call EPSSetStoppingTestFunction() first");
912: eps->stopping = eps->stoppinguser;
913: break;
914: default:
915: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
916: }
917: eps->stop = stop;
918: PetscFunctionReturn(PETSC_SUCCESS);
919: }
921: /*@
922: EPSGetStoppingTest - Gets the method used to decide the termination of the outer
923: loop of the eigensolver.
925: Not Collective
927: Input Parameter:
928: . eps - the linear eigensolver context
930: Output Parameter:
931: . stop - the type of stopping test
933: Level: advanced
935: .seealso: [](ch:eps), `EPSSetStoppingTest()`, `EPSStop`
936: @*/
937: PetscErrorCode EPSGetStoppingTest(EPS eps,EPSStop *stop)
938: {
939: PetscFunctionBegin;
941: PetscAssertPointer(stop,2);
942: *stop = eps->stop;
943: PetscFunctionReturn(PETSC_SUCCESS);
944: }
946: /*@
947: EPSSetProblemType - Specifies the type of the eigenvalue problem.
949: Logically Collective
951: Input Parameters:
952: + eps - the linear eigensolver context
953: - type - a known type of eigenvalue problem
955: Options Database Keys:
956: + -eps_hermitian - Hermitian eigenvalue problem
957: . -eps_gen_hermitian - generalized Hermitian eigenvalue problem
958: . -eps_non_hermitian - non-Hermitian eigenvalue problem
959: . -eps_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
960: . -eps_pos_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
961: with positive semi-definite $B$
962: . -eps_gen_indefinite - generalized Hermitian-indefinite eigenvalue problem
963: . -eps_bse - structured Bethe-Salpeter eigenvalue problem
964: . -eps_hamiltonian - structured Hamiltonian eigenvalue problem
965: - -eps_lrep - structured Linear Response eigenvalue problem
967: Notes:
968: This function must be used to instruct SLEPc to exploit symmetry or other
969: kind of structure. If no
970: problem type is specified, by default a non-Hermitian problem is assumed
971: (either standard or generalized). If the user knows that the problem is
972: Hermitian (i.e., $A=A^*$) or generalized Hermitian (i.e., $A=A^*$, $B=B^*$,
973: and $B$ positive definite) then it is recommended to set the problem type so
974: that the eigensolver can exploit these properties.
976: If the user does not call this function, the solver will use a reasonable
977: guess.
979: For structured problem types such as `EPS_BSE`, the matrices passed in via
980: `EPSSetOperators()` must have been created with the corresponding helper
981: function, i.e., `MatCreateBSE()`.
983: Level: intermediate
985: .seealso: [](ch:eps), `EPSSetOperators()`, `EPSSetType()`, `EPSGetProblemType()`, `EPSProblemType`
986: @*/
987: PetscErrorCode EPSSetProblemType(EPS eps,EPSProblemType type)
988: {
989: PetscFunctionBegin;
992: if (type == eps->problem_type) PetscFunctionReturn(PETSC_SUCCESS);
993: switch (type) {
994: case EPS_HEP:
995: eps->isgeneralized = PETSC_FALSE;
996: eps->ishermitian = PETSC_TRUE;
997: eps->ispositive = PETSC_FALSE;
998: eps->isstructured = PETSC_FALSE;
999: break;
1000: case EPS_NHEP:
1001: eps->isgeneralized = PETSC_FALSE;
1002: eps->ishermitian = PETSC_FALSE;
1003: eps->ispositive = PETSC_FALSE;
1004: eps->isstructured = PETSC_FALSE;
1005: break;
1006: case EPS_GHEP:
1007: eps->isgeneralized = PETSC_TRUE;
1008: eps->ishermitian = PETSC_TRUE;
1009: eps->ispositive = PETSC_TRUE;
1010: eps->isstructured = PETSC_FALSE;
1011: break;
1012: case EPS_GNHEP:
1013: eps->isgeneralized = PETSC_TRUE;
1014: eps->ishermitian = PETSC_FALSE;
1015: eps->ispositive = PETSC_FALSE;
1016: eps->isstructured = PETSC_FALSE;
1017: break;
1018: case EPS_PGNHEP:
1019: eps->isgeneralized = PETSC_TRUE;
1020: eps->ishermitian = PETSC_FALSE;
1021: eps->ispositive = PETSC_TRUE;
1022: eps->isstructured = PETSC_FALSE;
1023: break;
1024: case EPS_GHIEP:
1025: eps->isgeneralized = PETSC_TRUE;
1026: eps->ishermitian = PETSC_TRUE;
1027: eps->ispositive = PETSC_FALSE;
1028: eps->isstructured = PETSC_FALSE;
1029: break;
1030: case EPS_BSE:
1031: eps->isgeneralized = PETSC_FALSE;
1032: eps->ishermitian = PETSC_FALSE;
1033: eps->ispositive = PETSC_FALSE;
1034: eps->isstructured = PETSC_TRUE;
1035: break;
1036: case EPS_HAMILT:
1037: eps->isgeneralized = PETSC_FALSE;
1038: eps->ishermitian = PETSC_FALSE;
1039: eps->ispositive = PETSC_FALSE;
1040: eps->isstructured = PETSC_TRUE;
1041: break;
1042: case EPS_LREP:
1043: eps->isgeneralized = PETSC_FALSE;
1044: eps->ishermitian = PETSC_FALSE;
1045: eps->ispositive = PETSC_FALSE;
1046: eps->isstructured = PETSC_TRUE;
1047: break;
1048: default:
1049: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
1050: }
1051: eps->problem_type = type;
1052: eps->state = EPS_STATE_INITIAL;
1053: PetscFunctionReturn(PETSC_SUCCESS);
1054: }
1056: /*@
1057: EPSGetProblemType - Gets the problem type from the `EPS` object.
1059: Not Collective
1061: Input Parameter:
1062: . eps - the linear eigensolver context
1064: Output Parameter:
1065: . type - the problem type
1067: Level: intermediate
1069: .seealso: [](ch:eps), `EPSSetProblemType()`, `EPSProblemType`
1070: @*/
1071: PetscErrorCode EPSGetProblemType(EPS eps,EPSProblemType *type)
1072: {
1073: PetscFunctionBegin;
1075: PetscAssertPointer(type,2);
1076: *type = eps->problem_type;
1077: PetscFunctionReturn(PETSC_SUCCESS);
1078: }
1080: /*@
1081: EPSSetExtraction - Specifies the type of extraction technique to be employed
1082: by the eigensolver.
1084: Logically Collective
1086: Input Parameters:
1087: + eps - the linear eigensolver context
1088: - extr - a known type of extraction
1090: Options Database Keys:
1091: + -eps_ritz - Rayleigh-Ritz extraction
1092: . -eps_harmonic - harmonic Ritz extraction
1093: . -eps_harmonic_relative - harmonic Ritz extraction relative to the eigenvalue
1094: . -eps_harmonic_right - harmonic Ritz extraction for rightmost eigenvalues
1095: . -eps_harmonic_largest - harmonic Ritz extraction for largest magnitude (without target)
1096: . -eps_refined - refined Ritz extraction
1097: - -eps_refined_harmonic - refined harmonic Ritz extraction
1099: Notes:
1100: Not all eigensolvers support all types of extraction.
1102: By default, a standard Rayleigh-Ritz extraction is used. Other extractions
1103: may be useful when computing interior eigenvalues.
1105: Harmonic-type extractions are used in combination with a target, see `EPSSetTarget()`.
1107: Level: advanced
1109: .seealso: [](ch:eps), [](#sec:harmonic), `EPSSetTarget()`, `EPSGetExtraction()`, `EPSExtraction`
1110: @*/
1111: PetscErrorCode EPSSetExtraction(EPS eps,EPSExtraction extr)
1112: {
1113: PetscFunctionBegin;
1116: if (eps->extraction != extr) {
1117: eps->state = EPS_STATE_INITIAL;
1118: eps->extraction = extr;
1119: }
1120: PetscFunctionReturn(PETSC_SUCCESS);
1121: }
1123: /*@
1124: EPSGetExtraction - Gets the extraction type used by the `EPS` object.
1126: Not Collective
1128: Input Parameter:
1129: . eps - the linear eigensolver context
1131: Output Parameter:
1132: . extr - name of extraction type
1134: Level: advanced
1136: .seealso: [](ch:eps), `EPSSetExtraction()`, `EPSExtraction`
1137: @*/
1138: PetscErrorCode EPSGetExtraction(EPS eps,EPSExtraction *extr)
1139: {
1140: PetscFunctionBegin;
1142: PetscAssertPointer(extr,2);
1143: *extr = eps->extraction;
1144: PetscFunctionReturn(PETSC_SUCCESS);
1145: }
1147: /*@
1148: EPSSetBalance - Specifies the balancing technique to be employed by the
1149: eigensolver, and some parameters associated to it.
1151: Logically Collective
1153: Input Parameters:
1154: + eps - the linear eigensolver context
1155: . bal - the balancing method, see `EPSBalance` for possible values
1156: . its - number of iterations of the balancing algorithm
1157: - cutoff - cutoff value
1159: Options Database Keys:
1160: + -eps_balance (none|oneside|twoside|user) - the balancing method
1161: . -eps_balance_its its - number of iterations
1162: - -eps_balance_cutoff cutoff - cutoff value
1164: Notes:
1165: When balancing is enabled, the solver works implicitly with matrix $DAD^{-1}$,
1166: where $D$ is an appropriate diagonal matrix. This improves the accuracy of
1167: the computed results in some cases, see [](sec:balancing).
1169: Balancing makes sense only for non-Hermitian problems when the required
1170: precision is high (i.e., a small tolerance such as `1e-14`).
1172: By default, balancing is disabled. The two-sided method is much more
1173: effective than the one-sided counterpart, but it requires the system
1174: matrices to have the `MatMultTranspose()` operation defined. The methods
1175: are described in {cite:p}`Che00`.
1177: The parameter `its` is the number of iterations performed by the method. The
1178: `cutoff` value is used only in the two-side variant. Use `PETSC_DETERMINE` to assign
1179: a reasonably good value, or `PETSC_CURRENT` to leave the value unchanged.
1181: User-defined balancing is allowed provided that the corresponding matrix
1182: is set via `STSetBalanceMatrix()`.
1184: Level: intermediate
1186: .seealso: [](ch:eps), [](sec:balancing), `EPSGetBalance()`, `EPSBalance`, `STSetBalanceMatrix()`
1187: @*/
1188: PetscErrorCode EPSSetBalance(EPS eps,EPSBalance bal,PetscInt its,PetscReal cutoff)
1189: {
1190: PetscFunctionBegin;
1195: switch (bal) {
1196: case EPS_BALANCE_NONE:
1197: case EPS_BALANCE_ONESIDE:
1198: case EPS_BALANCE_TWOSIDE:
1199: case EPS_BALANCE_USER:
1200: if (eps->balance != bal) {
1201: eps->state = EPS_STATE_INITIAL;
1202: eps->balance = bal;
1203: }
1204: break;
1205: default:
1206: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'bal'");
1207: }
1208: if (its==PETSC_DETERMINE) eps->balance_its = 5;
1209: else if (its!=PETSC_CURRENT) {
1210: PetscCheck(its>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be > 0");
1211: eps->balance_its = its;
1212: }
1213: if (cutoff==(PetscReal)PETSC_DETERMINE) eps->balance_cutoff = 1e-8;
1214: else if (cutoff!=(PetscReal)PETSC_CURRENT) {
1215: PetscCheck(cutoff>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of cutoff. Must be > 0");
1216: eps->balance_cutoff = cutoff;
1217: }
1218: PetscFunctionReturn(PETSC_SUCCESS);
1219: }
1221: /*@
1222: EPSGetBalance - Gets the balancing type used by the `EPS` object, and the
1223: associated parameters.
1225: Not Collective
1227: Input Parameter:
1228: . eps - the linear eigensolver context
1230: Output Parameters:
1231: + bal - the balancing method
1232: . its - number of iterations of the balancing algorithm
1233: - cutoff - cutoff value
1235: Level: intermediate
1237: Note:
1238: The user can specify `NULL` for any parameter that is not needed.
1240: .seealso: [](ch:eps), `EPSSetBalance()`, `EPSBalance`
1241: @*/
1242: PetscErrorCode EPSGetBalance(EPS eps,EPSBalance *bal,PetscInt *its,PetscReal *cutoff)
1243: {
1244: PetscFunctionBegin;
1246: if (bal) *bal = eps->balance;
1247: if (its) *its = eps->balance_its;
1248: if (cutoff) *cutoff = eps->balance_cutoff;
1249: PetscFunctionReturn(PETSC_SUCCESS);
1250: }
1252: /*@
1253: EPSSetTwoSided - Sets the solver to use a two-sided variant so that left
1254: eigenvectors are also computed.
1256: Logically Collective
1258: Input Parameters:
1259: + eps - the linear eigensolver context
1260: - twosided - whether the two-sided variant is to be used or not
1262: Options Database Key:
1263: . -eps_two_sided (true|false) - toggles the twosided flag
1265: Notes:
1266: If the user sets `twosided`=`PETSC_TRUE` then the solver uses a variant of
1267: the algorithm that computes both right and left eigenvectors. This is
1268: usually much more costly. This option is not available in all solvers,
1269: see table [](#tab:support).
1271: When using two-sided solvers, the problem matrices must have both the
1272: `MATOP_MULT` and `MATOP_MULT_TRANSPOSE` operations defined.
1274: Level: advanced
1276: .seealso: [](ch:eps), `EPSGetTwoSided()`, `EPSGetLeftEigenvector()`
1277: @*/
1278: PetscErrorCode EPSSetTwoSided(EPS eps,PetscBool twosided)
1279: {
1280: PetscFunctionBegin;
1283: if (twosided!=eps->twosided) {
1284: eps->twosided = twosided;
1285: eps->state = EPS_STATE_INITIAL;
1286: }
1287: PetscFunctionReturn(PETSC_SUCCESS);
1288: }
1290: /*@
1291: EPSGetTwoSided - Returns the flag indicating whether a two-sided variant
1292: of the algorithm is being used or not.
1294: Not Collective
1296: Input Parameter:
1297: . eps - the linear eigensolver context
1299: Output Parameter:
1300: . twosided - the returned flag
1302: Level: advanced
1304: .seealso: [](ch:eps), `EPSSetTwoSided()`
1305: @*/
1306: PetscErrorCode EPSGetTwoSided(EPS eps,PetscBool *twosided)
1307: {
1308: PetscFunctionBegin;
1310: PetscAssertPointer(twosided,2);
1311: *twosided = eps->twosided;
1312: PetscFunctionReturn(PETSC_SUCCESS);
1313: }
1315: /*@
1316: EPSSetTrueResidual - Specifies if the solver must compute the true residual
1317: explicitly or not.
1319: Logically Collective
1321: Input Parameters:
1322: + eps - the linear eigensolver context
1323: - trueres - whether true residuals are required or not
1325: Options Database Key:
1326: . -eps_true_residual (true|false) - toggles the true residual
1328: Notes:
1329: If the user sets `trueres`=`PETSC_TRUE` then the solver explicitly computes
1330: the true residual for each eigenpair approximation, and uses it for
1331: convergence testing. Computing the residual is usually an expensive
1332: operation. Some solvers (e.g., Krylov solvers) can avoid this computation
1333: by using a cheap estimate of the residual norm, but this may sometimes
1334: give inaccurate results (especially if a spectral transform is being
1335: used). On the contrary, preconditioned eigensolvers (e.g., Davidson solvers)
1336: do rely on computing the true residual, so this option is irrelevant for them.
1338: Level: advanced
1340: .seealso: [](ch:eps), `EPSGetTrueResidual()`
1341: @*/
1342: PetscErrorCode EPSSetTrueResidual(EPS eps,PetscBool trueres)
1343: {
1344: PetscFunctionBegin;
1347: eps->trueres = trueres;
1348: PetscFunctionReturn(PETSC_SUCCESS);
1349: }
1351: /*@
1352: EPSGetTrueResidual - Returns the flag indicating whether true
1353: residuals must be computed explicitly or not.
1355: Not Collective
1357: Input Parameter:
1358: . eps - the linear eigensolver context
1360: Output Parameter:
1361: . trueres - the returned flag
1363: Level: advanced
1365: .seealso: [](ch:eps), `EPSSetTrueResidual()`
1366: @*/
1367: PetscErrorCode EPSGetTrueResidual(EPS eps,PetscBool *trueres)
1368: {
1369: PetscFunctionBegin;
1371: PetscAssertPointer(trueres,2);
1372: *trueres = eps->trueres;
1373: PetscFunctionReturn(PETSC_SUCCESS);
1374: }
1376: /*@
1377: EPSSetTrackAll - Specifies if the solver must compute the residual norm of all
1378: approximate eigenpairs or not.
1380: Logically Collective
1382: Input Parameters:
1383: + eps - the linear eigensolver context
1384: - trackall - whether to compute all residuals or not
1386: Notes:
1387: If the user sets `trackall`=`PETSC_TRUE` then the solver computes (or estimates)
1388: the residual norm for each eigenpair approximation. Computing the residual is
1389: usually an expensive operation and solvers commonly compute only the residual
1390: associated to the first unconverged eigenpair.
1392: The option `-eps_monitor_all` automatically activates this option.
1394: Level: developer
1396: .seealso: [](ch:eps), `EPSGetTrackAll()`
1397: @*/
1398: PetscErrorCode EPSSetTrackAll(EPS eps,PetscBool trackall)
1399: {
1400: PetscFunctionBegin;
1403: eps->trackall = trackall;
1404: PetscFunctionReturn(PETSC_SUCCESS);
1405: }
1407: /*@
1408: EPSGetTrackAll - Returns the flag indicating whether all residual norms must
1409: be computed or not.
1411: Not Collective
1413: Input Parameter:
1414: . eps - the linear eigensolver context
1416: Output Parameter:
1417: . trackall - the returned flag
1419: Level: developer
1421: .seealso: [](ch:eps), `EPSSetTrackAll()`
1422: @*/
1423: PetscErrorCode EPSGetTrackAll(EPS eps,PetscBool *trackall)
1424: {
1425: PetscFunctionBegin;
1427: PetscAssertPointer(trackall,2);
1428: *trackall = eps->trackall;
1429: PetscFunctionReturn(PETSC_SUCCESS);
1430: }
1432: /*@
1433: EPSSetPurify - Disable eigenvector purification (which is enabled by default).
1435: Logically Collective
1437: Input Parameters:
1438: + eps - the linear eigensolver context
1439: - purify - whether purification is done or not, use `PETSC_FALSE` to disable it
1441: Options Database Key:
1442: . -eps_purify (true|false) - toggles the purification flag
1444: Notes:
1445: By default, eigenvectors of generalized symmetric eigenproblems are purified
1446: in order to purge directions in the nullspace of matrix $B$. If the user knows
1447: that $B$ is non-singular, then purification can be safely deactivated and some
1448: computational cost is avoided (this is particularly important in interval computations).
1450: More details are given in section [](#sec:purif).
1452: Level: intermediate
1454: .seealso: [](ch:eps), [](#sec:purif), `EPSGetPurify()`, `EPSSetInterval()`
1455: @*/
1456: PetscErrorCode EPSSetPurify(EPS eps,PetscBool purify)
1457: {
1458: PetscFunctionBegin;
1461: if (purify!=eps->purify) {
1462: eps->purify = purify;
1463: eps->state = EPS_STATE_INITIAL;
1464: }
1465: PetscFunctionReturn(PETSC_SUCCESS);
1466: }
1468: /*@
1469: EPSGetPurify - Returns the flag indicating whether purification is activated
1470: or not.
1472: Not Collective
1474: Input Parameter:
1475: . eps - the linear eigensolver context
1477: Output Parameter:
1478: . purify - the returned flag
1480: Level: intermediate
1482: .seealso: [](ch:eps), `EPSSetPurify()`
1483: @*/
1484: PetscErrorCode EPSGetPurify(EPS eps,PetscBool *purify)
1485: {
1486: PetscFunctionBegin;
1488: PetscAssertPointer(purify,2);
1489: *purify = eps->purify;
1490: PetscFunctionReturn(PETSC_SUCCESS);
1491: }
1493: /*@
1494: EPSSetOptionsPrefix - Sets the prefix used for searching for all
1495: `EPS` options in the database.
1497: Logically Collective
1499: Input Parameters:
1500: + eps - the linear eigensolver context
1501: - prefix - the prefix string to prepend to all `EPS` option requests
1503: Notes:
1504: A hyphen (-) must NOT be given at the beginning of the prefix name.
1505: The first character of all runtime options is AUTOMATICALLY the
1506: hyphen.
1508: For example, to distinguish between the runtime options for two
1509: different `EPS` contexts, one could call
1510: .vb
1511: EPSSetOptionsPrefix(eps1,"eig1_")
1512: EPSSetOptionsPrefix(eps2,"eig2_")
1513: .ve
1515: Level: advanced
1517: .seealso: [](ch:eps), `EPSAppendOptionsPrefix()`, `EPSGetOptionsPrefix()`
1518: @*/
1519: PetscErrorCode EPSSetOptionsPrefix(EPS eps,const char prefix[])
1520: {
1521: PetscFunctionBegin;
1523: if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
1524: PetscCall(STSetOptionsPrefix(eps->st,prefix));
1525: if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
1526: PetscCall(BVSetOptionsPrefix(eps->V,prefix));
1527: if (!eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
1528: PetscCall(DSSetOptionsPrefix(eps->ds,prefix));
1529: if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
1530: PetscCall(RGSetOptionsPrefix(eps->rg,prefix));
1531: PetscCall(PetscObjectSetOptionsPrefix((PetscObject)eps,prefix));
1532: PetscFunctionReturn(PETSC_SUCCESS);
1533: }
1535: /*@
1536: EPSAppendOptionsPrefix - Appends to the prefix used for searching for all
1537: `EPS` options in the database.
1539: Logically Collective
1541: Input Parameters:
1542: + eps - the linear eigensolver context
1543: - prefix - the prefix string to prepend to all `EPS` option requests
1545: Notes:
1546: A hyphen (-) must NOT be given at the beginning of the prefix name.
1547: The first character of all runtime options is AUTOMATICALLY the hyphen.
1549: Level: advanced
1551: .seealso: [](ch:eps), `EPSSetOptionsPrefix()`, `EPSGetOptionsPrefix()`
1552: @*/
1553: PetscErrorCode EPSAppendOptionsPrefix(EPS eps,const char prefix[])
1554: {
1555: PetscFunctionBegin;
1557: if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
1558: PetscCall(STAppendOptionsPrefix(eps->st,prefix));
1559: if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
1560: PetscCall(BVAppendOptionsPrefix(eps->V,prefix));
1561: if (!eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
1562: PetscCall(DSAppendOptionsPrefix(eps->ds,prefix));
1563: if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
1564: PetscCall(RGAppendOptionsPrefix(eps->rg,prefix));
1565: PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)eps,prefix));
1566: PetscFunctionReturn(PETSC_SUCCESS);
1567: }
1569: /*@
1570: EPSGetOptionsPrefix - Gets the prefix used for searching for all
1571: `EPS` options in the database.
1573: Not Collective
1575: Input Parameter:
1576: . eps - the linear eigensolver context
1578: Output Parameter:
1579: . prefix - pointer to the prefix string used is returned
1581: Level: advanced
1583: .seealso: [](ch:eps), `EPSSetOptionsPrefix()`, `EPSAppendOptionsPrefix()`
1584: @*/
1585: PetscErrorCode EPSGetOptionsPrefix(EPS eps,const char *prefix[])
1586: {
1587: PetscFunctionBegin;
1589: PetscAssertPointer(prefix,2);
1590: PetscCall(PetscObjectGetOptionsPrefix((PetscObject)eps,prefix));
1591: PetscFunctionReturn(PETSC_SUCCESS);
1592: }