Actual source code: epsopts.c

slepc-main 2024-12-17
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,(PetscCtxDestroyFn*)dfunc));
 63:   if (trackall) PetscCall(EPSSetTrackAll(eps,PETSC_TRUE));
 64:   PetscFunctionReturn(PETSC_SUCCESS);
 65: }

 67: /*@
 68:    EPSSetFromOptions - Sets EPS options from the options database.
 69:    This routine must be called before EPSSetUp() if the user is to be
 70:    allowed to 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:     r = eps->thres;
146:     PetscCall(PetscOptionsReal("-eps_threshold_absolute","Absolute threshold","EPSSetThreshold",r,&r,&flg));
147:     if (flg) PetscCall(EPSSetThreshold(eps,r,PETSC_FALSE));
148:     PetscCall(PetscOptionsReal("-eps_threshold_relative","Relative threshold","EPSSetThreshold",r,&r,&flg));
149:     if (flg) PetscCall(EPSSetThreshold(eps,r,PETSC_TRUE));

151:     PetscCall(PetscOptionsBoolGroupBegin("-eps_conv_rel","Relative error convergence test","EPSSetConvergenceTest",&flg));
152:     if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_REL));
153:     PetscCall(PetscOptionsBoolGroup("-eps_conv_norm","Convergence test relative to the eigenvalue and the matrix norms","EPSSetConvergenceTest",&flg));
154:     if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_NORM));
155:     PetscCall(PetscOptionsBoolGroup("-eps_conv_abs","Absolute error convergence test","EPSSetConvergenceTest",&flg));
156:     if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_ABS));
157:     PetscCall(PetscOptionsBoolGroupEnd("-eps_conv_user","User-defined convergence test","EPSSetConvergenceTest",&flg));
158:     if (flg) PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_USER));

160:     PetscCall(PetscOptionsBoolGroupBegin("-eps_stop_basic","Stop iteration if all eigenvalues converged or max_it reached","EPSSetStoppingTest",&flg));
161:     if (flg) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_BASIC));
162:     PetscCall(PetscOptionsBoolGroup("-eps_stop_threshold","Stop iteration if a converged eigenvalue is below/above the threshold","EPSSetStoppingTest",&flg));
163:     if (flg) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_THRESHOLD));
164:     PetscCall(PetscOptionsBoolGroupEnd("-eps_stop_user","User-defined stopping test","EPSSetStoppingTest",&flg));
165:     if (flg) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_USER));

167:     i = eps->nev;
168:     PetscCall(PetscOptionsInt("-eps_nev","Number of eigenvalues to compute","EPSSetDimensions",eps->nev,&i,&flg1));
169:     if (!flg1) i = PETSC_CURRENT;
170:     j = eps->ncv;
171:     PetscCall(PetscOptionsInt("-eps_ncv","Number of basis vectors","EPSSetDimensions",eps->ncv,&j,&flg2));
172:     k = eps->mpd;
173:     PetscCall(PetscOptionsInt("-eps_mpd","Maximum dimension of projected problem","EPSSetDimensions",eps->mpd,&k,&flg3));
174:     if (flg1 || flg2 || flg3) PetscCall(EPSSetDimensions(eps,i,j,k));

176:     PetscCall(PetscOptionsBoolGroupBegin("-eps_largest_magnitude","Compute largest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg));
177:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE));
178:     PetscCall(PetscOptionsBoolGroup("-eps_smallest_magnitude","Compute smallest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg));
179:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_MAGNITUDE));
180:     PetscCall(PetscOptionsBoolGroup("-eps_largest_real","Compute eigenvalues with largest real parts","EPSSetWhichEigenpairs",&flg));
181:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL));
182:     PetscCall(PetscOptionsBoolGroup("-eps_smallest_real","Compute eigenvalues with smallest real parts","EPSSetWhichEigenpairs",&flg));
183:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
184:     PetscCall(PetscOptionsBoolGroup("-eps_largest_imaginary","Compute eigenvalues with largest imaginary parts","EPSSetWhichEigenpairs",&flg));
185:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_IMAGINARY));
186:     PetscCall(PetscOptionsBoolGroup("-eps_smallest_imaginary","Compute eigenvalues with smallest imaginary parts","EPSSetWhichEigenpairs",&flg));
187:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_IMAGINARY));
188:     PetscCall(PetscOptionsBoolGroup("-eps_target_magnitude","Compute eigenvalues closest to target","EPSSetWhichEigenpairs",&flg));
189:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
190:     PetscCall(PetscOptionsBoolGroup("-eps_target_real","Compute eigenvalues with real parts closest to target","EPSSetWhichEigenpairs",&flg));
191:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_REAL));
192:     PetscCall(PetscOptionsBoolGroup("-eps_target_imaginary","Compute eigenvalues with imaginary parts closest to target","EPSSetWhichEigenpairs",&flg));
193:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_IMAGINARY));
194:     PetscCall(PetscOptionsBoolGroupEnd("-eps_all","Compute all eigenvalues in an interval or a region","EPSSetWhichEigenpairs",&flg));
195:     if (flg) PetscCall(EPSSetWhichEigenpairs(eps,EPS_ALL));

197:     PetscCall(PetscOptionsScalar("-eps_target","Value of the target","EPSSetTarget",eps->target,&s,&flg));
198:     if (flg) {
199:       if (eps->which!=EPS_TARGET_REAL && eps->which!=EPS_TARGET_IMAGINARY) PetscCall(EPSSetWhichEigenpairs(eps,EPS_TARGET_MAGNITUDE));
200:       PetscCall(EPSSetTarget(eps,s));
201:     }

203:     k = 2;
204:     PetscCall(PetscOptionsRealArray("-eps_interval","Computational interval (two real values separated with a comma without spaces)","EPSSetInterval",array,&k,&flg));
205:     if (flg) {
206:       PetscCheck(k>1,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_SIZ,"Must pass two values in -eps_interval (comma-separated without spaces)");
207:       PetscCall(EPSSetWhichEigenpairs(eps,EPS_ALL));
208:       PetscCall(EPSSetInterval(eps,array[0],array[1]));
209:     }

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

217:     /* -----------------------------------------------------------------------*/
218:     /*
219:       Cancels all monitors hardwired into code before call to EPSSetFromOptions()
220:     */
221:     PetscCall(PetscOptionsBool("-eps_monitor_cancel","Remove any hardwired monitor routines","EPSMonitorCancel",PETSC_FALSE,&flg,&set));
222:     if (set && flg) PetscCall(EPSMonitorCancel(eps));
223:     PetscCall(EPSMonitorSetFromOptions(eps,"-eps_monitor","first_approximation",NULL,PETSC_FALSE));
224:     PetscCall(EPSMonitorSetFromOptions(eps,"-eps_monitor_all","all_approximations",NULL,PETSC_TRUE));
225:     PetscCall(EPSMonitorSetFromOptions(eps,"-eps_monitor_conv","convergence_history",NULL,PETSC_FALSE));

