Actual source code: epsopts.c

slepc-3.22.2 2024-12-02
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    EPS routines related to options that can be set via the command-line
 12:    or procedurally.
 13: */

 15: #include <slepc/private/epsimpl.h>
 16: #include <petscdraw.h>

 18: /*@C
 19:    EPSMonitorSetFromOptions - Sets a monitor function and viewer appropriate for the type
 20:    indicated by the user.

 22:    Collective

 24:    Input Parameters:
 25: +  eps      - the 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;

 47:   PetscFunctionBegin;
 48:   PetscCall(PetscOptionsCreateViewer(PetscObjectComm((PetscObject)eps),((PetscObject)eps)->options,((PetscObject)eps)->prefix,opt,&viewer,&format,&flg));
 49:   if (!flg) PetscFunctionReturn(PETSC_SUCCESS);

 51:   PetscCall(PetscViewerGetType(viewer,&vtype));
 52:   PetscCall(SlepcMonitorMakeKey_Internal(name,vtype,format,key));
 53:   PetscCall(PetscFunctionListFind(EPSMonitorList,key,&mfunc));
 54:   PetscCheck(mfunc,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Specified viewer and format not supported");
 55:   PetscCall(PetscFunctionListFind(EPSMonitorCreateList,key,&cfunc));
 56:   PetscCall(PetscFunctionListFind(EPSMonitorDestroyList,key,&dfunc));
 57:   if (!cfunc) cfunc = PetscViewerAndFormatCreate_Internal;
 58:   if (!dfunc) dfunc = PetscViewerAndFormatDestroy;

 60:   PetscCall((*cfunc)(viewer,format,ctx,&vf));
 61:   PetscCall(PetscViewerDestroy(&viewer));
 62:   PetscCall(EPSMonitorSet(eps,mfunc,vf,(PetscErrorCode(*)(void **))dfunc));
 63:   if (trackall) PetscCall(EPSSetTrackAll(eps,PETSC_TRUE));
 64:   PetscFunctionReturn(PETSC_SUCCESS);
 65: }

 67: /*@
 68:    EPSSetFromOptions - Sets EPS options from the options database.
 69:    This routine must be called before EPSSetUp() if the user is to be
 70:    allowed to set the solver type.

 72:    Collective

 74:    Input Parameters:
 75: .  eps - the eigensolver context

 77:    Notes:
 78:    To see all options, run your program with the -help option.

 80:    Level: beginner

 82: .seealso: EPSSetOptionsPrefix()
 83: @*/
 84: PetscErrorCode EPSSetFromOptions(EPS eps)
 85: {
 86:   char           type[256];
 87:   PetscBool      set,flg,flg1,flg2,flg3,bval;
 88:   PetscReal      r,array[2]={0,0};
 89:   PetscScalar    s;
 90:   PetscInt       i,j,k;
 91:   EPSBalance     bal;

 93:   PetscFunctionBegin;
 95:   PetscCall(EPSRegisterAll());
 96:   PetscObjectOptionsBegin((PetscObject)eps);
 97:     PetscCall(PetscOptionsFList("-eps_type","Eigensolver method","EPSSetType",EPSList,(char*)(((PetscObject)eps)->type_name?((PetscObject)eps)->type_name:EPSKRYLOVSCHUR),type,sizeof(type),&flg));
 98:     if (flg) PetscCall(EPSSetType(eps,type));
 99:     else if (!((PetscObject)eps)->type_name) PetscCall(EPSSetType(eps,EPSKRYLOVSCHUR));

101:     PetscCall(PetscOptionsBoolGroupBegin("-eps_hermitian","Hermitian eigenvalue problem","EPSSetProblemType",&flg));
102:     if (flg) PetscCall(EPSSetProblemType(eps,EPS_HEP));
103:     PetscCall(PetscOptionsBoolGroup("-eps_gen_hermitian","Generalized Hermitian eigenvalue problem","EPSSetProblemType",&flg));
104:     if (flg) PetscCall(EPSSetProblemType(eps,EPS_GHEP));
105:     PetscCall(PetscOptionsBoolGroup("-eps_non_hermitian","Non-Hermitian eigenvalue problem","EPSSetProblemType",&flg));
106:     if (flg) PetscCall(EPSSetProblemType(eps,EPS_NHEP));
107:     PetscCall(PetscOptionsBoolGroup("-eps_gen_non_hermitian","Generalized non-Hermitian eigenvalue problem","EPSSetProblemType",&flg));
108:     if (flg) PetscCall(EPSSetProblemType(eps,EPS_GNHEP));
109:     PetscCall(PetscOptionsBoolGroup("-eps_pos_gen_non_hermitian","Generalized non-Hermitian eigenvalue problem with positive semi-definite B","EPSSetProblemType",&flg));
110:     if (flg) PetscCall(EPSSetProblemType(eps,EPS_PGNHEP));
111:     PetscCall(PetscOptionsBoolGroup("-eps_gen_indefinite","Generalized Hermitian-indefinite eigenvalue problem","EPSSetProblemType",&flg));
112:     if (flg) PetscCall(EPSSetProblemType(eps,EPS_GHIEP));
113:     PetscCall(PetscOptionsBoolGroupEnd("-eps_bse","Structured Bethe-Salpeter eigenvalue problem","EPSSetProblemType",&flg));
114:     if (flg) PetscCall(EPSSetProblemType(eps,EPS_BSE));

116:     PetscCall(PetscOptionsBoolGroupBegin("-eps_ritz","Rayleigh-Ritz extraction","EPSSetExtraction",&flg));
117:     if (flg) PetscCall(EPSSetExtraction(eps,EPS_RITZ));
118:     PetscCall(PetscOptionsBoolGroup("-eps_harmonic","Harmonic Ritz extraction","EPSSetExtraction",&flg));
119:     if (flg) PetscCall(EPSSetExtraction(eps,EPS_HARMONIC));
120:     PetscCall(PetscOptionsBoolGroup("-eps_harmonic_relative","Relative harmonic Ritz extraction","EPSSetExtraction",&flg));
121:     if (flg) PetscCall(EPSSetExtraction(eps,EPS_HARMONIC_RELATIVE));
122:     PetscCall(PetscOptionsBoolGroup("-eps_harmonic_right","Right harmonic Ritz extraction","EPSSetExtraction",&flg));
123:     if (flg) PetscCall(EPSSetExtraction(eps,EPS_HARMONIC_RIGHT));
124:     PetscCall(PetscOptionsBoolGroup("-eps_harmonic_largest","Largest harmonic Ritz extraction","EPSSetExtraction",&flg));
125:     if (flg) PetscCall(EPSSetExtraction(eps,EPS_HARMONIC_LARGEST));
126:     PetscCall(PetscOptionsBoolGroup("-eps_refined","Refined Ritz extraction","EPSSetExtraction",&flg));
127:     if (flg) PetscCall(EPSSetExtraction(eps,EPS_REFINED));
128:     PetscCall(PetscOptionsBoolGroupEnd("-eps_refined_harmonic","Refined harmonic Ritz extraction","EPSSetExtraction",&flg));
129:     if (flg) PetscCall(EPSSetExtraction(eps,EPS_REFINED_HARMONIC));

131:     bal = eps->balance;
132:     PetscCall(PetscOptionsEnum("-eps_balance","Balancing method","EPSSetBalance",EPSBalanceTypes,(PetscEnum)bal,(PetscEnum*)&bal,&flg1));
133:     j = eps->balance_its;
134:     PetscCall(PetscOptionsInt("-eps_balance_its","Number of iterations in balancing","EPSSetBalance",eps->balance_its,&j,&flg2));
135:     r = eps->balance_cutoff;
136:     PetscCall(PetscOptionsReal("-eps_balance_cutoff","Cutoff value in balancing","EPSSetBalance",eps->balance_cutoff,&r,&flg3));
137:     if (flg1 || flg2 || flg3) PetscCall(EPSSetBalance(eps,bal,j,r));

139:     i = eps->max_it;
140:     PetscCall(PetscOptionsInt("-eps_max_it","Maximum number of iterations","EPSSetTolerances",eps->max_it,&i,&flg1));
141:     r = eps->tol;
142:     PetscCall(PetscOptionsReal("-eps_tol","Tolerance","EPSSetTolerances",SlepcDefaultTol(eps->tol),&r,&flg2));
143:     if (flg1 || flg2) PetscCall(EPSSetTolerances(eps,r,i));

145:     PetscCall(PetscOptionsBoolGroupBegin("-eps_conv_rel","Relative error convergence test","EPSSetConvergenceTest",&flg));
146:     if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_REL));
147:     PetscCall(PetscOptionsBoolGroup("-eps_conv_norm","Convergence test relative to the eigenvalue and the matrix norms","EPSSetConvergenceTest",&flg));
148:     if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_NORM));
149:     PetscCall(PetscOptionsBoolGroup("-eps_conv_abs","Absolute error convergence test","EPSSetConvergenceTest",&flg));
150:     if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_ABS));
151:     PetscCall(PetscOptionsBoolGroupEnd("-eps_conv_user","User-defined convergence test","EPSSetConvergenceTest",&flg));
152:     if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_USER));

