Actual source code: slepceps.h

  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:    User interface for the SLEPc linear eigenvalue solvers
 12: */

 14: #pragma once

 16: #include <slepcst.h>
 17: #include <slepcbv.h>
 18: #include <slepcds.h>
 19: #include <slepcrg.h>
 20: #include <slepclme.h>
 21: #include <petscsnes.h>

 23: /* SUBMANSEC = EPS */

 25: SLEPC_EXTERN PetscErrorCode EPSInitializePackage(void);
 26: SLEPC_EXTERN PetscErrorCode EPSFinalizePackage(void);

 28: /*S
 29:    EPS - SLEPc object that manages all the linear eigenvalue problem solvers.

 31:    Level: beginner

 33: .seealso: [](ch:eps), `EPSCreate()`, `ST`
 34: S*/
 35: typedef struct _p_EPS* EPS;

 37: /*J
 38:    EPSType - String with the name of a linear eigensolver.

 40:    Level: beginner

 42: .seealso: [](ch:eps), `EPSSetType()`, `EPS`
 43: J*/
 44: typedef const char *EPSType;
 45: #define EPSPOWER       "power"
 46: #define EPSSUBSPACE    "subspace"
 47: #define EPSARNOLDI     "arnoldi"
 48: #define EPSLANCZOS     "lanczos"
 49: #define EPSKRYLOVSCHUR "krylovschur"
 50: #define EPSGD          "gd"
 51: #define EPSJD          "jd"
 52: #define EPSRQCG        "rqcg"
 53: #define EPSLOBPCG      "lobpcg"
 54: #define EPSCISS        "ciss"
 55: #define EPSLYAPII      "lyapii"
 56: #define EPSLAPACK      "lapack"
 57: #define EPSARPACK      "arpack"
 58: #define EPSBLOPEX      "blopex"
 59: #define EPSPRIMME      "primme"
 60: #define EPSFEAST       "feast"
 61: #define EPSSCALAPACK   "scalapack"
 62: #define EPSELPA        "elpa"
 63: #define EPSELEMENTAL   "elemental"
 64: #define EPSEVSL        "evsl"
 65: #define EPSCHASE       "chase"

 67: /* Logging support */
 68: SLEPC_EXTERN PetscClassId EPS_CLASSID;

 70: /*E
 71:    EPSProblemType - Determines the type of eigenvalue problem.

 73:    Values:
 74: +  `EPS_HEP`    - Hermitian
 75: .  `EPS_GHEP`   - generalized Hermitian
 76: .  `EPS_NHEP`   - non-Hermitian
 77: .  `EPS_GNHEP`  - generalized non-Hermitian
 78: .  `EPS_PGNHEP` - generalized non-Hermitian with positive (semi-)definite $B$
 79: .  `EPS_GHIEP`  - generalized Hermitian-indefinite
 80: .  `EPS_BSE`    - structured Bethe-Salpeter
 81: .  `EPS_HAMILT` - structured Hamiltonian
 82: -  `EPS_LREP`   - structured Linear Response

 84:    Note:
 85:    In real scalars, one should read the term Hermitian as symmetric.

 87:    Level: intermediate

 89: .seealso: [](ch:eps), `EPSSetProblemType()`, `EPSGetProblemType()`
 90: E*/
 91: typedef enum { EPS_HEP    = 1,
 92:                EPS_GHEP   = 2,
 93:                EPS_NHEP   = 3,
 94:                EPS_GNHEP  = 4,
 95:                EPS_PGNHEP = 5,
 96:                EPS_GHIEP  = 6,
 97:                EPS_BSE    = 7,
 98:                EPS_HAMILT = 8,
 99:                EPS_LREP   = 9 } EPSProblemType;

101: /*MC
102:    EPS_HEP - A Hermitian eigenvalue problem.

104:    Note:
105:    The problem is formulated as $Ax=\lambda x$, where $A$ is real symmetric
106:    or complex Hermitian.

108:    Level: intermediate

110: .seealso: [](ch:eps), `EPSProblemType`, `EPSSetProblemType()`
111: M*/

113: /*MC
114:    EPS_GHEP - A generalized Hermitian eigenvalue problem.

116:    Note:
117:    The problem is formulated as $Ax=\lambda Bx$, where $A$ and $B$ are real
118:    symmetric or complex Hermitian, and $B$ is positive (semi-)definite.

120:    Level: intermediate

122: .seealso: [](ch:eps), `EPSProblemType`, `EPSSetProblemType()`
123: M*/

125: /*MC
126:    EPS_NHEP - A non-Hermitian eigenvalue problem.

128:    Note:
129:    The problem is formulated as $Ax=\lambda x$, where $A$ is non-symmetric
130:    (or non-Hermitian).

132:    Level: intermediate

134: .seealso: [](ch:eps), `EPSProblemType`, `EPSSetProblemType()`
135: M*/

137: /*MC
138:    EPS_GNHEP - A generalized non-Hermitian eigenvalue problem.

140:    Note:
141:    The problem is formulated as $Ax=\lambda Bx$, where $A$ or $B$ are
142:    non-symmetric (or non-Hermitian).

144:    Level: intermediate

146: .seealso: [](ch:eps), `EPSProblemType`, `EPSSetProblemType()`
147: M*/

149: /*MC
150:    EPS_PGNHEP - A generalized non-Hermitian eigenvalue problem with positive
151:    (semi-)definite $B$.

153:    Notes:
154:    The problem is formulated as $Ax=\lambda Bx$, where $A$ is non-symmetric
155:    (or non-Hermitian), but $B$ is symmetric (or Hermitian) and positive
156:    (semi-)definite.

158:    The problem will be solved with a non-Hermitian solver, but using an
159:    inner product induced by matrix $B$.

161:    Level: intermediate

163: .seealso: [](ch:eps), `EPSProblemType`, `EPSSetProblemType()`
164: M*/

166: /*MC
167:    EPS_GHIEP - A generalized Hermitian-indefinite eigenvalue problem.

169:    Notes:
170:    The problem is formulated as $Ax=\lambda Bx$, where both $A$ and $B$ are
171:    real symmetric or complex Hermitian, but $B$ is indefinite.

173:    The solver will try to exploit the symmetry by using an indefinite
174:    inner product, which may turn the computation numerically unstable.
175:    To avoid this, solve the problem as non-Hermitian.

177:    Level: intermediate

179: .seealso: [](ch:eps), `EPSProblemType`, `EPSSetProblemType()`
180: M*/

182: /*MC
183:    EPS_BSE - A structured Bethe-Salpeter eigenvalue problem.

185:    Notes:
186:    The problem is formulated as $Hx=\lambda x$, where $H$ has a Bethe-Salpeter
187:    structure,
188:      $$H = \begin{bmatrix}
189:         R & C \\
190:         -C^* & -R^T
191:         \end{bmatrix},$$
192:    where $R$ is Hermitian and $C$ is complex symmetric. Can also be used in
193:    the case of real matrices.

195:    A description of the properties of this problem can be found in {cite:p}`Alv25`
196:    and references therein.

198:    Level: intermediate

200: .seealso: [](ch:eps), [](sec:structured), `EPSProblemType`, `EPSSetProblemType()`
201: M*/

203: /*MC
204:    EPS_HAMILT - A structured Hamiltonian eigenvalue problem.

206:    Note:
207:    The problem is formulated as $Hx=\lambda x$, where $H$ has a Hamiltonian
208:    structure,
209:      $$H = \begin{bmatrix}
210:         A & B \\
211:         C & -A^*
212:         \end{bmatrix},$$
213:    where $A$, $B$ and $C$ are either real with $B=B^T$, $C=C^T$, or complex with
214:    $B=B^*$, $C=C^*$.

216:    Level: intermediate

218: .seealso: [](ch:eps), [](sec:structured), `EPSProblemType`, `EPSSetProblemType()`
219: M*/

221: /*MC
222:    EPS_LREP - A structured Linear Response eigenvalue problem.

224:    Note:
225:    The problem is formulated as $Hx=\lambda x$, where $H$ has the form
226:      $$H_1 = \begin{bmatrix}
227:         A & B \\
228:        -B & -A
229:         \end{bmatrix},$$
230:    where $A$ and $B$ are real symmetric. An alternative form is
231:      $$H_2 = \begin{bmatrix}
232:         0 & K \\
233:         M & 0
234:         \end{bmatrix},$$
235:    where $K=A-B$ and $M=A+B$.

237:    Level: intermediate

239: .seealso: [](ch:eps), [](sec:structured), `EPSProblemType`, `EPSSetProblemType()`, `MatCreateLREP()`
240: M*/

