Actual source code: slepcst.h
slepc-3.18.3 2023-03-24
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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: */
10: /*
11: Spectral transformation module for eigenvalue problems
12: */
14: #if !defined(SLEPCST_H)
15: #define SLEPCST_H
17: #include <slepcsys.h>
18: #include <slepcbv.h>
19: #include <petscksp.h>
21: /* SUBMANSEC = ST */
23: SLEPC_EXTERN PetscErrorCode STInitializePackage(void);
25: /*S
26: ST - Abstract SLEPc object that manages spectral transformations.
27: This object is accessed only in advanced applications.
29: Level: beginner
31: .seealso: STCreate(), EPS
32: S*/
33: typedef struct _p_ST* ST;
35: /*J
36: STType - String with the name of a SLEPc spectral transformation
38: Level: beginner
40: .seealso: STSetType(), ST
41: J*/
42: typedef const char* STType;
43: #define STSHELL "shell"
44: #define STSHIFT "shift"
45: #define STSINVERT "sinvert"
46: #define STCAYLEY "cayley"
47: #define STPRECOND "precond"
48: #define STFILTER "filter"
50: /* Logging support */
51: SLEPC_EXTERN PetscClassId ST_CLASSID;
53: SLEPC_EXTERN PetscErrorCode STCreate(MPI_Comm,ST*);
54: SLEPC_EXTERN PetscErrorCode STDestroy(ST*);
55: SLEPC_EXTERN PetscErrorCode STReset(ST);
56: SLEPC_EXTERN PetscErrorCode STSetType(ST,STType);
57: SLEPC_EXTERN PetscErrorCode STGetType(ST,STType*);
58: SLEPC_EXTERN PetscErrorCode STSetMatrices(ST,PetscInt,Mat*);
59: SLEPC_EXTERN PetscErrorCode STGetMatrix(ST,PetscInt,Mat*);
60: SLEPC_EXTERN PetscErrorCode STGetMatrixTransformed(ST,PetscInt,Mat*);
61: SLEPC_EXTERN PetscErrorCode STGetNumMatrices(ST,PetscInt*);
62: SLEPC_EXTERN PetscErrorCode STGetOperator(ST,Mat*);
63: SLEPC_EXTERN PetscErrorCode STRestoreOperator(ST,Mat*);
64: SLEPC_EXTERN PetscErrorCode STSetUp(ST);
65: SLEPC_EXTERN PetscErrorCode STSetFromOptions(ST);
66: SLEPC_EXTERN PetscErrorCode STView(ST,PetscViewer);
67: SLEPC_EXTERN PetscErrorCode STViewFromOptions(ST,PetscObject,const char[]);
69: PETSC_DEPRECATED_FUNCTION("Use STSetMatrices()") static inline PetscErrorCode STSetOperators(ST st,PetscInt n,Mat *A) {return STSetMatrices(st,n,A);}
70: PETSC_DEPRECATED_FUNCTION("Use STGetMatrix()") static inline PetscErrorCode STGetOperators(ST st,PetscInt k,Mat *A) {return STGetMatrix(st,k,A);}
71: PETSC_DEPRECATED_FUNCTION("Use STGetMatrixTransformed()") static inline PetscErrorCode STGetTOperators(ST st,PetscInt k,Mat *A) {return STGetMatrixTransformed(st,k,A);}
72: PETSC_DEPRECATED_FUNCTION("Use STGetOperator() followed by MatComputeOperator()") static inline PetscErrorCode STComputeExplicitOperator(ST st,Mat *A)
73: {
74: Mat Op;
75: STGetOperator(st,&Op);
76: MatComputeOperator(Op,MATAIJ,A);
77: STRestoreOperator(st,&Op);
78: return 0;
79: }
81: SLEPC_EXTERN PetscErrorCode STApply(ST,Vec,Vec);
82: SLEPC_EXTERN PetscErrorCode STApplyMat(ST,Mat,Mat);
83: SLEPC_EXTERN PetscErrorCode STApplyTranspose(ST,Vec,Vec);
84: SLEPC_EXTERN PetscErrorCode STApplyHermitianTranspose(ST,Vec,Vec);
85: SLEPC_EXTERN PetscErrorCode STMatMult(ST,PetscInt,Vec,Vec);
86: SLEPC_EXTERN PetscErrorCode STMatMultTranspose(ST,PetscInt,Vec,Vec);
87: SLEPC_EXTERN PetscErrorCode STMatSolve(ST,Vec,Vec);
88: SLEPC_EXTERN PetscErrorCode STMatSolveTranspose(ST,Vec,Vec);
89: SLEPC_EXTERN PetscErrorCode STMatMatSolve(ST,Mat,Mat);
90: SLEPC_EXTERN PetscErrorCode STGetBilinearForm(ST,Mat*);
91: SLEPC_EXTERN PetscErrorCode STMatSetUp(ST,PetscScalar,PetscScalar*);
92: SLEPC_EXTERN PetscErrorCode STPostSolve(ST);
93: SLEPC_EXTERN PetscErrorCode STResetMatrixState(ST);
94: SLEPC_EXTERN PetscErrorCode STSetWorkVecs(ST,PetscInt);
96: SLEPC_EXTERN PetscErrorCode STSetKSP(ST,KSP);
97: SLEPC_EXTERN PetscErrorCode STGetKSP(ST,KSP*);
98: SLEPC_EXTERN PetscErrorCode STSetShift(ST,PetscScalar);
99: SLEPC_EXTERN PetscErrorCode STGetShift(ST,PetscScalar*);
100: SLEPC_EXTERN PetscErrorCode STSetDefaultShift(ST,PetscScalar);
101: SLEPC_EXTERN PetscErrorCode STScaleShift(ST,PetscScalar);
102: SLEPC_EXTERN PetscErrorCode STSetBalanceMatrix(ST,Vec);
103: SLEPC_EXTERN PetscErrorCode STGetBalanceMatrix(ST,Vec*);
104: SLEPC_EXTERN PetscErrorCode STSetTransform(ST,PetscBool);
105: SLEPC_EXTERN PetscErrorCode STGetTransform(ST,PetscBool*);
107: SLEPC_EXTERN PetscErrorCode STSetOptionsPrefix(ST,const char*);
108: SLEPC_EXTERN PetscErrorCode STAppendOptionsPrefix(ST,const char*);
109: SLEPC_EXTERN PetscErrorCode STGetOptionsPrefix(ST,const char*[]);
111: SLEPC_EXTERN PetscErrorCode STBackTransform(ST,PetscInt,PetscScalar*,PetscScalar*);
112: SLEPC_EXTERN PetscErrorCode STIsInjective(ST,PetscBool*);
114: SLEPC_EXTERN PetscErrorCode STCheckNullSpace(ST,BV);
116: SLEPC_EXTERN PetscErrorCode STSetPreconditionerMat(ST,Mat);
117: SLEPC_EXTERN PetscErrorCode STGetPreconditionerMat(ST,Mat*);
118: SLEPC_EXTERN PetscErrorCode STSetSplitPreconditioner(ST,PetscInt,Mat[],MatStructure);
119: SLEPC_EXTERN PetscErrorCode STGetSplitPreconditionerTerm(ST,PetscInt,Mat*);
120: SLEPC_EXTERN PetscErrorCode STGetSplitPreconditionerInfo(ST,PetscInt*,MatStructure*);
122: SLEPC_EXTERN PetscErrorCode STMatCreateVecs(ST,Vec*,Vec*);
123: SLEPC_EXTERN PetscErrorCode STMatCreateVecsEmpty(ST,Vec*,Vec*);
124: SLEPC_EXTERN PetscErrorCode STMatGetSize(ST,PetscInt*,PetscInt*);
125: SLEPC_EXTERN PetscErrorCode STMatGetLocalSize(ST,PetscInt*,PetscInt*);
127: /*E
128: STMatMode - Determines how to handle the coefficient matrix associated
129: to the spectral transformation
131: Level: intermediate
133: .seealso: STSetMatMode(), STGetMatMode()
134: E*/
135: typedef enum { ST_MATMODE_COPY,
136: ST_MATMODE_INPLACE,
137: ST_MATMODE_SHELL } STMatMode;
138: SLEPC_EXTERN const char *STMatModes[];
140: SLEPC_EXTERN PetscErrorCode STSetMatMode(ST,STMatMode);
141: SLEPC_EXTERN PetscErrorCode STGetMatMode(ST,STMatMode*);
142: SLEPC_EXTERN PetscErrorCode STSetMatStructure(ST,MatStructure);
143: SLEPC_EXTERN PetscErrorCode STGetMatStructure(ST,MatStructure*);
145: SLEPC_EXTERN PetscFunctionList STList;
146: SLEPC_EXTERN PetscErrorCode STRegister(const char[],PetscErrorCode(*)(ST));
148: /* --------- options specific to particular spectral transformations-------- */
150: SLEPC_EXTERN PetscErrorCode STShellGetContext(ST,void*);
151: SLEPC_EXTERN PetscErrorCode STShellSetContext(ST,void*);
152: SLEPC_EXTERN PetscErrorCode STShellSetApply(ST,PetscErrorCode (*)(ST,Vec,Vec));
153: SLEPC_EXTERN PetscErrorCode STShellSetApplyTranspose(ST,PetscErrorCode (*)(ST,Vec,Vec));
154: SLEPC_EXTERN PetscErrorCode STShellSetBackTransform(ST,PetscErrorCode (*)(ST,PetscInt,PetscScalar*,PetscScalar*));
156: SLEPC_EXTERN PetscErrorCode STCayleyGetAntishift(ST,PetscScalar*);
157: SLEPC_EXTERN PetscErrorCode STCayleySetAntishift(ST,PetscScalar);
159: PETSC_DEPRECATED_FUNCTION("Use STGetPreconditionerMat()") static inline PetscErrorCode STPrecondGetMatForPC(ST st,Mat *A) {return STGetPreconditionerMat(st,A);}
160: PETSC_DEPRECATED_FUNCTION("Use STSetPreconditionerMat()") static inline PetscErrorCode STPrecondSetMatForPC(ST st,Mat A) {return STSetPreconditionerMat(st,A);}
161: SLEPC_EXTERN PetscErrorCode STPrecondGetKSPHasMat(ST,PetscBool*);
162: SLEPC_EXTERN PetscErrorCode STPrecondSetKSPHasMat(ST,PetscBool);
164: SLEPC_EXTERN PetscErrorCode STFilterSetInterval(ST,PetscReal,PetscReal);
165: SLEPC_EXTERN PetscErrorCode STFilterGetInterval(ST,PetscReal*,PetscReal*);
166: SLEPC_EXTERN PetscErrorCode STFilterSetRange(ST,PetscReal,PetscReal);
167: SLEPC_EXTERN PetscErrorCode STFilterGetRange(ST,PetscReal*,PetscReal*);
168: SLEPC_EXTERN PetscErrorCode STFilterSetDegree(ST,PetscInt);
169: SLEPC_EXTERN PetscErrorCode STFilterGetDegree(ST,PetscInt*);
170: SLEPC_EXTERN PetscErrorCode STFilterGetThreshold(ST,PetscReal*);
172: #endif