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

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

 86:    Level: intermediate

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

 99: /*MC
100:    EPS_HEP - A Hermitian eigenvalue problem.

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

106:    Level: intermediate

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

111: /*MC
112:    EPS_GHEP - A generalized Hermitian eigenvalue problem.

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

118:    Level: intermediate

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

123: /*MC
124:    EPS_NHEP - A non-Hermitian eigenvalue problem.

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

130:    Level: intermediate

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

135: /*MC
136:    EPS_GNHEP - A generalized non-Hermitian eigenvalue problem.

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

142:    Level: intermediate

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

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

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

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

159:    Level: intermediate

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

164: /*MC
165:    EPS_GHIEP - A generalized Hermitian-indefinite eigenvalue problem.

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

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

175:    Level: intermediate

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

180: /*MC
181:    EPS_BSE - A structured Bethe-Salpeter eigenvalue problem.

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

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

196:    Level: intermediate

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

201: /*MC
202:    EPS_HAMILT - A structured Hamiltonian eigenvalue problem.

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

214:    Level: intermediate

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

219: /*E
220:    EPSExtraction - Determines the type of extraction technique employed
221:    by the eigensolver.

223:    Values:
224: +  `EPS_RITZ`              - Rayleigh-Ritz extraction
225: .  `EPS_HARMONIC`          - harmonic Ritz extraction
226: .  `EPS_HARMONIC_RELATIVE` - harmonic Ritz extraction relative to the eigenvalue
227: .  `EPS_HARMONIC_RIGHT`    - harmonic Ritz extraction for rightmost eigenvalues
228: .  `EPS_HARMONIC_LARGEST`  - harmonic Ritz extraction for largest magnitude (without target)
229: .  `EPS_REFINED`           - refined Ritz extraction
230: -  `EPS_REFINED_HARMONIC`  - refined harmonic Ritz extraction

232:    Level: advanced

234: .seealso: [](ch:eps), `EPSSetExtraction()`, `EPSGetExtraction()`
235: E*/
236: typedef enum { EPS_RITZ,
237:                EPS_HARMONIC,
238:                EPS_HARMONIC_RELATIVE,
239:                EPS_HARMONIC_RIGHT,
240:                EPS_HARMONIC_LARGEST,
241:                EPS_REFINED,
242:                EPS_REFINED_HARMONIC } EPSExtraction;

244: /*MC
245:    EPS_RITZ - The standard Rayleigh-Ritz extraction.

247:    Note:
248:    This is the default way of computing eigenpair approximations from a
249:    given subspace.

251:    Level: advanced

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

256: /*MC
257:    EPS_HARMONIC - The harmonic Ritz extraction.

259:    Notes:
260:    This extraction method may provide better convergence when computing
261:    interior eigenvalues close to a given target.

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

266:    Level: advanced

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

271: /*MC
272:    EPS_HARMONIC_RELATIVE - The harmonic Ritz extraction relative to the eigenvalue.

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

277:    Level: advanced

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

282: /*MC
283:    EPS_HARMONIC_RIGHT - The harmonic Ritz extraction for rightmost eigenvalues.

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

288:    Level: advanced

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

293: /*MC
294:    EPS_HARMONIC_LARGEST - The harmonic Ritz extraction for largest magnitude
295:    eigenvalues (without target).

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

300:    Level: advanced

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

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

308:    Note:
309:    Currently implemented only in `EPSARNOLDI`.

311:    Level: advanced

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

316: /*MC
317:    EPS_REFINED_HARMONIC - The refined harmonic Ritz extraction.

319:    Note:
320:    This is a combination of `EPS_HARMONIC` and `EPS_REFINED`.

322:    Developer Note:
323:    Currently not implemented, reserved for future use.

325:    Level: advanced

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

330: /*E
331:    EPSWhich - Determines which part of the spectrum is requested.

333:    Values:
334: +  `EPS_LARGEST_MAGNITUDE`  - largest $|\lambda|$
335: .  `EPS_SMALLEST_MAGNITUDE` - smallest $|\lambda|$
336: .  `EPS_LARGEST_REAL`       - largest $\mathrm{Re}(\lambda)$
337: .  `EPS_SMALLEST_REAL`      - smallest $\mathrm{Re}(\lambda)$
338: .  `EPS_LARGEST_IMAGINARY`  - largest $\mathrm{Im}(\lambda)$
339: .  `EPS_SMALLEST_IMAGINARY` - smallest $\mathrm{Im}(\lambda)$
340: .  `EPS_TARGET_MAGNITUDE`   - smallest $|\lambda-\tau|$
341: .  `EPS_TARGET_REAL`        - smallest $|\mathrm{Re}(\lambda-\tau)|$
342: .  `EPS_TARGET_IMAGINARY`   - smallest $|\mathrm{Im}(\lambda-\tau)|$
343: .  `EPS_ALL`                - all $\lambda\in[a,b]$ or $\lambda\in\Omega$
344: -  `EPS_WHICH_USER`         - user-defined sorting criterion

346:    Notes:
347:    If SLEPc is compiled for real scalars `EPS_LARGEST_IMAGINARY` and
348:    `EPS_SMALLEST_IMAGINARY` use the absolute value of the imaginary part
349:    for eigenvalue selection.

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

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

356:    Level: intermediate

358: .seealso: [](ch:eps), `EPSSetWhichEigenpairs()`, `EPSSetTarget()`, `EPSSetInterval()`
359: E*/
360: typedef enum { EPS_LARGEST_MAGNITUDE  = 1,
361:                EPS_SMALLEST_MAGNITUDE = 2,
362:                EPS_LARGEST_REAL       = 3,
363:                EPS_SMALLEST_REAL      = 4,
364:                EPS_LARGEST_IMAGINARY  = 5,
365:                EPS_SMALLEST_IMAGINARY = 6,
366:                EPS_TARGET_MAGNITUDE   = 7,
367:                EPS_TARGET_REAL        = 8,
368:                EPS_TARGET_IMAGINARY   = 9,
369:                EPS_ALL                = 10,
370:                EPS_WHICH_USER         = 11 } EPSWhich;

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

375:    Values:
376: +  `EPS_BALANCE_NONE`    - no balancing matrix is used
377: .  `EPS_BALANCE_ONESIDE` - balancing matrix $D$ is computed with a one-sided Krylov method
378: .  `EPS_BALANCE_TWOSIDE` - balancing matrix $D$ is computed with a two-sided Krylov method
379: -  `EPS_BALANCE_USER`    - use a balancing matrix $D$ provided by the user

