Actual source code: stimpl.h

slepc-3.12.0 2019-09-30
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2019, 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_Apply,ST_ApplyTranspose,ST_MatSetUp,ST_MatMult,ST_MatMultTranspose,ST_MatSolve,ST_MatSolveTranspose;

 21: typedef struct _STOps *STOps;

 23: struct _STOps {
 24:   PetscErrorCode (*setup)(ST);
 25:   PetscErrorCode (*apply)(ST,Vec,Vec);
 26:   PetscErrorCode (*getbilinearform)(ST,Mat*);
 27:   PetscErrorCode (*applytrans)(ST,Vec,Vec);
 28:   PetscErrorCode (*setshift)(ST,PetscScalar);
 29:   PetscErrorCode (*setfromoptions)(PetscOptionItems*,ST);
 30:   PetscErrorCode (*postsolve)(ST);
 31:   PetscErrorCode (*backtransform)(ST,PetscInt,PetscScalar*,PetscScalar*);
 32:   PetscErrorCode (*destroy)(ST);
 33:   PetscErrorCode (*reset)(ST);
 34:   PetscErrorCode (*view)(ST,PetscViewer);
 35:   PetscErrorCode (*checknullspace)(ST,BV);
 36:   PetscErrorCode (*setdefaultksp)(ST);
 37: };

 39: /*
 40:      'Updated' state means STSetUp must be called because matrices have been
 41:      modified, but the pattern is the same (hence reuse symbolic factorization)
 42: */
 43: typedef enum { ST_STATE_INITIAL,
 44:                ST_STATE_SETUP,
 45:                ST_STATE_UPDATED } STStateType;

 47: struct _p_ST {
 48:   PETSCHEADER(struct _STOps);
 49:   /*------------------------- User parameters --------------------------*/
 50:   Mat              *A;               /* matrices that define the eigensystem */
 51:   PetscInt         nmat;             /* number of user-provided matrices */
 52:   PetscScalar      sigma;            /* value of the shift */
 53:   PetscScalar      defsigma;         /* default value of the shift */
 54:   STMatMode        matmode;          /* how the transformation matrix is handled */
 55:   MatStructure     str;              /* whether matrices have the same pattern or not */
 56:   PetscBool        transform;        /* whether transformed matrices are computed */
 57:   Vec              D;                /* diagonal matrix for balancing */

 59:   /*------------------------- Misc data --------------------------*/
 60:   KSP              ksp;              /* linear solver used in some ST's */
 61:   PetscBool        usesksp;          /* whether the KSP object is used or not */
 62:   PetscInt         nwork;            /* number of work vectors */
 63:   Vec              *work;            /* work vectors */
 64:   Vec              wb;               /* balancing requires an extra work vector */
 65:   STStateType      state;            /* initial -> setup -> with updated matrices */
 66:   PetscObjectState *Astate;          /* matrix state (to identify the original matrices) */
 67:   Mat              *T;               /* matrices resulting from transformation */
 68:   Mat              P;                /* matrix from which preconditioner is built */
 69:   PetscBool        sigma_set;        /* whether the user provided the shift or not */
 70:   void             *data;
 71: };

 73: /*
 74:     Macros to test valid ST arguments
 75: */
 76: #if !defined(PETSC_USE_DEBUG)

 78: #define STCheckMatrices(h,arg) do {} while (0)

 80: #else

 82: #define STCheckMatrices(h,arg) \
 83:   do { \
 84:     if (!h->A) SETERRQ1(PetscObjectComm((PetscObject)h),PETSC_ERR_ARG_WRONGSTATE,"ST matrices have not been set: Parameter #%d",arg); \
 85:   } while (0)

 87: #endif

 89: SLEPC_INTERN PetscErrorCode STGetBilinearForm_Default(ST,Mat*);
 90: SLEPC_INTERN PetscErrorCode STCheckNullSpace_Default(ST,BV);
 91: SLEPC_INTERN PetscErrorCode STMatShellCreate(ST,PetscScalar,PetscInt,PetscInt*,PetscScalar*,Mat*);
 92: SLEPC_INTERN PetscErrorCode STMatShellShift(Mat,PetscScalar);
 93: SLEPC_INTERN PetscErrorCode STMatSetHermitian(ST,Mat);
 94: SLEPC_INTERN PetscErrorCode STCheckFactorPackage(ST);
 95: SLEPC_INTERN PetscErrorCode STMatMAXPY_Private(ST,PetscScalar,PetscScalar,PetscInt,PetscScalar*,PetscBool,Mat*);
 96: SLEPC_INTERN PetscErrorCode STCoeffs_Monomial(ST,PetscScalar*);
 97: SLEPC_INTERN PetscErrorCode STSetDefaultKSP(ST);
 98: SLEPC_INTERN PetscErrorCode STSetDefaultKSP_Default(ST);
 99: SLEPC_INTERN PetscErrorCode STIsInjective_Shell(ST,PetscBool*);

101: #endif