242: /*E
243:    EPSExtraction - Determines the type of extraction technique employed
244:    by the eigensolver.

246:    Values:
247: +  `EPS_RITZ`              - Rayleigh-Ritz extraction
248: .  `EPS_HARMONIC`          - harmonic Ritz extraction
249: .  `EPS_HARMONIC_RELATIVE` - harmonic Ritz extraction relative to the eigenvalue
250: .  `EPS_HARMONIC_RIGHT`    - harmonic Ritz extraction for rightmost eigenvalues
251: .  `EPS_HARMONIC_LARGEST`  - harmonic Ritz extraction for largest magnitude (without target)
252: .  `EPS_REFINED`           - refined Ritz extraction
253: -  `EPS_REFINED_HARMONIC`  - refined harmonic Ritz extraction

255:    Level: advanced

257: .seealso: [](ch:eps), `EPSSetExtraction()`, `EPSGetExtraction()`
258: E*/
259: typedef enum { EPS_RITZ,
260:                EPS_HARMONIC,
261:                EPS_HARMONIC_RELATIVE,
262:                EPS_HARMONIC_RIGHT,
263:                EPS_HARMONIC_LARGEST,
264:                EPS_REFINED,
265:                EPS_REFINED_HARMONIC } EPSExtraction;

267: /*MC
268:    EPS_RITZ - The standard Rayleigh-Ritz extraction.

270:    Note:
271:    This is the default way of computing eigenpair approximations from a
272:    given subspace.

274:    Level: advanced

276: .seealso: [](ch:eps), `EPSExtraction`, `EPSSetExtraction()`
277: M*/

279: /*MC
280:    EPS_HARMONIC - The harmonic Ritz extraction.

282:    Notes:
283:    This extraction method may provide better convergence when computing
284:    interior eigenvalues close to a given target.

286:    For the particular case of Krylov-Schur, a detailed description can
287:    be found in {cite:p}`Rom09`.

289:    Level: advanced

291: .seealso: [](ch:eps), `EPSExtraction`, `EPSSetExtraction()`, `EPSSetTarget()`
292: M*/

294: /*MC
295:    EPS_HARMONIC_RELATIVE - The harmonic Ritz extraction relative to the eigenvalue.

297:    Note:
298:    This is a variation of `EPS_HARMONIC`, used in Davidson methods only.

300:    Level: advanced

302: .seealso: [](ch:eps), `EPSExtraction`, `EPSSetExtraction()`, `EPSSetTarget()`
303: M*/

305: /*MC
306:    EPS_HARMONIC_RIGHT - The harmonic Ritz extraction for rightmost eigenvalues.

308:    Note:
309:    This is a variation of `EPS_HARMONIC`, used in Davidson methods only.

311:    Level: advanced

313: .seealso: [](ch:eps), `EPSExtraction`, `EPSSetExtraction()`, `EPSSetTarget()`
314: M*/

316: /*MC
317:    EPS_HARMONIC_LARGEST - The harmonic Ritz extraction for largest magnitude
318:    eigenvalues (without target).

320:    Note:
321:    This is a variation of `EPS_HARMONIC`, used in Davidson methods only.

323:    Level: advanced

325: .seealso: [](ch:eps), `EPSExtraction`, `EPSSetExtraction()`
326: M*/

328: /*MC
329:    EPS_REFINED - The refined Ritz extraction method {cite:p}`Jia97`.

331:    Note:
332:    Currently implemented only in `EPSARNOLDI`.

334:    Level: advanced

336: .seealso: [](ch:eps), `EPSExtraction`, `EPSSetExtraction()`
337: M*/

339: /*MC
340:    EPS_REFINED_HARMONIC - The refined harmonic Ritz extraction.

342:    Note:
343:    This is a combination of `EPS_HARMONIC` and `EPS_REFINED`.

345:    Developer Note:
346:    Currently not implemented, reserved for future use.

348:    Level: advanced

350: .seealso: [](ch:eps), `EPSExtraction`, `EPSSetExtraction()`
351: M*/

353: /*E
354:    EPSWhich - Determines which part of the spectrum is requested.

356:    Values:
357: +  `EPS_LARGEST_MAGNITUDE`  - largest $|\lambda|$
358: .  `EPS_SMALLEST_MAGNITUDE` - smallest $|\lambda|$
359: .  `EPS_LARGEST_REAL`       - largest $\mathrm{Re}(\lambda)$
360: .  `EPS_SMALLEST_REAL`      - smallest $\mathrm{Re}(\lambda)$
361: .  `EPS_LARGEST_IMAGINARY`  - largest $\mathrm{Im}(\lambda)$
362: .  `EPS_SMALLEST_IMAGINARY` - smallest $\mathrm{Im}(\lambda)$
363: .  `EPS_TARGET_MAGNITUDE`   - smallest $|\lambda-\tau|$
364: .  `EPS_TARGET_REAL`        - smallest $|\mathrm{Re}(\lambda-\tau)|$
365: .  `EPS_TARGET_IMAGINARY`   - smallest $|\mathrm{Im}(\lambda-\tau)|$
366: .  `EPS_ALL`                - all $\lambda\in[a,b]$ or $\lambda\in\Omega$
367: -  `EPS_WHICH_USER`         - user-defined sorting criterion

369:    Notes:
370:    If SLEPc is compiled for real scalars `EPS_LARGEST_IMAGINARY` and
371:    `EPS_SMALLEST_IMAGINARY` use the absolute value of the imaginary part
372:    for eigenvalue selection.

374:    The target $\tau$ is a scalar value provided with `EPSSetTarget()`.

376:    The case `EPS_ALL` needs an interval $[a,b]$ given with `EPSSetInterval()`
377:    or a region $\Omega$ specified with an `RG` object.

379:    Level: intermediate

381: .seealso: [](ch:eps), `EPSSetWhichEigenpairs()`, `EPSSetTarget()`, `EPSSetInterval()`
382: E*/
383: typedef enum { EPS_LARGEST_MAGNITUDE  = 1,
384:                EPS_SMALLEST_MAGNITUDE = 2,
385:                EPS_LARGEST_REAL       = 3,
386:                EPS_SMALLEST_REAL      = 4,
387:                EPS_LARGEST_IMAGINARY  = 5,
388:                EPS_SMALLEST_IMAGINARY = 6,
389:                EPS_TARGET_MAGNITUDE   = 7,
390:                EPS_TARGET_REAL        = 8,
391:                EPS_TARGET_IMAGINARY   = 9,
392:                EPS_ALL                = 10,
393:                EPS_WHICH_USER         = 11 } EPSWhich;

395: /*E
396:    EPSBalance - The type of balancing used for non-Hermitian problems.

398:    Values:
399: +  `EPS_BALANCE_NONE`    - no balancing matrix is used
400: .  `EPS_BALANCE_ONESIDE` - balancing matrix $D$ is computed with a one-sided Krylov method
401: .  `EPS_BALANCE_TWOSIDE` - balancing matrix $D$ is computed with a two-sided Krylov method
402: -  `EPS_BALANCE_USER`    - use a balancing matrix $D$ provided by the user

404:    Level: intermediate

406: .seealso: [](ch:eps), [](sec:balancing), `EPSSetBalance()`
407: E*/
408: typedef enum { EPS_BALANCE_NONE,
409:                EPS_BALANCE_ONESIDE,
410:                EPS_BALANCE_TWOSIDE,
411:                EPS_BALANCE_USER } EPSBalance;
412: SLEPC_EXTERN const char *EPSBalanceTypes[];

414: /*E
415:    EPSErrorType - The error type used to assess the accuracy of computed solutions.

417:    Values:
418: +  `EPS_ERROR_ABSOLUTE` - compute error bound as $\|r\|$
419: .  `EPS_ERROR_RELATIVE` - compute error bound as $\|r\|/|\lambda|$
420: -  `EPS_ERROR_BACKWARD` - compute error bound as $\|r\|/(\|A\|+|\lambda|\|B\|)$

422:    Level: intermediate

424: .seealso: [](ch:eps), `EPSComputeError()`
425: E*/
426: typedef enum { EPS_ERROR_ABSOLUTE,
427:                EPS_ERROR_RELATIVE,
428:                EPS_ERROR_BACKWARD } EPSErrorType;
429: SLEPC_EXTERN const char *EPSErrorTypes[];