381:    Level: intermediate

383: .seealso: [](ch:eps), [](sec:balancing), `EPSSetBalance()`
384: E*/
385: typedef enum { EPS_BALANCE_NONE,
386:                EPS_BALANCE_ONESIDE,
387:                EPS_BALANCE_TWOSIDE,
388:                EPS_BALANCE_USER } EPSBalance;
389: SLEPC_EXTERN const char *EPSBalanceTypes[];

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

394:    Values:
395: +  `EPS_ERROR_ABSOLUTE` - compute error bound as $\|r\|$
396: .  `EPS_ERROR_RELATIVE` - compute error bound as $\|r\|/|\lambda|$
397: -  `EPS_ERROR_BACKWARD` - compute error bound as $\|r\|/(\|A\|+|\lambda|\|B\|)$

399:    Level: intermediate

401: .seealso: [](ch:eps), `EPSComputeError()`
402: E*/
403: typedef enum { EPS_ERROR_ABSOLUTE,
404:                EPS_ERROR_RELATIVE,
405:                EPS_ERROR_BACKWARD } EPSErrorType;
406: SLEPC_EXTERN const char *EPSErrorTypes[];

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

411:    Values:
412: +  `EPS_CONV_ABS`  - absolute convergence criterion, $\|r\|$
413: .  `EPS_CONV_REL`  - convergence criterion relative to eigenvalue, $\|r\|/|\lambda|$
414: .  `EPS_CONV_NORM` - convergence criterion relative to matrix norms, $\|r\|/(\|A\|+|\lambda|\|B\|)$
415: -  `EPS_CONV_USER` - convergence dictated by user-provided function

417:    Level: intermediate

419: .seealso: [](ch:eps), `EPSSetConvergenceTest()`, `EPSSetConvergenceTestFunction()`
420: E*/
421: typedef enum { EPS_CONV_ABS,
422:                EPS_CONV_REL,
423:                EPS_CONV_NORM,
424:                EPS_CONV_USER } EPSConv;

426: /*E
427:    EPSStop - The stopping test to decide the termination of the outer loop
428:    of the eigensolver.

430:    Values:
431: +  `EPS_STOP_BASIC`     - default stopping test
432: .  `EPS_STOP_USER`      - user-provided stopping test
433: -  `EPS_STOP_THRESHOLD` - threshold stopping test

435:    Level: advanced

437: .seealso: [](ch:eps), `EPSSetStoppingTest()`, `EPSSetStoppingTestFunction()`
438: E*/
439: typedef enum { EPS_STOP_BASIC,
440:                EPS_STOP_USER,
441:                EPS_STOP_THRESHOLD } EPSStop;

443: /*MC
444:    EPS_STOP_BASIC - The default stopping test.

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

451:    Level: advanced

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

456: /*MC
457:    EPS_STOP_USER - The user-provided stopping test.

459:    Note:
460:    Customized stopping test using the user-provided function given with
461:    `EPSSetStoppingTestFunction()`.

463:    Level: advanced

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

468: /*MC
469:    EPS_STOP_THRESHOLD - The threshold stopping test.

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

479:    Level: advanced

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

484: /*E
485:    EPSConvergedReason - Reason an eigensolver was determined to have converged
486:    or diverged.

488:    Values:
489: +  `EPS_CONVERGED_TOL`          - converged up to tolerance
490: .  `EPS_CONVERGED_USER`         - converged due to a user-defined condition
491: .  `EPS_DIVERGED_ITS`           - exceeded the maximum number of allowed iterations
492: .  `EPS_DIVERGED_BREAKDOWN`     - generic breakdown in method
493: .  `EPS_DIVERGED_SYMMETRY_LOST` - pseudo-Lanczos was not able to keep symmetry
494: -  `EPS_CONVERGED_ITERATING`    - the solver is still running

496:    Level: intermediate

498: .seealso: [](ch:eps), `EPSSolve()`, `EPSGetConvergedReason()`, `EPSSetTolerances()`
499: E*/
500: typedef enum {/* converged */
501:               EPS_CONVERGED_TOL                =  1,
502:               EPS_CONVERGED_USER               =  2,
503:               /* diverged */
504:               EPS_DIVERGED_ITS                 = -1,
505:               EPS_DIVERGED_BREAKDOWN           = -2,
506:               EPS_DIVERGED_SYMMETRY_LOST       = -3,
507:               EPS_CONVERGED_ITERATING          =  0} EPSConvergedReason;
508: SLEPC_EXTERN const char *const*EPSConvergedReasons;

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

514:    Level: intermediate

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

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

522:    Note:
523:    This happens only when a user-defined stopping test has been set with
524:    `EPSSetStoppingTestFunction()`.

526:    Level: intermediate

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

531: /*MC
532:    EPS_DIVERGED_ITS - Exceeded the maximum number of allowed iterations
533:    before the convergence criterion was satisfied.

535:    Level: intermediate

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

540: /*MC
541:    EPS_DIVERGED_BREAKDOWN - A breakdown in the solver was detected so the
542:    method could not continue.

544:    Level: intermediate

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

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

554:    Level: intermediate

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

559: /*MC
560:    EPS_CONVERGED_ITERATING - This value is returned if `EPSGetConvergedReason()` is called
561:    while `EPSSolve()` is still running.

563:    Level: intermediate

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

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

572:    Level: advanced

574: .seealso: [](ch:eps), `EPSSetStoppingTestFunction()`
575: S*/
576: struct _n_EPSStoppingCtx {
577:   PetscReal firstev;    /* the (absolute) value of the first converged eigenvalue */
578:   PetscReal lastev;     /* the (absolute) value of the last converged eigenvalue */
579:   PetscReal thres;      /* threshold set with EPSSetThreshold() */
580:   PetscBool threlative; /* threshold is relative */
581:   EPSWhich  which;      /* which eigenvalues are being computed */
582: };
583: typedef struct _n_EPSStoppingCtx* EPSStoppingCtx;

