Actual source code: slepcst.h

slepc-3.6.3 2016-03-29
Report Typos and Errors
  1: /*
  2:    Spectral transformation module for eigenvalue problems.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

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

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

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

 26: #include <slepcsys.h>
 27: #include <slepcbv.h>
 28: #include <petscksp.h>

 30: PETSC_EXTERN PetscErrorCode STInitializePackage(void);

 32: /*S
 33:     ST - Abstract SLEPc object that manages spectral transformations.
 34:     This object is accessed only in advanced applications.

 36:     Level: beginner

 38: .seealso:  STCreate(), EPS
 39: S*/
 40: typedef struct _p_ST* ST;

 42: /*J
 43:     STType - String with the name of a SLEPc spectral transformation

 45:     Level: beginner

 47: .seealso: STSetType(), ST
 48: J*/
 49: typedef const char* STType;
 50: #define STSHELL     "shell"
 51: #define STSHIFT     "shift"
 52: #define STSINVERT   "sinvert"
 53: #define STCAYLEY    "cayley"
 54: #define STPRECOND   "precond"

 56: /* Logging support */
 57: PETSC_EXTERN PetscClassId ST_CLASSID;

 59: PETSC_EXTERN PetscErrorCode STCreate(MPI_Comm,ST*);
 60: PETSC_EXTERN PetscErrorCode STDestroy(ST*);
 61: PETSC_EXTERN PetscErrorCode STReset(ST);
 62: PETSC_EXTERN PetscErrorCode STSetType(ST,STType);
 63: PETSC_EXTERN PetscErrorCode STGetType(ST,STType*);
 64: PETSC_EXTERN PetscErrorCode STSetOperators(ST,PetscInt,Mat*);
 65: PETSC_EXTERN PetscErrorCode STGetOperators(ST,PetscInt,Mat*);
 66: PETSC_EXTERN PetscErrorCode STGetTOperators(ST,PetscInt,Mat*);
 67: PETSC_EXTERN PetscErrorCode STGetNumMatrices(ST,PetscInt*);
 68: PETSC_EXTERN PetscErrorCode STSetUp(ST);
 69: PETSC_EXTERN PetscErrorCode STSetFromOptions(ST);
 70: PETSC_EXTERN PetscErrorCode STView(ST,PetscViewer);

 72: PETSC_EXTERN PetscErrorCode STApply(ST,Vec,Vec);
 73: PETSC_EXTERN PetscErrorCode STMatMult(ST,PetscInt,Vec,Vec);
 74: PETSC_EXTERN PetscErrorCode STMatMultTranspose(ST,PetscInt,Vec,Vec);
 75: PETSC_EXTERN PetscErrorCode STMatSolve(ST,Vec,Vec);
 76: PETSC_EXTERN PetscErrorCode STMatSolveTranspose(ST,Vec,Vec);
 77: PETSC_EXTERN PetscErrorCode STGetBilinearForm(ST,Mat*);
 78: PETSC_EXTERN PetscErrorCode STApplyTranspose(ST,Vec,Vec);
 79: PETSC_EXTERN PetscErrorCode STComputeExplicitOperator(ST,Mat*);
 80: PETSC_EXTERN PetscErrorCode STMatSetUp(ST,PetscScalar,PetscScalar*);
 81: PETSC_EXTERN PetscErrorCode STPostSolve(ST);

 83: PETSC_EXTERN PetscErrorCode STSetKSP(ST,KSP);
 84: PETSC_EXTERN PetscErrorCode STGetKSP(ST,KSP*);
 85: PETSC_EXTERN PetscErrorCode STSetShift(ST,PetscScalar);
 86: PETSC_EXTERN PetscErrorCode STGetShift(ST,PetscScalar*);
 87: PETSC_EXTERN PetscErrorCode STSetDefaultShift(ST,PetscScalar);
 88: PETSC_EXTERN PetscErrorCode STScaleShift(ST,PetscScalar);
 89: PETSC_EXTERN PetscErrorCode STSetBalanceMatrix(ST,Vec);
 90: PETSC_EXTERN PetscErrorCode STGetBalanceMatrix(ST,Vec*);
 91: PETSC_EXTERN PetscErrorCode STSetTransform(ST,PetscBool);
 92: PETSC_EXTERN PetscErrorCode STGetTransform(ST,PetscBool*);

 94: PETSC_EXTERN PetscErrorCode STSetOptionsPrefix(ST,const char*);
 95: PETSC_EXTERN PetscErrorCode STAppendOptionsPrefix(ST,const char*);
 96: PETSC_EXTERN PetscErrorCode STGetOptionsPrefix(ST,const char*[]);

 98: PETSC_EXTERN PetscErrorCode STBackTransform(ST,PetscInt,PetscScalar*,PetscScalar*);

100: PETSC_EXTERN PetscErrorCode STCheckNullSpace(ST,BV);

102: PETSC_EXTERN PetscErrorCode STMatCreateVecs(ST,Vec*,Vec*);
103: PETSC_EXTERN PetscErrorCode STMatGetSize(ST,PetscInt*,PetscInt*);
104: PETSC_EXTERN PetscErrorCode STMatGetLocalSize(ST,PetscInt*,PetscInt*);

106: /*E
107:     STMatMode - Determines how to handle the coefficient matrix associated
108:     to the spectral transformation

110:     Level: intermediate

112: .seealso: STSetMatMode(), STGetMatMode()
113: E*/
114: typedef enum { ST_MATMODE_COPY,
115:                ST_MATMODE_INPLACE,
116:                ST_MATMODE_SHELL } STMatMode;
117: PETSC_EXTERN PetscErrorCode STSetMatMode(ST,STMatMode);
118: PETSC_EXTERN PetscErrorCode STGetMatMode(ST,STMatMode*);
119: PETSC_EXTERN PetscErrorCode STSetMatStructure(ST,MatStructure);
120: PETSC_EXTERN PetscErrorCode STGetMatStructure(ST,MatStructure*);

122: PETSC_EXTERN PetscFunctionList STList;
123: PETSC_EXTERN PetscErrorCode STRegister(const char[],PetscErrorCode(*)(ST));

125: /* --------- options specific to particular spectral transformations-------- */

127: PETSC_EXTERN PetscErrorCode STShellGetContext(ST st,void **ctx);
128: PETSC_EXTERN PetscErrorCode STShellSetContext(ST st,void *ctx);
129: PETSC_EXTERN PetscErrorCode STShellSetApply(ST st,PetscErrorCode (*apply)(ST,Vec,Vec));
130: PETSC_EXTERN PetscErrorCode STShellSetApplyTranspose(ST st,PetscErrorCode (*applytrans)(ST,Vec,Vec));
131: PETSC_EXTERN PetscErrorCode STShellSetBackTransform(ST st,PetscErrorCode (*backtr)(ST,PetscInt,PetscScalar*,PetscScalar*));

133: PETSC_EXTERN PetscErrorCode STCayleyGetAntishift(ST,PetscScalar*);
134: PETSC_EXTERN PetscErrorCode STCayleySetAntishift(ST,PetscScalar);

136: PETSC_EXTERN PetscErrorCode STPrecondGetMatForPC(ST,Mat*);
137: PETSC_EXTERN PetscErrorCode STPrecondSetMatForPC(ST,Mat);
138: PETSC_EXTERN PetscErrorCode STPrecondGetKSPHasMat(ST,PetscBool*);
139: PETSC_EXTERN PetscErrorCode STPrecondSetKSPHasMat(ST,PetscBool);

141: #endif