227:     /* -----------------------------------------------------------------------*/
228:     PetscCall(PetscOptionsName("-eps_view","Print detailed information on solver used","EPSView",&set));
229:     PetscCall(PetscOptionsName("-eps_view_vectors","View computed eigenvectors","EPSVectorsView",&set));
230:     PetscCall(PetscOptionsName("-eps_view_values","View computed eigenvalues","EPSValuesView",&set));
231:     PetscCall(PetscOptionsName("-eps_converged_reason","Print reason for convergence, and number of iterations","EPSConvergedReasonView",&set));
232:     PetscCall(PetscOptionsName("-eps_error_absolute","Print absolute errors of each eigenpair","EPSErrorView",&set));
233:     PetscCall(PetscOptionsName("-eps_error_relative","Print relative errors of each eigenpair","EPSErrorView",&set));
234:     PetscCall(PetscOptionsName("-eps_error_backward","Print backward errors of each eigenpair","EPSErrorView",&set));

236:     PetscTryTypeMethod(eps,setfromoptions,PetscOptionsObject);
237:     PetscCall(PetscObjectProcessOptionsHandlers((PetscObject)eps,PetscOptionsObject));
238:   PetscOptionsEnd();

240:   if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
241:   PetscCall(BVSetFromOptions(eps->V));
242:   if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
243:   PetscCall(RGSetFromOptions(eps->rg));
244:   if (eps->useds) {
245:     if (!eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
246:     PetscCall(EPSSetDSType(eps));
247:     PetscCall(DSSetFromOptions(eps->ds));
248:   }
249:   if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
250:   PetscCall(EPSSetDefaultST(eps));
251:   PetscCall(STSetFromOptions(eps->st));
252:   PetscFunctionReturn(PETSC_SUCCESS);
253: }

255: /*@
256:    EPSGetTolerances - Gets the tolerance and maximum iteration count used
257:    by the EPS convergence tests.

259:    Not Collective

261:    Input Parameter:
262: .  eps - the eigensolver context

264:    Output Parameters:
265: +  tol - the convergence tolerance
266: -  maxits - maximum number of iterations

268:    Notes:
269:    The user can specify NULL for any parameter that is not needed.

271:    Level: intermediate

273: .seealso: EPSSetTolerances()
274: @*/
275: PetscErrorCode EPSGetTolerances(EPS eps,PetscReal *tol,PetscInt *maxits)
276: {
277:   PetscFunctionBegin;
279:   if (tol)    *tol    = eps->tol;
280:   if (maxits) *maxits = eps->max_it;
281:   PetscFunctionReturn(PETSC_SUCCESS);
282: }