585: SLEPC_EXTERN PetscErrorCode EPSCreate(MPI_Comm,EPS*);
586: SLEPC_EXTERN PetscErrorCode EPSDestroy(EPS*);
587: SLEPC_EXTERN PetscErrorCode EPSReset(EPS);
588: SLEPC_EXTERN PetscErrorCode EPSSetType(EPS,EPSType);
589: SLEPC_EXTERN PetscErrorCode EPSGetType(EPS,EPSType*);
590: SLEPC_EXTERN PetscErrorCode EPSSetProblemType(EPS,EPSProblemType);
591: SLEPC_EXTERN PetscErrorCode EPSGetProblemType(EPS,EPSProblemType*);
592: SLEPC_EXTERN PetscErrorCode EPSSetExtraction(EPS,EPSExtraction);
593: SLEPC_EXTERN PetscErrorCode EPSGetExtraction(EPS,EPSExtraction*);
594: SLEPC_EXTERN PetscErrorCode EPSSetBalance(EPS,EPSBalance,PetscInt,PetscReal);
595: SLEPC_EXTERN PetscErrorCode EPSGetBalance(EPS,EPSBalance*,PetscInt*,PetscReal*);
596: SLEPC_EXTERN PetscErrorCode EPSSetOperators(EPS,Mat,Mat);
597: SLEPC_EXTERN PetscErrorCode EPSGetOperators(EPS,Mat*,Mat*);
598: SLEPC_EXTERN PetscErrorCode EPSSetFromOptions(EPS);
599: SLEPC_EXTERN PetscErrorCode EPSSetDSType(EPS);
600: SLEPC_EXTERN PetscErrorCode EPSSetUp(EPS);
601: SLEPC_EXTERN PetscErrorCode EPSSolve(EPS);
602: SLEPC_EXTERN PetscErrorCode EPSView(EPS,PetscViewer);
603: SLEPC_EXTERN PetscErrorCode EPSViewFromOptions(EPS,PetscObject,const char[]);
604: SLEPC_EXTERN PetscErrorCode EPSErrorView(EPS,EPSErrorType,PetscViewer);
605: PETSC_DEPRECATED_FUNCTION(3, 6, 0, "EPSErrorView()", ) static inline PetscErrorCode EPSPrintSolution(EPS eps,PetscViewer v) {return EPSErrorView(eps,EPS_ERROR_RELATIVE,v);}
606: SLEPC_EXTERN PetscErrorCode EPSErrorViewFromOptions(EPS);
607: SLEPC_EXTERN PetscErrorCode EPSConvergedReasonView(EPS,PetscViewer);
608: SLEPC_EXTERN PetscErrorCode EPSConvergedReasonViewFromOptions(EPS);
609: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "EPSConvergedReasonView()", ) static inline PetscErrorCode EPSReasonView(EPS eps,PetscViewer v) {return EPSConvergedReasonView(eps,v);}
610: PETSC_DEPRECATED_FUNCTION(3, 14, 0, "EPSConvergedReasonViewFromOptions()", ) static inline PetscErrorCode EPSReasonViewFromOptions(EPS eps) {return EPSConvergedReasonViewFromOptions(eps);}
611: SLEPC_EXTERN PetscErrorCode EPSValuesView(EPS,PetscViewer);
612: SLEPC_EXTERN PetscErrorCode EPSValuesViewFromOptions(EPS);
613: SLEPC_EXTERN PetscErrorCode EPSVectorsView(EPS,PetscViewer);
614: SLEPC_EXTERN PetscErrorCode EPSVectorsViewFromOptions(EPS);

616: SLEPC_EXTERN PetscErrorCode EPSSetTarget(EPS,PetscScalar);
617: SLEPC_EXTERN PetscErrorCode EPSGetTarget(EPS,PetscScalar*);
618: SLEPC_EXTERN PetscErrorCode EPSSetInterval(EPS,PetscReal,PetscReal);
619: SLEPC_EXTERN PetscErrorCode EPSGetInterval(EPS,PetscReal*,PetscReal*);
620: SLEPC_EXTERN PetscErrorCode EPSSetST(EPS,ST);
621: SLEPC_EXTERN PetscErrorCode EPSGetST(EPS,ST*);
622: SLEPC_EXTERN PetscErrorCode EPSSetBV(EPS,BV);
623: SLEPC_EXTERN PetscErrorCode EPSGetBV(EPS,BV*);
624: SLEPC_EXTERN PetscErrorCode EPSSetRG(EPS,RG);
625: SLEPC_EXTERN PetscErrorCode EPSGetRG(EPS,RG*);
626: SLEPC_EXTERN PetscErrorCode EPSSetDS(EPS,DS);
627: SLEPC_EXTERN PetscErrorCode EPSGetDS(EPS,DS*);
628: SLEPC_EXTERN PetscErrorCode EPSSetTolerances(EPS,PetscReal,PetscInt);
629: SLEPC_EXTERN PetscErrorCode EPSGetTolerances(EPS,PetscReal*,PetscInt*);
630: SLEPC_EXTERN PetscErrorCode EPSSetDimensions(EPS,PetscInt,PetscInt,PetscInt);
631: SLEPC_EXTERN PetscErrorCode EPSGetDimensions(EPS,PetscInt*,PetscInt*,PetscInt*);

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

635: SLEPC_EXTERN PetscErrorCode EPSGetConverged(EPS,PetscInt*);
636: SLEPC_EXTERN PetscErrorCode EPSGetEigenpair(EPS,PetscInt,PetscScalar*,PetscScalar*,Vec,Vec);
637: SLEPC_EXTERN PetscErrorCode EPSGetEigenvalue(EPS,PetscInt,PetscScalar*,PetscScalar*);
638: SLEPC_EXTERN PetscErrorCode EPSGetEigenvector(EPS,PetscInt,Vec,Vec);
639: SLEPC_EXTERN PetscErrorCode EPSGetLeftEigenvector(EPS,PetscInt,Vec,Vec);

641: SLEPC_EXTERN PetscErrorCode EPSComputeError(EPS,PetscInt,EPSErrorType,PetscReal*);
642: 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);}
643: 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);}
644: SLEPC_EXTERN PetscErrorCode EPSGetInvariantSubspace(EPS,Vec[]);
645: SLEPC_EXTERN PetscErrorCode EPSGetErrorEstimate(EPS,PetscInt,PetscReal*);
646: SLEPC_EXTERN PetscErrorCode EPSGetIterationNumber(EPS,PetscInt*);

