Actual source code: stimpl.h

slepc-3.17.1 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */

 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 */
 60:   Mat              Pmat;             /* user-provided preconditioner matrix */
 61:   PetscBool        Pmat_set;         /* whether the user provided a preconditioner matrix or not  */
 62:   Mat              *Psplit;          /* matrices for the split preconditioner */
 63:   PetscInt         nsplit;           /* number of split preconditioner matrices */
 64:   MatStructure     strp;             /* pattern of split preconditioner matrices */

 66:   /*------------------------- Misc data --------------------------*/
 67:   KSP              ksp;              /* linear solver used in some ST's */
 68:   PetscBool        usesksp;          /* whether the KSP object is used or not */
 69:   PetscInt         nwork;            /* number of work vectors */
 70:   Vec              *work;            /* work vectors */
 71:   Vec              wb;               /* balancing requires an extra work vector */
 72:   Vec              wht;              /* extra work vector for hermitian transpose apply */
 73:   STStateType      state;            /* initial -> setup -> with updated matrices */
 74:   PetscObjectState *Astate;          /* matrix state (to identify the original matrices) */
 75:   Mat              *T;               /* matrices resulting from transformation */
 76:   Mat              Op;               /* shell matrix for operator = alpha*D*inv(P)*M*inv(D) */
 77:   PetscBool        opseized;         /* whether Op has been seized by user */
 78:   PetscBool        opready;          /* whether Op is up-to-date or need be computed  */
 79:   Mat              P;                /* matrix from which preconditioner is built */
 80:   Mat              M;                /* matrix corresponding to the non-inverted part of the operator */
 81:   PetscBool        sigma_set;        /* whether the user provided the shift or not */
 82:   PetscBool        asymm;            /* the user matrices are all symmetric */
 83:   PetscBool        aherm;            /* the user matrices are all hermitian */
 84:   void             *data;
 85: };

 87: /*
 88:     Macros to test valid ST arguments
 89: */
 90: #if !defined(PETSC_USE_DEBUG)

 92: #define STCheckMatrices(h,arg) do {(void)(h);} while (0)
 93: #define STCheckNotSeized(h,arg) do {(void)(h);} while (0)

 95: #else

 97: #define STCheckMatrices(h,arg) \
 98:   do { \
100:   } while (0)
101: #define STCheckNotSeized(h,arg) \
102:   do { \
104:   } while (0)

106: #endif

108: SLEPC_INTERN PetscErrorCode STGetBilinearForm_Default(ST,Mat*);
109: SLEPC_INTERN PetscErrorCode STCheckNullSpace_Default(ST,BV);
110: SLEPC_INTERN PetscErrorCode STMatShellCreate(ST,PetscScalar,PetscInt,PetscInt*,PetscScalar*,Mat*);
111: SLEPC_INTERN PetscErrorCode STMatShellShift(Mat,PetscScalar);
112: SLEPC_INTERN PetscErrorCode STCheckFactorPackage(ST);
113: SLEPC_INTERN PetscErrorCode STMatMAXPY_Private(ST,PetscScalar,PetscScalar,PetscInt,PetscScalar*,PetscBool,PetscBool,Mat*);
114: SLEPC_INTERN PetscErrorCode STCoeffs_Monomial(ST,PetscScalar*);
115: SLEPC_INTERN PetscErrorCode STSetDefaultKSP(ST);
116: SLEPC_INTERN PetscErrorCode STSetDefaultKSP_Default(ST);
117: SLEPC_INTERN PetscErrorCode STIsInjective_Shell(ST,PetscBool*);
118: SLEPC_INTERN PetscErrorCode STComputeOperator(ST);
119: SLEPC_INTERN PetscErrorCode STGetOperator_Private(ST,Mat*);
120: SLEPC_INTERN PetscErrorCode STApply_Generic(ST,Vec,Vec);
121: SLEPC_INTERN PetscErrorCode STApplyMat_Generic(ST,Mat,Mat);
122: SLEPC_INTERN PetscErrorCode STApplyTranspose_Generic(ST,Vec,Vec);

124: /*
125:   ST_KSPSetOperators - Sets the KSP matrices
126: */
127: static inline PetscErrorCode ST_KSPSetOperators(ST st,Mat A,Mat B)
128: {
129:   const char     *prefix;

131:   if (!st->ksp) STGetKSP(st,&st->ksp);
132:   STCheckFactorPackage(st);
133:   KSPSetOperators(st->ksp,A,B);
134:   MatGetOptionsPrefix(B,&prefix);
135:   if (!prefix) {
136:     /* set Mat prefix to be the same as KSP to enable setting command-line options (e.g. MUMPS)
137:        only applies if the Mat has no user-defined prefix */
138:     KSPGetOptionsPrefix(st->ksp,&prefix);
139:     MatSetOptionsPrefix(B,prefix);
140:   }
141:   PetscFunctionReturn(0);
142: }

144: #endif