1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2021, 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 on eps
24: Input Parameters:
25: + eps - the 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: EPSMonitorSet(), EPSSetTrackAll()
34: @*/
35: PetscErrorCode EPSMonitorSetFromOptions(EPS eps,const char opt[],const char name[],void *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;
46: PetscErrorCode ierr;
49: PetscOptionsGetViewer(PetscObjectComm((PetscObject)eps),((PetscObject)eps)->options,((PetscObject)eps)->prefix,opt,&viewer,&format,&flg);
50: if (!flg) return(0);
52: PetscViewerGetType(viewer,&vtype);
53: SlepcMonitorMakeKey_Internal(name,vtype,format,key);
54: PetscFunctionListFind(EPSMonitorList,key,&mfunc);
55: PetscFunctionListFind(EPSMonitorCreateList,key,&cfunc);
56: PetscFunctionListFind(EPSMonitorDestroyList,key,&dfunc);
57: if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
58: if (!dfunc) dfunc = PetscViewerAndFormatDestroy;
60: (*cfunc)(viewer,format,ctx,&vf);
61: PetscObjectDereference((PetscObject)viewer);
62: EPSMonitorSet(eps,mfunc,vf,(PetscErrorCode(*)(void **))dfunc);
63: if (trackall) {
64: EPSSetTrackAll(eps,PETSC_TRUE);
65: }
66: return(0);
67: }
69: /*@
70: EPSSetFromOptions - Sets EPS options from the options database.
71: This routine must be called before EPSSetUp() if the user is to be
72: allowed to set the solver type.
74: Collective on eps
76: Input Parameters:
77: . eps - the eigensolver context
79: Notes:
80: To see all options, run your program with the -help option.
82: Level: beginner
83: @*/
84: PetscErrorCode EPSSetFromOptions(EPS eps) 85: {
87: char type[256];
88: PetscBool set,flg,flg1,flg2,flg3,bval;
89: PetscReal r,array[2]={0,0};
90: PetscScalar s;
91: PetscInt i,j,k;
92: EPSBalance bal;
96: EPSRegisterAll();
97: PetscObjectOptionsBegin((PetscObject)eps);
98: PetscOptionsFList("-eps_type","Eigensolver method","EPSSetType",EPSList,(char*)(((PetscObject)eps)->type_name?((PetscObject)eps)->type_name:EPSKRYLOVSCHUR),type,sizeof(type),&flg);
99: if (flg) {
100: EPSSetType(eps,type);
101: } else if (!((PetscObject)eps)->type_name) {
102: EPSSetType(eps,EPSKRYLOVSCHUR);
103: }
105: PetscOptionsBoolGroupBegin("-eps_hermitian","Hermitian eigenvalue problem","EPSSetProblemType",&flg);
106: if (flg) { EPSSetProblemType(eps,EPS_HEP); }
107: PetscOptionsBoolGroup("-eps_gen_hermitian","Generalized Hermitian eigenvalue problem","EPSSetProblemType",&flg);
108: if (flg) { EPSSetProblemType(eps,EPS_GHEP); }
109: PetscOptionsBoolGroup("-eps_non_hermitian","Non-Hermitian eigenvalue problem","EPSSetProblemType",&flg);
110: if (flg) { EPSSetProblemType(eps,EPS_NHEP); }
111: PetscOptionsBoolGroup("-eps_gen_non_hermitian","Generalized non-Hermitian eigenvalue problem","EPSSetProblemType",&flg);
112: if (flg) { EPSSetProblemType(eps,EPS_GNHEP); }
113: PetscOptionsBoolGroup("-eps_pos_gen_non_hermitian","Generalized non-Hermitian eigenvalue problem with positive semi-definite B","EPSSetProblemType",&flg);
114: if (flg) { EPSSetProblemType(eps,EPS_PGNHEP); }
115: PetscOptionsBoolGroupEnd("-eps_gen_indefinite","Generalized Hermitian-indefinite eigenvalue problem","EPSSetProblemType",&flg);
116: if (flg) { EPSSetProblemType(eps,EPS_GHIEP); }
118: PetscOptionsBoolGroupBegin("-eps_ritz","Rayleigh-Ritz extraction","EPSSetExtraction",&flg);
119: if (flg) { EPSSetExtraction(eps,EPS_RITZ); }
120: PetscOptionsBoolGroup("-eps_harmonic","Harmonic Ritz extraction","EPSSetExtraction",&flg);
121: if (flg) { EPSSetExtraction(eps,EPS_HARMONIC); }
122: PetscOptionsBoolGroup("-eps_harmonic_relative","Relative harmonic Ritz extraction","EPSSetExtraction",&flg);
123: if (flg) { EPSSetExtraction(eps,EPS_HARMONIC_RELATIVE); }
124: PetscOptionsBoolGroup("-eps_harmonic_right","Right harmonic Ritz extraction","EPSSetExtraction",&flg);
125: if (flg) { EPSSetExtraction(eps,EPS_HARMONIC_RIGHT); }
126: PetscOptionsBoolGroup("-eps_harmonic_largest","Largest harmonic Ritz extraction","EPSSetExtraction",&flg);
127: if (flg) { EPSSetExtraction(eps,EPS_HARMONIC_LARGEST); }
128: PetscOptionsBoolGroup("-eps_refined","Refined Ritz extraction","EPSSetExtraction",&flg);
129: if (flg) { EPSSetExtraction(eps,EPS_REFINED); }
130: PetscOptionsBoolGroupEnd("-eps_refined_harmonic","Refined harmonic Ritz extraction","EPSSetExtraction",&flg);
131: if (flg) { EPSSetExtraction(eps,EPS_REFINED_HARMONIC); }
133: bal = eps->balance;
134: PetscOptionsEnum("-eps_balance","Balancing method","EPSSetBalance",EPSBalanceTypes,(PetscEnum)bal,(PetscEnum*)&bal,&flg1);
135: j = eps->balance_its;
136: PetscOptionsInt("-eps_balance_its","Number of iterations in balancing","EPSSetBalance",eps->balance_its,&j,&flg2);
137: r = eps->balance_cutoff;
138: PetscOptionsReal("-eps_balance_cutoff","Cutoff value in balancing","EPSSetBalance",eps->balance_cutoff,&r,&flg3);
139: if (flg1 || flg2 || flg3) { EPSSetBalance(eps,bal,j,r); }
141: i = eps->max_it;
142: PetscOptionsInt("-eps_max_it","Maximum number of iterations","EPSSetTolerances",eps->max_it,&i,&flg1);
143: r = eps->tol;
144: PetscOptionsReal("-eps_tol","Tolerance","EPSSetTolerances",eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol,&r,&flg2);
145: if (flg1 || flg2) { EPSSetTolerances(eps,r,i); }
147: PetscOptionsBoolGroupBegin("-eps_conv_rel","Relative error convergence test","EPSSetConvergenceTest",&flg);
148: if (flg) { EPSSetConvergenceTest(eps,EPS_CONV_REL); }
149: PetscOptionsBoolGroup("-eps_conv_norm","Convergence test relative to the eigenvalue and the matrix norms","EPSSetConvergenceTest",&flg);
150: if (flg) { EPSSetConvergenceTest(eps,EPS_CONV_NORM); }
151: PetscOptionsBoolGroup("-eps_conv_abs","Absolute error convergence test","EPSSetConvergenceTest",&flg);
152: if (flg) { EPSSetConvergenceTest(eps,EPS_CONV_ABS); }
153: PetscOptionsBoolGroupEnd("-eps_conv_user","User-defined convergence test","EPSSetConvergenceTest",&flg);
154: if (flg) { EPSSetConvergenceTest(eps,EPS_CONV_USER); }
156: PetscOptionsBoolGroupBegin("-eps_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","EPSSetStoppingTest",&flg);
157: if (flg) { EPSSetStoppingTest(eps,EPS_STOP_BASIC); }
158: PetscOptionsBoolGroupEnd("-eps_stop_user","User-defined stopping test","EPSSetStoppingTest",&flg);
159: if (flg) { EPSSetStoppingTest(eps,EPS_STOP_USER); }
161: i = eps->nev;
162: PetscOptionsInt("-eps_nev","Number of eigenvalues to compute","EPSSetDimensions",eps->nev,&i,&flg1);
163: j = eps->ncv;
164: PetscOptionsInt("-eps_ncv","Number of basis vectors","EPSSetDimensions",eps->ncv,&j,&flg2);
165: k = eps->mpd;
166: PetscOptionsInt("-eps_mpd","Maximum dimension of projected problem","EPSSetDimensions",eps->mpd,&k,&flg3);
167: if (flg1 || flg2 || flg3) { EPSSetDimensions(eps,i,j,k); }
169: PetscOptionsBoolGroupBegin("-eps_largest_magnitude","Compute largest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg);
170: if (flg) { EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE); }
171: PetscOptionsBoolGroup("-eps_smallest_magnitude","Compute smallest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg);
172: if (flg) { EPSSetWhichEigenpairs(eps,EPS_SMALLEST_MAGNITUDE); }
173: PetscOptionsBoolGroup("-eps_largest_real","Compute eigenvalues with largest real parts","EPSSetWhichEigenpairs",&flg);
174: if (flg) { EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL); }
175: PetscOptionsBoolGroup("-eps_smallest_real","Compute eigenvalues with smallest real parts","EPSSetWhichEigenpairs",&flg);
176: if (flg) { EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL); }
177: PetscOptionsBoolGroup("-eps_largest_imaginary","Compute eigenvalues with largest imaginary parts","EPSSetWhichEigenpairs",&flg);
178: if (flg) { EPSSetWhichEigenpairs(eps,EPS_LARGEST_IMAGINARY); }
179: PetscOptionsBoolGroup("-eps_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","EPSSetWhichEigenpairs",&flg);
180: if (flg) { EPSSetWhichEigenpairs(eps,EPS_SMALLEST_IMAGINARY); }
181: PetscOptionsBoolGroup("-eps_target_magnitude","Compute eigenvalues closest to target","EPSSetWhichEigenpairs",&flg);
182: if (flg) { EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE); }
183: PetscOptionsBoolGroup("-eps_target_real","Compute eigenvalues with real parts closest to target","EPSSetWhichEigenpairs",&flg);
184: if (flg) { EPSSetWhichEigenpairs(eps,EPS_TARGET_REAL); }
185: PetscOptionsBoolGroup("-eps_target_imaginary","Compute eigenvalues with imaginary parts closest to target","EPSSetWhichEigenpairs",&flg);
186: if (flg) { EPSSetWhichEigenpairs(eps,EPS_TARGET_IMAGINARY); }
187: PetscOptionsBoolGroupEnd("-eps_all","Compute all eigenvalues in an interval or a region","EPSSetWhichEigenpairs",&flg);
188: if (flg) { EPSSetWhichEigenpairs(eps,EPS_ALL); }
190: PetscOptionsScalar("-eps_target","Value of the target","EPSSetTarget",eps->target,&s,&flg);
191: if (flg) {
192: if (eps->which!=EPS_TARGET_REAL && eps->which!=EPS_TARGET_IMAGINARY) {
193: EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
194: }
195: EPSSetTarget(eps,s);
196: }
198: k = 2;
199: PetscOptionsRealArray("-eps_interval","Computational interval (two real values separated with a comma without spaces)","EPSSetInterval",array,&k,&flg);
200: if (flg) {
201: if (k<2) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_SIZ,"Must pass two values in -eps_interval (comma-separated without spaces)");
202: EPSSetWhichEigenpairs(eps,EPS_ALL);
203: EPSSetInterval(eps,array[0],array[1]);
204: }
206: PetscOptionsBool("-eps_true_residual","Compute true residuals explicitly","EPSSetTrueResidual",eps->trueres,&eps->trueres,NULL);
207: PetscOptionsBool("-eps_purify","Postprocess eigenvectors for purification","EPSSetPurify",eps->purify,&bval,&flg);
208: if (flg) { EPSSetPurify(eps,bval); }
209: PetscOptionsBool("-eps_two_sided","Use two-sided variant (to compute left eigenvectors)","EPSSetTwoSided",eps->twosided,&bval,&flg);
210: if (flg) { EPSSetTwoSided(eps,bval); }
212: /* -----------------------------------------------------------------------*/
213: /*
214: Cancels all monitors hardwired into code before call to EPSSetFromOptions()
215: */
216: PetscOptionsBool("-eps_monitor_cancel","Remove any hardwired monitor routines","EPSMonitorCancel",PETSC_FALSE,&flg,&set);
217: if (set && flg) { EPSMonitorCancel(eps); }
218: EPSMonitorSetFromOptions(eps,"-eps_monitor","first_approximation",NULL,PETSC_FALSE);
219: EPSMonitorSetFromOptions(eps,"-eps_monitor_all","all_approximations",NULL,PETSC_TRUE);
220: EPSMonitorSetFromOptions(eps,"-eps_monitor_conv","convergence_history",NULL,PETSC_FALSE);
222: /* -----------------------------------------------------------------------*/
223: PetscOptionsName("-eps_view","Print detailed information on solver used","EPSView",NULL);
224: PetscOptionsName("-eps_view_vectors","View computed eigenvectors","EPSVectorsView",NULL);
225: PetscOptionsName("-eps_view_values","View computed eigenvalues","EPSValuesView",NULL);
226: PetscOptionsName("-eps_converged_reason","Print reason for convergence, and number of iterations","EPSConvergedReasonView",NULL);
227: PetscOptionsName("-eps_error_absolute","Print absolute errors of each eigenpair","EPSErrorView",NULL);
228: PetscOptionsName("-eps_error_relative","Print relative errors of each eigenpair","EPSErrorView",NULL);
229: PetscOptionsName("-eps_error_backward","Print backward errors of each eigenpair","EPSErrorView",NULL);
231: if (eps->ops->setfromoptions) {
232: (*eps->ops->setfromoptions)(PetscOptionsObject,eps);
233: }
234: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)eps);
235: PetscOptionsEnd();
237: if (!eps->V) { EPSGetBV(eps,&eps->V); }
238: BVSetFromOptions(eps->V);
239: if (!eps->rg) { EPSGetRG(eps,&eps->rg); }
240: RGSetFromOptions(eps->rg);
241: if (eps->useds) {
242: if (!eps->ds) { EPSGetDS(eps,&eps->ds); }
243: DSSetFromOptions(eps->ds);
244: }
245: if (!eps->st) { EPSGetST(eps,&eps->st); }
246: EPSSetDefaultST(eps);
247: STSetFromOptions(eps->st);
248: return(0);
249: }
251: /*@C
252: EPSGetTolerances - Gets the tolerance and maximum iteration count used
253: by the EPS convergence tests.
255: Not Collective
257: Input Parameter:
258: . eps - the eigensolver context
260: Output Parameters:
261: + tol - the convergence tolerance
262: - maxits - maximum number of iterations
264: Notes:
265: The user can specify NULL for any parameter that is not needed.
267: Level: intermediate
269: .seealso: EPSSetTolerances()
270: @*/
271: PetscErrorCode EPSGetTolerances(EPS eps,PetscReal *tol,PetscInt *maxits)272: {
275: if (tol) *tol = eps->tol;
276: if (maxits) *maxits = eps->max_it;
277: return(0);
278: }
280: /*@
281: EPSSetTolerances - Sets the tolerance and maximum iteration count used
282: by the EPS convergence tests.
284: Logically Collective on eps
286: Input Parameters:
287: + eps - the eigensolver context
288: . tol - the convergence tolerance
289: - maxits - maximum number of iterations to use
291: Options Database Keys:
292: + -eps_tol <tol> - Sets the convergence tolerance
293: - -eps_max_it <maxits> - Sets the maximum number of iterations allowed
295: Notes:
296: Use PETSC_DEFAULT for either argument to assign a reasonably good value.
298: Level: intermediate
300: .seealso: EPSGetTolerances()
301: @*/
302: PetscErrorCode EPSSetTolerances(EPS eps,PetscReal tol,PetscInt maxits)303: {
308: if (tol == PETSC_DEFAULT) {
309: eps->tol = PETSC_DEFAULT;
310: eps->state = EPS_STATE_INITIAL;
311: } else {
312: if (tol <= 0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
313: eps->tol = tol;
314: }
315: if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
316: eps->max_it = PETSC_DEFAULT;
317: eps->state = EPS_STATE_INITIAL;
318: } else {
319: if (maxits <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
320: eps->max_it = maxits;
321: }
322: return(0);
323: }
325: /*@C
326: EPSGetDimensions - Gets the number of eigenvalues to compute
327: and the dimension of the subspace.
329: Not Collective
331: Input Parameter:
332: . eps - the eigensolver context
334: Output Parameters:
335: + nev - number of eigenvalues to compute
336: . ncv - the maximum dimension of the subspace to be used by the solver
337: - mpd - the maximum dimension allowed for the projected problem
339: Level: intermediate
341: .seealso: EPSSetDimensions()
342: @*/
343: PetscErrorCode EPSGetDimensions(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)344: {
347: if (nev) *nev = eps->nev;
348: if (ncv) *ncv = eps->ncv;
349: if (mpd) *mpd = eps->mpd;
350: return(0);
351: }
353: /*@
354: EPSSetDimensions - Sets the number of eigenvalues to compute
355: and the dimension of the subspace.
357: Logically Collective on eps
359: Input Parameters:
360: + eps - the eigensolver context
361: . nev - number of eigenvalues to compute
362: . ncv - the maximum dimension of the subspace to be used by the solver
363: - mpd - the maximum dimension allowed for the projected problem
365: Options Database Keys:
366: + -eps_nev <nev> - Sets the number of eigenvalues
367: . -eps_ncv <ncv> - Sets the dimension of the subspace
368: - -eps_mpd <mpd> - Sets the maximum projected dimension
370: Notes:
371: Use PETSC_DEFAULT for ncv and mpd to assign a reasonably good value, which is
372: dependent on the solution method.
374: The parameters ncv and mpd are intimately related, so that the user is advised
375: to set one of them at most. Normal usage is that
376: (a) in cases where nev is small, the user sets ncv (a reasonable default is 2*nev); and
377: (b) in cases where nev is large, the user sets mpd.
379: The value of ncv should always be between nev and (nev+mpd), typically
380: ncv=nev+mpd. If nev is not too large, mpd=nev is a reasonable choice, otherwise
381: a smaller value should be used.
383: When computing all eigenvalues in an interval, see EPSSetInterval(), these
384: parameters lose relevance, and tuning must be done with
385: EPSKrylovSchurSetDimensions().
387: Level: intermediate
389: .seealso: EPSGetDimensions(), EPSSetInterval(), EPSKrylovSchurSetDimensions()
390: @*/
391: PetscErrorCode EPSSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)392: {
398: if (nev<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
399: eps->nev = nev;
400: if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
401: eps->ncv = PETSC_DEFAULT;
402: } else {
403: if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
404: eps->ncv = ncv;
405: }
406: if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
407: eps->mpd = PETSC_DEFAULT;
408: } else {
409: if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
410: eps->mpd = mpd;
411: }
412: eps->state = EPS_STATE_INITIAL;
413: return(0);
414: }
416: /*@
417: EPSSetWhichEigenpairs - Specifies which portion of the spectrum is
418: to be sought.
420: Logically Collective on eps
422: Input Parameters:
423: + eps - eigensolver context obtained from EPSCreate()
424: - which - the portion of the spectrum to be sought
426: Possible values:
427: The parameter 'which' can have one of these values
429: + EPS_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
430: . EPS_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
431: . EPS_LARGEST_REAL - largest real parts
432: . EPS_SMALLEST_REAL - smallest real parts
433: . EPS_LARGEST_IMAGINARY - largest imaginary parts
434: . EPS_SMALLEST_IMAGINARY - smallest imaginary parts
435: . EPS_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
436: . EPS_TARGET_REAL - eigenvalues with real part closest to target
437: . EPS_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
438: . EPS_ALL - all eigenvalues contained in a given interval or region
439: - EPS_WHICH_USER - user defined ordering set with EPSSetEigenvalueComparison()
441: Options Database Keys:
442: + -eps_largest_magnitude - Sets largest eigenvalues in magnitude
443: . -eps_smallest_magnitude - Sets smallest eigenvalues in magnitude
444: . -eps_largest_real - Sets largest real parts
445: . -eps_smallest_real - Sets smallest real parts
446: . -eps_largest_imaginary - Sets largest imaginary parts
447: . -eps_smallest_imaginary - Sets smallest imaginary parts
448: . -eps_target_magnitude - Sets eigenvalues closest to target
449: . -eps_target_real - Sets real parts closest to target
450: . -eps_target_imaginary - Sets imaginary parts closest to target
451: - -eps_all - Sets all eigenvalues in an interval or region
453: Notes:
454: Not all eigensolvers implemented in EPS account for all the possible values
455: stated above. Also, some values make sense only for certain types of
456: problems. If SLEPc is compiled for real numbers EPS_LARGEST_IMAGINARY457: and EPS_SMALLEST_IMAGINARY use the absolute value of the imaginary part
458: for eigenvalue selection.
460: The target is a scalar value provided with EPSSetTarget().
462: The criterion EPS_TARGET_IMAGINARY is available only in case PETSc and
463: SLEPc have been built with complex scalars.
465: EPS_ALL is intended for use in combination with an interval (see
466: EPSSetInterval()), when all eigenvalues within the interval are requested,
467: or in the context of the CISS solver for computing all eigenvalues in a region.
468: In those cases, the number of eigenvalues is unknown, so the nev parameter
469: has a different sense, see EPSSetDimensions().
471: Level: intermediate
473: .seealso: EPSGetWhichEigenpairs(), EPSSetTarget(), EPSSetInterval(),
474: EPSSetDimensions(), EPSSetEigenvalueComparison(), EPSWhich475: @*/
476: PetscErrorCode EPSSetWhichEigenpairs(EPS eps,EPSWhich which)477: {
481: switch (which) {
482: case EPS_LARGEST_MAGNITUDE:
483: case EPS_SMALLEST_MAGNITUDE:
484: case EPS_LARGEST_REAL:
485: case EPS_SMALLEST_REAL:
486: case EPS_LARGEST_IMAGINARY:
487: case EPS_SMALLEST_IMAGINARY:
488: case EPS_TARGET_MAGNITUDE:
489: case EPS_TARGET_REAL:
490: #if defined(PETSC_USE_COMPLEX)
491: case EPS_TARGET_IMAGINARY:
492: #endif
493: case EPS_ALL:
494: case EPS_WHICH_USER:
495: if (eps->which != which) {
496: eps->state = EPS_STATE_INITIAL;
497: eps->which = which;
498: }
499: break;
500: #if !defined(PETSC_USE_COMPLEX)
501: case EPS_TARGET_IMAGINARY:
502: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"EPS_TARGET_IMAGINARY can be used only with complex scalars");
503: #endif
504: default:505: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
506: }
507: return(0);
508: }
510: /*@
511: EPSGetWhichEigenpairs - Returns which portion of the spectrum is to be
512: sought.
514: Not Collective
516: Input Parameter:
517: . eps - eigensolver context obtained from EPSCreate()
519: Output Parameter:
520: . which - the portion of the spectrum to be sought
522: Notes:
523: See EPSSetWhichEigenpairs() for possible values of 'which'.
525: Level: intermediate
527: .seealso: EPSSetWhichEigenpairs(), EPSWhich528: @*/
529: PetscErrorCode EPSGetWhichEigenpairs(EPS eps,EPSWhich *which)530: {
534: *which = eps->which;
535: return(0);
536: }
538: /*@C
539: EPSSetEigenvalueComparison - Specifies the eigenvalue comparison function
540: when EPSSetWhichEigenpairs() is set to EPS_WHICH_USER.
542: Logically Collective on eps
544: Input Parameters:
545: + eps - eigensolver context obtained from EPSCreate()
546: . func - a pointer to the comparison function
547: - ctx - a context pointer (the last parameter to the comparison function)
549: Calling Sequence of func:
550: $ func(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)
552: + ar - real part of the 1st eigenvalue
553: . ai - imaginary part of the 1st eigenvalue
554: . br - real part of the 2nd eigenvalue
555: . bi - imaginary part of the 2nd eigenvalue
556: . res - result of comparison
557: - ctx - optional context, as set by EPSSetEigenvalueComparison()
559: Note:
560: The returning parameter 'res' can be
561: + negative - if the 1st eigenvalue is preferred to the 2st one
562: . zero - if both eigenvalues are equally preferred
563: - positive - if the 2st eigenvalue is preferred to the 1st one
565: Level: advanced
567: .seealso: EPSSetWhichEigenpairs(), EPSWhich568: @*/
569: PetscErrorCode EPSSetEigenvalueComparison(EPS eps,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)570: {
573: eps->sc->comparison = func;
574: eps->sc->comparisonctx = ctx;
575: eps->which = EPS_WHICH_USER;
576: return(0);
577: }
579: /*@C
580: EPSSetArbitrarySelection - Specifies a function intended to look for
581: eigenvalues according to an arbitrary selection criterion. This criterion
582: can be based on a computation involving the current eigenvector approximation.
584: Logically Collective on eps
586: Input Parameters:
587: + eps - eigensolver context obtained from EPSCreate()
588: . func - a pointer to the evaluation function
589: - ctx - a context pointer (the last parameter to the evaluation function)
591: Calling Sequence of func:
592: $ func(PetscScalar er,PetscScalar ei,Vec xr,Vec xi,PetscScalar *rr,PetscScalar *ri,void *ctx)
594: + er - real part of the current eigenvalue approximation
595: . ei - imaginary part of the current eigenvalue approximation
596: . xr - real part of the current eigenvector approximation
597: . xi - imaginary part of the current eigenvector approximation
598: . rr - result of evaluation (real part)
599: . ri - result of evaluation (imaginary part)
600: - ctx - optional context, as set by EPSSetArbitrarySelection()
602: Notes:
603: This provides a mechanism to select eigenpairs by evaluating a user-defined
604: function. When a function has been provided, the default selection based on
605: sorting the eigenvalues is replaced by the sorting of the results of this
606: function (with the same sorting criterion given in EPSSetWhichEigenpairs()).
608: For instance, suppose you want to compute those eigenvectors that maximize
609: a certain computable expression. Then implement the computation using
610: the arguments xr and xi, and return the result in rr. Then set the standard
611: sorting by magnitude so that the eigenpair with largest value of rr is
612: selected.
614: This evaluation function is collective, that is, all processes call it and
615: it can use collective operations; furthermore, the computed result must
616: be the same in all processes.
618: The result of func is expressed as a complex number so that it is possible to
619: use the standard eigenvalue sorting functions, but normally only rr is used.
620: Set ri to zero unless it is meaningful in your application.
622: Level: advanced
624: .seealso: EPSSetWhichEigenpairs()
625: @*/
626: PetscErrorCode EPSSetArbitrarySelection(EPS eps,PetscErrorCode (*func)(PetscScalar,PetscScalar,Vec,Vec,PetscScalar*,PetscScalar*,void*),void* ctx)627: {
630: eps->arbitrary = func;
631: eps->arbitraryctx = ctx;
632: eps->state = EPS_STATE_INITIAL;
633: return(0);
634: }
636: /*@C
637: EPSSetConvergenceTestFunction - Sets a function to compute the error estimate
638: used in the convergence test.
640: Logically Collective on eps
642: Input Parameters:
643: + eps - eigensolver context obtained from EPSCreate()
644: . func - a pointer to the convergence test function
645: . ctx - context for private data for the convergence routine (may be null)
646: - destroy - a routine for destroying the context (may be null)
648: Calling Sequence of func:
649: $ func(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
651: + eps - eigensolver context obtained from EPSCreate()
652: . eigr - real part of the eigenvalue
653: . eigi - imaginary part of the eigenvalue
654: . res - residual norm associated to the eigenpair
655: . errest - (output) computed error estimate
656: - ctx - optional context, as set by EPSSetConvergenceTestFunction()
658: Note:
659: If the error estimate returned by the convergence test function is less than
660: the tolerance, then the eigenvalue is accepted as converged.
662: Level: advanced
664: .seealso: EPSSetConvergenceTest(), EPSSetTolerances()
665: @*/
666: PetscErrorCode EPSSetConvergenceTestFunction(EPS eps,PetscErrorCode (*func)(EPS,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))667: {
672: if (eps->convergeddestroy) {
673: (*eps->convergeddestroy)(eps->convergedctx);
674: }
675: eps->convergeduser = func;
676: eps->convergeddestroy = destroy;
677: eps->convergedctx = ctx;
678: if (func == EPSConvergedRelative) eps->conv = EPS_CONV_REL;
679: else if (func == EPSConvergedNorm) eps->conv = EPS_CONV_NORM;
680: else if (func == EPSConvergedAbsolute) eps->conv = EPS_CONV_ABS;
681: else {
682: eps->conv = EPS_CONV_USER;
683: eps->converged = eps->convergeduser;
684: }
685: return(0);
686: }
688: /*@
689: EPSSetConvergenceTest - Specifies how to compute the error estimate
690: used in the convergence test.
692: Logically Collective on eps
694: Input Parameters:
695: + eps - eigensolver context obtained from EPSCreate()
696: - conv - the type of convergence test
698: Options Database Keys:
699: + -eps_conv_abs - Sets the absolute convergence test
700: . -eps_conv_rel - Sets the convergence test relative to the eigenvalue
701: . -eps_conv_norm - Sets the convergence test relative to the matrix norms
702: - -eps_conv_user - Selects the user-defined convergence test
704: Note:
705: The parameter 'conv' can have one of these values
706: + EPS_CONV_ABS - absolute error ||r||
707: . EPS_CONV_REL - error relative to the eigenvalue l, ||r||/|l|
708: . EPS_CONV_NORM - error relative to the matrix norms, ||r||/(||A||+|l|*||B||)
709: - EPS_CONV_USER - function set by EPSSetConvergenceTestFunction()
711: Level: intermediate
713: .seealso: EPSGetConvergenceTest(), EPSSetConvergenceTestFunction(), EPSSetStoppingTest(), EPSConv714: @*/
715: PetscErrorCode EPSSetConvergenceTest(EPS eps,EPSConv conv)716: {
720: switch (conv) {
721: case EPS_CONV_ABS: eps->converged = EPSConvergedAbsolute; break;
722: case EPS_CONV_REL: eps->converged = EPSConvergedRelative; break;
723: case EPS_CONV_NORM: eps->converged = EPSConvergedNorm; break;
724: case EPS_CONV_USER:
725: if (!eps->convergeduser) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ORDER,"Must call EPSSetConvergenceTestFunction() first");
726: eps->converged = eps->convergeduser;
727: break;
728: default:729: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
730: }
731: eps->conv = conv;
732: return(0);
733: }
735: /*@
736: EPSGetConvergenceTest - Gets the method used to compute the error estimate
737: used in the convergence test.
739: Not Collective
741: Input Parameters:
742: . eps - eigensolver context obtained from EPSCreate()
744: Output Parameters:
745: . conv - the type of convergence test
747: Level: intermediate
749: .seealso: EPSSetConvergenceTest(), EPSConv750: @*/
751: PetscErrorCode EPSGetConvergenceTest(EPS eps,EPSConv *conv)752: {
756: *conv = eps->conv;
757: return(0);
758: }
760: /*@C
761: EPSSetStoppingTestFunction - Sets a function to decide when to stop the outer
762: iteration of the eigensolver.
764: Logically Collective on eps
766: Input Parameters:
767: + eps - eigensolver context obtained from EPSCreate()
768: . func - pointer to the stopping test function
769: . ctx - context for private data for the stopping routine (may be null)
770: - destroy - a routine for destroying the context (may be null)
772: Calling Sequence of func:
773: $ func(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
775: + eps - eigensolver context obtained from EPSCreate()
776: . its - current number of iterations
777: . max_it - maximum number of iterations
778: . nconv - number of currently converged eigenpairs
779: . nev - number of requested eigenpairs
780: . reason - (output) result of the stopping test
781: - ctx - optional context, as set by EPSSetStoppingTestFunction()
783: Note:
784: Normal usage is to first call the default routine EPSStoppingBasic() and then
785: set reason to EPS_CONVERGED_USER if some user-defined conditions have been
786: met. To let the eigensolver continue iterating, the result must be left as
787: EPS_CONVERGED_ITERATING.
789: Level: advanced
791: .seealso: EPSSetStoppingTest(), EPSStoppingBasic()
792: @*/
793: PetscErrorCode EPSSetStoppingTestFunction(EPS eps,PetscErrorCode (*func)(EPS,PetscInt,PetscInt,PetscInt,PetscInt,EPSConvergedReason*,void*),void* ctx,PetscErrorCode (*destroy)(void*))794: {
799: if (eps->stoppingdestroy) {
800: (*eps->stoppingdestroy)(eps->stoppingctx);
801: }
802: eps->stoppinguser = func;
803: eps->stoppingdestroy = destroy;
804: eps->stoppingctx = ctx;
805: if (func == EPSStoppingBasic) eps->stop = EPS_STOP_BASIC;
806: else {
807: eps->stop = EPS_STOP_USER;
808: eps->stopping = eps->stoppinguser;
809: }
810: return(0);
811: }
813: /*@
814: EPSSetStoppingTest - Specifies how to decide the termination of the outer
815: loop of the eigensolver.
817: Logically Collective on eps
819: Input Parameters:
820: + eps - eigensolver context obtained from EPSCreate()
821: - stop - the type of stopping test
823: Options Database Keys:
824: + -eps_stop_basic - Sets the default stopping test
825: - -eps_stop_user - Selects the user-defined stopping test
827: Note:
828: The parameter 'stop' can have one of these values
829: + EPS_STOP_BASIC - default stopping test
830: - EPS_STOP_USER - function set by EPSSetStoppingTestFunction()
832: Level: advanced
834: .seealso: EPSGetStoppingTest(), EPSSetStoppingTestFunction(), EPSSetConvergenceTest(), EPSStop835: @*/
836: PetscErrorCode EPSSetStoppingTest(EPS eps,EPSStop stop)837: {
841: switch (stop) {
842: case EPS_STOP_BASIC: eps->stopping = EPSStoppingBasic; break;
843: case EPS_STOP_USER:
844: if (!eps->stoppinguser) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ORDER,"Must call EPSSetStoppingTestFunction() first");
845: eps->stopping = eps->stoppinguser;
846: break;
847: default:848: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
849: }
850: eps->stop = stop;
851: return(0);
852: }
854: /*@
855: EPSGetStoppingTest - Gets the method used to decide the termination of the outer
856: loop of the eigensolver.
858: Not Collective
860: Input Parameters:
861: . eps - eigensolver context obtained from EPSCreate()
863: Output Parameters:
864: . stop - the type of stopping test
866: Level: advanced
868: .seealso: EPSSetStoppingTest(), EPSStop869: @*/
870: PetscErrorCode EPSGetStoppingTest(EPS eps,EPSStop *stop)871: {
875: *stop = eps->stop;
876: return(0);
877: }
879: /*@
880: EPSSetProblemType - Specifies the type of the eigenvalue problem.
882: Logically Collective on eps
884: Input Parameters:
885: + eps - the eigensolver context
886: - type - a known type of eigenvalue problem
888: Options Database Keys:
889: + -eps_hermitian - Hermitian eigenvalue problem
890: . -eps_gen_hermitian - generalized Hermitian eigenvalue problem
891: . -eps_non_hermitian - non-Hermitian eigenvalue problem
892: . -eps_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
893: - -eps_pos_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
894: with positive semi-definite B
896: Notes:
897: Allowed values for the problem type are: Hermitian (EPS_HEP), non-Hermitian
898: (EPS_NHEP), generalized Hermitian (EPS_GHEP), generalized non-Hermitian
899: (EPS_GNHEP), generalized non-Hermitian with positive semi-definite B
900: (EPS_PGNHEP), and generalized Hermitian-indefinite (EPS_GHIEP).
902: This function must be used to instruct SLEPc to exploit symmetry. If no
903: problem type is specified, by default a non-Hermitian problem is assumed
904: (either standard or generalized). If the user knows that the problem is
905: Hermitian (i.e. A=A^H) or generalized Hermitian (i.e. A=A^H, B=B^H, and
906: B positive definite) then it is recommended to set the problem type so
907: that eigensolver can exploit these properties.
909: Level: intermediate
911: .seealso: EPSSetOperators(), EPSSetType(), EPSGetProblemType(), EPSProblemType912: @*/
913: PetscErrorCode EPSSetProblemType(EPS eps,EPSProblemType type)914: {
918: if (type == eps->problem_type) return(0);
919: switch (type) {
920: case EPS_HEP:
921: eps->isgeneralized = PETSC_FALSE;
922: eps->ishermitian = PETSC_TRUE;
923: eps->ispositive = PETSC_FALSE;
924: break;
925: case EPS_NHEP:
926: eps->isgeneralized = PETSC_FALSE;
927: eps->ishermitian = PETSC_FALSE;
928: eps->ispositive = PETSC_FALSE;
929: break;
930: case EPS_GHEP:
931: eps->isgeneralized = PETSC_TRUE;
932: eps->ishermitian = PETSC_TRUE;
933: eps->ispositive = PETSC_TRUE;
934: break;
935: case EPS_GNHEP:
936: eps->isgeneralized = PETSC_TRUE;
937: eps->ishermitian = PETSC_FALSE;
938: eps->ispositive = PETSC_FALSE;
939: break;
940: case EPS_PGNHEP:
941: eps->isgeneralized = PETSC_TRUE;
942: eps->ishermitian = PETSC_FALSE;
943: eps->ispositive = PETSC_TRUE;
944: break;
945: case EPS_GHIEP:
946: eps->isgeneralized = PETSC_TRUE;
947: eps->ishermitian = PETSC_TRUE;
948: eps->ispositive = PETSC_FALSE;
949: break;
950: default:951: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
952: }
953: eps->problem_type = type;
954: eps->state = EPS_STATE_INITIAL;
955: return(0);
956: }
958: /*@
959: EPSGetProblemType - Gets the problem type from the EPS object.
961: Not Collective
963: Input Parameter:
964: . eps - the eigensolver context
966: Output Parameter:
967: . type - the problem type
969: Level: intermediate
971: .seealso: EPSSetProblemType(), EPSProblemType972: @*/
973: PetscErrorCode EPSGetProblemType(EPS eps,EPSProblemType *type)974: {
978: *type = eps->problem_type;
979: return(0);
980: }
982: /*@
983: EPSSetExtraction - Specifies the type of extraction technique to be employed
984: by the eigensolver.
986: Logically Collective on eps
988: Input Parameters:
989: + eps - the eigensolver context
990: - extr - a known type of extraction
992: Options Database Keys:
993: + -eps_ritz - Rayleigh-Ritz extraction
994: . -eps_harmonic - harmonic Ritz extraction
995: . -eps_harmonic_relative - harmonic Ritz extraction relative to the eigenvalue
996: . -eps_harmonic_right - harmonic Ritz extraction for rightmost eigenvalues
997: . -eps_harmonic_largest - harmonic Ritz extraction for largest magnitude
998: (without target)
999: . -eps_refined - refined Ritz extraction
1000: - -eps_refined_harmonic - refined harmonic Ritz extraction
1002: Notes:
1003: Not all eigensolvers support all types of extraction. See the SLEPc
1004: Users Manual for details.
1006: By default, a standard Rayleigh-Ritz extraction is used. Other extractions
1007: may be useful when computing interior eigenvalues.
1009: Harmonic-type extractions are used in combination with a 'target'.
1011: Level: advanced
1013: .seealso: EPSSetTarget(), EPSGetExtraction(), EPSExtraction1014: @*/
1015: PetscErrorCode EPSSetExtraction(EPS eps,EPSExtraction extr)1016: {
1020: eps->extraction = extr;
1021: return(0);
1022: }
1024: /*@
1025: EPSGetExtraction - Gets the extraction type used by the EPS object.
1027: Not Collective
1029: Input Parameter:
1030: . eps - the eigensolver context
1032: Output Parameter:
1033: . extr - name of extraction type
1035: Level: advanced
1037: .seealso: EPSSetExtraction(), EPSExtraction1038: @*/
1039: PetscErrorCode EPSGetExtraction(EPS eps,EPSExtraction *extr)1040: {
1044: *extr = eps->extraction;
1045: return(0);
1046: }
1048: /*@
1049: EPSSetBalance - Specifies the balancing technique to be employed by the
1050: eigensolver, and some parameters associated to it.
1052: Logically Collective on eps
1054: Input Parameters:
1055: + eps - the eigensolver context
1056: . bal - the balancing method, one of EPS_BALANCE_NONE, EPS_BALANCE_ONESIDE,
1057: EPS_BALANCE_TWOSIDE, or EPS_BALANCE_USER1058: . its - number of iterations of the balancing algorithm
1059: - cutoff - cutoff value
1061: Options Database Keys:
1062: + -eps_balance <method> - the balancing method, where <method> is one of
1063: 'none', 'oneside', 'twoside', or 'user'
1064: . -eps_balance_its <its> - number of iterations
1065: - -eps_balance_cutoff <cutoff> - cutoff value
1067: Notes:
1068: When balancing is enabled, the solver works implicitly with matrix DAD^-1,
1069: where D is an appropriate diagonal matrix. This improves the accuracy of
1070: the computed results in some cases. See the SLEPc Users Manual for details.
1072: Balancing makes sense only for non-Hermitian problems when the required
1073: precision is high (i.e. a small tolerance such as 1e-15).
1075: By default, balancing is disabled. The two-sided method is much more
1076: effective than the one-sided counterpart, but it requires the system
1077: matrices to have the MatMultTranspose operation defined.
1079: The parameter 'its' is the number of iterations performed by the method. The
1080: cutoff value is used only in the two-side variant. Use PETSC_DEFAULT to assign
1081: a reasonably good value.
1083: User-defined balancing is allowed provided that the corresponding matrix
1084: is set via STSetBalanceMatrix.
1086: Level: intermediate
1088: .seealso: EPSGetBalance(), EPSBalance, STSetBalanceMatrix()
1089: @*/
1090: PetscErrorCode EPSSetBalance(EPS eps,EPSBalance bal,PetscInt its,PetscReal cutoff)1091: {
1097: switch (bal) {
1098: case EPS_BALANCE_NONE:
1099: case EPS_BALANCE_ONESIDE:
1100: case EPS_BALANCE_TWOSIDE:
1101: case EPS_BALANCE_USER:
1102: if (eps->balance != bal) {
1103: eps->state = EPS_STATE_INITIAL;
1104: eps->balance = bal;
1105: }
1106: break;
1107: default:1108: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'bal'");
1109: }
1110: if (its==PETSC_DECIDE || its==PETSC_DEFAULT) eps->balance_its = 5;
1111: else {
1112: if (its <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be > 0");
1113: eps->balance_its = its;
1114: }
1115: if (cutoff==PETSC_DECIDE || cutoff==PETSC_DEFAULT) eps->balance_cutoff = 1e-8;
1116: else {
1117: if (cutoff <= 0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of cutoff. Must be > 0");
1118: eps->balance_cutoff = cutoff;
1119: }
1120: return(0);
1121: }
1123: /*@C
1124: EPSGetBalance - Gets the balancing type used by the EPS object, and the
1125: associated parameters.
1127: Not Collective
1129: Input Parameter:
1130: . eps - the eigensolver context
1132: Output Parameters:
1133: + bal - the balancing method
1134: . its - number of iterations of the balancing algorithm
1135: - cutoff - cutoff value
1137: Level: intermediate
1139: Note:
1140: The user can specify NULL for any parameter that is not needed.
1142: .seealso: EPSSetBalance(), EPSBalance1143: @*/
1144: PetscErrorCode EPSGetBalance(EPS eps,EPSBalance *bal,PetscInt *its,PetscReal *cutoff)1145: {
1148: if (bal) *bal = eps->balance;
1149: if (its) *its = eps->balance_its;
1150: if (cutoff) *cutoff = eps->balance_cutoff;
1151: return(0);
1152: }
1154: /*@
1155: EPSSetTwoSided - Sets the solver to use a two-sided variant so that left
1156: eigenvectors are also computed.
1158: Logically Collective on eps
1160: Input Parameters:
1161: + eps - the eigensolver context
1162: - twosided - whether the two-sided variant is to be used or not
1164: Options Database Keys:
1165: . -eps_two_sided <boolean> - Sets/resets the twosided flag
1167: Notes:
1168: If the user sets twosided=PETSC_TRUE then the solver uses a variant of
1169: the algorithm that computes both right and left eigenvectors. This is
1170: usually much more costly. This option is not available in all solvers.
1172: When using two-sided solvers, the problem matrices must have both the
1173: MatMult and MatMultTranspose operations defined.
1175: Level: advanced
1177: .seealso: EPSGetTwoSided(), EPSGetLeftEigenvector()
1178: @*/
1179: PetscErrorCode EPSSetTwoSided(EPS eps,PetscBool twosided)1180: {
1184: if (twosided!=eps->twosided) {
1185: eps->twosided = twosided;
1186: eps->state = EPS_STATE_INITIAL;
1187: }
1188: return(0);
1189: }
1191: /*@
1192: EPSGetTwoSided - Returns the flag indicating whether a two-sided variant
1193: of the algorithm is being used or not.
1195: Not Collective
1197: Input Parameter:
1198: . eps - the eigensolver context
1200: Output Parameter:
1201: . twosided - the returned flag
1203: Level: advanced
1205: .seealso: EPSSetTwoSided()
1206: @*/
1207: PetscErrorCode EPSGetTwoSided(EPS eps,PetscBool *twosided)1208: {
1212: *twosided = eps->twosided;
1213: return(0);
1214: }
1216: /*@
1217: EPSSetTrueResidual - Specifies if the solver must compute the true residual
1218: explicitly or not.
1220: Logically Collective on eps
1222: Input Parameters:
1223: + eps - the eigensolver context
1224: - trueres - whether true residuals are required or not
1226: Options Database Keys:
1227: . -eps_true_residual <boolean> - Sets/resets the boolean flag 'trueres'
1229: Notes:
1230: If the user sets trueres=PETSC_TRUE then the solver explicitly computes
1231: the true residual for each eigenpair approximation, and uses it for
1232: convergence testing. Computing the residual is usually an expensive
1233: operation. Some solvers (e.g., Krylov solvers) can avoid this computation
1234: by using a cheap estimate of the residual norm, but this may sometimes
1235: give inaccurate results (especially if a spectral transform is being
1236: used). On the contrary, preconditioned eigensolvers (e.g., Davidson solvers)
1237: do rely on computing the true residual, so this option is irrelevant for them.
1239: Level: advanced
1241: .seealso: EPSGetTrueResidual()
1242: @*/
1243: PetscErrorCode EPSSetTrueResidual(EPS eps,PetscBool trueres)1244: {
1248: eps->trueres = trueres;
1249: return(0);
1250: }
1252: /*@
1253: EPSGetTrueResidual - Returns the flag indicating whether true
1254: residuals must be computed explicitly or not.
1256: Not Collective
1258: Input Parameter:
1259: . eps - the eigensolver context
1261: Output Parameter:
1262: . trueres - the returned flag
1264: Level: advanced
1266: .seealso: EPSSetTrueResidual()
1267: @*/
1268: PetscErrorCode EPSGetTrueResidual(EPS eps,PetscBool *trueres)1269: {
1273: *trueres = eps->trueres;
1274: return(0);
1275: }
1277: /*@
1278: EPSSetTrackAll - Specifies if the solver must compute the residual norm of all
1279: approximate eigenpairs or not.
1281: Logically Collective on eps
1283: Input Parameters:
1284: + eps - the eigensolver context
1285: - trackall - whether to compute all residuals or not
1287: Notes:
1288: If the user sets trackall=PETSC_TRUE then the solver computes (or estimates)
1289: the residual norm for each eigenpair approximation. Computing the residual is
1290: usually an expensive operation and solvers commonly compute only the residual
1291: associated to the first unconverged eigenpair.
1293: The option '-eps_monitor_all' automatically activates this option.
1295: Level: developer
1297: .seealso: EPSGetTrackAll()
1298: @*/
1299: PetscErrorCode EPSSetTrackAll(EPS eps,PetscBool trackall)1300: {
1304: eps->trackall = trackall;
1305: return(0);
1306: }
1308: /*@
1309: EPSGetTrackAll - Returns the flag indicating whether all residual norms must
1310: be computed or not.
1312: Not Collective
1314: Input Parameter:
1315: . eps - the eigensolver context
1317: Output Parameter:
1318: . trackall - the returned flag
1320: Level: developer
1322: .seealso: EPSSetTrackAll()
1323: @*/
1324: PetscErrorCode EPSGetTrackAll(EPS eps,PetscBool *trackall)1325: {
1329: *trackall = eps->trackall;
1330: return(0);
1331: }
1333: /*@
1334: EPSSetPurify - Deactivate eigenvector purification (which is activated by default).
1336: Logically Collective on eps
1338: Input Parameters:
1339: + eps - the eigensolver context
1340: - purify - whether purification is required or not
1342: Options Database Keys:
1343: . -eps_purify <boolean> - Sets/resets the boolean flag 'purify'
1345: Notes:
1346: By default, eigenvectors of generalized symmetric eigenproblems are purified
1347: in order to purge directions in the nullspace of matrix B. If the user knows
1348: that B is non-singular, then purification can be safely deactivated and some
1349: computational cost is avoided (this is particularly important in interval computations).
1351: Level: intermediate
1353: .seealso: EPSGetPurify(), EPSSetInterval()
1354: @*/
1355: PetscErrorCode EPSSetPurify(EPS eps,PetscBool purify)1356: {
1360: if (purify!=eps->purify) {
1361: eps->purify = purify;
1362: eps->state = EPS_STATE_INITIAL;
1363: }
1364: return(0);
1365: }
1367: /*@
1368: EPSGetPurify - Returns the flag indicating whether purification is activated
1369: or not.
1371: Not Collective
1373: Input Parameter:
1374: . eps - the eigensolver context
1376: Output Parameter:
1377: . purify - the returned flag
1379: Level: intermediate
1381: .seealso: EPSSetPurify()
1382: @*/
1383: PetscErrorCode EPSGetPurify(EPS eps,PetscBool *purify)1384: {
1388: *purify = eps->purify;
1389: return(0);
1390: }
1392: /*@C
1393: EPSSetOptionsPrefix - Sets the prefix used for searching for all
1394: EPS options in the database.
1396: Logically Collective on eps
1398: Input Parameters:
1399: + eps - the eigensolver context
1400: - prefix - the prefix string to prepend to all EPS option requests
1402: Notes:
1403: A hyphen (-) must NOT be given at the beginning of the prefix name.
1404: The first character of all runtime options is AUTOMATICALLY the
1405: hyphen.
1407: For example, to distinguish between the runtime options for two
1408: different EPS contexts, one could call
1409: .vb
1410: EPSSetOptionsPrefix(eps1,"eig1_")
1411: EPSSetOptionsPrefix(eps2,"eig2_")
1412: .ve
1414: Level: advanced
1416: .seealso: EPSAppendOptionsPrefix(), EPSGetOptionsPrefix()
1417: @*/
1418: PetscErrorCode EPSSetOptionsPrefix(EPS eps,const char *prefix)1419: {
1424: if (!eps->st) { EPSGetST(eps,&eps->st); }
1425: STSetOptionsPrefix(eps->st,prefix);
1426: if (!eps->V) { EPSGetBV(eps,&eps->V); }
1427: BVSetOptionsPrefix(eps->V,prefix);
1428: if (!eps->ds) { EPSGetDS(eps,&eps->ds); }
1429: DSSetOptionsPrefix(eps->ds,prefix);
1430: if (!eps->rg) { EPSGetRG(eps,&eps->rg); }
1431: RGSetOptionsPrefix(eps->rg,prefix);
1432: PetscObjectSetOptionsPrefix((PetscObject)eps,prefix);
1433: return(0);
1434: }
1436: /*@C
1437: EPSAppendOptionsPrefix - Appends to the prefix used for searching for all
1438: EPS options in the database.
1440: Logically Collective on eps
1442: Input Parameters:
1443: + eps - the eigensolver context
1444: - prefix - the prefix string to prepend to all EPS option requests
1446: Notes:
1447: A hyphen (-) must NOT be given at the beginning of the prefix name.
1448: The first character of all runtime options is AUTOMATICALLY the hyphen.
1450: Level: advanced
1452: .seealso: EPSSetOptionsPrefix(), EPSGetOptionsPrefix()
1453: @*/
1454: PetscErrorCode EPSAppendOptionsPrefix(EPS eps,const char *prefix)1455: {
1460: if (!eps->st) { EPSGetST(eps,&eps->st); }
1461: STAppendOptionsPrefix(eps->st,prefix);
1462: if (!eps->V) { EPSGetBV(eps,&eps->V); }
1463: BVAppendOptionsPrefix(eps->V,prefix);
1464: if (!eps->ds) { EPSGetDS(eps,&eps->ds); }
1465: DSAppendOptionsPrefix(eps->ds,prefix);
1466: if (!eps->rg) { EPSGetRG(eps,&eps->rg); }
1467: RGAppendOptionsPrefix(eps->rg,prefix);
1468: PetscObjectAppendOptionsPrefix((PetscObject)eps,prefix);
1469: return(0);
1470: }
1472: /*@C
1473: EPSGetOptionsPrefix - Gets the prefix used for searching for all
1474: EPS options in the database.
1476: Not Collective
1478: Input Parameters:
1479: . eps - the eigensolver context
1481: Output Parameters:
1482: . prefix - pointer to the prefix string used is returned
1484: Note:
1485: On the Fortran side, the user should pass in a string 'prefix' of
1486: sufficient length to hold the prefix.
1488: Level: advanced
1490: .seealso: EPSSetOptionsPrefix(), EPSAppendOptionsPrefix()
1491: @*/
1492: PetscErrorCode EPSGetOptionsPrefix(EPS eps,const char *prefix[])1493: {
1499: PetscObjectGetOptionsPrefix((PetscObject)eps,prefix);
1500: return(0);
1501: }