648: SLEPC_EXTERN PetscErrorCode EPSSetWhichEigenpairs(EPS,EPSWhich);
649: SLEPC_EXTERN PetscErrorCode EPSGetWhichEigenpairs(EPS,EPSWhich*);
650: SLEPC_EXTERN PetscErrorCode EPSSetThreshold(EPS,PetscReal,PetscBool);
651: SLEPC_EXTERN PetscErrorCode EPSGetThreshold(EPS,PetscReal*,PetscBool*);
652: SLEPC_EXTERN PetscErrorCode EPSSetTwoSided(EPS,PetscBool);
653: SLEPC_EXTERN PetscErrorCode EPSGetTwoSided(EPS,PetscBool*);
654: SLEPC_EXTERN PetscErrorCode EPSSetTrueResidual(EPS,PetscBool);
655: SLEPC_EXTERN PetscErrorCode EPSGetTrueResidual(EPS,PetscBool*);
656: SLEPC_EXTERN PetscErrorCode EPSSetPurify(EPS,PetscBool);
657: SLEPC_EXTERN PetscErrorCode EPSGetPurify(EPS,PetscBool*);
658: SLEPC_EXTERN PetscErrorCode EPSIsGeneralized(EPS,PetscBool*);
659: SLEPC_EXTERN PetscErrorCode EPSIsHermitian(EPS,PetscBool*);
660: SLEPC_EXTERN PetscErrorCode EPSIsPositive(EPS,PetscBool*);
661: SLEPC_EXTERN PetscErrorCode EPSIsStructured(EPS,PetscBool*);

663: SLEPC_EXTERN PetscErrorCode EPSSetTrackAll(EPS,PetscBool);
664: SLEPC_EXTERN PetscErrorCode EPSGetTrackAll(EPS,PetscBool*);

666: SLEPC_EXTERN PetscErrorCode EPSSetDeflationSpace(EPS,PetscInt,Vec[]);
667: SLEPC_EXTERN PetscErrorCode EPSSetInitialSpace(EPS,PetscInt,Vec[]);
668: SLEPC_EXTERN PetscErrorCode EPSSetLeftInitialSpace(EPS,PetscInt,Vec[]);

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

673:    Calling Sequence:
674: +  eps    - the linear eigensolver context
675: .  its    - iteration number
676: .  nconv  - number of converged eigenpairs
677: .  eigr   - real part of the eigenvalues
678: .  eigi   - imaginary part of the eigenvalues
679: .  errest - relative error estimates for each eigenpair
680: .  nest   - number of error estimates
681: -  ctx    - optional monitoring context, as provided with `EPSMonitorSet()`

683:    Level: intermediate

685: .seealso: [](ch:eps), `EPSMonitorSet()`
686: S*/
687: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode EPSMonitorFn(EPS eps,PetscInt its,PetscInt nconv,PetscScalar eigr[],PetscScalar eigi[],PetscReal errest[],PetscInt nest,void *ctx);

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

692:    Calling Sequence:
693: +  eps    - the linear eigensolver context
694: .  its    - iteration number
695: .  nconv  - number of converged eigenpairs
696: .  eigr   - real part of the eigenvalues
697: .  eigi   - imaginary part of the eigenvalues
698: .  errest - relative error estimates for each eigenpair
699: .  nest   - number of error estimates
700: -  ctx    - `PetscViewerAndFormat` object

702:    Level: advanced

704:    Note:
705:    This is an `EPSMonitorFn` specialized for a context of `PetscViewerAndFormat`.

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

711: /*S
712:    EPSMonitorRegisterCreateFn - A function prototype for functions that do the
713:    creation when provided to `EPSMonitorRegister()`.

715:    Calling Sequence:
716: +  viewer - the viewer to be used with the `EPSMonitorRegisterFn`
717: .  format - the format of the viewer
718: .  ctx    - a context for the monitor
719: -  result - a `PetscViewerAndFormat` object

721:    Level: advanced

723: .seealso: [](ch:eps), `EPSMonitorRegisterFn`, `EPSMonitorSet()`, `EPSMonitorRegister()`, `EPSMonitorFn`, `EPSMonitorRegisterDestroyFn`
724: S*/
725: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode EPSMonitorRegisterCreateFn(PetscViewer viewer,PetscViewerFormat format,void *ctx,PetscViewerAndFormat **result);

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

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

734:    Level: advanced

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

740: SLEPC_EXTERN PetscErrorCode EPSMonitor(EPS,PetscInt,PetscInt,PetscScalar[],PetscScalar[],PetscReal[],PetscInt);
741: SLEPC_EXTERN PetscErrorCode EPSMonitorSet(EPS,EPSMonitorFn,void*,PetscCtxDestroyFn*);
742: SLEPC_EXTERN PetscErrorCode EPSMonitorCancel(EPS);
743: SLEPC_EXTERN PetscErrorCode EPSGetMonitorContext(EPS,void*);

745: SLEPC_EXTERN PetscErrorCode EPSMonitorSetFromOptions(EPS,const char[],const char[],void*,PetscBool);
746: SLEPC_EXTERN EPSMonitorRegisterFn        EPSMonitorFirst;
747: SLEPC_EXTERN EPSMonitorRegisterFn        EPSMonitorFirstDrawLG;
748: SLEPC_EXTERN EPSMonitorRegisterCreateFn  EPSMonitorFirstDrawLGCreate;
749: SLEPC_EXTERN EPSMonitorRegisterFn        EPSMonitorAll;
750: SLEPC_EXTERN EPSMonitorRegisterFn        EPSMonitorAllDrawLG;
751: SLEPC_EXTERN EPSMonitorRegisterCreateFn  EPSMonitorAllDrawLGCreate;
752: SLEPC_EXTERN EPSMonitorRegisterFn        EPSMonitorConverged;
753: SLEPC_EXTERN EPSMonitorRegisterCreateFn  EPSMonitorConvergedCreate;
754: SLEPC_EXTERN EPSMonitorRegisterFn        EPSMonitorConvergedDrawLG;
755: SLEPC_EXTERN EPSMonitorRegisterCreateFn  EPSMonitorConvergedDrawLGCreate;
756: SLEPC_EXTERN EPSMonitorRegisterDestroyFn EPSMonitorConvergedDestroy;

758: SLEPC_EXTERN PetscErrorCode EPSSetOptionsPrefix(EPS,const char[]);
759: SLEPC_EXTERN PetscErrorCode EPSAppendOptionsPrefix(EPS,const char[]);
760: SLEPC_EXTERN PetscErrorCode EPSGetOptionsPrefix(EPS,const char*[]);

