Actual source code: epsopts.c
slepc-3.11.2 2019-07-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2019, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: EPS routines related to options that can be set via the command-line
12: or procedurally.
13: */
15: #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
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: . name - the monitor option name
27: . help - message indicating what monitoring is done
28: . manual - manual page for the monitor
29: . monitor - the monitor function, whose context is a PetscViewerAndFormat
30: - trackall - whether this monitor tracks all eigenvalues or not
32: Level: developer
34: .seealso: EPSMonitorSet(), EPSSetTrackAll(), EPSConvMonitorSetFromOptions()
35: @*/
36: PetscErrorCode EPSMonitorSetFromOptions(EPS eps,const char name[],const char help[],const char manual[],PetscErrorCode (*monitor)(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,PetscViewerAndFormat*),PetscBool trackall)
37: {
38: PetscErrorCode ierr;
39: PetscBool flg;
40: PetscViewer viewer;
41: PetscViewerFormat format;
42: PetscViewerAndFormat *vf;
45: PetscOptionsGetViewer(PetscObjectComm((PetscObject)eps),((PetscObject)eps)->options,((PetscObject)eps)->prefix,name,&viewer,&format,&flg);
46: if (flg) {
47: PetscViewerAndFormatCreate(viewer,format,&vf);
48: PetscObjectDereference((PetscObject)viewer);
49: EPSMonitorSet(eps,(PetscErrorCode (*)(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))monitor,vf,(PetscErrorCode (*)(void**))PetscViewerAndFormatDestroy);
50: if (trackall) {
51: EPSSetTrackAll(eps,PETSC_TRUE);
52: }
53: }
54: return(0);
55: }
57: /*@C
58: EPSConvMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
59: indicated by the user (for monitors that only show iteration numbers of convergence).
61: Collective on EPS
63: Input Parameters:
64: + eps - the eigensolver context
65: . name - the monitor option name
66: . help - message indicating what monitoring is done
67: . manual - manual page for the monitor
68: - monitor - the monitor function, whose context is a SlepcConvMonitor
70: Level: developer
72: .seealso: EPSMonitorSet(), EPSMonitorSetFromOptions()
73: @*/
74: PetscErrorCode EPSConvMonitorSetFromOptions(EPS eps,const char name[],const char help[],const char manual[],PetscErrorCode (*monitor)(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,SlepcConvMonitor))
75: {
76: PetscErrorCode ierr;
77: PetscBool flg;
78: PetscViewer viewer;
79: PetscViewerFormat format;
80: SlepcConvMonitor ctx;
83: PetscOptionsGetViewer(PetscObjectComm((PetscObject)eps),((PetscObject)eps)->options,((PetscObject)eps)->prefix,name,&viewer,&format,&flg);
84: if (flg) {
85: SlepcConvMonitorCreate(viewer,format,&ctx);
86: PetscObjectDereference((PetscObject)viewer);
87: EPSMonitorSet(eps,(PetscErrorCode (*)(EPS,PetscInt,PetscInt,PetscScalar*,PetscScalar*,PetscReal*,PetscInt,void*))monitor,ctx,(PetscErrorCode (*)(void**))SlepcConvMonitorDestroy);
88: }
89: return(0);
90: }
92: /*@
93: EPSSetFromOptions - Sets EPS options from the options database.
94: This routine must be called before EPSSetUp() if the user is to be
95: allowed to set the solver type.
97: Collective on EPS
99: Input Parameters:
100: . eps - the eigensolver context
102: Notes:
103: To see all options, run your program with the -help option.
105: Level: beginner
106: @*/
107: PetscErrorCode EPSSetFromOptions(EPS eps)
108: {
110: char type[256];
111: PetscBool set,flg,flg1,flg2,flg3,bval;
112: PetscReal r,array[2]={0,0};
113: PetscScalar s;
114: PetscInt i,j,k;
115: PetscDrawLG lg;
116: EPSBalance bal;
120: EPSRegisterAll();
121: PetscObjectOptionsBegin((PetscObject)eps);
122: PetscOptionsFList("-eps_type","Eigensolver method","EPSSetType",EPSList,(char*)(((PetscObject)eps)->type_name?((PetscObject)eps)->type_name:EPSKRYLOVSCHUR),type,256,&flg);
123: if (flg) {
124: EPSSetType(eps,type);
125: } else if (!((PetscObject)eps)->type_name) {
126: EPSSetType(eps,EPSKRYLOVSCHUR);
127: }
129: PetscOptionsBoolGroupBegin("-eps_hermitian","Hermitian eigenvalue problem","EPSSetProblemType",&flg);
130: if (flg) { EPSSetProblemType(eps,EPS_HEP); }
131: PetscOptionsBoolGroup("-eps_gen_hermitian","Generalized Hermitian eigenvalue problem","EPSSetProblemType",&flg);
132: if (flg) { EPSSetProblemType(eps,EPS_GHEP); }
133: PetscOptionsBoolGroup("-eps_non_hermitian","Non-Hermitian eigenvalue problem","EPSSetProblemType",&flg);
134: if (flg) { EPSSetProblemType(eps,EPS_NHEP); }
135: PetscOptionsBoolGroup("-eps_gen_non_hermitian","Generalized non-Hermitian eigenvalue problem","EPSSetProblemType",&flg);
136: if (flg) { EPSSetProblemType(eps,EPS_GNHEP); }
137: PetscOptionsBoolGroup("-eps_pos_gen_non_hermitian","Generalized non-Hermitian eigenvalue problem with positive semi-definite B","EPSSetProblemType",&flg);
138: if (flg) { EPSSetProblemType(eps,EPS_PGNHEP); }
139: PetscOptionsBoolGroupEnd("-eps_gen_indefinite","Generalized Hermitian-indefinite eigenvalue problem","EPSSetProblemType",&flg);
140: if (flg) { EPSSetProblemType(eps,EPS_GHIEP); }
142: PetscOptionsBoolGroupBegin("-eps_ritz","Rayleigh-Ritz extraction","EPSSetExtraction",&flg);
143: if (flg) { EPSSetExtraction(eps,EPS_RITZ); }
144: PetscOptionsBoolGroup("-eps_harmonic","Harmonic Ritz extraction","EPSSetExtraction",&flg);
145: if (flg) { EPSSetExtraction(eps,EPS_HARMONIC); }
146: PetscOptionsBoolGroup("-eps_harmonic_relative","Relative harmonic Ritz extraction","EPSSetExtraction",&flg);
147: if (flg) { EPSSetExtraction(eps,EPS_HARMONIC_RELATIVE); }
148: PetscOptionsBoolGroup("-eps_harmonic_right","Right harmonic Ritz extraction","EPSSetExtraction",&flg);
149: if (flg) { EPSSetExtraction(eps,EPS_HARMONIC_RIGHT); }
150: PetscOptionsBoolGroup("-eps_harmonic_largest","Largest harmonic Ritz extraction","EPSSetExtraction",&flg);
151: if (flg) { EPSSetExtraction(eps,EPS_HARMONIC_LARGEST); }
152: PetscOptionsBoolGroup("-eps_refined","Refined Ritz extraction","EPSSetExtraction",&flg);
153: if (flg) { EPSSetExtraction(eps,EPS_REFINED); }
154: PetscOptionsBoolGroupEnd("-eps_refined_harmonic","Refined harmonic Ritz extraction","EPSSetExtraction",&flg);
155: if (flg) { EPSSetExtraction(eps,EPS_REFINED_HARMONIC); }
157: bal = eps->balance;
158: PetscOptionsEnum("-eps_balance","Balancing method","EPSSetBalance",EPSBalanceTypes,(PetscEnum)bal,(PetscEnum*)&bal,&flg1);
159: j = eps->balance_its;
160: PetscOptionsInt("-eps_balance_its","Number of iterations in balancing","EPSSetBalance",eps->balance_its,&j,&flg2);
161: r = eps->balance_cutoff;
162: PetscOptionsReal("-eps_balance_cutoff","Cutoff value in balancing","EPSSetBalance",eps->balance_cutoff,&r,&flg3);
163: if (flg1 || flg2 || flg3) { EPSSetBalance(eps,bal,j,r); }
165: i = eps->max_it? eps->max_it: PETSC_DEFAULT;
166: PetscOptionsInt("-eps_max_it","Maximum number of iterations","EPSSetTolerances",eps->max_it,&i,&flg1);
167: r = eps->tol;
168: PetscOptionsReal("-eps_tol","Tolerance","EPSSetTolerances",eps->tol==PETSC_DEFAULT?SLEPC_DEFAULT_TOL:eps->tol,&r,&flg2);
169: if (flg1 || flg2) { EPSSetTolerances(eps,r,i); }
171: PetscOptionsBoolGroupBegin("-eps_conv_rel","Relative error convergence test","EPSSetConvergenceTest",&flg);
172: if (flg) { EPSSetConvergenceTest(eps,EPS_CONV_REL); }
173: PetscOptionsBoolGroup("-eps_conv_norm","Convergence test relative to the eigenvalue and the matrix norms","EPSSetConvergenceTest",&flg);
174: if (flg) { EPSSetConvergenceTest(eps,EPS_CONV_NORM); }
175: PetscOptionsBoolGroup("-eps_conv_abs","Absolute error convergence test","EPSSetConvergenceTest",&flg);
176: if (flg) { EPSSetConvergenceTest(eps,EPS_CONV_ABS); }
177: PetscOptionsBoolGroupEnd("-eps_conv_user","User-defined convergence test","EPSSetConvergenceTest",&flg);
178: if (flg) { EPSSetConvergenceTest(eps,EPS_CONV_USER); }
180: PetscOptionsBoolGroupBegin("-eps_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","EPSSetStoppingTest",&flg);
181: if (flg) { EPSSetStoppingTest(eps,EPS_STOP_BASIC); }
182: PetscOptionsBoolGroupEnd("-eps_stop_user","User-defined stopping test","EPSSetStoppingTest",&flg);
183: if (flg) { EPSSetStoppingTest(eps,EPS_STOP_USER); }
185: i = eps->nev;
186: PetscOptionsInt("-eps_nev","Number of eigenvalues to compute","EPSSetDimensions",eps->nev,&i,&flg1);
187: j = eps->ncv? eps->ncv: PETSC_DEFAULT;
188: PetscOptionsInt("-eps_ncv","Number of basis vectors","EPSSetDimensions",eps->ncv,&j,&flg2);
189: k = eps->mpd? eps->mpd: PETSC_DEFAULT;
190: PetscOptionsInt("-eps_mpd","Maximum dimension of projected problem","EPSSetDimensions",eps->mpd,&k,&flg3);
191: if (flg1 || flg2 || flg3) { EPSSetDimensions(eps,i,j,k); }
193: PetscOptionsBoolGroupBegin("-eps_largest_magnitude","Compute largest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg);
194: if (flg) { EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE); }
195: PetscOptionsBoolGroup("-eps_smallest_magnitude","Compute smallest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg);
196: if (flg) { EPSSetWhichEigenpairs(eps,EPS_SMALLEST_MAGNITUDE); }
197: PetscOptionsBoolGroup("-eps_largest_real","Compute eigenvalues with largest real parts","EPSSetWhichEigenpairs",&flg);
198: if (flg) { EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL); }
199: PetscOptionsBoolGroup("-eps_smallest_real","Compute eigenvalues with smallest real parts","EPSSetWhichEigenpairs",&flg);
200: if (flg) { EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL); }
201: PetscOptionsBoolGroup("-eps_largest_imaginary","Compute eigenvalues with largest imaginary parts","EPSSetWhichEigenpairs",&flg);
202: if (flg) { EPSSetWhichEigenpairs(eps,EPS_LARGEST_IMAGINARY); }
203: PetscOptionsBoolGroup("-eps_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","EPSSetWhichEigenpairs",&flg);
204: if (flg) { EPSSetWhichEigenpairs(eps,EPS_SMALLEST_IMAGINARY); }
205: PetscOptionsBoolGroup("-eps_target_magnitude","Compute eigenvalues closest to target","EPSSetWhichEigenpairs",&flg);
206: if (flg) { EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE); }
207: PetscOptionsBoolGroup("-eps_target_real","Compute eigenvalues with real parts closest to target","EPSSetWhichEigenpairs",&flg);
208: if (flg) { EPSSetWhichEigenpairs(eps,EPS_TARGET_REAL); }
209: PetscOptionsBoolGroup("-eps_target_imaginary","Compute eigenvalues with imaginary parts closest to target","EPSSetWhichEigenpairs",&flg);
210: if (flg) { EPSSetWhichEigenpairs(eps,EPS_TARGET_IMAGINARY); }
211: PetscOptionsBoolGroupEnd("-eps_all","Compute all eigenvalues in an interval or a region","EPSSetWhichEigenpairs",&flg);
212: if (flg) { EPSSetWhichEigenpairs(eps,EPS_ALL); }
214: PetscOptionsScalar("-eps_target","Value of the target","EPSSetTarget",eps->target,&s,&flg);
215: if (flg) {
216: if (eps->which!=EPS_TARGET_REAL && eps->which!=EPS_TARGET_IMAGINARY) {
217: EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE);
218: }
219: EPSSetTarget(eps,s);
220: }
222: k = 2;
223: PetscOptionsRealArray("-eps_interval","Computational interval (two real values separated with a comma without spaces)","EPSSetInterval",array,&k,&flg);
224: if (flg) {
225: if (k<2) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_SIZ,"Must pass two values in -eps_interval (comma-separated without spaces)");
226: EPSSetWhichEigenpairs(eps,EPS_ALL);
227: EPSSetInterval(eps,array[0],array[1]);
228: }
230: PetscOptionsBool("-eps_true_residual","Compute true residuals explicitly","EPSSetTrueResidual",eps->trueres,&eps->trueres,NULL);
231: PetscOptionsBool("-eps_purify","Postprocess eigenvectors for purification","EPSSetPurify",eps->purify,&bval,&flg);
232: if (flg) { EPSSetPurify(eps,bval); }
233: PetscOptionsBool("-eps_two_sided","Use two-sided variant (to compute left eigenvectors)","EPSSetTwoSided",eps->twosided,&bval,&flg);
234: if (flg) { EPSSetTwoSided(eps,bval); }
236: /* -----------------------------------------------------------------------*/
237: /*
238: Cancels all monitors hardwired into code before call to EPSSetFromOptions()
239: */
240: PetscOptionsBool("-eps_monitor_cancel","Remove any hardwired monitor routines","EPSMonitorCancel",PETSC_FALSE,&flg,&set);
241: if (set && flg) {
242: EPSMonitorCancel(eps);
243: }
244: /*
245: Text monitors
246: */
247: EPSMonitorSetFromOptions(eps,"-eps_monitor","Monitor first unconverged approximate eigenvalue and error estimate","EPSMonitorFirst",EPSMonitorFirst,PETSC_FALSE);
248: EPSConvMonitorSetFromOptions(eps,"-eps_monitor_conv","Monitor approximate eigenvalues and error estimates as they converge","EPSMonitorConverged",EPSMonitorConverged);
249: EPSMonitorSetFromOptions(eps,"-eps_monitor_all","Monitor approximate eigenvalues and error estimates","EPSMonitorAll",EPSMonitorAll,PETSC_TRUE);
250: /*
251: Line graph monitors
252: */
253: PetscOptionsBool("-eps_monitor_lg","Monitor first unconverged approximate eigenvalue and error estimate graphically","EPSMonitorSet",PETSC_FALSE,&flg,&set);
254: if (set && flg) {
255: EPSMonitorLGCreate(PetscObjectComm((PetscObject)eps),NULL,"Error estimates",PETSC_DECIDE,PETSC_DECIDE,300,300,&lg);
256: EPSMonitorSet(eps,EPSMonitorLG,lg,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
257: }
258: PetscOptionsBool("-eps_monitor_lg_all","Monitor error estimates graphically","EPSMonitorSet",PETSC_FALSE,&flg,&set);
259: if (set && flg) {
260: EPSMonitorLGCreate(PetscObjectComm((PetscObject)eps),NULL,"Error estimates",PETSC_DECIDE,PETSC_DECIDE,300,300,&lg);
261: EPSMonitorSet(eps,EPSMonitorLGAll,lg,(PetscErrorCode (*)(void**))PetscDrawLGDestroy);
262: EPSSetTrackAll(eps,PETSC_TRUE);
263: }
265: /* -----------------------------------------------------------------------*/
266: PetscOptionsName("-eps_view","Print detailed information on solver used","EPSView",NULL);
267: PetscOptionsName("-eps_view_vectors","View computed eigenvectors","EPSVectorsView",NULL);
268: PetscOptionsName("-eps_view_values","View computed eigenvalues","EPSValuesView",NULL);
269: PetscOptionsName("-eps_converged_reason","Print reason for convergence, and number of iterations","EPSReasonView",NULL);
270: PetscOptionsName("-eps_error_absolute","Print absolute errors of each eigenpair","EPSErrorView",NULL);
271: PetscOptionsName("-eps_error_relative","Print relative errors of each eigenpair","EPSErrorView",NULL);
272: PetscOptionsName("-eps_error_backward","Print backward errors of each eigenpair","EPSErrorView",NULL);
274: if (eps->ops->setfromoptions) {
275: (*eps->ops->setfromoptions)(PetscOptionsObject,eps);
276: }
277: PetscObjectProcessOptionsHandlers(PetscOptionsObject,(PetscObject)eps);
278: PetscOptionsEnd();
280: if (!eps->V) { EPSGetBV(eps,&eps->V); }
281: BVSetFromOptions(eps->V);
282: if (!eps->rg) { EPSGetRG(eps,&eps->rg); }
283: RGSetFromOptions(eps->rg);
284: if (eps->useds) {
285: if (!eps->ds) { EPSGetDS(eps,&eps->ds); }
286: DSSetFromOptions(eps->ds);
287: }
288: if (!eps->st) { EPSGetST(eps,&eps->st); }
289: EPSSetDefaultST(eps);
290: STSetFromOptions(eps->st);
291: return(0);
292: }
294: /*@C
295: EPSGetTolerances - Gets the tolerance and maximum iteration count used
296: by the EPS convergence tests.
298: Not Collective
300: Input Parameter:
301: . eps - the eigensolver context
303: Output Parameters:
304: + tol - the convergence tolerance
305: - maxits - maximum number of iterations
307: Notes:
308: The user can specify NULL for any parameter that is not needed.
310: Level: intermediate
312: .seealso: EPSSetTolerances()
313: @*/
314: PetscErrorCode EPSGetTolerances(EPS eps,PetscReal *tol,PetscInt *maxits)
315: {
318: if (tol) *tol = eps->tol;
319: if (maxits) *maxits = eps->max_it;
320: return(0);
321: }
323: /*@
324: EPSSetTolerances - Sets the tolerance and maximum iteration count used
325: by the EPS convergence tests.
327: Logically Collective on EPS
329: Input Parameters:
330: + eps - the eigensolver context
331: . tol - the convergence tolerance
332: - maxits - maximum number of iterations to use
334: Options Database Keys:
335: + -eps_tol <tol> - Sets the convergence tolerance
336: - -eps_max_it <maxits> - Sets the maximum number of iterations allowed
338: Notes:
339: Use PETSC_DEFAULT for either argument to assign a reasonably good value.
341: Level: intermediate
343: .seealso: EPSGetTolerances()
344: @*/
345: PetscErrorCode EPSSetTolerances(EPS eps,PetscReal tol,PetscInt maxits)
346: {
351: if (tol == PETSC_DEFAULT) {
352: eps->tol = PETSC_DEFAULT;
353: eps->state = EPS_STATE_INITIAL;
354: } else {
355: if (tol <= 0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
356: eps->tol = tol;
357: }
358: if (maxits == PETSC_DEFAULT || maxits == PETSC_DECIDE) {
359: eps->max_it = 0;
360: eps->state = EPS_STATE_INITIAL;
361: } else {
362: if (maxits <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
363: eps->max_it = maxits;
364: }
365: return(0);
366: }
368: /*@C
369: EPSGetDimensions - Gets the number of eigenvalues to compute
370: and the dimension of the subspace.
372: Not Collective
374: Input Parameter:
375: . eps - the eigensolver context
377: Output Parameters:
378: + nev - number of eigenvalues to compute
379: . ncv - the maximum dimension of the subspace to be used by the solver
380: - mpd - the maximum dimension allowed for the projected problem
382: Level: intermediate
384: .seealso: EPSSetDimensions()
385: @*/
386: PetscErrorCode EPSGetDimensions(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
387: {
390: if (nev) *nev = eps->nev;
391: if (ncv) *ncv = eps->ncv;
392: if (mpd) *mpd = eps->mpd;
393: return(0);
394: }
396: /*@
397: EPSSetDimensions - Sets the number of eigenvalues to compute
398: and the dimension of the subspace.
400: Logically Collective on EPS
402: Input Parameters:
403: + eps - the eigensolver context
404: . nev - number of eigenvalues to compute
405: . ncv - the maximum dimension of the subspace to be used by the solver
406: - mpd - the maximum dimension allowed for the projected problem
408: Options Database Keys:
409: + -eps_nev <nev> - Sets the number of eigenvalues
410: . -eps_ncv <ncv> - Sets the dimension of the subspace
411: - -eps_mpd <mpd> - Sets the maximum projected dimension
413: Notes:
414: Use PETSC_DEFAULT for ncv and mpd to assign a reasonably good value, which is
415: dependent on the solution method.
417: The parameters ncv and mpd are intimately related, so that the user is advised
418: to set one of them at most. Normal usage is that
419: (a) in cases where nev is small, the user sets ncv (a reasonable default is 2*nev); and
420: (b) in cases where nev is large, the user sets mpd.
422: The value of ncv should always be between nev and (nev+mpd), typically
423: ncv=nev+mpd. If nev is not too large, mpd=nev is a reasonable choice, otherwise
424: a smaller value should be used.
426: When computing all eigenvalues in an interval, see EPSSetInterval(), these
427: parameters lose relevance, and tuning must be done with
428: EPSKrylovSchurSetDimensions().
430: Level: intermediate
432: .seealso: EPSGetDimensions(), EPSSetInterval(), EPSKrylovSchurSetDimensions()
433: @*/
434: PetscErrorCode EPSSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
435: {
441: if (nev<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
442: eps->nev = nev;
443: if (ncv == PETSC_DECIDE || ncv == PETSC_DEFAULT) {
444: eps->ncv = 0;
445: } else {
446: if (ncv<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
447: eps->ncv = ncv;
448: }
449: if (mpd == PETSC_DECIDE || mpd == PETSC_DEFAULT) {
450: eps->mpd = 0;
451: } else {
452: if (mpd<1) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
453: eps->mpd = mpd;
454: }
455: eps->state = EPS_STATE_INITIAL;
456: return(0);
457: }
459: /*@
460: EPSSetWhichEigenpairs - Specifies which portion of the spectrum is
461: to be sought.
463: Logically Collective on EPS
465: Input Parameters:
466: + eps - eigensolver context obtained from EPSCreate()
467: - which - the portion of the spectrum to be sought
469: Possible values:
470: The parameter 'which' can have one of these values
472: + EPS_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
473: . EPS_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
474: . EPS_LARGEST_REAL - largest real parts
475: . EPS_SMALLEST_REAL - smallest real parts
476: . EPS_LARGEST_IMAGINARY - largest imaginary parts
477: . EPS_SMALLEST_IMAGINARY - smallest imaginary parts
478: . EPS_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
479: . EPS_TARGET_REAL - eigenvalues with real part closest to target
480: . EPS_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
481: . EPS_ALL - all eigenvalues contained in a given interval or region
482: - EPS_WHICH_USER - user defined ordering set with EPSSetEigenvalueComparison()
484: Options Database Keys:
485: + -eps_largest_magnitude - Sets largest eigenvalues in magnitude
486: . -eps_smallest_magnitude - Sets smallest eigenvalues in magnitude
487: . -eps_largest_real - Sets largest real parts
488: . -eps_smallest_real - Sets smallest real parts
489: . -eps_largest_imaginary - Sets largest imaginary parts
490: . -eps_smallest_imaginary - Sets smallest imaginary parts
491: . -eps_target_magnitude - Sets eigenvalues closest to target
492: . -eps_target_real - Sets real parts closest to target
493: . -eps_target_imaginary - Sets imaginary parts closest to target
494: - -eps_all - Sets all eigenvalues in an interval or region
496: Notes:
497: Not all eigensolvers implemented in EPS account for all the possible values
498: stated above. Also, some values make sense only for certain types of
499: problems. If SLEPc is compiled for real numbers EPS_LARGEST_IMAGINARY
500: and EPS_SMALLEST_IMAGINARY use the absolute value of the imaginary part
501: for eigenvalue selection.
503: The target is a scalar value provided with EPSSetTarget().
505: The criterion EPS_TARGET_IMAGINARY is available only in case PETSc and
506: SLEPc have been built with complex scalars.
508: EPS_ALL is intended for use in combination with an interval (see
509: EPSSetInterval()), when all eigenvalues within the interval are requested,
510: or in the context of the CISS solver for computing all eigenvalues in a region.
511: In those cases, the number of eigenvalues is unknown, so the nev parameter
512: has a different sense, see EPSSetDimensions().
514: Level: intermediate
516: .seealso: EPSGetWhichEigenpairs(), EPSSetTarget(), EPSSetInterval(),
517: EPSSetDimensions(), EPSSetEigenvalueComparison(), EPSWhich
518: @*/
519: PetscErrorCode EPSSetWhichEigenpairs(EPS eps,EPSWhich which)
520: {
524: switch (which) {
525: case EPS_LARGEST_MAGNITUDE:
526: case EPS_SMALLEST_MAGNITUDE:
527: case EPS_LARGEST_REAL:
528: case EPS_SMALLEST_REAL:
529: case EPS_LARGEST_IMAGINARY:
530: case EPS_SMALLEST_IMAGINARY:
531: case EPS_TARGET_MAGNITUDE:
532: case EPS_TARGET_REAL:
533: #if defined(PETSC_USE_COMPLEX)
534: case EPS_TARGET_IMAGINARY:
535: #endif
536: case EPS_ALL:
537: case EPS_WHICH_USER:
538: if (eps->which != which) {
539: eps->state = EPS_STATE_INITIAL;
540: eps->which = which;
541: }
542: break;
543: default:
544: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
545: }
546: return(0);
547: }
549: /*@
550: EPSGetWhichEigenpairs - Returns which portion of the spectrum is to be
551: sought.
553: Not Collective
555: Input Parameter:
556: . eps - eigensolver context obtained from EPSCreate()
558: Output Parameter:
559: . which - the portion of the spectrum to be sought
561: Notes:
562: See EPSSetWhichEigenpairs() for possible values of 'which'.
564: Level: intermediate
566: .seealso: EPSSetWhichEigenpairs(), EPSWhich
567: @*/
568: PetscErrorCode EPSGetWhichEigenpairs(EPS eps,EPSWhich *which)
569: {
573: *which = eps->which;
574: return(0);
575: }
577: /*@C
578: EPSSetEigenvalueComparison - Specifies the eigenvalue comparison function
579: when EPSSetWhichEigenpairs() is set to EPS_WHICH_USER.
581: Logically Collective on EPS
583: Input Parameters:
584: + eps - eigensolver context obtained from EPSCreate()
585: . func - a pointer to the comparison function
586: - ctx - a context pointer (the last parameter to the comparison function)
588: Calling Sequence of func:
589: $ func(PetscScalar ar,PetscScalar ai,PetscScalar br,PetscScalar bi,PetscInt *res,void *ctx)
591: + ar - real part of the 1st eigenvalue
592: . ai - imaginary part of the 1st eigenvalue
593: . br - real part of the 2nd eigenvalue
594: . bi - imaginary part of the 2nd eigenvalue
595: . res - result of comparison
596: - ctx - optional context, as set by EPSSetEigenvalueComparison()
598: Note:
599: The returning parameter 'res' can be
600: + negative - if the 1st eigenvalue is preferred to the 2st one
601: . zero - if both eigenvalues are equally preferred
602: - positive - if the 2st eigenvalue is preferred to the 1st one
604: Level: advanced
606: .seealso: EPSSetWhichEigenpairs(), EPSWhich
607: @*/
608: PetscErrorCode EPSSetEigenvalueComparison(EPS eps,PetscErrorCode (*func)(PetscScalar,PetscScalar,PetscScalar,PetscScalar,PetscInt*,void*),void* ctx)
609: {
612: eps->sc->comparison = func;
613: eps->sc->comparisonctx = ctx;
614: eps->which = EPS_WHICH_USER;
615: return(0);
616: }
618: /*@C
619: EPSSetArbitrarySelection - Specifies a function intended to look for
620: eigenvalues according to an arbitrary selection criterion. This criterion
621: can be based on a computation involving the current eigenvector approximation.
623: Logically Collective on EPS
625: Input Parameters:
626: + eps - eigensolver context obtained from EPSCreate()
627: . func - a pointer to the evaluation function
628: - ctx - a context pointer (the last parameter to the evaluation function)
630: Calling Sequence of func:
631: $ func(PetscScalar er,PetscScalar ei,Vec xr,Vec xi,PetscScalar *rr,PetscScalar *ri,void *ctx)
633: + er - real part of the current eigenvalue approximation
634: . ei - imaginary part of the current eigenvalue approximation
635: . xr - real part of the current eigenvector approximation
636: . xi - imaginary part of the current eigenvector approximation
637: . rr - result of evaluation (real part)
638: . ri - result of evaluation (imaginary part)
639: - ctx - optional context, as set by EPSSetArbitrarySelection()
641: Notes:
642: This provides a mechanism to select eigenpairs by evaluating a user-defined
643: function. When a function has been provided, the default selection based on
644: sorting the eigenvalues is replaced by the sorting of the results of this
645: function (with the same sorting criterion given in EPSSetWhichEigenpairs()).
647: For instance, suppose you want to compute those eigenvectors that maximize
648: a certain computable expression. Then implement the computation using
649: the arguments xr and xi, and return the result in rr. Then set the standard
650: sorting by magnitude so that the eigenpair with largest value of rr is
651: selected.
653: This evaluation function is collective, that is, all processes call it and
654: it can use collective operations; furthermore, the computed result must
655: be the same in all processes.
657: The result of func is expressed as a complex number so that it is possible to
658: use the standard eigenvalue sorting functions, but normally only rr is used.
659: Set ri to zero unless it is meaningful in your application.
661: Level: advanced
663: .seealso: EPSSetWhichEigenpairs()
664: @*/
665: PetscErrorCode EPSSetArbitrarySelection(EPS eps,PetscErrorCode (*func)(PetscScalar,PetscScalar,Vec,Vec,PetscScalar*,PetscScalar*,void*),void* ctx)
666: {
669: eps->arbitrary = func;
670: eps->arbitraryctx = ctx;
671: eps->state = EPS_STATE_INITIAL;
672: return(0);
673: }
675: /*@C
676: EPSSetConvergenceTestFunction - Sets a function to compute the error estimate
677: used in the convergence test.
679: Logically Collective on EPS
681: Input Parameters:
682: + eps - eigensolver context obtained from EPSCreate()
683: . func - a pointer to the convergence test function
684: . ctx - context for private data for the convergence routine (may be null)
685: - destroy - a routine for destroying the context (may be null)
687: Calling Sequence of func:
688: $ func(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx)
690: + eps - eigensolver context obtained from EPSCreate()
691: . eigr - real part of the eigenvalue
692: . eigi - imaginary part of the eigenvalue
693: . res - residual norm associated to the eigenpair
694: . errest - (output) computed error estimate
695: - ctx - optional context, as set by EPSSetConvergenceTestFunction()
697: Note:
698: If the error estimate returned by the convergence test function is less than
699: the tolerance, then the eigenvalue is accepted as converged.
701: Level: advanced
703: .seealso: EPSSetConvergenceTest(), EPSSetTolerances()
704: @*/
705: PetscErrorCode EPSSetConvergenceTestFunction(EPS eps,PetscErrorCode (*func)(EPS,PetscScalar,PetscScalar,PetscReal,PetscReal*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
706: {
711: if (eps->convergeddestroy) {
712: (*eps->convergeddestroy)(eps->convergedctx);
713: }
714: eps->convergeduser = func;
715: eps->convergeddestroy = destroy;
716: eps->convergedctx = ctx;
717: if (func == EPSConvergedRelative) eps->conv = EPS_CONV_REL;
718: else if (func == EPSConvergedNorm) eps->conv = EPS_CONV_NORM;
719: else if (func == EPSConvergedAbsolute) eps->conv = EPS_CONV_ABS;
720: else {
721: eps->conv = EPS_CONV_USER;
722: eps->converged = eps->convergeduser;
723: }
724: return(0);
725: }
727: /*@
728: EPSSetConvergenceTest - Specifies how to compute the error estimate
729: used in the convergence test.
731: Logically Collective on EPS
733: Input Parameters:
734: + eps - eigensolver context obtained from EPSCreate()
735: - conv - the type of convergence test
737: Options Database Keys:
738: + -eps_conv_abs - Sets the absolute convergence test
739: . -eps_conv_rel - Sets the convergence test relative to the eigenvalue
740: . -eps_conv_norm - Sets the convergence test relative to the matrix norms
741: - -eps_conv_user - Selects the user-defined convergence test
743: Note:
744: The parameter 'conv' can have one of these values
745: + EPS_CONV_ABS - absolute error ||r||
746: . EPS_CONV_REL - error relative to the eigenvalue l, ||r||/|l|
747: . EPS_CONV_NORM - error relative to the matrix norms, ||r||/(||A||+|l|*||B||)
748: - EPS_CONV_USER - function set by EPSSetConvergenceTestFunction()
750: Level: intermediate
752: .seealso: EPSGetConvergenceTest(), EPSSetConvergenceTestFunction(), EPSSetStoppingTest(), EPSConv
753: @*/
754: PetscErrorCode EPSSetConvergenceTest(EPS eps,EPSConv conv)
755: {
759: switch (conv) {
760: case EPS_CONV_ABS: eps->converged = EPSConvergedAbsolute; break;
761: case EPS_CONV_REL: eps->converged = EPSConvergedRelative; break;
762: case EPS_CONV_NORM: eps->converged = EPSConvergedNorm; break;
763: case EPS_CONV_USER:
764: if (!eps->convergeduser) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ORDER,"Must call EPSSetConvergenceTestFunction() first");
765: eps->converged = eps->convergeduser;
766: break;
767: default:
768: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
769: }
770: eps->conv = conv;
771: return(0);
772: }
774: /*@
775: EPSGetConvergenceTest - Gets the method used to compute the error estimate
776: used in the convergence test.
778: Not Collective
780: Input Parameters:
781: . eps - eigensolver context obtained from EPSCreate()
783: Output Parameters:
784: . conv - the type of convergence test
786: Level: intermediate
788: .seealso: EPSSetConvergenceTest(), EPSConv
789: @*/
790: PetscErrorCode EPSGetConvergenceTest(EPS eps,EPSConv *conv)
791: {
795: *conv = eps->conv;
796: return(0);
797: }
799: /*@C
800: EPSSetStoppingTestFunction - Sets a function to decide when to stop the outer
801: iteration of the eigensolver.
803: Logically Collective on EPS
805: Input Parameters:
806: + eps - eigensolver context obtained from EPSCreate()
807: . func - pointer to the stopping test function
808: . ctx - context for private data for the stopping routine (may be null)
809: - destroy - a routine for destroying the context (may be null)
811: Calling Sequence of func:
812: $ func(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx)
814: + eps - eigensolver context obtained from EPSCreate()
815: . its - current number of iterations
816: . max_it - maximum number of iterations
817: . nconv - number of currently converged eigenpairs
818: . nev - number of requested eigenpairs
819: . reason - (output) result of the stopping test
820: - ctx - optional context, as set by EPSSetStoppingTestFunction()
822: Note:
823: Normal usage is to first call the default routine EPSStoppingBasic() and then
824: set reason to EPS_CONVERGED_USER if some user-defined conditions have been
825: met. To let the eigensolver continue iterating, the result must be left as
826: EPS_CONVERGED_ITERATING.
828: Level: advanced
830: .seealso: EPSSetStoppingTest(), EPSStoppingBasic()
831: @*/
832: PetscErrorCode EPSSetStoppingTestFunction(EPS eps,PetscErrorCode (*func)(EPS,PetscInt,PetscInt,PetscInt,PetscInt,EPSConvergedReason*,void*),void* ctx,PetscErrorCode (*destroy)(void*))
833: {
838: if (eps->stoppingdestroy) {
839: (*eps->stoppingdestroy)(eps->stoppingctx);
840: }
841: eps->stoppinguser = func;
842: eps->stoppingdestroy = destroy;
843: eps->stoppingctx = ctx;
844: if (func == EPSStoppingBasic) eps->stop = EPS_STOP_BASIC;
845: else {
846: eps->stop = EPS_STOP_USER;
847: eps->stopping = eps->stoppinguser;
848: }
849: return(0);
850: }
852: /*@
853: EPSSetStoppingTest - Specifies how to decide the termination of the outer
854: loop of the eigensolver.
856: Logically Collective on EPS
858: Input Parameters:
859: + eps - eigensolver context obtained from EPSCreate()
860: - stop - the type of stopping test
862: Options Database Keys:
863: + -eps_stop_basic - Sets the default stopping test
864: - -eps_stop_user - Selects the user-defined stopping test
866: Note:
867: The parameter 'stop' can have one of these values
868: + EPS_STOP_BASIC - default stopping test
869: - EPS_STOP_USER - function set by EPSSetStoppingTestFunction()
871: Level: advanced
873: .seealso: EPSGetStoppingTest(), EPSSetStoppingTestFunction(), EPSSetConvergenceTest(), EPSStop
874: @*/
875: PetscErrorCode EPSSetStoppingTest(EPS eps,EPSStop stop)
876: {
880: switch (stop) {
881: case EPS_STOP_BASIC: eps->stopping = EPSStoppingBasic; break;
882: case EPS_STOP_USER:
883: if (!eps->stoppinguser) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ORDER,"Must call EPSSetStoppingTestFunction() first");
884: eps->stopping = eps->stoppinguser;
885: break;
886: default:
887: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
888: }
889: eps->stop = stop;
890: return(0);
891: }
893: /*@
894: EPSGetStoppingTest - Gets the method used to decide the termination of the outer
895: loop of the eigensolver.
897: Not Collective
899: Input Parameters:
900: . eps - eigensolver context obtained from EPSCreate()
902: Output Parameters:
903: . stop - the type of stopping test
905: Level: advanced
907: .seealso: EPSSetStoppingTest(), EPSStop
908: @*/
909: PetscErrorCode EPSGetStoppingTest(EPS eps,EPSStop *stop)
910: {
914: *stop = eps->stop;
915: return(0);
916: }
918: /*@
919: EPSSetProblemType - Specifies the type of the eigenvalue problem.
921: Logically Collective on EPS
923: Input Parameters:
924: + eps - the eigensolver context
925: - type - a known type of eigenvalue problem
927: Options Database Keys:
928: + -eps_hermitian - Hermitian eigenvalue problem
929: . -eps_gen_hermitian - generalized Hermitian eigenvalue problem
930: . -eps_non_hermitian - non-Hermitian eigenvalue problem
931: . -eps_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
932: - -eps_pos_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
933: with positive semi-definite B
935: Notes:
936: Allowed values for the problem type are: Hermitian (EPS_HEP), non-Hermitian
937: (EPS_NHEP), generalized Hermitian (EPS_GHEP), generalized non-Hermitian
938: (EPS_GNHEP), generalized non-Hermitian with positive semi-definite B
939: (EPS_PGNHEP), and generalized Hermitian-indefinite (EPS_GHIEP).
941: This function must be used to instruct SLEPc to exploit symmetry. If no
942: problem type is specified, by default a non-Hermitian problem is assumed
943: (either standard or generalized). If the user knows that the problem is
944: Hermitian (i.e. A=A^H) or generalized Hermitian (i.e. A=A^H, B=B^H, and
945: B positive definite) then it is recommended to set the problem type so
946: that eigensolver can exploit these properties.
948: Level: intermediate
950: .seealso: EPSSetOperators(), EPSSetType(), EPSGetProblemType(), EPSProblemType
951: @*/
952: PetscErrorCode EPSSetProblemType(EPS eps,EPSProblemType type)
953: {
957: if (type == eps->problem_type) return(0);
958: switch (type) {
959: case EPS_HEP:
960: eps->isgeneralized = PETSC_FALSE;
961: eps->ishermitian = PETSC_TRUE;
962: eps->ispositive = PETSC_FALSE;
963: break;
964: case EPS_NHEP:
965: eps->isgeneralized = PETSC_FALSE;
966: eps->ishermitian = PETSC_FALSE;
967: eps->ispositive = PETSC_FALSE;
968: break;
969: case EPS_GHEP:
970: eps->isgeneralized = PETSC_TRUE;
971: eps->ishermitian = PETSC_TRUE;
972: eps->ispositive = PETSC_TRUE;
973: break;
974: case EPS_GNHEP:
975: eps->isgeneralized = PETSC_TRUE;
976: eps->ishermitian = PETSC_FALSE;
977: eps->ispositive = PETSC_FALSE;
978: break;
979: case EPS_PGNHEP:
980: eps->isgeneralized = PETSC_TRUE;
981: eps->ishermitian = PETSC_FALSE;
982: eps->ispositive = PETSC_TRUE;
983: break;
984: case EPS_GHIEP:
985: eps->isgeneralized = PETSC_TRUE;
986: eps->ishermitian = PETSC_TRUE;
987: eps->ispositive = PETSC_FALSE;
988: break;
989: default:
990: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
991: }
992: eps->problem_type = type;
993: eps->state = EPS_STATE_INITIAL;
994: return(0);
995: }
997: /*@
998: EPSGetProblemType - Gets the problem type from the EPS object.
1000: Not Collective
1002: Input Parameter:
1003: . eps - the eigensolver context
1005: Output Parameter:
1006: . type - the problem type
1008: Level: intermediate
1010: .seealso: EPSSetProblemType(), EPSProblemType
1011: @*/
1012: PetscErrorCode EPSGetProblemType(EPS eps,EPSProblemType *type)
1013: {
1017: *type = eps->problem_type;
1018: return(0);
1019: }
1021: /*@
1022: EPSSetExtraction - Specifies the type of extraction technique to be employed
1023: by the eigensolver.
1025: Logically Collective on EPS
1027: Input Parameters:
1028: + eps - the eigensolver context
1029: - extr - a known type of extraction
1031: Options Database Keys:
1032: + -eps_ritz - Rayleigh-Ritz extraction
1033: . -eps_harmonic - harmonic Ritz extraction
1034: . -eps_harmonic_relative - harmonic Ritz extraction relative to the eigenvalue
1035: . -eps_harmonic_right - harmonic Ritz extraction for rightmost eigenvalues
1036: . -eps_harmonic_largest - harmonic Ritz extraction for largest magnitude
1037: (without target)
1038: . -eps_refined - refined Ritz extraction
1039: - -eps_refined_harmonic - refined harmonic Ritz extraction
1041: Notes:
1042: Not all eigensolvers support all types of extraction. See the SLEPc
1043: Users Manual for details.
1045: By default, a standard Rayleigh-Ritz extraction is used. Other extractions
1046: may be useful when computing interior eigenvalues.
1048: Harmonic-type extractions are used in combination with a 'target'.
1050: Level: advanced
1052: .seealso: EPSSetTarget(), EPSGetExtraction(), EPSExtraction
1053: @*/
1054: PetscErrorCode EPSSetExtraction(EPS eps,EPSExtraction extr)
1055: {
1059: eps->extraction = extr;
1060: return(0);
1061: }
1063: /*@
1064: EPSGetExtraction - Gets the extraction type used by the EPS object.
1066: Not Collective
1068: Input Parameter:
1069: . eps - the eigensolver context
1071: Output Parameter:
1072: . extr - name of extraction type
1074: Level: advanced
1076: .seealso: EPSSetExtraction(), EPSExtraction
1077: @*/
1078: PetscErrorCode EPSGetExtraction(EPS eps,EPSExtraction *extr)
1079: {
1083: *extr = eps->extraction;
1084: return(0);
1085: }
1087: /*@
1088: EPSSetBalance - Specifies the balancing technique to be employed by the
1089: eigensolver, and some parameters associated to it.
1091: Logically Collective on EPS
1093: Input Parameters:
1094: + eps - the eigensolver context
1095: . bal - the balancing method, one of EPS_BALANCE_NONE, EPS_BALANCE_ONESIDE,
1096: EPS_BALANCE_TWOSIDE, or EPS_BALANCE_USER
1097: . its - number of iterations of the balancing algorithm
1098: - cutoff - cutoff value
1100: Options Database Keys:
1101: + -eps_balance <method> - the balancing method, where <method> is one of
1102: 'none', 'oneside', 'twoside', or 'user'
1103: . -eps_balance_its <its> - number of iterations
1104: - -eps_balance_cutoff <cutoff> - cutoff value
1106: Notes:
1107: When balancing is enabled, the solver works implicitly with matrix DAD^-1,
1108: where D is an appropriate diagonal matrix. This improves the accuracy of
1109: the computed results in some cases. See the SLEPc Users Manual for details.
1111: Balancing makes sense only for non-Hermitian problems when the required
1112: precision is high (i.e. a small tolerance such as 1e-15).
1114: By default, balancing is disabled. The two-sided method is much more
1115: effective than the one-sided counterpart, but it requires the system
1116: matrices to have the MatMultTranspose operation defined.
1118: The parameter 'its' is the number of iterations performed by the method. The
1119: cutoff value is used only in the two-side variant. Use PETSC_DEFAULT to assign
1120: a reasonably good value.
1122: User-defined balancing is allowed provided that the corresponding matrix
1123: is set via STSetBalanceMatrix.
1125: Level: intermediate
1127: .seealso: EPSGetBalance(), EPSBalance, STSetBalanceMatrix()
1128: @*/
1129: PetscErrorCode EPSSetBalance(EPS eps,EPSBalance bal,PetscInt its,PetscReal cutoff)
1130: {
1136: switch (bal) {
1137: case EPS_BALANCE_NONE:
1138: case EPS_BALANCE_ONESIDE:
1139: case EPS_BALANCE_TWOSIDE:
1140: case EPS_BALANCE_USER:
1141: if (eps->balance != bal) {
1142: eps->state = EPS_STATE_INITIAL;
1143: eps->balance = bal;
1144: }
1145: break;
1146: default:
1147: SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'bal'");
1148: }
1149: if (its==PETSC_DECIDE || its==PETSC_DEFAULT) eps->balance_its = 5;
1150: else {
1151: if (its <= 0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be > 0");
1152: eps->balance_its = its;
1153: }
1154: if (cutoff==PETSC_DECIDE || cutoff==PETSC_DEFAULT) eps->balance_cutoff = 1e-8;
1155: else {
1156: if (cutoff <= 0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of cutoff. Must be > 0");
1157: eps->balance_cutoff = cutoff;
1158: }
1159: return(0);
1160: }
1162: /*@C
1163: EPSGetBalance - Gets the balancing type used by the EPS object, and the
1164: associated parameters.
1166: Not Collective
1168: Input Parameter:
1169: . eps - the eigensolver context
1171: Output Parameters:
1172: + bal - the balancing method
1173: . its - number of iterations of the balancing algorithm
1174: - cutoff - cutoff value
1176: Level: intermediate
1178: Note:
1179: The user can specify NULL for any parameter that is not needed.
1181: .seealso: EPSSetBalance(), EPSBalance
1182: @*/
1183: PetscErrorCode EPSGetBalance(EPS eps,EPSBalance *bal,PetscInt *its,PetscReal *cutoff)
1184: {
1187: if (bal) *bal = eps->balance;
1188: if (its) *its = eps->balance_its;
1189: if (cutoff) *cutoff = eps->balance_cutoff;
1190: return(0);
1191: }
1193: /*@
1194: EPSSetTwoSided - Sets the solver to use a two-sided variant so that left
1195: eigenvectors are also computed.
1197: Logically Collective on EPS
1199: Input Parameters:
1200: + eps - the eigensolver context
1201: - twosided - whether the two-sided variant is to be used or not
1203: Options Database Keys:
1204: . -eps_two_sided <boolean> - Sets/resets the twosided flag
1206: Notes:
1207: If the user sets twosided=PETSC_TRUE then the solver uses a variant of
1208: the algorithm that computes both right and left eigenvectors. This is
1209: usually much more costly. This option is not available in all solvers.
1211: When using two-sided solvers, the problem matrices must have both the
1212: MatMult and MatMultTranspose operations defined.
1214: Level: advanced
1216: .seealso: EPSGetTwoSided(), EPSGetLeftEigenvector()
1217: @*/
1218: PetscErrorCode EPSSetTwoSided(EPS eps,PetscBool twosided)
1219: {
1223: if (twosided!=eps->twosided) {
1224: eps->twosided = twosided;
1225: eps->state = EPS_STATE_INITIAL;
1226: }
1227: return(0);
1228: }
1230: /*@
1231: EPSGetTwoSided - Returns the flag indicating whether a two-sided variant
1232: of the algorithm is being used or not.
1234: Not Collective
1236: Input Parameter:
1237: . eps - the eigensolver context
1239: Output Parameter:
1240: . twosided - the returned flag
1242: Level: advanced
1244: .seealso: EPSSetTwoSided()
1245: @*/
1246: PetscErrorCode EPSGetTwoSided(EPS eps,PetscBool *twosided)
1247: {
1251: *twosided = eps->twosided;
1252: return(0);
1253: }
1255: /*@
1256: EPSSetTrueResidual - Specifies if the solver must compute the true residual
1257: explicitly or not.
1259: Logically Collective on EPS
1261: Input Parameters:
1262: + eps - the eigensolver context
1263: - trueres - whether true residuals are required or not
1265: Options Database Keys:
1266: . -eps_true_residual <boolean> - Sets/resets the boolean flag 'trueres'
1268: Notes:
1269: If the user sets trueres=PETSC_TRUE then the solver explicitly computes
1270: the true residual for each eigenpair approximation, and uses it for
1271: convergence testing. Computing the residual is usually an expensive
1272: operation. Some solvers (e.g., Krylov solvers) can avoid this computation
1273: by using a cheap estimate of the residual norm, but this may sometimes
1274: give inaccurate results (especially if a spectral transform is being
1275: used). On the contrary, preconditioned eigensolvers (e.g., Davidson solvers)
1276: do rely on computing the true residual, so this option is irrelevant for them.
1278: Level: advanced
1280: .seealso: EPSGetTrueResidual()
1281: @*/
1282: PetscErrorCode EPSSetTrueResidual(EPS eps,PetscBool trueres)
1283: {
1287: eps->trueres = trueres;
1288: return(0);
1289: }
1291: /*@
1292: EPSGetTrueResidual - Returns the flag indicating whether true
1293: residuals must be computed explicitly or not.
1295: Not Collective
1297: Input Parameter:
1298: . eps - the eigensolver context
1300: Output Parameter:
1301: . trueres - the returned flag
1303: Level: advanced
1305: .seealso: EPSSetTrueResidual()
1306: @*/
1307: PetscErrorCode EPSGetTrueResidual(EPS eps,PetscBool *trueres)
1308: {
1312: *trueres = eps->trueres;
1313: return(0);
1314: }
1316: /*@
1317: EPSSetTrackAll - Specifies if the solver must compute the residual norm of all
1318: approximate eigenpairs or not.
1320: Logically Collective on EPS
1322: Input Parameters:
1323: + eps - the eigensolver context
1324: - trackall - whether to compute all residuals or not
1326: Notes:
1327: If the user sets trackall=PETSC_TRUE then the solver computes (or estimates)
1328: the residual norm for each eigenpair approximation. Computing the residual is
1329: usually an expensive operation and solvers commonly compute only the residual
1330: associated to the first unconverged eigenpair.
1332: The options '-eps_monitor_all' and '-eps_monitor_lg_all' automatically
1333: activate this option.
1335: Level: developer
1337: .seealso: EPSGetTrackAll()
1338: @*/
1339: PetscErrorCode EPSSetTrackAll(EPS eps,PetscBool trackall)
1340: {
1344: eps->trackall = trackall;
1345: return(0);
1346: }
1348: /*@
1349: EPSGetTrackAll - Returns the flag indicating whether all residual norms must
1350: be computed or not.
1352: Not Collective
1354: Input Parameter:
1355: . eps - the eigensolver context
1357: Output Parameter:
1358: . trackall - the returned flag
1360: Level: developer
1362: .seealso: EPSSetTrackAll()
1363: @*/
1364: PetscErrorCode EPSGetTrackAll(EPS eps,PetscBool *trackall)
1365: {
1369: *trackall = eps->trackall;
1370: return(0);
1371: }
1373: /*@
1374: EPSSetPurify - Deactivate eigenvector purification (which is activated by default).
1376: Logically Collective on EPS
1378: Input Parameters:
1379: + eps - the eigensolver context
1380: - purify - whether purification is required or not
1382: Options Database Keys:
1383: . -eps_purify <boolean> - Sets/resets the boolean flag 'purify'
1385: Notes:
1386: By default, eigenvectors of generalized symmetric eigenproblems are purified
1387: in order to purge directions in the nullspace of matrix B. If the user knows
1388: that B is non-singular, then purification can be safely deactivated and some
1389: computational cost is avoided (this is particularly important in interval computations).
1391: Level: intermediate
1393: .seealso: EPSGetPurify(), EPSSetInterval()
1394: @*/
1395: PetscErrorCode EPSSetPurify(EPS eps,PetscBool purify)
1396: {
1400: if (purify!=eps->purify) {
1401: eps->purify = purify;
1402: eps->state = EPS_STATE_INITIAL;
1403: }
1404: return(0);
1405: }
1407: /*@
1408: EPSGetPurify - Returns the flag indicating whether purification is activated
1409: or not.
1411: Not Collective
1413: Input Parameter:
1414: . eps - the eigensolver context
1416: Output Parameter:
1417: . purify - the returned flag
1419: Level: intermediate
1421: .seealso: EPSSetPurify()
1422: @*/
1423: PetscErrorCode EPSGetPurify(EPS eps,PetscBool *purify)
1424: {
1428: *purify = eps->purify;
1429: return(0);
1430: }
1432: /*@C
1433: EPSSetOptionsPrefix - Sets the prefix used for searching for all
1434: EPS options in the database.
1436: Logically Collective on EPS
1438: Input Parameters:
1439: + eps - the eigensolver context
1440: - prefix - the prefix string to prepend to all EPS option requests
1442: Notes:
1443: A hyphen (-) must NOT be given at the beginning of the prefix name.
1444: The first character of all runtime options is AUTOMATICALLY the
1445: hyphen.
1447: For example, to distinguish between the runtime options for two
1448: different EPS contexts, one could call
1449: .vb
1450: EPSSetOptionsPrefix(eps1,"eig1_")
1451: EPSSetOptionsPrefix(eps2,"eig2_")
1452: .ve
1454: Level: advanced
1456: .seealso: EPSAppendOptionsPrefix(), EPSGetOptionsPrefix()
1457: @*/
1458: PetscErrorCode EPSSetOptionsPrefix(EPS eps,const char *prefix)
1459: {
1464: if (!eps->st) { EPSGetST(eps,&eps->st); }
1465: STSetOptionsPrefix(eps->st,prefix);
1466: if (!eps->V) { EPSGetBV(eps,&eps->V); }
1467: BVSetOptionsPrefix(eps->V,prefix);
1468: if (!eps->ds) { EPSGetDS(eps,&eps->ds); }
1469: DSSetOptionsPrefix(eps->ds,prefix);
1470: if (!eps->rg) { EPSGetRG(eps,&eps->rg); }
1471: RGSetOptionsPrefix(eps->rg,prefix);
1472: PetscObjectSetOptionsPrefix((PetscObject)eps,prefix);
1473: return(0);
1474: }
1476: /*@C
1477: EPSAppendOptionsPrefix - Appends to the prefix used for searching for all
1478: EPS options in the database.
1480: Logically Collective on EPS
1482: Input Parameters:
1483: + eps - the eigensolver context
1484: - prefix - the prefix string to prepend to all EPS option requests
1486: Notes:
1487: A hyphen (-) must NOT be given at the beginning of the prefix name.
1488: The first character of all runtime options is AUTOMATICALLY the hyphen.
1490: Level: advanced
1492: .seealso: EPSSetOptionsPrefix(), EPSGetOptionsPrefix()
1493: @*/
1494: PetscErrorCode EPSAppendOptionsPrefix(EPS eps,const char *prefix)
1495: {
1500: if (!eps->st) { EPSGetST(eps,&eps->st); }
1501: STAppendOptionsPrefix(eps->st,prefix);
1502: if (!eps->V) { EPSGetBV(eps,&eps->V); }
1503: BVAppendOptionsPrefix(eps->V,prefix);
1504: if (!eps->ds) { EPSGetDS(eps,&eps->ds); }
1505: DSAppendOptionsPrefix(eps->ds,prefix);
1506: if (!eps->rg) { EPSGetRG(eps,&eps->rg); }
1507: RGAppendOptionsPrefix(eps->rg,prefix);
1508: PetscObjectAppendOptionsPrefix((PetscObject)eps,prefix);
1509: return(0);
1510: }
1512: /*@C
1513: EPSGetOptionsPrefix - Gets the prefix used for searching for all
1514: EPS options in the database.
1516: Not Collective
1518: Input Parameters:
1519: . eps - the eigensolver context
1521: Output Parameters:
1522: . prefix - pointer to the prefix string used is returned
1524: Note:
1525: On the Fortran side, the user should pass in a string 'prefix' of
1526: sufficient length to hold the prefix.
1528: Level: advanced
1530: .seealso: EPSSetOptionsPrefix(), EPSAppendOptionsPrefix()
1531: @*/
1532: PetscErrorCode EPSGetOptionsPrefix(EPS eps,const char *prefix[])
1533: {
1539: PetscObjectGetOptionsPrefix((PetscObject)eps,prefix);
1540: return(0);
1541: }