431: /*E
432:    EPSConv - The convergence criterion to be used by the solver.

434:    Values:
435: +  `EPS_CONV_ABS`  - absolute convergence criterion, $\|r\|$
436: .  `EPS_CONV_REL`  - convergence criterion relative to eigenvalue, $\|r\|/|\lambda|$
437: .  `EPS_CONV_NORM` - convergence criterion relative to matrix norms, $\|r\|/(\|A\|+|\lambda|\|B\|)$
438: -  `EPS_CONV_USER` - convergence dictated by user-provided function

440:    Level: intermediate

442: .seealso: [](ch:eps), `EPSSetConvergenceTest()`, `EPSSetConvergenceTestFunction()`
443: E*/
444: typedef enum { EPS_CONV_ABS,
445:                EPS_CONV_REL,
446:                EPS_CONV_NORM,
447:                EPS_CONV_USER } EPSConv;

449: /*E
450:    EPSStop - The stopping test to decide the termination of the outer loop
451:    of the eigensolver.

453:    Values:
454: +  `EPS_STOP_BASIC`     - default stopping test
455: .  `EPS_STOP_USER`      - user-provided stopping test
456: -  `EPS_STOP_THRESHOLD` - threshold stopping test

458:    Level: advanced

460: .seealso: [](ch:eps), `EPSSetStoppingTest()`, `EPSSetStoppingTestFunction()`
461: E*/
462: typedef enum { EPS_STOP_BASIC,
463:                EPS_STOP_USER,
464:                EPS_STOP_THRESHOLD } EPSStop;

466: /*MC
467:    EPS_STOP_BASIC - The default stopping test.

469:    Note:
470:    By default, the termination of the outer loop is decided by calling
471:    `EPSStoppingBasic()`, which will stop if all requested eigenvalues are converged,
472:    or if the maximum number of iterations has been reached.

474:    Level: advanced

476: .seealso: [](ch:eps), `EPSStop`, `EPSSetStoppingTest()`, `EPSStoppingBasic()`
477: M*/

479: /*MC
480:    EPS_STOP_USER - The user-provided stopping test.

482:    Note:
483:    Customized stopping test using the user-provided function given with
484:    `EPSSetStoppingTestFunction()`.

486:    Level: advanced

488: .seealso: [](ch:eps), `EPSStop`, `EPSSetStoppingTest()`, `EPSSetStoppingTestFunction()`
489: M*/

491: /*MC
492:    EPS_STOP_THRESHOLD - The threshold stopping test.

494:    Note:
495:    When a threshold has been provided with `EPSSetThreshold()`, the termination
496:    of the outer loop is decided by calling `EPSStoppingThreshold()`, which will
497:    stop when one of the computed eigenvalues is not above/below the threshold.
498:    If a number of wanted eigenvalues has been specified via `EPSSetDimensions()`
499:    then it is also taken into account, and the solver will stop when one of the
500:    two conditions (threshold or number of converged values) is met.

502:    Level: advanced

504: .seealso: [](ch:eps), `EPSStop`, `EPSSetStoppingTest()`, `EPSStoppingThreshold()`, `EPSSetThreshold()`, `EPSSetDimensions()`
505: M*/

507: /*E
508:    EPSConvergedReason - Reason an eigensolver was determined to have converged
509:    or diverged.

511:    Values:
512: +  `EPS_CONVERGED_TOL`          - converged up to tolerance
513: .  `EPS_CONVERGED_USER`         - converged due to a user-defined condition
514: .  `EPS_DIVERGED_ITS`           - exceeded the maximum number of allowed iterations
515: .  `EPS_DIVERGED_BREAKDOWN`     - generic breakdown in method
516: .  `EPS_DIVERGED_SYMMETRY_LOST` - pseudo-Lanczos was not able to keep symmetry
517: -  `EPS_CONVERGED_ITERATING`    - the solver is still running

519:    Level: intermediate

521: .seealso: [](ch:eps), `EPSSolve()`, `EPSGetConvergedReason()`, `EPSSetTolerances()`
522: E*/
523: typedef enum {/* converged */
524:               EPS_CONVERGED_TOL                =  1,
525:               EPS_CONVERGED_USER               =  2,
526:               /* diverged */
527:               EPS_DIVERGED_ITS                 = -1,
528:               EPS_DIVERGED_BREAKDOWN           = -2,
529:               EPS_DIVERGED_SYMMETRY_LOST       = -3,
530:               EPS_CONVERGED_ITERATING          =  0} EPSConvergedReason;
531: SLEPC_EXTERN const char *const*EPSConvergedReasons;

533: /*MC
534:    EPS_CONVERGED_TOL - The computed error estimates, based on residual norms,
535:    for all requested eigenvalues are below the tolerance.

537:    Level: intermediate

539: .seealso: [](ch:eps), `EPSSolve()`, `EPSGetConvergedReason()`, `EPSConvergedReason`
540: M*/

542: /*MC
543:    EPS_CONVERGED_USER - The solver was declared converged due to a user-defined condition.

545:    Note:
546:    This happens only when a user-defined stopping test has been set with
547:    `EPSSetStoppingTestFunction()`.

549:    Level: intermediate

551: .seealso: [](ch:eps), `EPSSolve()`, `EPSGetConvergedReason()`, `EPSConvergedReason`, `EPSSetStoppingTestFunction()`
552: M*/

554: /*MC
555:    EPS_DIVERGED_ITS - Exceeded the maximum number of allowed iterations
556:    before the convergence criterion was satisfied.

558:    Level: intermediate

560: .seealso: [](ch:eps), `EPSSolve()`, `EPSGetConvergedReason()`, `EPSConvergedReason`
561: M*/

563: /*MC
564:    EPS_DIVERGED_BREAKDOWN - A breakdown in the solver was detected so the
565:    method could not continue.

567:    Level: intermediate

569: .seealso: [](ch:eps), `EPSSolve()`, `EPSGetConvergedReason()`, `EPSConvergedReason`
570: M*/

572: /*MC
573:    EPS_DIVERGED_SYMMETRY_LOST - The selected solver uses a pseudo-Lanczos recurrence,
574:    which is numerically unstable, and a symmetry test revealed that instability
575:    had appeared so the solver could not continue.

577:    Level: intermediate

579: .seealso: [](ch:eps), `EPSSolve()`, `EPSGetConvergedReason()`, `EPSConvergedReason`
580: M*/

582: /*MC
583:    EPS_CONVERGED_ITERATING - This value is returned if `EPSGetConvergedReason()` is called
584:    while `EPSSolve()` is still running.

586:    Level: intermediate

588: .seealso: [](ch:eps), `EPSSolve()`, `EPSGetConvergedReason()`, `EPSConvergedReason`
589: M*/

591: /*S
592:    EPSStoppingCtx - Data structure (C struct) to hold additional information to
593:    be used in some stopping test functions.

595:    Level: advanced

597: .seealso: [](ch:eps), `EPSSetStoppingTestFunction()`
598: S*/
599: struct _n_EPSStoppingCtx {
600:   PetscReal firstev;    /* the (absolute) value of the first converged eigenvalue */
601:   PetscReal lastev;     /* the (absolute) value of the last converged eigenvalue */
602:   PetscReal firstnc;    /* the (absolute) value of the first non-converged eigenvalue */
603:   PetscReal errest;     /* the error estimate of the first non-converged eigenvalue */
604:   PetscReal thres;      /* threshold set with EPSSetThreshold() */
605:   PetscBool threlative; /* threshold is relative */
606:   PetscInt  its;        /* iteration at which condition is met */
607:   PetscInt  napprox;    /* number of available approximations */
608:   EPSWhich  which;      /* which eigenvalues are being computed */
609: };
610: typedef struct _n_EPSStoppingCtx* EPSStoppingCtx;