762: SLEPC_EXTERN PetscFunctionList EPSList;
763: SLEPC_EXTERN PetscFunctionList EPSMonitorList;
764: SLEPC_EXTERN PetscFunctionList EPSMonitorCreateList;
765: SLEPC_EXTERN PetscFunctionList EPSMonitorDestroyList;
766: SLEPC_EXTERN PetscErrorCode EPSRegister(const char[],PetscErrorCode(*)(EPS));
767: SLEPC_EXTERN PetscErrorCode EPSMonitorRegister(const char[],PetscViewerType,PetscViewerFormat,EPSMonitorRegisterFn*,EPSMonitorRegisterCreateFn*,EPSMonitorRegisterDestroyFn*);

769: SLEPC_EXTERN PetscErrorCode EPSSetWorkVecs(EPS,PetscInt);
770: SLEPC_EXTERN PetscErrorCode EPSAllocateSolution(EPS,PetscInt);
771: SLEPC_EXTERN PetscErrorCode EPSReallocateSolution(EPS,PetscInt);

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

777:    Calling Sequence:
778: +  eps    - the linear eigensolver context
779: .  eigr   - real part of the eigenvalue
780: .  eigi   - imaginary part of the eigenvalue
781: .  res    - residual norm associated to the eigenpair
782: .  errest - [output] computed error estimate
783: -  ctx    - optional convergence context, as set by `EPSSetConvergenceTestFunction()`

785:    Level: advanced

787: .seealso: [](ch:eps), `EPSSetConvergenceTest()`, `EPSSetConvergenceTestFunction()`
788: S*/
789: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode EPSConvergenceTestFn(EPS eps,PetscScalar eigr,PetscScalar eigi,PetscReal res,PetscReal *errest,void *ctx);

791: SLEPC_EXTERN PetscErrorCode EPSSetConvergenceTest(EPS,EPSConv);
792: SLEPC_EXTERN PetscErrorCode EPSGetConvergenceTest(EPS,EPSConv*);
793: SLEPC_EXTERN EPSConvergenceTestFn EPSConvergedAbsolute;
794: SLEPC_EXTERN EPSConvergenceTestFn EPSConvergedRelative;
795: SLEPC_EXTERN EPSConvergenceTestFn EPSConvergedNorm;
796: SLEPC_EXTERN PetscErrorCode EPSSetConvergenceTestFunction(EPS,EPSConvergenceTestFn*,void*,PetscCtxDestroyFn*);

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

802:    Calling Sequence:
803: +  eps    - the linear eigensolver context
804: .  its    - current number of iterations
805: .  max_it - maximum number of iterations
806: .  nconv  - number of currently converged eigenpairs
807: .  nev    - number of requested eigenpairs
808: .  reason - [output] result of the stopping test
809: -  ctx    - optional stopping context, as set by `EPSSetStoppingTestFunction()`

811:    Note:
812:    A positive value of `reason` indicates that the iteration has finished successfully
813:    (converged), and a negative value indicates an error condition (diverged). If
814:    the iteration needs to be continued, `reason` must be set to `EPS_CONVERGED_ITERATING`
815:    (zero).

817:    Level: advanced

819: .seealso: [](ch:eps), `EPSSetStoppingTest()`, `EPSSetStoppingTestFunction()`
820: S*/
821: PETSC_EXTERN_TYPEDEF typedef PetscErrorCode EPSStoppingTestFn(EPS eps,PetscInt its,PetscInt max_it,PetscInt nconv,PetscInt nev,EPSConvergedReason *reason,void *ctx);

823: SLEPC_EXTERN PetscErrorCode EPSSetStoppingTest(EPS,EPSStop);
824: SLEPC_EXTERN PetscErrorCode EPSGetStoppingTest(EPS,EPSStop*);
825: SLEPC_EXTERN EPSStoppingTestFn EPSStoppingBasic;
826: SLEPC_EXTERN EPSStoppingTestFn EPSStoppingThreshold;
827: SLEPC_EXTERN PetscErrorCode EPSSetStoppingTestFunction(EPS,EPSStoppingTestFn*,void*,PetscCtxDestroyFn*);

829: SLEPC_EXTERN PetscErrorCode EPSSetEigenvalueComparison(EPS,SlepcEigenvalueComparisonFn*,void*);
830: SLEPC_EXTERN PetscErrorCode EPSSetArbitrarySelection(EPS,SlepcArbitrarySelectionFn*,void*);

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

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

837:    Values:
838: +  `EPS_POWER_SHIFT_CONSTANT`  - constant shift
839: .  `EPS_POWER_SHIFT_RAYLEIGH`  - variable shift using Rayleigh quotient
840: -  `EPS_POWER_SHIFT_WILKINSON` - variable shift using Wilkinson's approach

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

845:    Level: advanced

847: .seealso: [](ch:eps), `EPSPowerSetShiftType()`, `EPSPowerGetShiftType()`
848: E*/
849: typedef enum { EPS_POWER_SHIFT_CONSTANT,
850:                EPS_POWER_SHIFT_RAYLEIGH,
851:                EPS_POWER_SHIFT_WILKINSON } EPSPowerShiftType;
852: SLEPC_EXTERN const char *EPSPowerShiftTypes[];