154:     PetscCall(PetscOptionsBoolGroupBegin("-eps_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","EPSSetStoppingTest",&flg));
155:     if (flg) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_BASIC));
156:     PetscCall(PetscOptionsBoolGroupEnd("-eps_stop_user","User-defined stopping test","EPSSetStoppingTest",&flg));
157:     if (flg) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_USER));

159:     i = eps->nev;
160:     PetscCall(PetscOptionsInt("-eps_nev","Number of eigenvalues to compute","EPSSetDimensions",eps->nev,&i,&flg1));
161:     j = eps->ncv;
162:     PetscCall(PetscOptionsInt("-eps_ncv","Number of basis vectors","EPSSetDimensions",eps->ncv,&j,&flg2));
163:     k = eps->mpd;
164:     PetscCall(PetscOptionsInt("-eps_mpd","Maximum dimension of projected problem","EPSSetDimensions",eps->mpd,&k,&flg3));
165:     if (flg1 || flg2 || flg3) PetscCall(EPSSetDimensions(eps,i,j,k));

167:     PetscCall(PetscOptionsBoolGroupBegin("-eps_largest_magnitude","Compute largest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg));
168:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE));
169:     PetscCall(PetscOptionsBoolGroup("-eps_smallest_magnitude","Compute smallest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg));
170:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_MAGNITUDE));
171:     PetscCall(PetscOptionsBoolGroup("-eps_largest_real","Compute eigenvalues with largest real parts","EPSSetWhichEigenpairs",&flg));
172:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL));
173:     PetscCall(PetscOptionsBoolGroup("-eps_smallest_real","Compute eigenvalues with smallest real parts","EPSSetWhichEigenpairs",&flg));
174:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
175:     PetscCall(PetscOptionsBoolGroup("-eps_largest_imaginary","Compute eigenvalues with largest imaginary parts","EPSSetWhichEigenpairs",&flg));
176:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_IMAGINARY));
177:     PetscCall(PetscOptionsBoolGroup("-eps_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","EPSSetWhichEigenpairs",&flg));
178:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_IMAGINARY));
179:     PetscCall(PetscOptionsBoolGroup("-eps_target_magnitude","Compute eigenvalues closest to target","EPSSetWhichEigenpairs",&flg));
180:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
181:     PetscCall(PetscOptionsBoolGroup("-eps_target_real","Compute eigenvalues with real parts closest to target","EPSSetWhichEigenpairs",&flg));
182:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_REAL));
183:     PetscCall(PetscOptionsBoolGroup("-eps_target_imaginary","Compute eigenvalues with imaginary parts closest to target","EPSSetWhichEigenpairs",&flg));
184:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_IMAGINARY));
185:     PetscCall(PetscOptionsBoolGroupEnd("-eps_all","Compute all eigenvalues in an interval or a region","EPSSetWhichEigenpairs",&flg));
186:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_ALL));

188:     PetscCall(PetscOptionsScalar("-eps_target","Value of the target","EPSSetTarget",eps->target,&s,&flg));
189:     if (flg) {
190:       if (eps->which!=EPS_TARGET_REAL && eps->which!=EPS_TARGET_IMAGINARY) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
191:       PetscCall(EPSSetTarget(eps,s));
192:     }

194:     k = 2;
195:     PetscCall(PetscOptionsRealArray("-eps_interval","Computational interval (two real values separated with a comma without spaces)","EPSSetInterval",array,&k,&flg));
196:     if (flg) {
197:       PetscCheck(k>1,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_SIZ,"Must pass two values in -eps_interval (comma-separated without spaces)");
198:       PetscCall(EPSSetWhichEigenpairs(eps,EPS_ALL));
199:       PetscCall(EPSSetInterval(eps,array[0],array[1]));
200:     }

202:     PetscCall(PetscOptionsBool("-eps_true_residual","Compute true residuals explicitly","EPSSetTrueResidual",eps->trueres,&eps->trueres,NULL));
203:     PetscCall(PetscOptionsBool("-eps_purify","Postprocess eigenvectors for purification","EPSSetPurify",eps->purify,&bval,&flg));
204:     if (flg) PetscCall(EPSSetPurify(eps,bval));
205:     PetscCall(PetscOptionsBool("-eps_two_sided","Use two-sided variant (to compute left eigenvectors)","EPSSetTwoSided",eps->twosided,&bval,&flg));
206:     if (flg) PetscCall(EPSSetTwoSided(eps,bval));

208:     /* -----------------------------------------------------------------------*/
209:     /*
210:       Cancels all monitors hardwired into code before call to EPSSetFromOptions()
211:     */
212:     PetscCall(PetscOptionsBool("-eps_monitor_cancel","Remove any hardwired monitor routines","EPSMonitorCancel",PETSC_FALSE,&flg,&set));
213:     if (set && flg) PetscCall(EPSMonitorCancel(eps));
214:     PetscCall(EPSMonitorSetFromOptions(eps,"-eps_monitor","first_approximation",NULL,PETSC_FALSE));
215:     PetscCall(EPSMonitorSetFromOptions(eps,"-eps_monitor_all","all_approximations",NULL,PETSC_TRUE));
216:     PetscCall(EPSMonitorSetFromOptions(eps,"-eps_monitor_conv","convergence_history",NULL,PETSC_FALSE));

218:     /* -----------------------------------------------------------------------*/
219:     PetscCall(PetscOptionsName("-eps_view","Print detailed information on solver used","EPSView",&set));
220:     PetscCall(PetscOptionsName("-eps_view_vectors","View computed eigenvectors","EPSVectorsView",&set));
221:     PetscCall(PetscOptionsName("-eps_view_values","View computed eigenvalues","EPSValuesView",&set));
222:     PetscCall(PetscOptionsName("-eps_converged_reason","Print reason for convergence, and number of iterations","EPSConvergedReasonView",&set));
223:     PetscCall(PetscOptionsName("-eps_error_absolute","Print absolute errors of each eigenpair","EPSErrorView",&set));
224:     PetscCall(PetscOptionsName("-eps_error_relative","Print relative errors of each eigenpair","EPSErrorView",&set));
225:     PetscCall(PetscOptionsName("-eps_error_backward","Print backward errors of each eigenpair","EPSErrorView",&set));