612: SLEPC_EXTERN PetscErrorCode EPSCreate(MPI_Comm,EPS*);
613: SLEPC_EXTERN PetscErrorCode EPSDestroy(EPS*);
614: SLEPC_EXTERN PetscErrorCode EPSReset(EPS);
615: SLEPC_EXTERN PetscErrorCode EPSSetType(EPS,EPSType);
616: SLEPC_EXTERN PetscErrorCode EPSGetType(EPS,EPSType*);
617: SLEPC_EXTERN PetscErrorCode EPSSetProblemType(EPS,EPSProblemType);
618: SLEPC_EXTERN PetscErrorCode EPSGetProblemType(EPS,EPSProblemType*);
619: SLEPC_EXTERN PetscErrorCode EPSSetExtraction(EPS,EPSExtraction);
620: SLEPC_EXTERN PetscErrorCode EPSGetExtraction(EPS,EPSExtraction*);
621: SLEPC_EXTERN PetscErrorCode EPSSetBalance(EPS,EPSBalance,PetscInt,PetscReal);
622: SLEPC_EXTERN PetscErrorCode EPSGetBalance(EPS,EPSBalance*,PetscInt*,PetscReal*);
623: SLEPC_EXTERN PetscErrorCode EPSSetOperators(EPS,Mat,Mat);
624: SLEPC_EXTERN PetscErrorCode EPSGetOperators(EPS,Mat*,Mat*);
625: SLEPC_EXTERN PetscErrorCode EPSSetFromOptions(EPS);
626: SLEPC_EXTERN PetscErrorCode EPSSetDSType(EPS);
627: SLEPC_EXTERN PetscErrorCode EPSSetUp(EPS);
628: SLEPC_EXTERN PetscErrorCode EPSSolve(EPS);
629: SLEPC_EXTERN PetscErrorCode EPSView(EPS,PetscViewer);
630: SLEPC_EXTERN PetscErrorCode EPSViewFromOptions(EPS,PetscObject,const char[]);
631: SLEPC_EXTERN PetscErrorCode EPSErrorView(EPS,EPSErrorType,PetscViewer);
632: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "EPSErrorView()", ) static inline PetscErrorCode EPSPrintSolution(EPS eps,PetscViewer v) {return EPSErrorView(eps,EPS_ERROR_RELATIVE,v);}
633: SLEPC_EXTERN PetscErrorCode EPSErrorViewFromOptions(EPS);
634: SLEPC_EXTERN PetscErrorCode EPSConvergedReasonView(EPS,PetscViewer);
635: SLEPC_EXTERN PetscErrorCode EPSConvergedReasonViewFromOptions(EPS);
636: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "EPSConvergedReasonView()", ) static inline PetscErrorCode EPSReasonView(EPS eps,PetscViewer v) {return EPSConvergedReasonView(eps,v);}
637: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "EPSConvergedReasonViewFromOptions()", ) static inline PetscErrorCode EPSReasonViewFromOptions(EPS eps) {return EPSConvergedReasonViewFromOptions(eps);}
638: SLEPC_EXTERN PetscErrorCode EPSValuesView(EPS,PetscViewer);
639: SLEPC_EXTERN PetscErrorCode EPSValuesViewFromOptions(EPS);
640: SLEPC_EXTERN PetscErrorCode EPSVectorsView(EPS,PetscViewer);
641: SLEPC_EXTERN PetscErrorCode EPSVectorsViewFromOptions(EPS);

643: SLEPC_EXTERN PetscErrorCode EPSSetTarget(EPS,PetscScalar);
644: SLEPC_EXTERN PetscErrorCode EPSGetTarget(EPS,PetscScalar*);
645: SLEPC_EXTERN PetscErrorCode EPSSetInterval(EPS,PetscReal,PetscReal);
646: SLEPC_EXTERN PetscErrorCode EPSGetInterval(EPS,PetscReal*,PetscReal*);
647: SLEPC_EXTERN PetscErrorCode EPSSetST(EPS,ST);
648: SLEPC_EXTERN PetscErrorCode EPSGetST(EPS,ST*);
649: SLEPC_EXTERN PetscErrorCode EPSSetBV(EPS,BV);
650: SLEPC_EXTERN PetscErrorCode EPSGetBV(EPS,BV*);
651: SLEPC_EXTERN PetscErrorCode EPSSetRG(EPS,RG);
652: SLEPC_EXTERN PetscErrorCode EPSGetRG(EPS,RG*);
653: SLEPC_EXTERN PetscErrorCode EPSSetDS(EPS,DS);
654: SLEPC_EXTERN PetscErrorCode EPSGetDS(EPS,DS*);
655: SLEPC_EXTERN PetscErrorCode EPSSetTolerances(EPS,PetscReal,PetscInt);
656: SLEPC_EXTERN PetscErrorCode EPSGetTolerances(EPS,PetscReal*,PetscInt*);
657: SLEPC_EXTERN PetscErrorCode EPSSetDimensions(EPS,PetscInt,PetscInt,PetscInt);
658: SLEPC_EXTERN PetscErrorCode EPSGetDimensions(EPS,PetscInt*,PetscInt*,PetscInt*);

660: SLEPC_EXTERN PetscErrorCode EPSGetConvergedReason(EPS,EPSConvergedReason*);

662: SLEPC_EXTERN PetscErrorCode EPSGetConverged(EPS,PetscInt*);
663: SLEPC_EXTERN PetscErrorCode EPSGetEigenpair(EPS,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);
664: SLEPC_EXTERN PetscErrorCode EPSGetEigenvalue(EPS,PetscInt,PetscScalar*,PetscScalar*);
665: SLEPC_EXTERN PetscErrorCode EPSGetEigenvector(EPS,PetscInt,Vec,Vec);
666: SLEPC_EXTERN PetscErrorCode EPSGetLeftEigenvector(EPS,PetscInt,Vec,Vec);

668: SLEPC_EXTERN PetscErrorCode EPSComputeError(EPS,PetscInt,EPSErrorType,PetscReal*);
669: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "EPSComputeError()", ) static inline PetscErrorCode EPSComputeRelativeError(EPS eps,PetscInt i,PetscReal *r) {return EPSComputeError(eps,i,EPS_ERROR_RELATIVE,r);}
670: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "EPSComputeError() with EPS_ERROR_ABSOLUTE", ) static inline PetscErrorCode EPSComputeResidualNorm(EPS eps,PetscInt i,PetscReal *r) {return EPSComputeError(eps,i,EPS_ERROR_ABSOLUTE,r);}
671: SLEPC_EXTERN PetscErrorCode EPSGetInvariantSubspace(EPS,Vec[]);
672: SLEPC_EXTERN PetscErrorCode EPSGetErrorEstimate(EPS,PetscInt,PetscReal*);
673: SLEPC_EXTERN PetscErrorCode EPSGetIterationNumber(EPS,PetscInt*);

675: SLEPC_EXTERN PetscErrorCode EPSSetWhichEigenpairs(EPS,EPSWhich);
676: SLEPC_EXTERN PetscErrorCode EPSGetWhichEigenpairs(EPS,EPSWhich*);
677: SLEPC_EXTERN PetscErrorCode EPSSetThreshold(EPS,PetscReal,PetscBool);
678: SLEPC_EXTERN PetscErrorCode EPSGetThreshold(EPS,PetscReal*,PetscBool*);
679: SLEPC_EXTERN PetscErrorCode EPSSetTwoSided(EPS,PetscBool);
680: SLEPC_EXTERN PetscErrorCode EPSGetTwoSided(EPS,PetscBool*);
681: SLEPC_EXTERN PetscErrorCode EPSSetTrueResidual(EPS,PetscBool);
682: SLEPC_EXTERN PetscErrorCode EPSGetTrueResidual(EPS,PetscBool*);
683: SLEPC_EXTERN PetscErrorCode EPSSetPurify(EPS,PetscBool);
684: SLEPC_EXTERN PetscErrorCode EPSGetPurify(EPS,PetscBool*);
685: SLEPC_EXTERN PetscErrorCode EPSIsGeneralized(EPS,PetscBool*);
686: SLEPC_EXTERN PetscErrorCode EPSIsHermitian(EPS,PetscBool*);
687: SLEPC_EXTERN PetscErrorCode EPSIsPositive(EPS,PetscBool*);
688: SLEPC_EXTERN PetscErrorCode EPSIsStructured(EPS,PetscBool*);

690: SLEPC_EXTERN PetscErrorCode EPSSetTrackAll(EPS,PetscBool);
691: SLEPC_EXTERN PetscErrorCode EPSGetTrackAll(EPS,PetscBool*);

693: SLEPC_EXTERN PetscErrorCode EPSSetDeflationSpace(EPS,PetscInt,Vec[]);
694: SLEPC_EXTERN PetscErrorCode EPSSetInitialSpace(EPS,PetscInt,Vec[]);
695: SLEPC_EXTERN PetscErrorCode EPSSetLeftInitialSpace(EPS,PetscInt,Vec[]);

