Actual source code: stimpl.h

slepc-3.9.2 2018-07-02
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2018, 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(_STIMPL)
 12: #define _STIMPL

 14: #include <slepcst.h>
 15: #include <slepc/private/slepcimpl.h>

 17: PETSC_EXTERN PetscBool STRegisterAllCalled;
 18: PETSC_EXTERN PetscErrorCode STRegisterAll(void);
 19: PETSC_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:   PetscObjectState *Astate;          /* state (to identify the original matrices) */
 52:   Mat              *T;               /* matrices resulting from transformation */
 53:   Mat              P;                /* matrix from which preconditioner is built */
 54:   PetscInt         nmat;             /* number of matrices */
 55:   PetscScalar      sigma;            /* value of the shift */
 56:   PetscBool        sigma_set;        /* whether the user provided the shift or not */
 57:   PetscScalar      defsigma;         /* default value of the shift */
 58:   STMatMode        shift_matrix;
 59:   MatStructure     str;              /* whether matrices have the same pattern or not */
 60:   PetscBool        transform;        /* whether transformed matrices are computed */

 62:   /*------------------------- Misc data --------------------------*/
 63:   KSP              ksp;
 64:   PetscInt         nwork;            /* number of work vectors */
 65:   Vec              *work;            /* work vectors */
 66:   Vec              D;                /* diagonal matrix for balancing */
 67:   Vec              wb;               /* balancing requires an extra work vector */
 68:   void             *data;
 69:   STStateType      state;            /* initial -> setup -> with updated matrices */
 70: };

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

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

 79: #else

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

 86: #endif

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

100: #endif