227:     PetscTryTypeMethod(eps,setfromoptions,PetscOptionsObject);
228:     PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)eps,PetscOptionsObject));
229:   PetscOptionsEnd();

231:   if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
232:   PetscCall(BVSetFromOptions(eps->V));
233:   if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
234:   PetscCall(RGSetFromOptions(eps->rg));
235:   if (eps->useds) {
236:     if (!eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
237:     PetscCall(EPSSetDSType(eps));
238:     PetscCall(DSSetFromOptions(eps->ds));
239:   }
240:   if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
241:   PetscCall(EPSSetDefaultST(eps));
242:   PetscCall(STSetFromOptions(eps->st));
243:   PetscFunctionReturn(PETSC_SUCCESS);
244: }

246: /*@
247:    EPSGetTolerances - Gets the tolerance and maximum iteration count used
248:    by the EPS convergence tests.

250:    Not Collective

252:    Input Parameter:
253: .  eps - the eigensolver context

255:    Output Parameters:
256: +  tol - the convergence tolerance
257: -  maxits - maximum number of iterations

259:    Notes:
260:    The user can specify NULL for any parameter that is not needed.

262:    Level: intermediate

264: .seealso: EPSSetTolerances()
265: @*/
266: PetscErrorCode EPSGetTolerances(EPS eps,PetscReal *tol,PetscInt *maxits)
267: {
268:   PetscFunctionBegin;
270:   if (tol)    *tol    = eps->tol;
271:   if (maxits) *maxits = eps->max_it;
272:   PetscFunctionReturn(PETSC_SUCCESS);
273: }