697: /*S
698:    EPSMonitorFn - A function prototype for functions provided to `EPSMonitorSet()`.

700:    Calling Sequence:
701: +  eps    - the linear eigensolver context
702: .  its    - iteration number
703: .  nconv  - number of converged eigenpairs
704: .  eigr   - real part of the eigenvalues
705: .  eigi   - imaginary part of the eigenvalues
706: .  errest - relative error estimates for each eigenpair
707: .  nest   - number of error estimates
708: -  ctx    - optional monitoring context, as provided with `EPSMonitorSet()`

710:    Level: intermediate

712: .seealso: [](ch:eps), `EPSMonitorSet()`
713: S*/
714: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode EPSMonitorFn(EPS eps,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,PetscCtx ctx);

716: /*S
717:    EPSMonitorRegisterFn - A function prototype for functions provided to `EPSMonitorRegister()`.

719:    Calling Sequence:
720: +  eps    - the linear eigensolver context
721: .  its    - iteration number
722: .  nconv  - number of converged eigenpairs
723: .  eigr   - real part of the eigenvalues
724: .  eigi   - imaginary part of the eigenvalues
725: .  errest - relative error estimates for each eigenpair
726: .  nest   - number of error estimates
727: -  ctx    - `PetscViewerAndFormat` object

729:    Level: advanced

731:    Note:
732:    This is an `EPSMonitorFn` specialized for a context of `PetscViewerAndFormat`.

734: .seealso: [](ch:eps), `EPSMonitorSet()`, `EPSMonitorRegister()`, `EPSMonitorFn`, `EPSMonitorRegisterCreateFn`, `EPSMonitorRegisterDestroyFn`
735: S*/
736: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode EPSMonitorRegisterFn(EPS eps,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,PetscViewerAndFormat *ctx);

738: /*S
739:    EPSMonitorRegisterCreateFn - A function prototype for functions that do the
740:    creation when provided to `EPSMonitorRegister()`.

742:    Calling Sequence:
743: +  viewer - the viewer to be used with the `EPSMonitorRegisterFn`
744: .  format - the format of the viewer
745: .  ctx    - a context for the monitor
746: -  result - a `PetscViewerAndFormat` object

748:    Level: advanced

750: .seealso: [](ch:eps), `EPSMonitorRegisterFn`, `EPSMonitorSet()`, `EPSMonitorRegister()`, `EPSMonitorFn`, `EPSMonitorRegisterDestroyFn`
751: S*/
752: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode EPSMonitorRegisterCreateFn(PetscViewer viewer,PetscViewerFormat format,PetscCtx ctx,PetscViewerAndFormat **result);

754: /*S
755:    EPSMonitorRegisterDestroyFn - A function prototype for functions that do the after
756:    use destruction when provided to `EPSMonitorRegister()`.

758:    Calling Sequence:
759: .  vf - a `PetscViewerAndFormat` object to be destroyed, including any context

761:    Level: advanced

763: .seealso: [](ch:eps), `EPSMonitorRegisterFn`, `EPSMonitorSet()`, `EPSMonitorRegister()`, `EPSMonitorFn`, `EPSMonitorRegisterCreateFn`
764: S*/
765: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode EPSMonitorRegisterDestroyFn(PetscViewerAndFormat **result);

767: SLEPC_EXTERN PetscErrorCode EPSMonitor(EPS,PetscInt,PetscInt,PetscScalar[],PetscScalar[],PetscReal[],PetscInt);
768: SLEPC_EXTERN PetscErrorCode EPSMonitorSet(EPS,EPSMonitorFn,PetscCtx,PetscCtxDestroyFn*);
769: SLEPC_EXTERN PetscErrorCode EPSMonitorCancel(EPS);
770: SLEPC_EXTERN PetscErrorCode EPSGetMonitorContext(EPS,PetscCtxRt);

772: SLEPC_EXTERN PetscErrorCode EPSMonitorSetFromOptions(EPS,const char[],const char[],PetscCtx,PetscBool);
773: SLEPC_EXTERN EPSMonitorRegisterFn        EPSMonitorFirst;
774: SLEPC_EXTERN EPSMonitorRegisterFn        EPSMonitorFirstDrawLG;
775: SLEPC_EXTERN EPSMonitorRegisterCreateFn  EPSMonitorFirstDrawLGCreate;
776: SLEPC_EXTERN EPSMonitorRegisterFn        EPSMonitorAll;
777: SLEPC_EXTERN EPSMonitorRegisterFn        EPSMonitorAllDrawLG;
778: SLEPC_EXTERN EPSMonitorRegisterCreateFn  EPSMonitorAllDrawLGCreate;
779: SLEPC_EXTERN EPSMonitorRegisterFn        EPSMonitorConverged;
780: SLEPC_EXTERN EPSMonitorRegisterCreateFn  EPSMonitorConvergedCreate;
781: SLEPC_EXTERN EPSMonitorRegisterFn        EPSMonitorConvergedDrawLG;
782: SLEPC_EXTERN EPSMonitorRegisterCreateFn  EPSMonitorConvergedDrawLGCreate;

784: SLEPC_EXTERN PetscErrorCode EPSSetOptionsPrefix(EPS,const char[]);
785: SLEPC_EXTERN PetscErrorCode EPSAppendOptionsPrefix(EPS,const char[]);
786: SLEPC_EXTERN PetscErrorCode EPSGetOptionsPrefix(EPS,const char*[]);

788: SLEPC_EXTERN PetscFunctionList EPSList;
789: SLEPC_EXTERN PetscFunctionList EPSMonitorList;
790: SLEPC_EXTERN PetscFunctionList EPSMonitorCreateList;
791: SLEPC_EXTERN PetscFunctionList EPSMonitorDestroyList;
792: SLEPC_EXTERN PetscErrorCode EPSRegister(const char[],PetscErrorCode(*)(EPS));
793: SLEPC_EXTERN PetscErrorCode EPSMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,EPSMonitorRegisterFn*,EPSMonitorRegisterCreateFn*,EPSMonitorRegisterDestroyFn*);

795: SLEPC_EXTERN PetscErrorCode EPSSetWorkVecs(EPS,PetscInt);
796: SLEPC_EXTERN PetscErrorCode EPSAllocateSolution(EPS,PetscInt);
797: SLEPC_EXTERN PetscErrorCode EPSReallocateSolution(EPS,PetscInt);

799: /*S
800:    EPSConvergenceTestFn - A prototype of an `EPS` convergence test function that
801:    would be passed to `EPSSetConvergenceTestFunction()`.

803:    Calling Sequence:
804: +  eps    - the linear eigensolver context
805: .  eigr   - real part of the eigenvalue
806: .  eigi   - imaginary part of the eigenvalue
807: .  res    - residual norm associated to the eigenpair
808: .  errest - [output] computed error estimate
809: -  ctx    - optional convergence context, as set by `EPSSetConvergenceTestFunction()`

811:    Level: advanced

813: .seealso: [](ch:eps), `EPSSetConvergenceTest()`, `EPSSetConvergenceTestFunction()`
814: S*/
815: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode EPSConvergenceTestFn(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,PetscCtx ctx);

817: SLEPC_EXTERN PetscErrorCode EPSSetConvergenceTest(EPS,EPSConv);
818: SLEPC_EXTERN PetscErrorCode EPSGetConvergenceTest(EPS,EPSConv*);
819: SLEPC_EXTERN EPSConvergenceTestFn EPSConvergedAbsolute;
820: SLEPC_EXTERN EPSConvergenceTestFn EPSConvergedRelative;
821: SLEPC_EXTERN EPSConvergenceTestFn EPSConvergedNorm;
822: SLEPC_EXTERN PetscErrorCode EPSSetConvergenceTestFunction(EPS,EPSConvergenceTestFn*,PetscCtx,PetscCtxDestroyFn*);

824: /*S
825:    EPSStoppingTestFn - A prototype of an `EPS` stopping test function that would
826:    be passed to `EPSSetStoppingTestFunction()`.

828:    Calling Sequence:
829: +  eps    - the linear eigensolver context
830: .  its    - current number of iterations
831: .  max_it - maximum number of iterations
832: .  nconv  - number of currently converged eigenpairs
833: .  nev    - number of requested eigenpairs
834: .  reason - [output] result of the stopping test
835: -  ctx    - optional stopping context, as set by `EPSSetStoppingTestFunction()`

837:    Note:
838:    A positive value of `reason` indicates that the iteration has finished successfully
839:    (converged), and a negative value indicates an error condition (diverged). If
840:    the iteration needs to be continued, `reason` must be set to `EPS_CONVERGED_ITERATING`
841:    (zero).

843:    Level: advanced

845: .seealso: [](ch:eps), `EPSSetStoppingTest()`, `EPSSetStoppingTestFunction()`
846: S*/
847: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode EPSStoppingTestFn(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,PetscCtx ctx);

