Actual source code: slepcst.h
1: /*
2: Spectral transformation module for eigenvalue problems.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2009, Universidad Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
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 "petscksp.h"
31: /*S
32: ST - Abstract SLEPc object that manages spectral transformations.
33: This object is accessed only in advanced applications.
35: Level: beginner
37: .seealso: STCreate(), EPS
38: S*/
39: typedef struct _p_ST* ST;
41: /*E
42: STType - String with the name of a SLEPc spectral transformation
44: Level: beginner
46: .seealso: STSetType(), ST
47: E*/
48: #define STType char*
49: #define STSHELL "shell"
50: #define STSHIFT "shift"
51: #define STSINV "sinvert"
52: #define STCAYLEY "cayley"
53: #define STFOLD "fold"
55: EXTERN PetscErrorCode STCreate(MPI_Comm,ST*);
56: EXTERN PetscErrorCode STDestroy(ST);
57: EXTERN PetscErrorCode STSetType(ST,const STType);
58: EXTERN PetscErrorCode STGetType(ST,const STType*);
59: EXTERN PetscErrorCode STSetOperators(ST,Mat,Mat);
60: EXTERN PetscErrorCode STGetOperators(ST,Mat*,Mat*);
61: EXTERN PetscErrorCode STSetUp(ST);
62: EXTERN PetscErrorCode STSetFromOptions(ST);
63: EXTERN PetscErrorCode STView(ST,PetscViewer);
65: EXTERN PetscErrorCode STApply(ST,Vec,Vec);
66: EXTERN PetscErrorCode STGetBilinearForm(ST,Mat*);
67: EXTERN PetscErrorCode STApplyTranspose(ST,Vec,Vec);
68: EXTERN PetscErrorCode STComputeExplicitOperator(ST,Mat*);
69: EXTERN PetscErrorCode STPostSolve(ST);
71: EXTERN PetscErrorCode STInitializePackage(char*);
73: EXTERN PetscErrorCode STSetKSP(ST,KSP);
74: EXTERN PetscErrorCode STGetKSP(ST,KSP*);
75: EXTERN PetscErrorCode STAssociatedKSPSolve(ST,Vec,Vec);
76: EXTERN PetscErrorCode STSetShift(ST,PetscScalar);
77: EXTERN PetscErrorCode STGetShift(ST,PetscScalar*);
79: EXTERN PetscErrorCode STSetOptionsPrefix(ST,const char*);
80: EXTERN PetscErrorCode STAppendOptionsPrefix(ST,const char*);
81: EXTERN PetscErrorCode STGetOptionsPrefix(ST,const char*[]);
83: EXTERN PetscErrorCode STBackTransform(ST,PetscScalar*,PetscScalar*);
85: EXTERN PetscErrorCode STCheckNullSpace(ST,PetscInt,const Vec[]);
87: EXTERN PetscErrorCode STGetOperationCounters(ST,PetscInt*,PetscInt*);
88: EXTERN PetscErrorCode STResetOperationCounters(ST);
90: /*E
91: STMatMode - determines how to handle the coefficient matrix associated
92: to the spectral transformation
94: Level: intermediate
96: .seealso: STSetMatMode(), STGetMatMode()
97: E*/
98: typedef enum { STMATMODE_COPY, STMATMODE_INPLACE,
99: STMATMODE_SHELL } STMatMode;
100: EXTERN PetscErrorCode STSetMatMode(ST,STMatMode);
101: EXTERN PetscErrorCode STGetMatMode(ST,STMatMode*);
102: EXTERN PetscErrorCode STSetMatStructure(ST,MatStructure);
104: EXTERN PetscErrorCode STRegister(const char*,const char*,const char*,PetscErrorCode(*)(ST));
105: #if defined(PETSC_USE_DYNAMIC_LIBRARIES)
106: #define STRegisterDynamic(a,b,c,d) STRegister(a,b,c,0)
107: #else
108: #define STRegisterDynamic(a,b,c,d) STRegister(a,b,c,d)
109: #endif
110: EXTERN PetscErrorCode STRegisterDestroy(void);
112: /* --------- options specific to particular spectral transformations-------- */
114: EXTERN PetscErrorCode STShellGetContext(ST st,void **ctx);
115: EXTERN PetscErrorCode STShellSetContext(ST st,void *ctx);
116: EXTERN PetscErrorCode STShellSetApply(ST st,PetscErrorCode (*apply)(void*,Vec,Vec));
117: EXTERN PetscErrorCode STShellSetApplyTranspose(ST st,PetscErrorCode (*applytrans)(void*,Vec,Vec));
118: EXTERN PetscErrorCode STShellSetBackTransform(ST st,PetscErrorCode (*backtr)(void*,PetscScalar*,PetscScalar*));
119: EXTERN PetscErrorCode STShellSetName(ST,const char[]);
120: EXTERN PetscErrorCode STShellGetName(ST,char*[]);
122: EXTERN PetscErrorCode STCayleySetAntishift(ST,PetscScalar);
124: EXTERN PetscErrorCode STFoldSetLeftSide(ST st,PetscTruth left);
127: #endif