284: /*@
285:    EPSSetTolerances - Sets the tolerance and maximum iteration count used
286:    by the EPS convergence tests.

288:    Logically Collective

290:    Input Parameters:
291: +  eps - the eigensolver context
292: .  tol - the convergence tolerance
293: -  maxits - maximum number of iterations to use

295:    Options Database Keys:
296: +  -eps_tol <tol> - Sets the convergence tolerance
297: -  -eps_max_it <maxits> - Sets the maximum number of iterations allowed

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

305:    Level: intermediate

307: .seealso: EPSGetTolerances()
308: @*/
309: PetscErrorCode EPSSetTolerances(EPS eps,PetscReal tol,PetscInt maxits)
310: {
311:   PetscFunctionBegin;
315:   if (tol == (PetscReal)PETSC_DETERMINE) {
316:     eps->tol   = PETSC_DETERMINE;
317:     eps->state = EPS_STATE_INITIAL;
318:   } else if (tol != (PetscReal)PETSC_CURRENT) {
319:     PetscCheck(tol>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of tol. Must be > 0");
320:     eps->tol = tol;
321:   }
322:   if (maxits == PETSC_DETERMINE) {
323:     eps->max_it = PETSC_DETERMINE;
324:     eps->state  = EPS_STATE_INITIAL;
325:   } else if (maxits == PETSC_UNLIMITED) {
326:     eps->max_it = PETSC_INT_MAX;
327:   } else if (maxits != PETSC_CURRENT) {
328:     PetscCheck(maxits>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of maxits. Must be > 0");
329:     eps->max_it = maxits;
330:   }
331:   PetscFunctionReturn(PETSC_SUCCESS);
332: }

334: /*@
335:    EPSGetDimensions - Gets the number of eigenvalues to compute
336:    and the dimension of the subspace.

338:    Not Collective

340:    Input Parameter:
341: .  eps - the eigensolver context

343:    Output Parameters:
344: +  nev - number of eigenvalues to compute
345: .  ncv - the maximum dimension of the subspace to be used by the solver
346: -  mpd - the maximum dimension allowed for the projected problem

348:    Level: intermediate

350: .seealso: EPSSetDimensions()
351: @*/
352: PetscErrorCode EPSGetDimensions(EPS eps,PetscInt *nev,PetscInt *ncv,PetscInt *mpd)
353: {
354:   PetscFunctionBegin;
356:   if (nev) *nev = eps->nev? eps->nev: 1;
357:   if (ncv) *ncv = eps->ncv;
358:   if (mpd) *mpd = eps->mpd;
359:   PetscFunctionReturn(PETSC_SUCCESS);
360: }

362: /*@
363:    EPSSetDimensions - Sets the number of eigenvalues to compute
364:    and the dimension of the subspace.

366:    Logically Collective

368:    Input Parameters:
369: +  eps - the eigensolver context
370: .  nev - number of eigenvalues to compute
371: .  ncv - the maximum dimension of the subspace to be used by the solver
372: -  mpd - the maximum dimension allowed for the projected problem

374:    Options Database Keys:
375: +  -eps_nev <nev> - Sets the number of eigenvalues
376: .  -eps_ncv <ncv> - Sets the dimension of the subspace
377: -  -eps_mpd <mpd> - Sets the maximum projected dimension

379:    Notes:
380:    Use PETSC_DETERMINE for ncv and mpd to assign a reasonably good value, which is
381:    dependent on the solution method. For any of the arguments, use PETSC_CURRENT
382:    to preserve the current value.

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

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

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

397:    Level: intermediate

399: .seealso: EPSGetDimensions(), EPSSetInterval(), EPSKrylovSchurSetDimensions()
400: @*/
401: PetscErrorCode EPSSetDimensions(EPS eps,PetscInt nev,PetscInt ncv,PetscInt mpd)
402: {
403:   PetscFunctionBegin;
408:   if (nev != PETSC_CURRENT) {
409:     PetscCheck(nev>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of nev. Must be > 0");
410:     eps->nev = nev;
411:   }
412:   if (ncv == PETSC_DETERMINE) {
413:     eps->ncv = PETSC_DETERMINE;
414:   } else if (ncv != PETSC_CURRENT) {
415:     PetscCheck(ncv>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of ncv. Must be > 0");
416:     eps->ncv = ncv;
417:   }
418:   if (mpd == PETSC_DETERMINE) {
419:     eps->mpd = PETSC_DETERMINE;
420:   } else if (mpd != PETSC_CURRENT) {
421:     PetscCheck(mpd>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of mpd. Must be > 0");
422:     eps->mpd = mpd;
423:   }
424:   eps->state = EPS_STATE_INITIAL;
425:   PetscFunctionReturn(PETSC_SUCCESS);
426: }

428: /*@
429:    EPSSetWhichEigenpairs - Specifies which portion of the spectrum is
430:    to be sought.

432:    Logically Collective

434:    Input Parameters:
435: +  eps   - eigensolver context obtained from EPSCreate()
436: -  which - the portion of the spectrum to be sought

438:    Options Database Keys:
439: +   -eps_largest_magnitude - Sets largest eigenvalues in magnitude
440: .   -eps_smallest_magnitude - Sets smallest eigenvalues in magnitude
441: .   -eps_largest_real - Sets largest real parts
442: .   -eps_smallest_real - Sets smallest real parts
443: .   -eps_largest_imaginary - Sets largest imaginary parts
444: .   -eps_smallest_imaginary - Sets smallest imaginary parts
445: .   -eps_target_magnitude - Sets eigenvalues closest to target
446: .   -eps_target_real - Sets real parts closest to target
447: .   -eps_target_imaginary - Sets imaginary parts closest to target
448: -   -eps_all - Sets all eigenvalues in an interval or region

450:    Notes:
451:    The parameter 'which' can have one of these values

453: +     EPS_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
454: .     EPS_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
455: .     EPS_LARGEST_REAL - largest real parts
456: .     EPS_SMALLEST_REAL - smallest real parts
457: .     EPS_LARGEST_IMAGINARY - largest imaginary parts
458: .     EPS_SMALLEST_IMAGINARY - smallest imaginary parts
459: .     EPS_TARGET_MAGNITUDE - eigenvalues closest to the target (in magnitude)
460: .     EPS_TARGET_REAL - eigenvalues with real part closest to target
461: .     EPS_TARGET_IMAGINARY - eigenvalues with imaginary part closest to target
462: .     EPS_ALL - all eigenvalues contained in a given interval or region
463: -     EPS_WHICH_USER - user defined ordering set with EPSSetEigenvalueComparison()

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

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

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

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

482:    Level: intermediate

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

521: /*@
522:    EPSGetWhichEigenpairs - Returns which portion of the spectrum is to be
523:    sought.

525:    Not Collective

527:    Input Parameter:
528: .  eps - eigensolver context obtained from EPSCreate()

530:    Output Parameter:
531: .  which - the portion of the spectrum to be sought

533:    Notes:
534:    See EPSSetWhichEigenpairs() for possible values of 'which'.

536:    Level: intermediate

538: .seealso: EPSSetWhichEigenpairs(), EPSWhich
539: @*/
540: PetscErrorCode EPSGetWhichEigenpairs(EPS eps,EPSWhich *which)
541: {
542:   PetscFunctionBegin;
544:   PetscAssertPointer(which,2);
545:   *which = eps->which;
546:   PetscFunctionReturn(PETSC_SUCCESS);
547: }

549: /*@
550:    EPSSetThreshold - Sets the threshold used in the threshold stopping test.

552:    Logically Collective

554:    Input Parameters:
555: +  eps   - the eigenvalue solver context
556: .  thres - the threshold value
557: -  rel   - whether the threshold is relative or not

559:    Options Database Keys:
560: +  -eps_threshold_absolute <thres> - Sets an absolute threshold
561: -  -eps_threshold_relative <thres> - Sets a relative threshold

563:    Notes:
564:    This function internally calls EPSSetStoppingTest() to set a special stopping
565:    test based on the threshold, where eigenvalues are computed in sequence until
566:    one of the computed eigenvalues is below the threshold (in magnitude). This is
567:    the interpretation in case of searching for largest eigenvalues in magnitude,
568:    see EPSSetWhichEigenpairs().

570:    If the solver is configured to compute smallest magnitude eigenvalues, then the
571:    threshold must be interpreted in the opposite direction, i.e., the computation
572:    will stop when one of the computed values is above the threshold (in magnitude).

574:    The threshold can also be used when computing largest/smallest real eigenvalues
575:    (i.e, rightmost or leftmost), in which case the threshold is allowed to be
576:    negative. The solver will stop when one of the computed eigenvalues is above
577:    or below the threshold (considering the real part of the eigenvalue). This mode
578:    is allowed only in problem types whose eigenvalues are always real (e.g., HEP).

580:    In the case of largest magnitude eigenvalues, the threshold can be made relative
581:    with respect to the dominant eigenvalue. Otherwise, the argument rel should be
582:    PETSC_FALSE.

584:    An additional use case is with target magnitude selection of eigenvalues (e.g.,
585:    with shift-and-invert), but this must be used with caution to avoid unexpected
586:    behaviour. With an absolute threshold, the solver will assume that leftmost
587:    eigenvalues are being computed (e.g., with target=0 for a problem with real
588:    positive eigenvalues). In case of a relative threshold, a value of threshold<1
589:    implies that the wanted eigenvalues are the largest ones, and otherwise the
590:    solver assumes that smallest eigenvalues are being computed.

592:    The test against the threshold is done for converged eigenvalues, which
593:    implies that the final number of converged eigenvalues will be at least
594:    one more than the actual number of values below/above the threshold.

596:    Since the number of computed eigenvalues is not known a priori, the solver
597:    will need to reallocate the basis of vectors internally, to have enough room
598:    to accommodate all the eigenvectors. Hence, this option must be used with
599:    caution to avoid out-of-memory problems. The recommendation is to set the value
600:    of ncv to be larger than the estimated number of eigenvalues, to minimize the
601:    number of reallocations.

603:    If a number of wanted eigenvalues has been set with EPSSetDimensions()
604:    it is also taken into account and the solver will stop when one of the two
605:    conditions (threshold or number of converged values) is met.

607:    Use EPSSetStoppingTest() to return to the usual computation of a fixed number
608:    of eigenvalues.

610:    Level: advanced

612: .seealso: EPSGetThreshold(), EPSSetStoppingTest(), EPSSetDimensions(), EPSSetWhichEigenpairs(), EPSSetProblemType()
613: @*/
614: PetscErrorCode EPSSetThreshold(EPS eps,PetscReal thres,PetscBool rel)
615: {
616:   PetscFunctionBegin;
620:   if (eps->thres != thres || eps->threlative != rel) {
621:     eps->thres = thres;
622:     eps->threlative = rel;
623:     eps->state = EPS_STATE_INITIAL;
624:     PetscCall(EPSSetStoppingTest(eps,EPS_STOP_THRESHOLD));
625:   }
626:   PetscFunctionReturn(PETSC_SUCCESS);
627: }

629: /*@
630:    EPSGetThreshold - Gets the threshold used by the threshold stopping test.

632:    Not Collective

634:    Input Parameter:
635: .  eps - the eigenvalue solver context

637:    Output Parameters:
638: +  thres - the threshold
639: -  rel   - whether the threshold is relative or not

641:    Level: advanced

643: .seealso: EPSSetThreshold()
644: @*/
645: PetscErrorCode EPSGetThreshold(EPS eps,PetscReal *thres,PetscBool *rel)
646: {
647:   PetscFunctionBegin;
649:   if (thres) *thres = eps->thres;
650:   if (rel)   *rel   = eps->threlative;
651:   PetscFunctionReturn(PETSC_SUCCESS);
652: }

654: /*@C
655:    EPSSetEigenvalueComparison - Specifies the eigenvalue comparison function
656:    when EPSSetWhichEigenpairs() is set to EPS_WHICH_USER.

658:    Logically Collective

660:    Input Parameters:
661: +  eps  - eigensolver context obtained from EPSCreate()
662: .  func - the comparison function, see EPSEigenvalueComparisonFn for the calling sequence
663: -  ctx  - a context pointer (the last parameter to the comparison function)

665:    Note:
666:    The returning parameter 'res' can be
667: +  negative - if the 1st eigenvalue is preferred to the 2st one
668: .  zero     - if both eigenvalues are equally preferred
669: -  positive - if the 2st eigenvalue is preferred to the 1st one

671:    Level: advanced

673: .seealso: EPSSetWhichEigenpairs(), EPSWhich
674: @*/
675: PetscErrorCode EPSSetEigenvalueComparison(EPS eps,SlepcEigenvalueComparisonFn *func,void* ctx)
676: {
677:   PetscFunctionBegin;
679:   eps->sc->comparison    = func;
680:   eps->sc->comparisonctx = ctx;
681:   eps->which             = EPS_WHICH_USER;
682:   PetscFunctionReturn(PETSC_SUCCESS);
683: }

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

690:    Logically Collective

692:    Input Parameters:
693: +  eps  - eigensolver context obtained from EPSCreate()
694: .  func - the arbitrary selection function, see SlepcArbitrarySelectionFn for a calling sequence
695: -  ctx  - a context pointer (the last parameter to the arbitrary selection function)

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

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

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

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

717:    Level: advanced

719: .seealso: EPSSetWhichEigenpairs()
720: @*/
721: PetscErrorCode EPSSetArbitrarySelection(EPS eps,SlepcArbitrarySelectionFn *func,void* ctx)
722: {
723:   PetscFunctionBegin;
725:   eps->arbitrary    = func;
726:   eps->arbitraryctx = ctx;
727:   eps->state        = EPS_STATE_INITIAL;
728:   PetscFunctionReturn(PETSC_SUCCESS);
729: }

731: /*@C
732:    EPSSetConvergenceTestFunction - Sets a function to compute the error estimate
733:    used in the convergence test.

735:    Logically Collective

737:    Input Parameters:
738: +  eps     - eigensolver context obtained from EPSCreate()
739: .  func    - convergence test function, see EPSConvergenceTestFn for the calling sequence
740: .  ctx     - context for private data for the convergence routine (may be NULL)
741: -  destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence

743:    Note:
744:    If the error estimate returned by the convergence test function is less than
745:    the tolerance, then the eigenvalue is accepted as converged.

747:    Level: advanced

749: .seealso: EPSSetConvergenceTest(), EPSSetTolerances()
750: @*/
751: PetscErrorCode EPSSetConvergenceTestFunction(EPS eps,EPSConvergenceTestFn *func,void* ctx,PetscCtxDestroyFn *destroy)
752: {
753:   PetscFunctionBegin;
755:   if (eps->convergeddestroy) PetscCall((*eps->convergeddestroy)(&eps->convergedctx));
756:   eps->convergeduser    = func;
757:   eps->convergeddestroy = destroy;
758:   eps->convergedctx     = ctx;
759:   if (func == EPSConvergedRelative) eps->conv = EPS_CONV_REL;
760:   else if (func == EPSConvergedNorm) eps->conv = EPS_CONV_NORM;
761:   else if (func == EPSConvergedAbsolute) eps->conv = EPS_CONV_ABS;
762:   else {
763:     eps->conv      = EPS_CONV_USER;
764:     eps->converged = eps->convergeduser;
765:   }
766:   PetscFunctionReturn(PETSC_SUCCESS);
767: }

769: /*@
770:    EPSSetConvergenceTest - Specifies how to compute the error estimate
771:    used in the convergence test.

773:    Logically Collective

775:    Input Parameters:
776: +  eps  - eigensolver context obtained from EPSCreate()
777: -  conv - the type of convergence test

779:    Options Database Keys:
780: +  -eps_conv_abs  - Sets the absolute convergence test
781: .  -eps_conv_rel  - Sets the convergence test relative to the eigenvalue
782: .  -eps_conv_norm - Sets the convergence test relative to the matrix norms
783: -  -eps_conv_user - Selects the user-defined convergence test

785:    Note:
786:    The parameter 'conv' can have one of these values
787: +     EPS_CONV_ABS  - absolute error ||r||
788: .     EPS_CONV_REL  - error relative to the eigenvalue l, ||r||/|l|
789: .     EPS_CONV_NORM - error relative to the matrix norms, ||r||/(||A||+|l|*||B||)
790: -     EPS_CONV_USER - function set by EPSSetConvergenceTestFunction()

792:    Level: intermediate

794: .seealso: EPSGetConvergenceTest(), EPSSetConvergenceTestFunction(), EPSSetStoppingTest(), EPSConv
795: @*/
796: PetscErrorCode EPSSetConvergenceTest(EPS eps,EPSConv conv)
797: {
798:   PetscFunctionBegin;
801:   switch (conv) {
802:     case EPS_CONV_ABS:  eps->converged = EPSConvergedAbsolute; break;
803:     case EPS_CONV_REL:  eps->converged = EPSConvergedRelative; break;
804:     case EPS_CONV_NORM: eps->converged = EPSConvergedNorm; break;
805:     case EPS_CONV_USER:
806:       PetscCheck(eps->convergeduser,PetscObjectComm((PetscObject)eps),PETSC_ERR_ORDER,"Must call EPSSetConvergenceTestFunction() first");
807:       eps->converged = eps->convergeduser;
808:       break;
809:     default:
810:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'conv' value");
811:   }
812:   eps->conv = conv;
813:   PetscFunctionReturn(PETSC_SUCCESS);
814: }

816: /*@
817:    EPSGetConvergenceTest - Gets the method used to compute the error estimate
818:    used in the convergence test.

820:    Not Collective

822:    Input Parameters:
823: .  eps   - eigensolver context obtained from EPSCreate()

825:    Output Parameters:
826: .  conv  - the type of convergence test

828:    Level: intermediate

830: .seealso: EPSSetConvergenceTest(), EPSConv
831: @*/
832: PetscErrorCode EPSGetConvergenceTest(EPS eps,EPSConv *conv)
833: {
834:   PetscFunctionBegin;
836:   PetscAssertPointer(conv,2);
837:   *conv = eps->conv;
838:   PetscFunctionReturn(PETSC_SUCCESS);
839: }

841: /*@C
842:    EPSSetStoppingTestFunction - Sets a function to decide when to stop the outer
843:    iteration of the eigensolver.

845:    Logically Collective

847:    Input Parameters:
848: +  eps     - eigensolver context obtained from EPSCreate()
849: .  func    - stopping test function, see EPSStoppingTestFn for the calling sequence
850: .  ctx     - context for private data for the stopping routine (may be NULL)
851: -  destroy - a routine for destroying the context (may be NULL), see PetscCtxDestroyFn for the calling sequence

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

859:    Level: advanced

861: .seealso: EPSSetStoppingTest(), EPSStoppingBasic()
862: @*/
863: PetscErrorCode EPSSetStoppingTestFunction(EPS eps,EPSStoppingTestFn *func,void* ctx,PetscCtxDestroyFn *destroy)
864: {
865:   PetscFunctionBegin;
867:   if (eps->stoppingdestroy) PetscCall((*eps->stoppingdestroy)(&eps->stoppingctx));
868:   eps->stoppinguser    = func;
869:   eps->stoppingdestroy = destroy;
870:   eps->stoppingctx     = ctx;
871:   if (func == EPSStoppingBasic) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_BASIC));
872:   else if (func == EPSStoppingThreshold) PetscCall(EPSSetStoppingTest(eps,EPS_STOP_THRESHOLD));
873:   else {
874:     eps->stop     = EPS_STOP_USER;
875:     eps->stopping = eps->stoppinguser;
876:   }
877:   PetscFunctionReturn(PETSC_SUCCESS);
878: }

880: /*@
881:    EPSSetStoppingTest - Specifies how to decide the termination of the outer
882:    loop of the eigensolver.

884:    Logically Collective

886:    Input Parameters:
887: +  eps  - eigensolver context obtained from EPSCreate()
888: -  stop - the type of stopping test

890:    Options Database Keys:
891: +  -eps_stop_basic     - Sets the default stopping test
892: .  -eps_stop_threshold - Sets the threshold stopping test
893: -  -eps_stop_user      - Selects the user-defined stopping test

895:    Note:
896:    The parameter 'stop' can have one of these values
897: +     EPS_STOP_BASIC     - default stopping test
898: .     EPS_STOP_THRESHOLD - threshold stopping test)
899: -     EPS_STOP_USER      - function set by EPSSetStoppingTestFunction()

901:    Level: advanced

903: .seealso: EPSGetStoppingTest(), EPSSetStoppingTestFunction(), EPSSetConvergenceTest(), EPSStop
904: @*/
905: PetscErrorCode EPSSetStoppingTest(EPS eps,EPSStop stop)
906: {
907:   PetscFunctionBegin;
910:   switch (stop) {
911:     case EPS_STOP_BASIC: eps->stopping = EPSStoppingBasic; break;
912:     case EPS_STOP_THRESHOLD: eps->stopping = EPSStoppingThreshold; break;
913:     case EPS_STOP_USER:
914:       PetscCheck(eps->stoppinguser,PetscObjectComm((PetscObject)eps),PETSC_ERR_ORDER,"Must call EPSSetStoppingTestFunction() first");
915:       eps->stopping = eps->stoppinguser;
916:       break;
917:     default:
918:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid 'stop' value");
919:   }
920:   eps->stop = stop;
921:   PetscFunctionReturn(PETSC_SUCCESS);
922: }

924: /*@
925:    EPSGetStoppingTest - Gets the method used to decide the termination of the outer
926:    loop of the eigensolver.

928:    Not Collective

930:    Input Parameters:
931: .  eps   - eigensolver context obtained from EPSCreate()

933:    Output Parameters:
934: .  stop  - the type of stopping test

936:    Level: advanced

938: .seealso: EPSSetStoppingTest(), EPSStop
939: @*/
940: PetscErrorCode EPSGetStoppingTest(EPS eps,EPSStop *stop)
941: {
942:   PetscFunctionBegin;
944:   PetscAssertPointer(stop,2);
945:   *stop = eps->stop;
946:   PetscFunctionReturn(PETSC_SUCCESS);
947: }

949: /*@
950:    EPSSetProblemType - Specifies the type of the eigenvalue problem.

952:    Logically Collective

954:    Input Parameters:
955: +  eps      - the eigensolver context
956: -  type     - a known type of eigenvalue problem

958:    Options Database Keys:
959: +  -eps_hermitian - Hermitian eigenvalue problem
960: .  -eps_gen_hermitian - generalized Hermitian eigenvalue problem
961: .  -eps_non_hermitian - non-Hermitian eigenvalue problem
962: .  -eps_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
963: .  -eps_pos_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
964:    with positive semi-definite B
965: .  -eps_gen_indefinite - generalized Hermitian-indefinite eigenvalue problem
966: -  -eps_bse - structured Bethe-Salpeter eigenvalue problem

968:    Notes:
969:    This function must be used to instruct SLEPc to exploit symmetry or other
970:    kind of structure. If no
971:    problem type is specified, by default a non-Hermitian problem is assumed
972:    (either standard or generalized). If the user knows that the problem is
973:    Hermitian (i.e. A=A^H) or generalized Hermitian (i.e. A=A^H, B=B^H, and
974:    B positive definite) then it is recommended to set the problem type so
975:    that eigensolver can exploit these properties.

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

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

984:    Level: intermediate

986: .seealso: EPSSetOperators(), EPSSetType(), EPSGetProblemType(), EPSProblemType
987: @*/
988: PetscErrorCode EPSSetProblemType(EPS eps,EPSProblemType type)
989: {
990:   PetscFunctionBegin;
993:   if (type == eps->problem_type) PetscFunctionReturn(PETSC_SUCCESS);
994:   switch (type) {
995:     case EPS_HEP:
996:       eps->isgeneralized = PETSC_FALSE;
997:       eps->ishermitian = PETSC_TRUE;
998:       eps->ispositive = PETSC_FALSE;
999:       eps->isstructured = PETSC_FALSE;
1000:       break;
1001:     case EPS_NHEP:
1002:       eps->isgeneralized = PETSC_FALSE;
1003:       eps->ishermitian = PETSC_FALSE;
1004:       eps->ispositive = PETSC_FALSE;
1005:       eps->isstructured = PETSC_FALSE;
1006:       break;
1007:     case EPS_GHEP:
1008:       eps->isgeneralized = PETSC_TRUE;
1009:       eps->ishermitian = PETSC_TRUE;
1010:       eps->ispositive = PETSC_TRUE;
1011:       eps->isstructured = PETSC_FALSE;
1012:       break;
1013:     case EPS_GNHEP:
1014:       eps->isgeneralized = PETSC_TRUE;
1015:       eps->ishermitian = PETSC_FALSE;
1016:       eps->ispositive = PETSC_FALSE;
1017:       eps->isstructured = PETSC_FALSE;
1018:       break;
1019:     case EPS_PGNHEP:
1020:       eps->isgeneralized = PETSC_TRUE;
1021:       eps->ishermitian = PETSC_FALSE;
1022:       eps->ispositive = PETSC_TRUE;
1023:       eps->isstructured = PETSC_FALSE;
1024:       break;
1025:     case EPS_GHIEP:
1026:       eps->isgeneralized = PETSC_TRUE;
1027:       eps->ishermitian = PETSC_TRUE;
1028:       eps->ispositive = PETSC_FALSE;
1029:       eps->isstructured = PETSC_FALSE;
1030:       break;
1031:     case EPS_BSE:
1032:       eps->isgeneralized = PETSC_FALSE;
1033:       eps->ishermitian = PETSC_FALSE;
1034:       eps->ispositive = PETSC_FALSE;
1035:       eps->isstructured = PETSC_TRUE;
1036:       break;
1037:     default:
1038:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type");
1039:   }
1040:   eps->problem_type = type;
1041:   eps->state = EPS_STATE_INITIAL;
1042:   PetscFunctionReturn(PETSC_SUCCESS);
1043: }

1045: /*@
1046:    EPSGetProblemType - Gets the problem type from the EPS object.

1048:    Not Collective

1050:    Input Parameter:
1051: .  eps - the eigensolver context

1053:    Output Parameter:
1054: .  type - the problem type

1056:    Level: intermediate

1058: .seealso: EPSSetProblemType(), EPSProblemType
1059: @*/
1060: PetscErrorCode EPSGetProblemType(EPS eps,EPSProblemType *type)
1061: {
1062:   PetscFunctionBegin;
1064:   PetscAssertPointer(type,2);
1065:   *type = eps->problem_type;
1066:   PetscFunctionReturn(PETSC_SUCCESS);
1067: }

1069: /*@
1070:    EPSSetExtraction - Specifies the type of extraction technique to be employed
1071:    by the eigensolver.

1073:    Logically Collective

1075:    Input Parameters:
1076: +  eps  - the eigensolver context
1077: -  extr - a known type of extraction

1079:    Options Database Keys:
1080: +  -eps_ritz - Rayleigh-Ritz extraction
1081: .  -eps_harmonic - harmonic Ritz extraction
1082: .  -eps_harmonic_relative - harmonic Ritz extraction relative to the eigenvalue
1083: .  -eps_harmonic_right - harmonic Ritz extraction for rightmost eigenvalues
1084: .  -eps_harmonic_largest - harmonic Ritz extraction for largest magnitude
1085:    (without target)
1086: .  -eps_refined - refined Ritz extraction
1087: -  -eps_refined_harmonic - refined harmonic Ritz extraction

1089:    Notes:
1090:    Not all eigensolvers support all types of extraction. See the SLEPc
1091:    Users Manual for details.

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

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

1098:    Level: advanced

1100: .seealso: EPSSetTarget(), EPSGetExtraction(), EPSExtraction
1101: @*/
1102: PetscErrorCode EPSSetExtraction(EPS eps,EPSExtraction extr)
1103: {
1104:   PetscFunctionBegin;
1107:   if (eps->extraction != extr) {
1108:     eps->state      = EPS_STATE_INITIAL;
1109:     eps->extraction = extr;
1110:   }
1111:   PetscFunctionReturn(PETSC_SUCCESS);
1112: }

1114: /*@
1115:    EPSGetExtraction - Gets the extraction type used by the EPS object.

1117:    Not Collective

1119:    Input Parameter:
1120: .  eps - the eigensolver context

1122:    Output Parameter:
1123: .  extr - name of extraction type

1125:    Level: advanced

1127: .seealso: EPSSetExtraction(), EPSExtraction
1128: @*/
1129: PetscErrorCode EPSGetExtraction(EPS eps,EPSExtraction *extr)
1130: {
1131:   PetscFunctionBegin;
1133:   PetscAssertPointer(extr,2);
1134:   *extr = eps->extraction;
1135:   PetscFunctionReturn(PETSC_SUCCESS);
1136: }

1138: /*@
1139:    EPSSetBalance - Specifies the balancing technique to be employed by the
1140:    eigensolver, and some parameters associated to it.

1142:    Logically Collective

1144:    Input Parameters:
1145: +  eps    - the eigensolver context
1146: .  bal    - the balancing method, one of EPS_BALANCE_NONE, EPS_BALANCE_ONESIDE,
1147:             EPS_BALANCE_TWOSIDE, or EPS_BALANCE_USER
1148: .  its    - number of iterations of the balancing algorithm
1149: -  cutoff - cutoff value

1151:    Options Database Keys:
1152: +  -eps_balance <method> - the balancing method, where <method> is one of
1153:                            'none', 'oneside', 'twoside', or 'user'
1154: .  -eps_balance_its <its> - number of iterations
1155: -  -eps_balance_cutoff <cutoff> - cutoff value

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

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

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

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

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

1176:    Level: intermediate

1178: .seealso: EPSGetBalance(), EPSBalance, STSetBalanceMatrix()
1179: @*/
1180: PetscErrorCode EPSSetBalance(EPS eps,EPSBalance bal,PetscInt its,PetscReal cutoff)
1181: {
1182:   PetscFunctionBegin;
1187:   switch (bal) {
1188:     case EPS_BALANCE_NONE:
1189:     case EPS_BALANCE_ONESIDE:
1190:     case EPS_BALANCE_TWOSIDE:
1191:     case EPS_BALANCE_USER:
1192:       if (eps->balance != bal) {
1193:         eps->state = EPS_STATE_INITIAL;
1194:         eps->balance = bal;
1195:       }
1196:       break;
1197:     default:
1198:       SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Invalid value of argument 'bal'");
1199:   }
1200:   if (its==PETSC_DETERMINE) eps->balance_its = 5;
1201:   else if (its!=PETSC_CURRENT) {
1202:     PetscCheck(its>0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of its. Must be > 0");
1203:     eps->balance_its = its;
1204:   }
1205:   if (cutoff==(PetscReal)PETSC_DETERMINE) eps->balance_cutoff = 1e-8;
1206:   else if (cutoff!=(PetscReal)PETSC_CURRENT) {
1207:     PetscCheck(cutoff>0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Illegal value of cutoff. Must be > 0");
1208:     eps->balance_cutoff = cutoff;
1209:   }
1210:   PetscFunctionReturn(PETSC_SUCCESS);
1211: }

1213: /*@
1214:    EPSGetBalance - Gets the balancing type used by the EPS object, and the
1215:    associated parameters.

1217:    Not Collective

1219:    Input Parameter:
1220: .  eps - the eigensolver context

1222:    Output Parameters:
1223: +  bal    - the balancing method
1224: .  its    - number of iterations of the balancing algorithm
1225: -  cutoff - cutoff value

1227:    Level: intermediate

1229:    Note:
1230:    The user can specify NULL for any parameter that is not needed.

1232: .seealso: EPSSetBalance(), EPSBalance
1233: @*/
1234: PetscErrorCode EPSGetBalance(EPS eps,EPSBalance *bal,PetscInt *its,PetscReal *cutoff)
1235: {
1236:   PetscFunctionBegin;
1238:   if (bal)    *bal = eps->balance;
1239:   if (its)    *its = eps->balance_its;
1240:   if (cutoff) *cutoff = eps->balance_cutoff;
1241:   PetscFunctionReturn(PETSC_SUCCESS);
1242: }

1244: /*@
1245:    EPSSetTwoSided - Sets the solver to use a two-sided variant so that left
1246:    eigenvectors are also computed.

1248:    Logically Collective

1250:    Input Parameters:
1251: +  eps      - the eigensolver context
1252: -  twosided - whether the two-sided variant is to be used or not

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

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

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

1265:    Level: advanced

1267: .seealso: EPSGetTwoSided(), EPSGetLeftEigenvector()
1268: @*/
1269: PetscErrorCode EPSSetTwoSided(EPS eps,PetscBool twosided)
1270: {
1271:   PetscFunctionBegin;
1274:   if (twosided!=eps->twosided) {
1275:     eps->twosided = twosided;
1276:     eps->state    = EPS_STATE_INITIAL;
1277:   }
1278:   PetscFunctionReturn(PETSC_SUCCESS);
1279: }

1281: /*@
1282:    EPSGetTwoSided - Returns the flag indicating whether a two-sided variant
1283:    of the algorithm is being used or not.

1285:    Not Collective

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

1290:    Output Parameter:
1291: .  twosided - the returned flag

1293:    Level: advanced

1295: .seealso: EPSSetTwoSided()
1296: @*/
1297: PetscErrorCode EPSGetTwoSided(EPS eps,PetscBool *twosided)
1298: {
1299:   PetscFunctionBegin;
1301:   PetscAssertPointer(twosided,2);
1302:   *twosided = eps->twosided;
1303:   PetscFunctionReturn(PETSC_SUCCESS);
1304: }

1306: /*@
1307:    EPSSetTrueResidual - Specifies if the solver must compute the true residual
1308:    explicitly or not.

1310:    Logically Collective

1312:    Input Parameters:
1313: +  eps     - the eigensolver context
1314: -  trueres - whether true residuals are required or not

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

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

1329:    Level: advanced

1331: .seealso: EPSGetTrueResidual()
1332: @*/
1333: PetscErrorCode EPSSetTrueResidual(EPS eps,PetscBool trueres)
1334: {
1335:   PetscFunctionBegin;
1338:   eps->trueres = trueres;
1339:   PetscFunctionReturn(PETSC_SUCCESS);
1340: }

1342: /*@
1343:    EPSGetTrueResidual - Returns the flag indicating whether true
1344:    residuals must be computed explicitly or not.

1346:    Not Collective

1348:    Input Parameter:
1349: .  eps - the eigensolver context

1351:    Output Parameter:
1352: .  trueres - the returned flag

1354:    Level: advanced

1356: .seealso: EPSSetTrueResidual()
1357: @*/
1358: PetscErrorCode EPSGetTrueResidual(EPS eps,PetscBool *trueres)
1359: {
1360:   PetscFunctionBegin;
1362:   PetscAssertPointer(trueres,2);
1363:   *trueres = eps->trueres;
1364:   PetscFunctionReturn(PETSC_SUCCESS);
1365: }

1367: /*@
1368:    EPSSetTrackAll - Specifies if the solver must compute the residual norm of all
1369:    approximate eigenpairs or not.

1371:    Logically Collective

1373:    Input Parameters:
1374: +  eps      - the eigensolver context
1375: -  trackall - whether to compute all residuals or not

1377:    Notes:
1378:    If the user sets trackall=PETSC_TRUE then the solver computes (or estimates)
1379:    the residual norm for each eigenpair approximation. Computing the residual is
1380:    usually an expensive operation and solvers commonly compute only the residual
1381:    associated to the first unconverged eigenpair.

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

1385:    Level: developer

1387: .seealso: EPSGetTrackAll()
1388: @*/
1389: PetscErrorCode EPSSetTrackAll(EPS eps,PetscBool trackall)
1390: {
1391:   PetscFunctionBegin;
1394:   eps->trackall = trackall;
1395:   PetscFunctionReturn(PETSC_SUCCESS);
1396: }

1398: /*@
1399:    EPSGetTrackAll - Returns the flag indicating whether all residual norms must
1400:    be computed or not.

1402:    Not Collective

1404:    Input Parameter:
1405: .  eps - the eigensolver context

1407:    Output Parameter:
1408: .  trackall - the returned flag

1410:    Level: developer

1412: .seealso: EPSSetTrackAll()
1413: @*/
1414: PetscErrorCode EPSGetTrackAll(EPS eps,PetscBool *trackall)
1415: {
1416:   PetscFunctionBegin;
1418:   PetscAssertPointer(trackall,2);
1419:   *trackall = eps->trackall;
1420:   PetscFunctionReturn(PETSC_SUCCESS);
1421: }

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

1426:    Logically Collective

1428:    Input Parameters:
1429: +  eps    - the eigensolver context
1430: -  purify - whether purification is required or not

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

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

1441:    Level: intermediate

1443: .seealso: EPSGetPurify(), EPSSetInterval()
1444: @*/
1445: PetscErrorCode EPSSetPurify(EPS eps,PetscBool purify)
1446: {
1447:   PetscFunctionBegin;
1450:   if (purify!=eps->purify) {
1451:     eps->purify = purify;
1452:     eps->state  = EPS_STATE_INITIAL;
1453:   }
1454:   PetscFunctionReturn(PETSC_SUCCESS);
1455: }

1457: /*@
1458:    EPSGetPurify - Returns the flag indicating whether purification is activated
1459:    or not.

1461:    Not Collective

1463:    Input Parameter:
1464: .  eps - the eigensolver context

1466:    Output Parameter:
1467: .  purify - the returned flag

1469:    Level: intermediate

1471: .seealso: EPSSetPurify()
1472: @*/
1473: PetscErrorCode EPSGetPurify(EPS eps,PetscBool *purify)
1474: {
1475:   PetscFunctionBegin;
1477:   PetscAssertPointer(purify,2);
1478:   *purify = eps->purify;
1479:   PetscFunctionReturn(PETSC_SUCCESS);
1480: }

1482: /*@
1483:    EPSSetOptionsPrefix - Sets the prefix used for searching for all
1484:    EPS options in the database.

1486:    Logically Collective

1488:    Input Parameters:
1489: +  eps - the eigensolver context
1490: -  prefix - the prefix string to prepend to all EPS option requests

1492:    Notes:
1493:    A hyphen (-) must NOT be given at the beginning of the prefix name.
1494:    The first character of all runtime options is AUTOMATICALLY the
1495:    hyphen.

1497:    For example, to distinguish between the runtime options for two
1498:    different EPS contexts, one could call
1499: .vb
1500:       EPSSetOptionsPrefix(eps1,"eig1_")
1501:       EPSSetOptionsPrefix(eps2,"eig2_")
1502: .ve

1504:    Level: advanced

1506: .seealso: EPSAppendOptionsPrefix(), EPSGetOptionsPrefix()
1507: @*/
1508: PetscErrorCode EPSSetOptionsPrefix(EPS eps,const char *prefix)
1509: {
1510:   PetscFunctionBegin;
1512:   if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
1513:   PetscCall(STSetOptionsPrefix(eps->st,prefix));
1514:   if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
1515:   PetscCall(BVSetOptionsPrefix(eps->V,prefix));
1516:   if (!eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
1517:   PetscCall(DSSetOptionsPrefix(eps->ds,prefix));
1518:   if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
1519:   PetscCall(RGSetOptionsPrefix(eps->rg,prefix));
1520:   PetscCall(PetscObjectSetOptionsPrefix((PetscObject)eps,prefix));
1521:   PetscFunctionReturn(PETSC_SUCCESS);
1522: }

1524: /*@
1525:    EPSAppendOptionsPrefix - Appends to the prefix used for searching for all
1526:    EPS options in the database.

1528:    Logically Collective

1530:    Input Parameters:
1531: +  eps - the eigensolver context
1532: -  prefix - the prefix string to prepend to all EPS option requests

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

1538:    Level: advanced

1540: .seealso: EPSSetOptionsPrefix(), EPSGetOptionsPrefix()
1541: @*/
1542: PetscErrorCode EPSAppendOptionsPrefix(EPS eps,const char *prefix)
1543: {
1544:   PetscFunctionBegin;
1546:   if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
1547:   PetscCall(STAppendOptionsPrefix(eps->st,prefix));
1548:   if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
1549:   PetscCall(BVAppendOptionsPrefix(eps->V,prefix));
1550:   if (!eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
1551:   PetscCall(DSAppendOptionsPrefix(eps->ds,prefix));
1552:   if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
1553:   PetscCall(RGAppendOptionsPrefix(eps->rg,prefix));
1554:   PetscCall(PetscObjectAppendOptionsPrefix((PetscObject)eps,prefix));
1555:   PetscFunctionReturn(PETSC_SUCCESS);
1556: }

1558: /*@
1559:    EPSGetOptionsPrefix - Gets the prefix used for searching for all
1560:    EPS options in the database.

1562:    Not Collective

1564:    Input Parameters:
1565: .  eps - the eigensolver context

1567:    Output Parameters:
1568: .  prefix - pointer to the prefix string used is returned

1570:    Note:
1571:    On the Fortran side, the user should pass in a string 'prefix' of
1572:    sufficient length to hold the prefix.

1574:    Level: advanced

1576: .seealso: EPSSetOptionsPrefix(), EPSAppendOptionsPrefix()
1577: @*/
1578: PetscErrorCode EPSGetOptionsPrefix(EPS eps,const char *prefix[])
1579: {
1580:   PetscFunctionBegin;
1582:   PetscAssertPointer(prefix,2);
1583:   PetscCall(PetscObjectGetOptionsPrefix((PetscObject)eps,prefix));
1584:   PetscFunctionReturn(PETSC_SUCCESS);
1585: }