849: SLEPC_EXTERN PetscErrorCode EPSSetStoppingTest(EPS,EPSStop);
850: SLEPC_EXTERN PetscErrorCode EPSGetStoppingTest(EPS,EPSStop*);
851: SLEPC_EXTERN EPSStoppingTestFn EPSStoppingBasic;
852: SLEPC_EXTERN EPSStoppingTestFn EPSStoppingThreshold;
853: SLEPC_EXTERN PetscErrorCode EPSSetStoppingTestFunction(EPS,EPSStoppingTestFn*,PetscCtx,PetscCtxDestroyFn*);

855: SLEPC_EXTERN PetscErrorCode EPSSetEigenvalueComparison(EPS,SlepcEigenvalueComparisonFn*,PetscCtx);
856: SLEPC_EXTERN PetscErrorCode EPSSetArbitrarySelection(EPS,SlepcArbitrarySelectionFn*,PetscCtx);
857: SLEPC_EXTERN PetscErrorCode EPSSetArbitrarySelectionContextDestroy(EPS,PetscCtxDestroyFn*);

859: /* --------- options specific to particular eigensolvers -------- */

861: /*E
862:    EPSPowerShiftType - The type of shift used in the Power iteration solver.

864:    Values:
865: +  `EPS_POWER_SHIFT_CONSTANT`  - constant shift
866: .  `EPS_POWER_SHIFT_RAYLEIGH`  - variable shift using Rayleigh quotient
867: -  `EPS_POWER_SHIFT_WILKINSON` - variable shift using Wilkinson's approach

869:    Note:
870:    Details of the three variants can be found in {cite:p}`Her05`.

872:    Level: advanced

874: .seealso: [](ch:eps), `EPSPowerSetShiftType()`, `EPSPowerGetShiftType()`
875: E*/
876: typedef enum { EPS_POWER_SHIFT_CONSTANT,
877:                EPS_POWER_SHIFT_RAYLEIGH,
878:                EPS_POWER_SHIFT_WILKINSON } EPSPowerShiftType;
879: SLEPC_EXTERN const char *EPSPowerShiftTypes[];

881: /*MC
882:    EPS_POWER_SHIFT_CONSTANT - The power iteration will use a constant shift.

884:    Note:
885:    Together with `STSINVERT`, the `EPSPOWER` solver implements the inverse iteration
886:    method, i.e., it will apply $(A-\sigma I)^{-1}$ at each iteration, by solving
887:    a linear system. By default, the shift $\sigma$ is constant and given by the
888:    user with `EPSSetTarget()`.

890:    Details of the three variants can be found in {cite:p}`Her05`.

892:    Level: advanced

894: .seealso: [](ch:eps), `EPSPowerShiftType`, `EPSPowerSetShiftType()`, `STSetShift()`, `EPSSetTarget()`, `EPS_POWER_SHIFT_RAYLEIGH`, `EPS_POWER_SHIFT_WILKINSON`
895: M*/

897: /*MC
898:    EPS_POWER_SHIFT_RAYLEIGH - The power iteration will use a variable shift
899:    computed with the Rayleigh quotient.

901:    Notes:
902:    Together with `STSINVERT`, the `EPSPOWER` solver implements the inverse iteration
903:    method, i.e., it will apply $(A-\sigma I)^{-1}$ at each iteration, by solving
904:    a linear system. With this strategy, the value of the shift will be updated at
905:    each iteration as $\sigma=\frac{x^*Ax}{x^*x}$, where $x$ is the current eigenvector
906:    approximation. The resulting iteration is the RQI method.

908:    Updating the shift may involve a high computational cost if the linear solve
909:    is done via a factorization.

911:    Details of the three variants can be found in {cite:p}`Her05`.

913:    Level: advanced

915: .seealso: [](ch:eps), `EPSPowerShiftType`, `EPSPowerSetShiftType()`, `STSetShift()`, `EPSSetTarget()`, `EPS_POWER_SHIFT_CONSTANT`, `EPS_POWER_SHIFT_WILKINSON`
916: M*/

918: /*MC
919:    EPS_POWER_SHIFT_WILKINSON - The power iteration will use a variable shift
920:    computed with Wilkinson's approach.

922:    Note:
923:    Together with `STSINVERT`, the `EPSPOWER` solver implements the inverse iteration
924:    method, i.e., it will apply $(A-\sigma I)^{-1}$ at each iteration, by solving
925:    a linear system. With this strategy, the value of the shift will be updated at
926:    each iteration as proposed by Wilkinson, see {cite:p}`Par80{8.10}`.

928:    Updating the shift may involve a high computational cost if the linear solve
929:    is done via a factorization.

931:    Details of the three variants can be found in {cite:p}`Her05`.

933:    Level: advanced

935: .seealso: [](ch:eps), `EPSPowerShiftType`, `EPSPowerSetShiftType()`, `STSetShift()`, `EPSSetTarget()`, `EPS_POWER_SHIFT_CONSTANT`, `EPS_POWER_SHIFT_RAYLEIGH`
936: M*/

938: SLEPC_EXTERN PetscErrorCode EPSPowerSetShiftType(EPS,EPSPowerShiftType);
939: SLEPC_EXTERN PetscErrorCode EPSPowerGetShiftType(EPS,EPSPowerShiftType*);
940: SLEPC_EXTERN PetscErrorCode EPSPowerSetNonlinear(EPS,PetscBool);
941: SLEPC_EXTERN PetscErrorCode EPSPowerGetNonlinear(EPS,PetscBool*);
942: SLEPC_EXTERN PetscErrorCode EPSPowerSetUpdate(EPS,PetscBool);
943: SLEPC_EXTERN PetscErrorCode EPSPowerGetUpdate(EPS,PetscBool*);
944: SLEPC_EXTERN PetscErrorCode EPSPowerSetSignNormalization(EPS,PetscBool);
945: SLEPC_EXTERN PetscErrorCode EPSPowerGetSignNormalization(EPS,PetscBool*);
946: SLEPC_EXTERN PetscErrorCode EPSPowerSetSNES(EPS,SNES);
947: SLEPC_EXTERN PetscErrorCode EPSPowerGetSNES(EPS,SNES*);

949: SLEPC_EXTERN PetscErrorCode EPSArnoldiSetDelayed(EPS,PetscBool);
950: SLEPC_EXTERN PetscErrorCode EPSArnoldiGetDelayed(EPS,PetscBool*);

952: /*E
953:    EPSKrylovSchurBSEType - The method to be used in the Krylov-Schur solver
954:    for the case of BSE structured eigenproblems.

956:    Values:
957: +  `EPS_KRYLOVSCHUR_BSE_SHAO`         - a Lanczos method proposed by Shao and coauthors
958: .  `EPS_KRYLOVSCHUR_BSE_GRUNING`      - a Lanczos method proposed by Gruning and coauthors
959: -  `EPS_KRYLOVSCHUR_BSE_PROJECTEDBSE` - a Lanczos method resulting is a projected problem with BSE structure

961:    Note:
962:    All variants are implemented in combination with a thick restart, see
963:    the details in {cite:p}`Alv25`.

965:    Level: advanced

967: .seealso: [](ch:eps), `EPSKrylovSchurSetBSEType()`, `EPSKrylovSchurGetBSEType()`
968: E*/
969: typedef enum { EPS_KRYLOVSCHUR_BSE_SHAO,
970:                EPS_KRYLOVSCHUR_BSE_GRUNING,
971:                EPS_KRYLOVSCHUR_BSE_PROJECTEDBSE } EPSKrylovSchurBSEType;
972: SLEPC_EXTERN const char *EPSKrylovSchurBSETypes[];