854: /*MC
855:    EPS_POWER_SHIFT_CONSTANT - The power iteration will use a constant shift.

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

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

865:    Level: advanced

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

870: /*MC
871:    EPS_POWER_SHIFT_RAYLEIGH - The power iteration will use a variable shift
872:    computed with the Rayleigh quotient.

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

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

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

886:    Level: advanced

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

891: /*MC
892:    EPS_POWER_SHIFT_WILKINSON - The power iteration will use a variable shift
893:    computed with Wilkinson's approach.

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

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

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

906:    Level: advanced

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

911: SLEPC_EXTERN PetscErrorCode EPSPowerSetShiftType(EPS,EPSPowerShiftType);
912: SLEPC_EXTERN PetscErrorCode EPSPowerGetShiftType(EPS,EPSPowerShiftType*);
913: SLEPC_EXTERN PetscErrorCode EPSPowerSetNonlinear(EPS,PetscBool);
914: SLEPC_EXTERN PetscErrorCode EPSPowerGetNonlinear(EPS,PetscBool*);
915: SLEPC_EXTERN PetscErrorCode EPSPowerSetUpdate(EPS,PetscBool);
916: SLEPC_EXTERN PetscErrorCode EPSPowerGetUpdate(EPS,PetscBool*);
917: SLEPC_EXTERN PetscErrorCode EPSPowerSetSignNormalization(EPS,PetscBool);
918: SLEPC_EXTERN PetscErrorCode EPSPowerGetSignNormalization(EPS,PetscBool*);
919: SLEPC_EXTERN PetscErrorCode EPSPowerSetSNES(EPS,SNES);
920: SLEPC_EXTERN PetscErrorCode EPSPowerGetSNES(EPS,SNES*);

922: SLEPC_EXTERN PetscErrorCode EPSArnoldiSetDelayed(EPS,PetscBool);
923: SLEPC_EXTERN PetscErrorCode EPSArnoldiGetDelayed(EPS,PetscBool*);

925: /*E
926:    EPSKrylovSchurBSEType - The method to be used in the Krylov-Schur solver
927:    for the case of BSE structured eigenproblems.

929:    Values:
930: +  `EPS_KRYLOVSCHUR_BSE_SHAO`         - a Lanczos method proposed by Shao and coauthors
931: .  `EPS_KRYLOVSCHUR_BSE_GRUNING`      - a Lanczos method proposed by Gruning and coauthors
932: -  `EPS_KRYLOVSCHUR_BSE_PROJECTEDBSE` - a Lanczos method resulting is a projected problem with BSE structure

934:    Note:
935:    All variants are implemented in combination with a thick restart, see
936:    the details in {cite:p}`Alv25`.

938:    Level: advanced

940: .seealso: [](ch:eps), `EPSKrylovSchurSetBSEType()`, `EPSKrylovSchurGetBSEType()`
941: E*/
942: typedef enum { EPS_KRYLOVSCHUR_BSE_SHAO,
943:                EPS_KRYLOVSCHUR_BSE_GRUNING,
944:                EPS_KRYLOVSCHUR_BSE_PROJECTEDBSE } EPSKrylovSchurBSEType;
945: SLEPC_EXTERN const char *EPSKrylovSchurBSETypes[];

947: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurSetBSEType(EPS,EPSKrylovSchurBSEType);
948: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetBSEType(EPS,EPSKrylovSchurBSEType*);
949: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurSetRestart(EPS,PetscReal);
950: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetRestart(EPS,PetscReal*);
951: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurSetLocking(EPS,PetscBool);
952: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetLocking(EPS,PetscBool*);
953: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurSetPartitions(EPS,PetscInt);
954: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetPartitions(EPS,PetscInt*);
955: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurSetDetectZeros(EPS,PetscBool);
956: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetDetectZeros(EPS,PetscBool*);
957: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurSetDimensions(EPS,PetscInt,PetscInt,PetscInt);
958: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetDimensions(EPS,PetscInt*,PetscInt*,PetscInt*);
959: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurSetSubintervals(EPS,PetscReal[]);
960: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetSubintervals(EPS,PetscReal*[]);
961: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetInertias(EPS,PetscInt*,PetscReal*[],PetscInt*[]);
962: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetSubcommInfo(EPS,PetscInt*,PetscInt*,Vec*);
963: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetSubcommPairs(EPS,PetscInt,PetscScalar*,Vec);
964: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetSubcommMats(EPS,Mat*,Mat*);
965: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurUpdateSubcommMats(EPS,PetscScalar,PetscScalar,Mat,PetscScalar,PetscScalar, Mat,MatStructure,PetscBool);
966: SLEPC_EXTERN PetscErrorCode EPSKrylovSchurGetKSP(EPS,KSP*);

968: /*E
969:    EPSLanczosReorthogType - The type of reorthogonalization used in `EPSLANCZOS`.

971:    Values:
972: +  `EPS_LANCZOS_REORTHOG_LOCAL`     - local orthogonalization, only involves previous two vectors
973: .  `EPS_LANCZOS_REORTHOG_FULL`      - full orthogonalization against all previous vectors
974: .  `EPS_LANCZOS_REORTHOG_SELECTIVE` - selective reorthogonalization against nearly converged Ritz vectors
975: .  `EPS_LANCZOS_REORTHOG_PERIODIC`  - periodic reorthogonalization against all previous vectors
976: .  `EPS_LANCZOS_REORTHOG_PARTIAL`   - partial reorthogonalization against as subset of previous vectors
977: -  `EPS_LANCZOS_REORTHOG_DELAYED`   - full orthogonalization with delayed reorthogonalization

979:    Note:
980:    Details of the different reorthogonalization strategies can be found in
981:    {cite:p}`Her06`.

983:    Level: advanced

985: .seealso: [](ch:eps), `EPSLanczosSetReorthog()`, `EPSLanczosGetReorthog()`
986: E*/
987: typedef enum { EPS_LANCZOS_REORTHOG_LOCAL,
988:                EPS_LANCZOS_REORTHOG_FULL,
989:                EPS_LANCZOS_REORTHOG_SELECTIVE,
990:                EPS_LANCZOS_REORTHOG_PERIODIC,
991:                EPS_LANCZOS_REORTHOG_PARTIAL,
992:                EPS_LANCZOS_REORTHOG_DELAYED } EPSLanczosReorthogType;
993: SLEPC_EXTERN const char *EPSLanczosReorthogTypes[];

995: SLEPC_EXTERN PetscErrorCode EPSLanczosSetReorthog(EPS,EPSLanczosReorthogType);
996: SLEPC_EXTERN PetscErrorCode EPSLanczosGetReorthog(EPS,EPSLanczosReorthogType*);

998: /*E
999:    EPSPRIMMEMethod - The method selected in the PRIMME library.

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

1004:    Level: advanced

1006: .seealso: [](ch:eps), `EPSPRIMMESetMethod()`, `EPSPRIMMEGetMethod()`
1007: E*/
1008: typedef enum { EPS_PRIMME_DYNAMIC             = 1,
1009:                EPS_PRIMME_DEFAULT_MIN_TIME    = 2,
1010:                EPS_PRIMME_DEFAULT_MIN_MATVECS = 3,
1011:                EPS_PRIMME_ARNOLDI             = 4,
1012:                EPS_PRIMME_GD                  = 5,
1013:                EPS_PRIMME_GD_PLUSK            = 6,
1014:                EPS_PRIMME_GD_OLSEN_PLUSK      = 7,
1015:                EPS_PRIMME_JD_OLSEN_PLUSK      = 8,
1016:                EPS_PRIMME_RQI                 = 9,
1017:                EPS_PRIMME_JDQR                = 10,
1018:                EPS_PRIMME_JDQMR               = 11,
1019:                EPS_PRIMME_JDQMR_ETOL          = 12,
1020:                EPS_PRIMME_SUBSPACE_ITERATION  = 13,
1021:                EPS_PRIMME_LOBPCG_ORTHOBASIS   = 14,
1022:                EPS_PRIMME_LOBPCG_ORTHOBASISW  = 15 } EPSPRIMMEMethod;
1023: SLEPC_EXTERN const char *EPSPRIMMEMethods[];