275: /*@
276:    EPSSetTolerances - Sets the tolerance and maximum iteration count used
277:    by the EPS convergence tests.

279:    Logically Collective

281:    Input Parameters:
282: +  eps - the eigensolver context
283: .  tol - the convergence tolerance
284: -  maxits - maximum number of iterations to use

286:    Options Database Keys:
287: +  -eps_tol <tol> - Sets the convergence tolerance
288: -  -eps_max_it <maxits> - Sets the maximum number of iterations allowed

290:    Notes:
291:    Use PETSC_CURRENT to retain the current value of any of the parameters.
292:    Use PETSC_DETERMINE for either argument to assign a default value computed
293:    internally (may be different in each solver).
294:    For maxits use PETSC_UMLIMITED to indicate there is no upper bound on this value.

296:    Level: intermediate

298: .seealso: EPSGetTolerances()
299: @*/
300: PetscErrorCode EPSSetTolerances(EPS eps,PetscReal tol,PetscInt maxits)
301: {
302:   PetscFunctionBegin;
306:   if (tol == (PetscReal)PETSC_DETERMINE) {
307:     eps->tol   = PETSC_DETERMINE;
308:     eps->state = EPS_STATE_INITIAL;
309:   } else if (tol != (PetscReal)PETSC_CURRENT) {
310:     PetscCheck(tol>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
311:     eps->tol = tol;
312:   }
313:   if (maxits == PETSC_DETERMINE) {
314:     eps->max_it = PETSC_DETERMINE;
315:     eps->state  = EPS_STATE_INITIAL;
316:   } else if (maxits == PETSC_UNLIMITED) {
317:     eps->max_it = PETSC_INT_MAX;
318:   } else if (maxits != PETSC_CURRENT) {
319:     PetscCheck(maxits>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
320:     eps->max_it = maxits;
321:   }
322:   PetscFunctionReturn(PETSC_SUCCESS);
323: }

325: /*@
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: {
345:   PetscFunctionBegin;
347:   if (nev) *nev = eps->nev;
348:   if (ncv) *ncv = eps->ncv;
349:   if (mpd) *mpd = eps->mpd;
350:   PetscFunctionReturn(PETSC_SUCCESS);
351: }

353: /*@
354:    EPSSetDimensions - Sets the number of eigenvalues to compute
355:    and the dimension of the subspace.

357:    Logically Collective

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_DETERMINE for ncv and mpd to assign a reasonably good value, which is
372:    dependent on the solution method. For any of the arguments, use PETSC_CURRENT
373:    to preserve the current value.

375:    The parameters ncv and mpd are intimately related, so that the user is advised
376:    to set one of them at most. Normal usage is that
377:    (a) in cases where nev is small, the user sets ncv (a reasonable default is 2*nev); and
378:    (b) in cases where nev is large, the user sets mpd.

380:    The value of ncv should always be between nev and (nev+mpd), typically
381:    ncv=nev+mpd. If nev is not too large, mpd=nev is a reasonable choice, otherwise
382:    a smaller value should be used.

384:    When computing all eigenvalues in an interval, see EPSSetInterval(), these
385:    parameters lose relevance, and tuning must be done with
386:    EPSKrylovSchurSetDimensions().

388:    Level: intermediate

390: .seealso: EPSGetDimensions(), EPSSetInterval(), EPSKrylovSchurSetDimensions()
391: @*/
392: PetscErrorCode EPSSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
393: {
394:   PetscFunctionBegin;
399:   if (nev != PETSC_CURRENT) {
400:     PetscCheck(nev>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
401:     eps->nev = nev;
402:   }
403:   if (ncv == PETSC_DETERMINE) {
404:     eps->ncv = PETSC_DETERMINE;
405:   } else if (ncv != PETSC_CURRENT) {
406:     PetscCheck(ncv>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
407:     eps->ncv = ncv;
408:   }
409:   if (mpd == PETSC_DETERMINE) {
410:     eps->mpd = PETSC_DETERMINE;
411:   } else if (mpd != PETSC_CURRENT) {
412:     PetscCheck(mpd>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
413:     eps->mpd = mpd;
414:   }
415:   eps->state = EPS_STATE_INITIAL;
416:   PetscFunctionReturn(PETSC_SUCCESS);
417: }

419: /*@
420:    EPSSetWhichEigenpairs - Specifies which portion of the spectrum is
421:    to be sought.

423:    Logically Collective

425:    Input Parameters:
426: +  eps   - eigensolver context obtained from EPSCreate()
427: -  which - the portion of the spectrum to be sought

429:    Options Database Keys:
430: +   -eps_largest_magnitude - Sets largest eigenvalues in magnitude
431: .   -eps_smallest_magnitude - Sets smallest eigenvalues in magnitude
432: .   -eps_largest_real - Sets largest real parts
433: .   -eps_smallest_real - Sets smallest real parts
434: .   -eps_largest_imaginary - Sets largest imaginary parts
435: .   -eps_smallest_imaginary - Sets smallest imaginary parts
436: .   -eps_target_magnitude - Sets eigenvalues closest to target
437: .   -eps_target_real - Sets real parts closest to target
438: .   -eps_target_imaginary - Sets imaginary parts closest to target
439: -   -eps_all - Sets all eigenvalues in an interval or region

441:    Notes:
442:    The parameter 'which' can have one of these values

444: +     EPS_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
445: .     EPS_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
446: .     EPS_LARGEST_REAL - largest real parts
447: .     EPS_SMALLEST_REAL - smallest real parts
448: .     EPS_LARGEST_IMAGINARY - largest imaginary parts
449: .     EPS_SMALLEST_IMAGINARY - smallest imaginary parts
450: .     EPS_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
451: .     EPS_TARGET_REAL - eigenvalues with real part closest to target
452: .     EPS_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
453: .     EPS_ALL - all eigenvalues contained in a given interval or region
454: -     EPS_WHICH_USER - user defined ordering set with EPSSetEigenvalueComparison()

456:    Not all eigensolvers implemented in EPS account for all the possible values
457:    stated above. Also, some values make sense only for certain types of
458:    problems. If SLEPc is compiled for real numbers EPS_LARGEST_IMAGINARY
459:    and EPS_SMALLEST_IMAGINARY use the absolute value of the imaginary part
460:    for eigenvalue selection.

462:    The target is a scalar value provided with EPSSetTarget().

464:    The criterion EPS_TARGET_IMAGINARY is available only in case PETSc and
465:    SLEPc have been built with complex scalars.

467:    EPS_ALL is intended for use in combination with an interval (see
468:    EPSSetInterval()), when all eigenvalues within the interval are requested,
469:    or in the context of the CISS solver for computing all eigenvalues in a region.
470:    In those cases, the number of eigenvalues is unknown, so the nev parameter
471:    has a different sense, see EPSSetDimensions().

473:    Level: intermediate

475: .seealso: EPSGetWhichEigenpairs(), EPSSetTarget(), EPSSetInterval(),
476:           EPSSetDimensions(), EPSSetEigenvalueComparison(), EPSWhich
477: @*/
478: PetscErrorCode EPSSetWhichEigenpairs(EPS eps,EPSWhich which)
479: {
480:   PetscFunctionBegin;
483:   switch (which) {
484:     case EPS_LARGEST_MAGNITUDE:
485:     case EPS_SMALLEST_MAGNITUDE:
486:     case EPS_LARGEST_REAL:
487:     case EPS_SMALLEST_REAL:
488:     case EPS_LARGEST_IMAGINARY:
489:     case EPS_SMALLEST_IMAGINARY:
490:     case EPS_TARGET_MAGNITUDE:
491:     case EPS_TARGET_REAL:
492: #if defined(PETSC_USE_COMPLEX)
493:     case EPS_TARGET_IMAGINARY:
494: #endif
495:     case EPS_ALL:
496:     case EPS_WHICH_USER:
497:       if (eps->which != which) {
498:         eps->state = EPS_STATE_INITIAL;
499:         eps->which = which;
500:       }
501:       break;
502: #if !defined(PETSC_USE_COMPLEX)
503:     case EPS_TARGET_IMAGINARY:
504:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"EPS_TARGET_IMAGINARY can be used only with complex scalars");
505: #endif
506:     default:
507:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'which' value");
508:   }
509:   PetscFunctionReturn(PETSC_SUCCESS);
510: }

512: /*@
513:    EPSGetWhichEigenpairs - Returns which portion of the spectrum is to be
514:    sought.

516:    Not Collective

518:    Input Parameter:
519: .  eps - eigensolver context obtained from EPSCreate()

521:    Output Parameter:
522: .  which - the portion of the spectrum to be sought

524:    Notes:
525:    See EPSSetWhichEigenpairs() for possible values of 'which'.

527:    Level: intermediate

529: .seealso: EPSSetWhichEigenpairs(), EPSWhich
530: @*/
531: PetscErrorCode EPSGetWhichEigenpairs(EPS eps,EPSWhich *which)
532: {
533:   PetscFunctionBegin;
535:   PetscAssertPointer(which,2);
536:   *which = eps->which;
537:   PetscFunctionReturn(PETSC_SUCCESS);
538: }

540: /*@C
541:    EPSSetEigenvalueComparison - Specifies the eigenvalue comparison function
542:    when EPSSetWhichEigenpairs() is set to EPS_WHICH_USER.

544:    Logically Collective

546:    Input Parameters:
547: +  eps  - eigensolver context obtained from EPSCreate()
548: .  func - the comparison function, see EPSEigenvalueComparisonFn for the calling sequence
549: -  ctx  - a context pointer (the last parameter to the comparison function)

551:    Note:
552:    The returning parameter 'res' can be
553: +  negative - if the 1st eigenvalue is preferred to the 2st one
554: .  zero     - if both eigenvalues are equally preferred
555: -  positive - if the 2st eigenvalue is preferred to the 1st one

557:    Level: advanced

559: .seealso: EPSSetWhichEigenpairs(), EPSWhich
560: @*/
561: PetscErrorCode EPSSetEigenvalueComparison(EPS eps,SlepcEigenvalueComparisonFn *func,void* ctx)
562: {
563:   PetscFunctionBegin;
565:   eps->sc->comparison    = func;
566:   eps->sc->comparisonctx = ctx;
567:   eps->which             = EPS_WHICH_USER;
568:   PetscFunctionReturn(PETSC_SUCCESS);
569: }

571: /*@C
572:    EPSSetArbitrarySelection - Specifies a function intended to look for
573:    eigenvalues according to an arbitrary selection criterion. This criterion
574:    can be based on a computation involving the current eigenvector approximation.

576:    Logically Collective

578:    Input Parameters:
579: +  eps  - eigensolver context obtained from EPSCreate()
580: .  func - the arbitrary selection function, see SlepcArbitrarySelectionFn for a calling sequence
581: -  ctx  - a context pointer (the last parameter to the arbitrary selection function)

583:    Notes:
584:    This provides a mechanism to select eigenpairs by evaluating a user-defined
585:    function. When a function has been provided, the default selection based on
586:    sorting the eigenvalues is replaced by the sorting of the results of this
587:    function (with the same sorting criterion given in EPSSetWhichEigenpairs()).

589:    For instance, suppose you want to compute those eigenvectors that maximize
590:    a certain computable expression. Then implement the computation using
591:    the arguments xr and xi, and return the result in rr. Then set the standard
592:    sorting by magnitude so that the eigenpair with largest value of rr is
593:    selected.

595:    This evaluation function is collective, that is, all processes call it and
596:    it can use collective operations; furthermore, the computed result must
597:    be the same in all processes.

599:    The result of func is expressed as a complex number so that it is possible to
600:    use the standard eigenvalue sorting functions, but normally only rr is used.
601:    Set ri to zero unless it is meaningful in your application.

603:    Level: advanced

605: .seealso: EPSSetWhichEigenpairs()
606: @*/
607: PetscErrorCode EPSSetArbitrarySelection(EPS eps,SlepcArbitrarySelectionFn *func,void* ctx)
608: {
609:   PetscFunctionBegin;
611:   eps->arbitrary    = func;
612:   eps->arbitraryctx = ctx;
613:   eps->state        = EPS_STATE_INITIAL;
614:   PetscFunctionReturn(PETSC_SUCCESS);
615: }

617: /*@C
618:    EPSSetConvergenceTestFunction - Sets a function to compute the error estimate
619:    used in the convergence test.

621:    Logically Collective

623:    Input Parameters:
624: +  eps     - eigensolver context obtained from EPSCreate()
625: .  func    - convergence test function, see EPSConvergenceTestFn for the calling sequence
626: .  ctx     - context for private data for the convergence routine (may be NULL)
627: -  destroy - a routine for destroying the context (may be NULL)

629:    Note:
630:    If the error estimate returned by the convergence test function is less than
631:    the tolerance, then the eigenvalue is accepted as converged.

633:    Level: advanced

635: .seealso: EPSSetConvergenceTest(), EPSSetTolerances()
636: @*/
637: PetscErrorCode EPSSetConvergenceTestFunction(EPS eps,EPSConvergenceTestFn *func,void* ctx,PetscErrorCode (*destroy)(void*))
638: {
639:   PetscFunctionBegin;
641:   if (eps->convergeddestroy) PetscCall((*eps->convergeddestroy)(eps->convergedctx));
642:   eps->convergeduser    = func;
643:   eps->convergeddestroy = destroy;
644:   eps->convergedctx     = ctx;
645:   if (func == EPSConvergedRelative) eps->conv = EPS_CONV_REL;
646:   else if (func == EPSConvergedNorm) eps->conv = EPS_CONV_NORM;
647:   else if (func == EPSConvergedAbsolute) eps->conv = EPS_CONV_ABS;
648:   else {
649:     eps->conv      = EPS_CONV_USER;
650:     eps->converged = eps->convergeduser;
651:   }
652:   PetscFunctionReturn(PETSC_SUCCESS);
653: }

655: /*@
656:    EPSSetConvergenceTest - Specifies how to compute the error estimate
657:    used in the convergence test.

659:    Logically Collective

661:    Input Parameters:
662: +  eps  - eigensolver context obtained from EPSCreate()
663: -  conv - the type of convergence test

665:    Options Database Keys:
666: +  -eps_conv_abs  - Sets the absolute convergence test
667: .  -eps_conv_rel  - Sets the convergence test relative to the eigenvalue
668: .  -eps_conv_norm - Sets the convergence test relative to the matrix norms
669: -  -eps_conv_user - Selects the user-defined convergence test

671:    Note:
672:    The parameter 'conv' can have one of these values
673: +     EPS_CONV_ABS  - absolute error ||r||
674: .     EPS_CONV_REL  - error relative to the eigenvalue l, ||r||/|l|
675: .     EPS_CONV_NORM - error relative to the matrix norms, ||r||/(||A||+|l|*||B||)
676: -     EPS_CONV_USER - function set by EPSSetConvergenceTestFunction()

678:    Level: intermediate

680: .seealso: EPSGetConvergenceTest(), EPSSetConvergenceTestFunction(), EPSSetStoppingTest(), EPSConv
681: @*/
682: PetscErrorCode EPSSetConvergenceTest(EPS eps,EPSConv conv)
683: {
684:   PetscFunctionBegin;
687:   switch (conv) {
688:     case EPS_CONV_ABS:  eps->converged = EPSConvergedAbsolute; break;
689:     case EPS_CONV_REL:  eps->converged = EPSConvergedRelative; break;
690:     case EPS_CONV_NORM: eps->converged = EPSConvergedNorm; break;
691:     case EPS_CONV_USER:
692:       PetscCheck(eps->convergeduser,PetscObjectComm((PetscObject)eps),PETSC_ERR_ORDER,"Must call EPSSetConvergenceTestFunction() first");
693:       eps->converged = eps->convergeduser;
694:       break;
695:     default:
696:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
697:   }
698:   eps->conv = conv;
699:   PetscFunctionReturn(PETSC_SUCCESS);
700: }

702: /*@
703:    EPSGetConvergenceTest - Gets the method used to compute the error estimate
704:    used in the convergence test.

706:    Not Collective

708:    Input Parameters:
709: .  eps   - eigensolver context obtained from EPSCreate()

711:    Output Parameters:
712: .  conv  - the type of convergence test

714:    Level: intermediate

716: .seealso: EPSSetConvergenceTest(), EPSConv
717: @*/
718: PetscErrorCode EPSGetConvergenceTest(EPS eps,EPSConv *conv)
719: {
720:   PetscFunctionBegin;
722:   PetscAssertPointer(conv,2);
723:   *conv = eps->conv;
724:   PetscFunctionReturn(PETSC_SUCCESS);
725: }

727: /*@C
728:    EPSSetStoppingTestFunction - Sets a function to decide when to stop the outer
729:    iteration of the eigensolver.

731:    Logically Collective

733:    Input Parameters:
734: +  eps     - eigensolver context obtained from EPSCreate()
735: .  func    - stopping test function, see EPSStoppingTestFn for the calling sequence
736: .  ctx     - context for private data for the stopping routine (may be NULL)
737: -  destroy - a routine for destroying the context (may be NULL)

739:    Note:
740:    Normal usage is to first call the default routine EPSStoppingBasic() and then
741:    set reason to EPS_CONVERGED_USER if some user-defined conditions have been
742:    met. To let the eigensolver continue iterating, the result must be left as
743:    EPS_CONVERGED_ITERATING.

745:    Level: advanced

747: .seealso: EPSSetStoppingTest(), EPSStoppingBasic()
748: @*/
749: PetscErrorCode EPSSetStoppingTestFunction(EPS eps,EPSStoppingTestFn *func,void* ctx,PetscErrorCode (*destroy)(void*))
750: {
751:   PetscFunctionBegin;
753:   if (eps->stoppingdestroy) PetscCall((*eps->stoppingdestroy)(eps->stoppingctx));
754:   eps->stoppinguser    = func;
755:   eps->stoppingdestroy = destroy;
756:   eps->stoppingctx     = ctx;
757:   if (func == EPSStoppingBasic) eps->stop = EPS_STOP_BASIC;
758:   else {
759:     eps->stop     = EPS_STOP_USER;
760:     eps->stopping = eps->stoppinguser;
761:   }
762:   PetscFunctionReturn(PETSC_SUCCESS);
763: }

765: /*@
766:    EPSSetStoppingTest - Specifies how to decide the termination of the outer
767:    loop of the eigensolver.

769:    Logically Collective

771:    Input Parameters:
772: +  eps  - eigensolver context obtained from EPSCreate()
773: -  stop - the type of stopping test

775:    Options Database Keys:
776: +  -eps_stop_basic - Sets the default stopping test
777: -  -eps_stop_user  - Selects the user-defined stopping test

779:    Note:
780:    The parameter 'stop' can have one of these values
781: +     EPS_STOP_BASIC - default stopping test
782: -     EPS_STOP_USER  - function set by EPSSetStoppingTestFunction()

784:    Level: advanced

786: .seealso: EPSGetStoppingTest(), EPSSetStoppingTestFunction(), EPSSetConvergenceTest(), EPSStop
787: @*/
788: PetscErrorCode EPSSetStoppingTest(EPS eps,EPSStop stop)
789: {
790:   PetscFunctionBegin;
793:   switch (stop) {
794:     case EPS_STOP_BASIC: eps->stopping = EPSStoppingBasic; break;
795:     case EPS_STOP_USER:
796:       PetscCheck(eps->stoppinguser,PetscObjectComm((PetscObject)eps),PETSC_ERR_ORDER,"Must call EPSSetStoppingTestFunction() first");
797:       eps->stopping = eps->stoppinguser;
798:       break;
799:     default:
800:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
801:   }
802:   eps->stop = stop;
803:   PetscFunctionReturn(PETSC_SUCCESS);
804: }

806: /*@
807:    EPSGetStoppingTest - Gets the method used to decide the termination of the outer
808:    loop of the eigensolver.

810:    Not Collective

812:    Input Parameters:
813: .  eps   - eigensolver context obtained from EPSCreate()

815:    Output Parameters:
816: .  stop  - the type of stopping test

818:    Level: advanced

820: .seealso: EPSSetStoppingTest(), EPSStop
821: @*/
822: PetscErrorCode EPSGetStoppingTest(EPS eps,EPSStop *stop)
823: {
824:   PetscFunctionBegin;
826:   PetscAssertPointer(stop,2);
827:   *stop = eps->stop;
828:   PetscFunctionReturn(PETSC_SUCCESS);
829: }

831: /*@
832:    EPSSetProblemType - Specifies the type of the eigenvalue problem.

834:    Logically Collective

836:    Input Parameters:
837: +  eps      - the eigensolver context
838: -  type     - a known type of eigenvalue problem

840:    Options Database Keys:
841: +  -eps_hermitian - Hermitian eigenvalue problem
842: .  -eps_gen_hermitian - generalized Hermitian eigenvalue problem
843: .  -eps_non_hermitian - non-Hermitian eigenvalue problem
844: .  -eps_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
845: .  -eps_pos_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
846:    with positive semi-definite B
847: .  -eps_gen_indefinite - generalized Hermitian-indefinite eigenvalue problem
848: -  -eps_bse - structured Bethe-Salpeter eigenvalue problem

850:    Notes:
851:    This function must be used to instruct SLEPc to exploit symmetry or other
852:    kind of structure. If no
853:    problem type is specified, by default a non-Hermitian problem is assumed
854:    (either standard or generalized). If the user knows that the problem is
855:    Hermitian (i.e. A=A^H) or generalized Hermitian (i.e. A=A^H, B=B^H, and
856:    B positive definite) then it is recommended to set the problem type so
857:    that eigensolver can exploit these properties.

859:    If the user does not call this function, the solver will use a reasonable
860:    guess.

862:    For structured problem types such as EPS_BSE, the matrices passed in via
863:    EPSSetOperators() must have been created with the corresponding helper
864:    function, i.e., MatCreateBSE().

866:    Level: intermediate

868: .seealso: EPSSetOperators(), EPSSetType(), EPSGetProblemType(), EPSProblemType
869: @*/
870: PetscErrorCode EPSSetProblemType(EPS eps,EPSProblemType type)
871: {
872:   PetscFunctionBegin;
875:   if (type == eps->problem_type) PetscFunctionReturn(PETSC_SUCCESS);
876:   switch (type) {
877:     case EPS_HEP:
878:       eps->isgeneralized = PETSC_FALSE;
879:       eps->ishermitian = PETSC_TRUE;
880:       eps->ispositive = PETSC_FALSE;
881:       eps->isstructured = PETSC_FALSE;
882:       break;
883:     case EPS_NHEP:
884:       eps->isgeneralized = PETSC_FALSE;
885:       eps->ishermitian = PETSC_FALSE;
886:       eps->ispositive = PETSC_FALSE;
887:       eps->isstructured = PETSC_FALSE;
888:       break;
889:     case EPS_GHEP:
890:       eps->isgeneralized = PETSC_TRUE;
891:       eps->ishermitian = PETSC_TRUE;
892:       eps->ispositive = PETSC_TRUE;
893:       eps->isstructured = PETSC_FALSE;
894:       break;
895:     case EPS_GNHEP:
896:       eps->isgeneralized = PETSC_TRUE;
897:       eps->ishermitian = PETSC_FALSE;
898:       eps->ispositive = PETSC_FALSE;
899:       eps->isstructured = PETSC_FALSE;
900:       break;
901:     case EPS_PGNHEP:
902:       eps->isgeneralized = PETSC_TRUE;
903:       eps->ishermitian = PETSC_FALSE;
904:       eps->ispositive = PETSC_TRUE;
905:       eps->isstructured = PETSC_FALSE;
906:       break;
907:     case EPS_GHIEP:
908:       eps->isgeneralized = PETSC_TRUE;
909:       eps->ishermitian = PETSC_TRUE;
910:       eps->ispositive = PETSC_FALSE;
911:       eps->isstructured = PETSC_FALSE;
912:       break;
913:     case EPS_BSE:
914:       eps->isgeneralized = PETSC_FALSE;
915:       eps->ishermitian = PETSC_FALSE;
916:       eps->ispositive = PETSC_FALSE;
917:       eps->isstructured = PETSC_TRUE;
918:       break;
919:     default:
920:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
921:   }
922:   eps->problem_type = type;
923:   eps->state = EPS_STATE_INITIAL;
924:   PetscFunctionReturn(PETSC_SUCCESS);
925: }

927: /*@
928:    EPSGetProblemType - Gets the problem type from the EPS object.

930:    Not Collective

932:    Input Parameter:
933: .  eps - the eigensolver context

935:    Output Parameter:
936: .  type - the problem type

938:    Level: intermediate

940: .seealso: EPSSetProblemType(), EPSProblemType
941: @*/
942: PetscErrorCode EPSGetProblemType(EPS eps,EPSProblemType *type)
943: {
944:   PetscFunctionBegin;
946:   PetscAssertPointer(type,2);
947:   *type = eps->problem_type;
948:   PetscFunctionReturn(PETSC_SUCCESS);
949: }

951: /*@
952:    EPSSetExtraction - Specifies the type of extraction technique to be employed
953:    by the eigensolver.

955:    Logically Collective

957:    Input Parameters:
958: +  eps  - the eigensolver context
959: -  extr - a known type of extraction

961:    Options Database Keys:
962: +  -eps_ritz - Rayleigh-Ritz extraction
963: .  -eps_harmonic - harmonic Ritz extraction
964: .  -eps_harmonic_relative - harmonic Ritz extraction relative to the eigenvalue
965: .  -eps_harmonic_right - harmonic Ritz extraction for rightmost eigenvalues
966: .  -eps_harmonic_largest - harmonic Ritz extraction for largest magnitude
967:    (without target)
968: .  -eps_refined - refined Ritz extraction
969: -  -eps_refined_harmonic - refined harmonic Ritz extraction

971:    Notes:
972:    Not all eigensolvers support all types of extraction. See the SLEPc
973:    Users Manual for details.

975:    By default, a standard Rayleigh-Ritz extraction is used. Other extractions
976:    may be useful when computing interior eigenvalues.

978:    Harmonic-type extractions are used in combination with a 'target'.

980:    Level: advanced

982: .seealso: EPSSetTarget(), EPSGetExtraction(), EPSExtraction
983: @*/
984: PetscErrorCode EPSSetExtraction(EPS eps,EPSExtraction extr)
985: {
986:   PetscFunctionBegin;
989:   if (eps->extraction != extr) {
990:     eps->state      = EPS_STATE_INITIAL;
991:     eps->extraction = extr;
992:   }
993:   PetscFunctionReturn(PETSC_SUCCESS);
994: }

996: /*@
997:    EPSGetExtraction - Gets the extraction type used by the EPS object.

999:    Not Collective

1001:    Input Parameter:
1002: .  eps - the eigensolver context

1004:    Output Parameter:
1005: .  extr - name of extraction type

1007:    Level: advanced

1009: .seealso: EPSSetExtraction(), EPSExtraction
1010: @*/
1011: PetscErrorCode EPSGetExtraction(EPS eps,EPSExtraction *extr)
1012: {
1013:   PetscFunctionBegin;
1015:   PetscAssertPointer(extr,2);
1016:   *extr = eps->extraction;
1017:   PetscFunctionReturn(PETSC_SUCCESS);
1018: }

1020: /*@
1021:    EPSSetBalance - Specifies the balancing technique to be employed by the
1022:    eigensolver, and some parameters associated to it.

1024:    Logically Collective

1026:    Input Parameters:
1027: +  eps    - the eigensolver context
1028: .  bal    - the balancing method, one of EPS_BALANCE_NONE, EPS_BALANCE_ONESIDE,
1029:             EPS_BALANCE_TWOSIDE, or EPS_BALANCE_USER
1030: .  its    - number of iterations of the balancing algorithm
1031: -  cutoff - cutoff value

1033:    Options Database Keys:
1034: +  -eps_balance <method> - the balancing method, where <method> is one of
1035:                            'none', 'oneside', 'twoside', or 'user'
1036: .  -eps_balance_its <its> - number of iterations
1037: -  -eps_balance_cutoff <cutoff> - cutoff value

1039:    Notes:
1040:    When balancing is enabled, the solver works implicitly with matrix DAD^-1,
1041:    where D is an appropriate diagonal matrix. This improves the accuracy of
1042:    the computed results in some cases. See the SLEPc Users Manual for details.

1044:    Balancing makes sense only for non-Hermitian problems when the required
1045:    precision is high (i.e. a small tolerance such as 1e-15).

1047:    By default, balancing is disabled. The two-sided method is much more
1048:    effective than the one-sided counterpart, but it requires the system
1049:    matrices to have the MatMultTranspose operation defined.

1051:    The parameter 'its' is the number of iterations performed by the method. The
1052:    cutoff value is used only in the two-side variant. Use PETSC_DETERMINE to assign
1053:    a reasonably good value, or PETSC_CURRENT to leave the value unchanged.

1055:    User-defined balancing is allowed provided that the corresponding matrix
1056:    is set via STSetBalanceMatrix.

1058:    Level: intermediate

1060: .seealso: EPSGetBalance(), EPSBalance, STSetBalanceMatrix()
1061: @*/
1062: PetscErrorCode EPSSetBalance(EPS eps,EPSBalance bal,PetscInt its,PetscReal cutoff)
1063: {
1064:   PetscFunctionBegin;
1069:   switch (bal) {
1070:     case EPS_BALANCE_NONE:
1071:     case EPS_BALANCE_ONESIDE:
1072:     case EPS_BALANCE_TWOSIDE:
1073:     case EPS_BALANCE_USER:
1074:       if (eps->balance != bal) {
1075:         eps->state = EPS_STATE_INITIAL;
1076:         eps->balance = bal;
1077:       }
1078:       break;
1079:     default:
1080:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'bal'");
1081:   }
1082:   if (its==PETSC_DETERMINE) eps->balance_its = 5;
1083:   else if (its!=PETSC_CURRENT) {
1084:     PetscCheck(its>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be > 0");
1085:     eps->balance_its = its;
1086:   }
1087:   if (cutoff==(PetscReal)PETSC_DETERMINE) eps->balance_cutoff = 1e-8;
1088:   else if (cutoff!=(PetscReal)PETSC_CURRENT) {
1089:     PetscCheck(cutoff>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of cutoff. Must be > 0");
1090:     eps->balance_cutoff = cutoff;
1091:   }
1092:   PetscFunctionReturn(PETSC_SUCCESS);
1093: }

1095: /*@
1096:    EPSGetBalance - Gets the balancing type used by the EPS object, and the
1097:    associated parameters.

1099:    Not Collective

1101:    Input Parameter:
1102: .  eps - the eigensolver context

1104:    Output Parameters:
1105: +  bal    - the balancing method
1106: .  its    - number of iterations of the balancing algorithm
1107: -  cutoff - cutoff value

1109:    Level: intermediate

1111:    Note:
1112:    The user can specify NULL for any parameter that is not needed.

1114: .seealso: EPSSetBalance(), EPSBalance
1115: @*/
1116: PetscErrorCode EPSGetBalance(EPS eps,EPSBalance *bal,PetscInt *its,PetscReal *cutoff)
1117: {
1118:   PetscFunctionBegin;
1120:   if (bal)    *bal = eps->balance;
1121:   if (its)    *its = eps->balance_its;
1122:   if (cutoff) *cutoff = eps->balance_cutoff;
1123:   PetscFunctionReturn(PETSC_SUCCESS);
1124: }

1126: /*@
1127:    EPSSetTwoSided - Sets the solver to use a two-sided variant so that left
1128:    eigenvectors are also computed.

1130:    Logically Collective

1132:    Input Parameters:
1133: +  eps      - the eigensolver context
1134: -  twosided - whether the two-sided variant is to be used or not

1136:    Options Database Keys:
1137: .  -eps_two_sided <boolean> - Sets/resets the twosided flag

1139:    Notes:
1140:    If the user sets twosided=PETSC_TRUE then the solver uses a variant of
1141:    the algorithm that computes both right and left eigenvectors. This is
1142:    usually much more costly. This option is not available in all solvers.

1144:    When using two-sided solvers, the problem matrices must have both the
1145:    MatMult and MatMultTranspose operations defined.

1147:    Level: advanced

1149: .seealso: EPSGetTwoSided(), EPSGetLeftEigenvector()
1150: @*/
1151: PetscErrorCode EPSSetTwoSided(EPS eps,PetscBool twosided)
1152: {
1153:   PetscFunctionBegin;
1156:   if (twosided!=eps->twosided) {
1157:     eps->twosided = twosided;
1158:     eps->state    = EPS_STATE_INITIAL;
1159:   }
1160:   PetscFunctionReturn(PETSC_SUCCESS);
1161: }

1163: /*@
1164:    EPSGetTwoSided - Returns the flag indicating whether a two-sided variant
1165:    of the algorithm is being used or not.

1167:    Not Collective

1169:    Input Parameter:
1170: .  eps - the eigensolver context

1172:    Output Parameter:
1173: .  twosided - the returned flag

1175:    Level: advanced

1177: .seealso: EPSSetTwoSided()
1178: @*/
1179: PetscErrorCode EPSGetTwoSided(EPS eps,PetscBool *twosided)
1180: {
1181:   PetscFunctionBegin;
1183:   PetscAssertPointer(twosided,2);
1184:   *twosided = eps->twosided;
1185:   PetscFunctionReturn(PETSC_SUCCESS);
1186: }

1188: /*@
1189:    EPSSetTrueResidual - Specifies if the solver must compute the true residual
1190:    explicitly or not.

1192:    Logically Collective

1194:    Input Parameters:
1195: +  eps     - the eigensolver context
1196: -  trueres - whether true residuals are required or not

1198:    Options Database Keys:
1199: .  -eps_true_residual <boolean> - Sets/resets the boolean flag 'trueres'

1201:    Notes:
1202:    If the user sets trueres=PETSC_TRUE then the solver explicitly computes
1203:    the true residual for each eigenpair approximation, and uses it for
1204:    convergence testing. Computing the residual is usually an expensive
1205:    operation. Some solvers (e.g., Krylov solvers) can avoid this computation
1206:    by using a cheap estimate of the residual norm, but this may sometimes
1207:    give inaccurate results (especially if a spectral transform is being
1208:    used). On the contrary, preconditioned eigensolvers (e.g., Davidson solvers)
1209:    do rely on computing the true residual, so this option is irrelevant for them.

1211:    Level: advanced

1213: .seealso: EPSGetTrueResidual()
1214: @*/
1215: PetscErrorCode EPSSetTrueResidual(EPS eps,PetscBool trueres)
1216: {
1217:   PetscFunctionBegin;
1220:   eps->trueres = trueres;
1221:   PetscFunctionReturn(PETSC_SUCCESS);
1222: }

1224: /*@
1225:    EPSGetTrueResidual - Returns the flag indicating whether true
1226:    residuals must be computed explicitly or not.

1228:    Not Collective

1230:    Input Parameter:
1231: .  eps - the eigensolver context

1233:    Output Parameter:
1234: .  trueres - the returned flag

1236:    Level: advanced

1238: .seealso: EPSSetTrueResidual()
1239: @*/
1240: PetscErrorCode EPSGetTrueResidual(EPS eps,PetscBool *trueres)
1241: {
1242:   PetscFunctionBegin;
1244:   PetscAssertPointer(trueres,2);
1245:   *trueres = eps->trueres;
1246:   PetscFunctionReturn(PETSC_SUCCESS);
1247: }

1249: /*@
1250:    EPSSetTrackAll - Specifies if the solver must compute the residual norm of all
1251:    approximate eigenpairs or not.

1253:    Logically Collective

1255:    Input Parameters:
1256: +  eps      - the eigensolver context
1257: -  trackall - whether to compute all residuals or not

1259:    Notes:
1260:    If the user sets trackall=PETSC_TRUE then the solver computes (or estimates)
1261:    the residual norm for each eigenpair approximation. Computing the residual is
1262:    usually an expensive operation and solvers commonly compute only the residual
1263:    associated to the first unconverged eigenpair.

1265:    The option '-eps_monitor_all' automatically activates this option.

1267:    Level: developer

1269: .seealso: EPSGetTrackAll()
1270: @*/
1271: PetscErrorCode EPSSetTrackAll(EPS eps,PetscBool trackall)
1272: {
1273:   PetscFunctionBegin;
1276:   eps->trackall = trackall;
1277:   PetscFunctionReturn(PETSC_SUCCESS);
1278: }

1280: /*@
1281:    EPSGetTrackAll - Returns the flag indicating whether all residual norms must
1282:    be computed or not.

1284:    Not Collective

1286:    Input Parameter:
1287: .  eps - the eigensolver context

1289:    Output Parameter:
1290: .  trackall - the returned flag

1292:    Level: developer

1294: .seealso: EPSSetTrackAll()
1295: @*/
1296: PetscErrorCode EPSGetTrackAll(EPS eps,PetscBool *trackall)
1297: {
1298:   PetscFunctionBegin;
1300:   PetscAssertPointer(trackall,2);
1301:   *trackall = eps->trackall;
1302:   PetscFunctionReturn(PETSC_SUCCESS);
1303: }

1305: /*@
1306:    EPSSetPurify - Deactivate eigenvector purification (which is activated by default).

1308:    Logically Collective

1310:    Input Parameters:
1311: +  eps    - the eigensolver context
1312: -  purify - whether purification is required or not

1314:    Options Database Keys:
1315: .  -eps_purify <boolean> - Sets/resets the boolean flag 'purify'

1317:    Notes:
1318:    By default, eigenvectors of generalized symmetric eigenproblems are purified
1319:    in order to purge directions in the nullspace of matrix B. If the user knows
1320:    that B is non-singular, then purification can be safely deactivated and some
1321:    computational cost is avoided (this is particularly important in interval computations).

1323:    Level: intermediate

1325: .seealso: EPSGetPurify(), EPSSetInterval()
1326: @*/
1327: PetscErrorCode EPSSetPurify(EPS eps,PetscBool purify)
1328: {
1329:   PetscFunctionBegin;
1332:   if (purify!=eps->purify) {
1333:     eps->purify = purify;
1334:     eps->state  = EPS_STATE_INITIAL;
1335:   }
1336:   PetscFunctionReturn(PETSC_SUCCESS);
1337: }

1339: /*@
1340:    EPSGetPurify - Returns the flag indicating whether purification is activated
1341:    or not.

1343:    Not Collective

1345:    Input Parameter:
1346: .  eps - the eigensolver context

1348:    Output Parameter:
1349: .  purify - the returned flag

1351:    Level: intermediate

1353: .seealso: EPSSetPurify()
1354: @*/
1355: PetscErrorCode EPSGetPurify(EPS eps,PetscBool *purify)
1356: {
1357:   PetscFunctionBegin;
1359:   PetscAssertPointer(purify,2);
1360:   *purify = eps->purify;
1361:   PetscFunctionReturn(PETSC_SUCCESS);
1362: }

1364: /*@
1365:    EPSSetOptionsPrefix - Sets the prefix used for searching for all
1366:    EPS options in the database.

1368:    Logically Collective

1370:    Input Parameters:
1371: +  eps - the eigensolver context
1372: -  prefix - the prefix string to prepend to all EPS option requests

1374:    Notes:
1375:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1376:    The first character of all runtime options is AUTOMATICALLY the
1377:    hyphen.

1379:    For example, to distinguish between the runtime options for two
1380:    different EPS contexts, one could call
1381: .vb
1382:       EPSSetOptionsPrefix(eps1,"eig1_")
1383:       EPSSetOptionsPrefix(eps2,"eig2_")
1384: .ve

1386:    Level: advanced

1388: .seealso: EPSAppendOptionsPrefix(), EPSGetOptionsPrefix()
1389: @*/
1390: PetscErrorCode EPSSetOptionsPrefix(EPS eps,const char *prefix)
1391: {
1392:   PetscFunctionBegin;
1394:   if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
1395:   PetscCall(STSetOptionsPrefix(eps->st,prefix));
1396:   if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
1397:   PetscCall(BVSetOptionsPrefix(eps->V,prefix));
1398:   if (!eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
1399:   PetscCall(DSSetOptionsPrefix(eps->ds,prefix));
1400:   if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
1401:   PetscCall(RGSetOptionsPrefix(eps->rg,prefix));
1402:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)eps,prefix));
1403:   PetscFunctionReturn(PETSC_SUCCESS);
1404: }

1406: /*@
1407:    EPSAppendOptionsPrefix - Appends to the prefix used for searching for all
1408:    EPS options in the database.

1410:    Logically Collective

1412:    Input Parameters:
1413: +  eps - the eigensolver context
1414: -  prefix - the prefix string to prepend to all EPS option requests

1416:    Notes:
1417:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1418:    The first character of all runtime options is AUTOMATICALLY the hyphen.

1420:    Level: advanced

1422: .seealso: EPSSetOptionsPrefix(), EPSGetOptionsPrefix()
1423: @*/
1424: PetscErrorCode EPSAppendOptionsPrefix(EPS eps,const char *prefix)
1425: {
1426:   PetscFunctionBegin;
1428:   if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
1429:   PetscCall(STAppendOptionsPrefix(eps->st,prefix));
1430:   if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
1431:   PetscCall(BVAppendOptionsPrefix(eps->V,prefix));
1432:   if (!eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
1433:   PetscCall(DSAppendOptionsPrefix(eps->ds,prefix));
1434:   if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
1435:   PetscCall(RGAppendOptionsPrefix(eps->rg,prefix));
1436:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)eps,prefix));
1437:   PetscFunctionReturn(PETSC_SUCCESS);
1438: }

1440: /*@
1441:    EPSGetOptionsPrefix - Gets the prefix used for searching for all
1442:    EPS options in the database.

1444:    Not Collective

1446:    Input Parameters:
1447: .  eps - the eigensolver context

1449:    Output Parameters:
1450: .  prefix - pointer to the prefix string used is returned

1452:    Note:
1453:    On the Fortran side, the user should pass in a string 'prefix' of
1454:    sufficient length to hold the prefix.

1456:    Level: advanced

1458: .seealso: EPSSetOptionsPrefix(), EPSAppendOptionsPrefix()
1459: @*/
1460: PetscErrorCode EPSGetOptionsPrefix(EPS eps,const char *prefix[])
1461: {
1462:   PetscFunctionBegin;
1464:   PetscAssertPointer(prefix,2);
1465:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)eps,prefix));
1466:   PetscFunctionReturn(PETSC_SUCCESS);
1467: }