974: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurSetBSEType(EPS,EPSKrylovSchurBSEType);
975: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetBSEType(EPS,EPSKrylovSchurBSEType*);
976: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurSetRestart(EPS,PetscReal);
977: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetRestart(EPS,PetscReal*);
978: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurSetLocking(EPS,PetscBool);
979: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetLocking(EPS,PetscBool*);
980: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurSetPartitions(EPS,PetscInt);
981: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetPartitions(EPS,PetscInt*);
982: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurSetDetectZeros(EPS,PetscBool);
983: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetDetectZeros(EPS,PetscBool*);
984: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurSetDimensions(EPS,PetscInt,PetscInt,PetscInt);
985: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetDimensions(EPS,PetscInt*,PetscInt*,PetscInt*);
986: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurSetSubintervals(EPS,PetscReal[]);
987: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetSubintervals(EPS,PetscReal*[]);
988: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetInertias(EPS,PetscInt*,PetscReal*[],PetscInt*[]);
989: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetSubcommInfo(EPS,PetscInt*,PetscInt*,Vec*);
990: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetSubcommPairs(EPS,PetscInt,PetscScalar*,Vec);
991: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetSubcommMats(EPS,Mat*,Mat*);
992: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurUpdateSubcommMats(EPS,PetscScalar,PetscScalar,Mat,PetscScalar,PetscScalar, Mat,MatStructure,PetscBool);
993: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetKSP(EPS,KSP*);

995: /*E
996:    EPSLanczosReorthogType - The type of reorthogonalization used in `EPSLANCZOS`.

998:    Values:
999: +  `EPS_LANCZOS_REORTHOG_LOCAL`     - local orthogonalization, only involves previous two vectors
1000: .  `EPS_LANCZOS_REORTHOG_FULL`      - full orthogonalization against all previous vectors
1001: .  `EPS_LANCZOS_REORTHOG_SELECTIVE` - selective reorthogonalization against nearly converged Ritz vectors
1002: .  `EPS_LANCZOS_REORTHOG_PERIODIC`  - periodic reorthogonalization against all previous vectors
1003: .  `EPS_LANCZOS_REORTHOG_PARTIAL`   - partial reorthogonalization against as subset of previous vectors
1004: -  `EPS_LANCZOS_REORTHOG_DELAYED`   - full orthogonalization with delayed reorthogonalization

1006:    Note:
1007:    Details of the different reorthogonalization strategies can be found in
1008:    {cite:p}`Her06`.

1010:    Level: advanced

1012: .seealso: [](ch:eps), `EPSLanczosSetReorthog()`, `EPSLanczosGetReorthog()`
1013: E*/
1014: typedef enum { EPS_LANCZOS_REORTHOG_LOCAL,
1015:                EPS_LANCZOS_REORTHOG_FULL,
1016:                EPS_LANCZOS_REORTHOG_SELECTIVE,
1017:                EPS_LANCZOS_REORTHOG_PERIODIC,
1018:                EPS_LANCZOS_REORTHOG_PARTIAL,
1019:                EPS_LANCZOS_REORTHOG_DELAYED } EPSLanczosReorthogType;
1020: SLEPC_EXTERN const char *EPSLanczosReorthogTypes[];

1022: SLEPC_EXTERN PetscErrorCode EPSLanczosSetReorthog(EPS,EPSLanczosReorthogType);
1023: SLEPC_EXTERN PetscErrorCode EPSLanczosGetReorthog(EPS,EPSLanczosReorthogType*);

1025: /*E
1026:    EPSPRIMMEMethod - The method selected in the PRIMME library.

1028:    Note:
1029:    See the documentation of PRIMME {cite:p}`Sta10` for a description of the methods.

1031:    Level: advanced

1033: .seealso: [](ch:eps), `EPSPRIMMESetMethod()`, `EPSPRIMMEGetMethod()`
1034: E*/
1035: typedef enum { EPS_PRIMME_DYNAMIC             = 1,
1036:                EPS_PRIMME_DEFAULT_MIN_TIME    = 2,
1037:                EPS_PRIMME_DEFAULT_MIN_MATVECS = 3,
1038:                EPS_PRIMME_ARNOLDI             = 4,
1039:                EPS_PRIMME_GD                  = 5,
1040:                EPS_PRIMME_GD_PLUSK            = 6,
1041:                EPS_PRIMME_GD_OLSEN_PLUSK      = 7,
1042:                EPS_PRIMME_JD_OLSEN_PLUSK      = 8,
1043:                EPS_PRIMME_RQI                 = 9,
1044:                EPS_PRIMME_JDQR                = 10,
1045:                EPS_PRIMME_JDQMR               = 11,
1046:                EPS_PRIMME_JDQMR_ETOL          = 12,
1047:                EPS_PRIMME_SUBSPACE_ITERATION  = 13,
1048:                EPS_PRIMME_LOBPCG_ORTHOBASIS   = 14,
1049:                EPS_PRIMME_LOBPCG_ORTHOBASISW  = 15 } EPSPRIMMEMethod;
1050: SLEPC_EXTERN const char *EPSPRIMMEMethods[];

1052: SLEPC_EXTERN PetscErrorCode EPSPRIMMESetBlockSize(EPS,PetscInt);
1053: SLEPC_EXTERN PetscErrorCode EPSPRIMMEGetBlockSize(EPS,PetscInt*);
1054: SLEPC_EXTERN PetscErrorCode EPSPRIMMESetMethod(EPS,EPSPRIMMEMethod);
1055: SLEPC_EXTERN PetscErrorCode EPSPRIMMEGetMethod(EPS,EPSPRIMMEMethod*);

1057: SLEPC_EXTERN PetscErrorCode EPSGDSetKrylovStart(EPS,PetscBool);
1058: SLEPC_EXTERN PetscErrorCode EPSGDGetKrylovStart(EPS,PetscBool*);
1059: SLEPC_EXTERN PetscErrorCode EPSGDSetBlockSize(EPS,PetscInt);
1060: SLEPC_EXTERN PetscErrorCode EPSGDGetBlockSize(EPS,PetscInt*);
1061: SLEPC_EXTERN PetscErrorCode EPSGDSetRestart(EPS,PetscInt,PetscInt);
1062: SLEPC_EXTERN PetscErrorCode EPSGDGetRestart(EPS,PetscInt*,PetscInt*);
1063: SLEPC_EXTERN PetscErrorCode EPSGDSetInitialSize(EPS,PetscInt);
1064: SLEPC_EXTERN PetscErrorCode EPSGDGetInitialSize(EPS,PetscInt*);
1065: SLEPC_EXTERN PetscErrorCode EPSGDSetBOrth(EPS,PetscBool);
1066: SLEPC_EXTERN PetscErrorCode EPSGDGetBOrth(EPS,PetscBool*);
1067: SLEPC_EXTERN PetscErrorCode EPSGDSetDoubleExpansion(EPS,PetscBool);
1068: SLEPC_EXTERN PetscErrorCode EPSGDGetDoubleExpansion(EPS,PetscBool*);

1070: SLEPC_EXTERN PetscErrorCode EPSJDSetKrylovStart(EPS,PetscBool);
1071: SLEPC_EXTERN PetscErrorCode EPSJDGetKrylovStart(EPS,PetscBool*);
1072: SLEPC_EXTERN PetscErrorCode EPSJDSetBlockSize(EPS,PetscInt);
1073: SLEPC_EXTERN PetscErrorCode EPSJDGetBlockSize(EPS,PetscInt*);
1074: SLEPC_EXTERN PetscErrorCode EPSJDSetRestart(EPS,PetscInt,PetscInt);
1075: SLEPC_EXTERN PetscErrorCode EPSJDGetRestart(EPS,PetscInt*,PetscInt*);
1076: SLEPC_EXTERN PetscErrorCode EPSJDSetInitialSize(EPS,PetscInt);
1077: SLEPC_EXTERN PetscErrorCode EPSJDGetInitialSize(EPS,PetscInt*);
1078: SLEPC_EXTERN PetscErrorCode EPSJDSetFix(EPS,PetscReal);
1079: SLEPC_EXTERN PetscErrorCode EPSJDGetFix(EPS,PetscReal*);
1080: SLEPC_EXTERN PetscErrorCode EPSJDSetConstCorrectionTol(EPS,PetscBool);
1081: SLEPC_EXTERN PetscErrorCode EPSJDGetConstCorrectionTol(EPS,PetscBool*);
1082: SLEPC_EXTERN PetscErrorCode EPSJDSetBOrth(EPS,PetscBool);
1083: SLEPC_EXTERN PetscErrorCode EPSJDGetBOrth(EPS,PetscBool*);

1085: SLEPC_EXTERN PetscErrorCode EPSRQCGSetReset(EPS,PetscInt);
1086: SLEPC_EXTERN PetscErrorCode EPSRQCGGetReset(EPS,PetscInt*);