1025: SLEPC_EXTERN PetscErrorCode EPSPRIMMESetBlockSize(EPS,PetscInt);
1026: SLEPC_EXTERN PetscErrorCode EPSPRIMMEGetBlockSize(EPS,PetscInt*);
1027: SLEPC_EXTERN PetscErrorCode EPSPRIMMESetMethod(EPS,EPSPRIMMEMethod);
1028: SLEPC_EXTERN PetscErrorCode EPSPRIMMEGetMethod(EPS,EPSPRIMMEMethod*);

1030: SLEPC_EXTERN PetscErrorCode EPSGDSetKrylovStart(EPS,PetscBool);
1031: SLEPC_EXTERN PetscErrorCode EPSGDGetKrylovStart(EPS,PetscBool*);
1032: SLEPC_EXTERN PetscErrorCode EPSGDSetBlockSize(EPS,PetscInt);
1033: SLEPC_EXTERN PetscErrorCode EPSGDGetBlockSize(EPS,PetscInt*);
1034: SLEPC_EXTERN PetscErrorCode EPSGDSetRestart(EPS,PetscInt,PetscInt);
1035: SLEPC_EXTERN PetscErrorCode EPSGDGetRestart(EPS,PetscInt*,PetscInt*);
1036: SLEPC_EXTERN PetscErrorCode EPSGDSetInitialSize(EPS,PetscInt);
1037: SLEPC_EXTERN PetscErrorCode EPSGDGetInitialSize(EPS,PetscInt*);
1038: SLEPC_EXTERN PetscErrorCode EPSGDSetBOrth(EPS,PetscBool);
1039: SLEPC_EXTERN PetscErrorCode EPSGDGetBOrth(EPS,PetscBool*);
1040: SLEPC_EXTERN PetscErrorCode EPSGDSetDoubleExpansion(EPS,PetscBool);
1041: SLEPC_EXTERN PetscErrorCode EPSGDGetDoubleExpansion(EPS,PetscBool*);

1043: SLEPC_EXTERN PetscErrorCode EPSJDSetKrylovStart(EPS,PetscBool);
1044: SLEPC_EXTERN PetscErrorCode EPSJDGetKrylovStart(EPS,PetscBool*);
1045: SLEPC_EXTERN PetscErrorCode EPSJDSetBlockSize(EPS,PetscInt);
1046: SLEPC_EXTERN PetscErrorCode EPSJDGetBlockSize(EPS,PetscInt*);
1047: SLEPC_EXTERN PetscErrorCode EPSJDSetRestart(EPS,PetscInt,PetscInt);
1048: SLEPC_EXTERN PetscErrorCode EPSJDGetRestart(EPS,PetscInt*,PetscInt*);
1049: SLEPC_EXTERN PetscErrorCode EPSJDSetInitialSize(EPS,PetscInt);
1050: SLEPC_EXTERN PetscErrorCode EPSJDGetInitialSize(EPS,PetscInt*);
1051: SLEPC_EXTERN PetscErrorCode EPSJDSetFix(EPS,PetscReal);
1052: SLEPC_EXTERN PetscErrorCode EPSJDGetFix(EPS,PetscReal*);
1053: SLEPC_EXTERN PetscErrorCode EPSJDSetConstCorrectionTol(EPS,PetscBool);
1054: SLEPC_EXTERN PetscErrorCode EPSJDGetConstCorrectionTol(EPS,PetscBool*);
1055: SLEPC_EXTERN PetscErrorCode EPSJDSetBOrth(EPS,PetscBool);
1056: SLEPC_EXTERN PetscErrorCode EPSJDGetBOrth(EPS,PetscBool*);

1058: SLEPC_EXTERN PetscErrorCode EPSRQCGSetReset(EPS,PetscInt);
1059: SLEPC_EXTERN PetscErrorCode EPSRQCGGetReset(EPS,PetscInt*);

1061: SLEPC_EXTERN PetscErrorCode EPSLOBPCGSetBlockSize(EPS,PetscInt);
1062: SLEPC_EXTERN PetscErrorCode EPSLOBPCGGetBlockSize(EPS,PetscInt*);
1063: SLEPC_EXTERN PetscErrorCode EPSLOBPCGSetRestart(EPS,PetscReal);
1064: SLEPC_EXTERN PetscErrorCode EPSLOBPCGGetRestart(EPS,PetscReal*);
1065: SLEPC_EXTERN PetscErrorCode EPSLOBPCGSetLocking(EPS,PetscBool);
1066: SLEPC_EXTERN PetscErrorCode EPSLOBPCGGetLocking(EPS,PetscBool*);

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

1071:    Values:
1072: +  `EPS_CISS_QUADRULE_TRAPEZOIDAL` - trapezoidal rule
1073: -  `EPS_CISS_QUADRULE_CHEBYSHEV`   - Gauss quadrature on Chebyshev points

1075:    Note:
1076:    For a detailed description see {cite:p}`Mae16`.

1078:    Level: advanced

1080: .seealso: [](ch:eps), `EPSCISSSetQuadRule()`, `EPSCISSGetQuadRule()`
1081: E*/
1082: typedef enum { EPS_CISS_QUADRULE_TRAPEZOIDAL = 1,
1083:                EPS_CISS_QUADRULE_CHEBYSHEV   = 2 } EPSCISSQuadRule;
1084: SLEPC_EXTERN const char *EPSCISSQuadRules[];

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

1089:    Values:
1090: +  `EPS_CISS_EXTRACTION_RITZ`   - Ritz approximations from Rayleigh-Ritz projection
1091: -  `EPS_CISS_EXTRACTION_HANKEL` - use Hankel pencil as in the original Sakurai-Sugiura method

1093:    Note:
1094:    For a detailed description see {cite:p}`Mae16`.

