Actual source code: stimpl.h
slepc-3.14.0 2020-09-30
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: */
11: #if !defined(SLEPCSTIMPL_H)
12: #define SLEPCSTIMPL_H
14: #include <slepcst.h>
15: #include <slepc/private/slepcimpl.h>
17: SLEPC_EXTERN PetscBool STRegisterAllCalled;
18: SLEPC_EXTERN PetscErrorCode STRegisterAll(void);
19: SLEPC_EXTERN PetscLogEvent ST_SetUp,ST_ComputeOperator,ST_Apply,ST_ApplyTranspose,ST_MatSetUp,ST_MatMult,ST_MatMultTranspose,ST_MatSolve,ST_MatSolveTranspose;
21: typedef struct _STOps *STOps;
23: struct _STOps {
24: PetscErrorCode (*apply)(ST,Vec,Vec);
25: PetscErrorCode (*applymat)(ST,Mat,Mat);
26: PetscErrorCode (*applytrans)(ST,Vec,Vec);
27: PetscErrorCode (*backtransform)(ST,PetscInt,PetscScalar*,PetscScalar*);
28: PetscErrorCode (*setshift)(ST,PetscScalar);
29: PetscErrorCode (*getbilinearform)(ST,Mat*);
30: PetscErrorCode (*setup)(ST);
31: PetscErrorCode (*computeoperator)(ST);
32: PetscErrorCode (*setfromoptions)(PetscOptionItems*,ST);
33: PetscErrorCode (*postsolve)(ST);
34: PetscErrorCode (*destroy)(ST);
35: PetscErrorCode (*reset)(ST);
36: PetscErrorCode (*view)(ST,PetscViewer);
37: PetscErrorCode (*checknullspace)(ST,BV);
38: PetscErrorCode (*setdefaultksp)(ST);
39: };
41: /*
42: 'Updated' state means STSetUp must be called because matrices have been
43: modified, but the pattern is the same (hence reuse symbolic factorization)
44: */
45: typedef enum { ST_STATE_INITIAL,
46: ST_STATE_SETUP,
47: ST_STATE_UPDATED } STStateType;
49: struct _p_ST {
50: PETSCHEADER(struct _STOps);
51: /*------------------------- User parameters --------------------------*/
52: Mat *A; /* matrices that define the eigensystem */
53: PetscInt nmat; /* number of user-provided matrices */
54: PetscScalar sigma; /* value of the shift */
55: PetscScalar defsigma; /* default value of the shift */
56: STMatMode matmode; /* how the transformation matrix is handled */
57: MatStructure str; /* whether matrices have the same pattern or not */
58: PetscBool transform; /* whether transformed matrices are computed */
59: Vec D; /* diagonal matrix for balancing */
61: /*------------------------- Misc data --------------------------*/
62: KSP ksp; /* linear solver used in some ST's */
63: PetscBool usesksp; /* whether the KSP object is used or not */
64: PetscInt nwork; /* number of work vectors */
65: Vec *work; /* work vectors */
66: Vec wb; /* balancing requires an extra work vector */
67: Vec wht; /* extra work vector for hermitian transpose apply */
68: STStateType state; /* initial -> setup -> with updated matrices */
69: PetscObjectState *Astate; /* matrix state (to identify the original matrices) */
70: Mat *T; /* matrices resulting from transformation */
71: Mat Op; /* shell matrix for operator = alpha*D*inv(P)*M*inv(D) */
72: PetscBool opseized; /* whether Op has been seized by user */
73: PetscBool opready; /* whether Op is up-to-date or need be computed */
74: Mat P; /* matrix from which preconditioner is built */
75: Mat M; /* matrix corresponding to the non-inverted part of the operator */
76: PetscBool sigma_set; /* whether the user provided the shift or not */
77: PetscBool asymm; /* the user matrices are all symmetric */
78: PetscBool aherm; /* the user matrices are all hermitian */
79: void *data;
80: };
82: /*
83: Macros to test valid ST arguments
84: */
85: #if !defined(PETSC_USE_DEBUG)
87: #define STCheckMatrices(h,arg) do {} while (0)
88: #define STCheckNotSeized(h,arg) do {} while (0)
90: #else
92: #define STCheckMatrices(h,arg) \
93: do { \
94: if (!(h)->A) SETERRQ1(PetscObjectComm((PetscObject)(h)),PETSC_ERR_ARG_WRONGSTATE,"ST matrices have not been set: Parameter #%d",arg); \
95: } while (0)
96: #define STCheckNotSeized(h,arg) \
97: do { \
98: if (h->opseized) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"Must call STRestoreOperator() first: Parameter #%d",arg); \
99: } while (0)
101: #endif
103: SLEPC_INTERN PetscErrorCode STGetBilinearForm_Default(ST,Mat*);
104: SLEPC_INTERN PetscErrorCode STCheckNullSpace_Default(ST,BV);
105: SLEPC_INTERN PetscErrorCode STMatShellCreate(ST,PetscScalar,PetscInt,PetscInt*,PetscScalar*,Mat*);
106: SLEPC_INTERN PetscErrorCode STMatShellShift(Mat,PetscScalar);
107: SLEPC_INTERN PetscErrorCode STCheckFactorPackage(ST);
108: SLEPC_INTERN PetscErrorCode STMatMAXPY_Private(ST,PetscScalar,PetscScalar,PetscInt,PetscScalar*,PetscBool,Mat*);
109: SLEPC_INTERN PetscErrorCode STCoeffs_Monomial(ST,PetscScalar*);
110: SLEPC_INTERN PetscErrorCode STSetDefaultKSP(ST);
111: SLEPC_INTERN PetscErrorCode STSetDefaultKSP_Default(ST);
112: SLEPC_INTERN PetscErrorCode STIsInjective_Shell(ST,PetscBool*);
113: SLEPC_INTERN PetscErrorCode STComputeOperator(ST);
114: SLEPC_INTERN PetscErrorCode STGetOperator_Private(ST,Mat*);
115: SLEPC_INTERN PetscErrorCode STApply_Generic(ST,Vec,Vec);
116: SLEPC_INTERN PetscErrorCode STApplyMat_Generic(ST,Mat,Mat);
117: SLEPC_INTERN PetscErrorCode STApplyTranspose_Generic(ST,Vec,Vec);
119: #endif