1088: SLEPC_EXTERN PetscErrorCode EPSLOBPCGSetBlockSize(EPS,PetscInt);
1089: SLEPC_EXTERN PetscErrorCode EPSLOBPCGGetBlockSize(EPS,PetscInt*);
1090: SLEPC_EXTERN PetscErrorCode EPSLOBPCGSetRestart(EPS,PetscReal);
1091: SLEPC_EXTERN PetscErrorCode EPSLOBPCGGetRestart(EPS,PetscReal*);
1092: SLEPC_EXTERN PetscErrorCode EPSLOBPCGSetLocking(EPS,PetscBool);
1093: SLEPC_EXTERN PetscErrorCode EPSLOBPCGGetLocking(EPS,PetscBool*);

1095: /*E
1096:    EPSCISSQuadRule - The quadrature rule used in the `EPSCISS` solver.

1098:    Values:
1099: +  `EPS_CISS_QUADRULE_TRAPEZOIDAL` - trapezoidal rule
1100: -  `EPS_CISS_QUADRULE_CHEBYSHEV`   - Gauss quadrature on Chebyshev points

1102:    Note:
1103:    For a detailed description see {cite:p}`Mae16`.

1105:    Level: advanced

1107: .seealso: [](ch:eps), `EPSCISSSetQuadRule()`, `EPSCISSGetQuadRule()`
1108: E*/
1109: typedef enum { EPS_CISS_QUADRULE_TRAPEZOIDAL = 1,
1110:                EPS_CISS_QUADRULE_CHEBYSHEV   = 2 } EPSCISSQuadRule;
1111: SLEPC_EXTERN const char *EPSCISSQuadRules[];

1113: /*E
1114:    EPSCISSExtraction - The extraction technique used in the `EPSCISS` solver.

1116:    Values:
1117: +  `EPS_CISS_EXTRACTION_RITZ`   - Ritz approximations from Rayleigh-Ritz projection
1118: -  `EPS_CISS_EXTRACTION_HANKEL` - use Hankel pencil as in the original Sakurai-Sugiura method

1120:    Note:
1121:    For a detailed description see {cite:p}`Mae16`.

1123:    Level: advanced

1125: .seealso: [](ch:eps), `EPSCISSSetExtraction()`, `EPSCISSGetExtraction()`
1126: E*/
1127: typedef enum { EPS_CISS_EXTRACTION_RITZ,
1128:                EPS_CISS_EXTRACTION_HANKEL } EPSCISSExtraction;
1129: SLEPC_EXTERN const char *EPSCISSExtractions[];

1131: SLEPC_EXTERN PetscErrorCode EPSCISSSetExtraction(EPS,EPSCISSExtraction);
1132: SLEPC_EXTERN PetscErrorCode EPSCISSGetExtraction(EPS,EPSCISSExtraction*);
1133: SLEPC_EXTERN PetscErrorCode EPSCISSSetQuadRule(EPS,EPSCISSQuadRule);
1134: SLEPC_EXTERN PetscErrorCode EPSCISSGetQuadRule(EPS,EPSCISSQuadRule*);
1135: SLEPC_EXTERN PetscErrorCode EPSCISSSetSizes(EPS,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool);
1136: SLEPC_EXTERN PetscErrorCode EPSCISSGetSizes(EPS,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*);
1137: SLEPC_EXTERN PetscErrorCode EPSCISSSetThreshold(EPS,PetscReal,PetscReal);
1138: SLEPC_EXTERN PetscErrorCode EPSCISSGetThreshold(EPS,PetscReal*,PetscReal*);
1139: SLEPC_EXTERN PetscErrorCode EPSCISSSetRefinement(EPS,PetscInt,PetscInt);
1140: SLEPC_EXTERN PetscErrorCode EPSCISSGetRefinement(EPS,PetscInt*,PetscInt*);
1141: SLEPC_EXTERN PetscErrorCode EPSCISSSetUseST(EPS,PetscBool);
1142: SLEPC_EXTERN PetscErrorCode EPSCISSGetUseST(EPS,PetscBool*);
1143: SLEPC_EXTERN PetscErrorCode EPSCISSGetKSPs(EPS,PetscInt*,KSP*[]);

1145: SLEPC_EXTERN PetscErrorCode EPSLyapIISetLME(EPS,LME);
1146: SLEPC_EXTERN PetscErrorCode EPSLyapIIGetLME(EPS,LME*);
1147: SLEPC_EXTERN PetscErrorCode EPSLyapIISetRanks(EPS,PetscInt,PetscInt);
1148: SLEPC_EXTERN PetscErrorCode EPSLyapIIGetRanks(EPS,PetscInt*,PetscInt*);

1150: SLEPC_EXTERN PetscErrorCode EPSBLOPEXSetBlockSize(EPS,PetscInt);
1151: SLEPC_EXTERN PetscErrorCode EPSBLOPEXGetBlockSize(EPS,PetscInt*);

1153: /*E
1154:    EPSEVSLDOSMethod - The method to approximate the density of states (DOS)
1155:    in the `EPSEVSL` solver.

1157:    Values:
1158: +  `EPS_EVSL_DOS_KPM`     - Kernel Polynomial Method
1159: -  `EPS_EVSL_DOS_LANCZOS` - Lanczos method

1161:    Note:
1162:    See the documentation of EVSL {cite:p}`Li19` for an explanation.

1164:    Level: advanced

1166: .seealso: [](ch:eps), `EPSEVSLSetDOSParameters()`, `EPSEVSLGetDOSParameters()`
1167: E*/
1168: typedef enum { EPS_EVSL_DOS_KPM,
1169:                EPS_EVSL_DOS_LANCZOS } EPSEVSLDOSMethod;
1170: SLEPC_EXTERN const char *EPSEVSLDOSMethods[];

1172: /*E
1173:    EPSEVSLDamping - The damping type used in the `EPSEVSL` solver.

1175:    Values:
1176: +  `EPS_EVSL_DAMPING_NONE`    - no damping
1177: .  `EPS_EVSL_DAMPING_JACKSON` - Jackson damping
1178: -  `EPS_EVSL_DAMPING_SIGMA`   - Lanczos damping

1180:    Note:
1181:    See the documentation of EVSL {cite:p}`Li19` for an explanation.

1183:    Level: advanced

1185: .seealso: [](ch:eps), `EPSEVSLSetDOSParameters()`, `EPSEVSLGetDOSParameters()`
1186: E*/
1187: typedef enum { EPS_EVSL_DAMPING_NONE,
1188:                EPS_EVSL_DAMPING_JACKSON,
1189:                EPS_EVSL_DAMPING_SIGMA } EPSEVSLDamping;
1190: SLEPC_EXTERN const char *EPSEVSLDampings[];

1192: SLEPC_EXTERN PetscErrorCode EPSEVSLSetRange(EPS,PetscReal,PetscReal);
1193: SLEPC_EXTERN PetscErrorCode EPSEVSLGetRange(EPS,PetscReal*,PetscReal*);
1194: SLEPC_EXTERN PetscErrorCode EPSEVSLSetSlices(EPS,PetscInt);
1195: SLEPC_EXTERN PetscErrorCode EPSEVSLGetSlices(EPS,PetscInt*);
1196: SLEPC_EXTERN PetscErrorCode EPSEVSLSetDOSParameters(EPS,EPSEVSLDOSMethod,PetscInt,PetscInt,PetscInt,PetscInt);
1197: SLEPC_EXTERN PetscErrorCode EPSEVSLGetDOSParameters(EPS,EPSEVSLDOSMethod*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
1198: SLEPC_EXTERN PetscErrorCode EPSEVSLSetPolParameters(EPS,PetscInt,PetscReal);
1199: SLEPC_EXTERN PetscErrorCode EPSEVSLGetPolParameters(EPS,PetscInt*,PetscReal*);
1200: SLEPC_EXTERN PetscErrorCode EPSEVSLSetDamping(EPS,EPSEVSLDamping);
1201: SLEPC_EXTERN PetscErrorCode EPSEVSLGetDamping(EPS,EPSEVSLDamping*);

1203: SLEPC_EXTERN PetscErrorCode EPSFEASTSetNumPoints(EPS,PetscInt);
1204: SLEPC_EXTERN PetscErrorCode EPSFEASTGetNumPoints(EPS,PetscInt*);

1206: SLEPC_EXTERN PetscErrorCode EPSCHASESetDegree(EPS,PetscInt,PetscBool);
1207: SLEPC_EXTERN PetscErrorCode EPSCHASEGetDegree(EPS,PetscInt*,PetscBool*);