1096:    Level: advanced

1098: .seealso: [](ch:eps), `EPSCISSSetExtraction()`, `EPSCISSGetExtraction()`
1099: E*/
1100: typedef enum { EPS_CISS_EXTRACTION_RITZ,
1101:                EPS_CISS_EXTRACTION_HANKEL } EPSCISSExtraction;
1102: SLEPC_EXTERN const char *EPSCISSExtractions[];

1104: SLEPC_EXTERN PetscErrorCode EPSCISSSetExtraction(EPS,EPSCISSExtraction);
1105: SLEPC_EXTERN PetscErrorCode EPSCISSGetExtraction(EPS,EPSCISSExtraction*);
1106: SLEPC_EXTERN PetscErrorCode EPSCISSSetQuadRule(EPS,EPSCISSQuadRule);
1107: SLEPC_EXTERN PetscErrorCode EPSCISSGetQuadRule(EPS,EPSCISSQuadRule*);
1108: SLEPC_EXTERN PetscErrorCode EPSCISSSetSizes(EPS,PetscInt,PetscInt,PetscInt,PetscInt,PetscInt,PetscBool);
1109: SLEPC_EXTERN PetscErrorCode EPSCISSGetSizes(EPS,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscInt*,PetscBool*);
1110: SLEPC_EXTERN PetscErrorCode EPSCISSSetThreshold(EPS,PetscReal,PetscReal);
1111: SLEPC_EXTERN PetscErrorCode EPSCISSGetThreshold(EPS,PetscReal*,PetscReal*);
1112: SLEPC_EXTERN PetscErrorCode EPSCISSSetRefinement(EPS,PetscInt,PetscInt);
1113: SLEPC_EXTERN PetscErrorCode EPSCISSGetRefinement(EPS,PetscInt*,PetscInt*);
1114: SLEPC_EXTERN PetscErrorCode EPSCISSSetUseST(EPS,PetscBool);
1115: SLEPC_EXTERN PetscErrorCode EPSCISSGetUseST(EPS,PetscBool*);
1116: SLEPC_EXTERN PetscErrorCode EPSCISSGetKSPs(EPS,PetscInt*,KSP*[]);

1118: SLEPC_EXTERN PetscErrorCode EPSLyapIISetLME(EPS,LME);
1119: SLEPC_EXTERN PetscErrorCode EPSLyapIIGetLME(EPS,LME*);
1120: SLEPC_EXTERN PetscErrorCode EPSLyapIISetRanks(EPS,PetscInt,PetscInt);
1121: SLEPC_EXTERN PetscErrorCode EPSLyapIIGetRanks(EPS,PetscInt*,PetscInt*);

1123: SLEPC_EXTERN PetscErrorCode EPSBLOPEXSetBlockSize(EPS,PetscInt);
1124: SLEPC_EXTERN PetscErrorCode EPSBLOPEXGetBlockSize(EPS,PetscInt*);

1126: /*E
1127:    EPSEVSLDOSMethod - The method to approximate the density of states (DOS)
1128:    in the `EPSEVSL` solver.

1130:    Values:
1131: +  `EPS_EVSL_DOS_KPM`     - Kernel Polynomial Method
1132: -  `EPS_EVSL_DOS_LANCZOS` - Lanczos method

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

1137:    Level: advanced

1139: .seealso: [](ch:eps), `EPSEVSLSetDOSParameters()`, `EPSEVSLGetDOSParameters()`
1140: E*/
1141: typedef enum { EPS_EVSL_DOS_KPM,
1142:                EPS_EVSL_DOS_LANCZOS } EPSEVSLDOSMethod;
1143: SLEPC_EXTERN const char *EPSEVSLDOSMethods[];

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

1148:    Values:
1149: +  `EPS_EVSL_DAMPING_NONE`    - no damping
1150: .  `EPS_EVSL_DAMPING_JACKSON` - Jackson damping
1151: -  `EPS_EVSL_DAMPING_SIGMA`   - Lanczos damping

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

1156:    Level: advanced

1158: .seealso: [](ch:eps), `EPSEVSLSetDOSParameters()`, `EPSEVSLGetDOSParameters()`
1159: E*/
1160: typedef enum { EPS_EVSL_DAMPING_NONE,
1161:                EPS_EVSL_DAMPING_JACKSON,
1162:                EPS_EVSL_DAMPING_SIGMA } EPSEVSLDamping;
1163: SLEPC_EXTERN const char *EPSEVSLDampings[];

1165: SLEPC_EXTERN PetscErrorCode EPSEVSLSetRange(EPS,PetscReal,PetscReal);
1166: SLEPC_EXTERN PetscErrorCode EPSEVSLGetRange(EPS,PetscReal*,PetscReal*);
1167: SLEPC_EXTERN PetscErrorCode EPSEVSLSetSlices(EPS,PetscInt);
1168: SLEPC_EXTERN PetscErrorCode EPSEVSLGetSlices(EPS,PetscInt*);
1169: SLEPC_EXTERN PetscErrorCode EPSEVSLSetDOSParameters(EPS,EPSEVSLDOSMethod,PetscInt,PetscInt,PetscInt,PetscInt);
1170: SLEPC_EXTERN PetscErrorCode EPSEVSLGetDOSParameters(EPS,EPSEVSLDOSMethod*,PetscInt*,PetscInt*,PetscInt*,PetscInt*);
1171: SLEPC_EXTERN PetscErrorCode EPSEVSLSetPolParameters(EPS,PetscInt,PetscReal);
1172: SLEPC_EXTERN PetscErrorCode EPSEVSLGetPolParameters(EPS,PetscInt*,PetscReal*);
1173: SLEPC_EXTERN PetscErrorCode EPSEVSLSetDamping(EPS,EPSEVSLDamping);
1174: SLEPC_EXTERN PetscErrorCode EPSEVSLGetDamping(EPS,EPSEVSLDamping*);

1176: SLEPC_EXTERN PetscErrorCode EPSFEASTSetNumPoints(EPS,PetscInt);
1177: SLEPC_EXTERN PetscErrorCode EPSFEASTGetNumPoints(EPS,PetscInt*);

1179: SLEPC_EXTERN PetscErrorCode EPSCHASESetDegree(EPS,PetscInt,PetscBool);
1180: SLEPC_EXTERN PetscErrorCode EPSCHASEGetDegree(EPS,PetscInt*,PetscBool*);