Actual source code: stimpl.h

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) , Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.

  8:    SLEPc is free software: you can redistribute it and/or modify it under  the
  9:    terms of version 3 of the GNU Lesser General Public License as published by
 10:    the Free Software Foundation.

 12:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 13:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 14:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 15:    more details.

 17:    You  should have received a copy of the GNU Lesser General  Public  License
 18:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 19:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 20: */

 22: #if !defined(_STIMPL)
 23: #define _STIMPL

 25: #include <slepcst.h>
 26: #include <slepc-private/slepcimpl.h>

 28: PETSC_EXTERN PetscLogEvent ST_SetUp,ST_Apply,ST_ApplyTranspose,ST_MatSetUp,ST_MatMult,ST_MatMultTranspose,ST_MatSolve,ST_MatSolveTranspose;

 30: typedef struct _STOps *STOps;

 32: struct _STOps {
 33:   PetscErrorCode (*setup)(ST);
 34:   PetscErrorCode (*apply)(ST,Vec,Vec);
 35:   PetscErrorCode (*getbilinearform)(ST,Mat*);
 36:   PetscErrorCode (*applytrans)(ST,Vec,Vec);
 37:   PetscErrorCode (*setshift)(ST,PetscScalar);
 38:   PetscErrorCode (*setfromoptions)(ST);
 39:   PetscErrorCode (*postsolve)(ST);
 40:   PetscErrorCode (*backtransform)(ST,PetscInt,PetscScalar*,PetscScalar*);
 41:   PetscErrorCode (*destroy)(ST);
 42:   PetscErrorCode (*reset)(ST);
 43:   PetscErrorCode (*view)(ST,PetscViewer);
 44:   PetscErrorCode (*checknullspace)(ST,BV);
 45: };

 47: struct _p_ST {
 48:   PETSCHEADER(struct _STOps);
 49:   /*------------------------- User parameters --------------------------*/
 50:   Mat          *A;               /* Matrices that define the eigensystem */
 51:   PetscInt     *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:   Vec          w;
 65:   Vec          D;                /* diagonal matrix for balancing */
 66:   Vec          wb;               /* balancing requires an extra work vector */
 67:   PetscInt     linearits;        /* number of linear iterations */
 68:   PetscInt     applys;           /* number of operator applies */
 69:   void         *data;
 70:   PetscInt     setupcalled;
 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: PETSC_INTERN PetscErrorCode STGetBilinearForm_Default(ST,Mat*);
 90: PETSC_INTERN PetscErrorCode STCheckNullSpace_Default(ST,BV);
 91: PETSC_INTERN PetscErrorCode STMatShellCreate(ST,PetscScalar,PetscInt,PetscInt*,PetscScalar*,Mat*);
 92: PETSC_INTERN PetscErrorCode STMatShellShift(Mat,PetscScalar);
 93: PETSC_INTERN PetscErrorCode STMatSetHermitian(ST,Mat);
 94: PETSC_INTERN PetscErrorCode STMatMAXPY_Private(ST,PetscScalar,PetscScalar,PetscInt,PetscScalar*,PetscBool,Mat*);
 95: PETSC_INTERN PetscErrorCode STCoeffs_Monomial(ST,PetscScalar*);

 97: #endif