Actual source code: slepcst.h

slepc-main 2024-12-17
Report Typos and Errors
  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: #pragma once

 16: #include <slepcsys.h>
 17: #include <slepcbv.h>
 18: #include <petscksp.h>

 20: /* SUBMANSEC = ST */

 22: SLEPC_EXTERN PetscErrorCode STInitializePackage(void);
 23: SLEPC_EXTERN PetscErrorCode STFinalizePackage(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(3, 15, 0, "STSetMatrices()", ) static inline PetscErrorCode STSetOperators(ST st,PetscInt n,Mat *A) {return STSetMatrices(st,n,A);}
 70: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "STGetMatrix()", ) static inline PetscErrorCode STGetOperators(ST st,PetscInt k,Mat *A) {return STGetMatrix(st,k,A);}
 71: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "STGetMatrixTransformed()", ) static inline PetscErrorCode STGetTOperators(ST st,PetscInt k,Mat *A) {return STGetMatrixTransformed(st,k,A);}
 72: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "STGetOperator() followed by MatComputeOperator()", ) static inline PetscErrorCode STComputeExplicitOperator(ST st,Mat *A)
 73: {
 74:   Mat Op;

 76:   PetscFunctionBegin;
 77:   PetscCall(STGetOperator(st,&Op));
 78:   PetscCall(MatComputeOperator(Op,MATAIJ,A));
 79:   PetscCall(STRestoreOperator(st,&Op));
 80:   PetscFunctionReturn(PETSC_SUCCESS);
 81: }

 83: SLEPC_EXTERN PetscErrorCode STApply(ST,Vec,Vec);
 84: SLEPC_EXTERN PetscErrorCode STApplyMat(ST,Mat,Mat);
 85: SLEPC_EXTERN PetscErrorCode STApplyTranspose(ST,Vec,Vec);
 86: SLEPC_EXTERN PetscErrorCode STApplyHermitianTranspose(ST,Vec,Vec);
 87: SLEPC_EXTERN PetscErrorCode STMatMult(ST,PetscInt,Vec,Vec);
 88: SLEPC_EXTERN PetscErrorCode STMatMultTranspose(ST,PetscInt,Vec,Vec);
 89: SLEPC_EXTERN PetscErrorCode STMatMultHermitianTranspose(ST,PetscInt,Vec,Vec);
 90: SLEPC_EXTERN PetscErrorCode STMatSolve(ST,Vec,Vec);
 91: SLEPC_EXTERN PetscErrorCode STMatSolveTranspose(ST,Vec,Vec);
 92: SLEPC_EXTERN PetscErrorCode STMatSolveHermitianTranspose(ST,Vec,Vec);
 93: SLEPC_EXTERN PetscErrorCode STMatMatSolve(ST,Mat,Mat);
 94: SLEPC_EXTERN PetscErrorCode STGetBilinearForm(ST,Mat*);
 95: SLEPC_EXTERN PetscErrorCode STMatSetUp(ST,PetscScalar,PetscScalar*);
 96: SLEPC_EXTERN PetscErrorCode STPostSolve(ST);
 97: SLEPC_EXTERN PetscErrorCode STResetMatrixState(ST);
 98: SLEPC_EXTERN PetscErrorCode STSetWorkVecs(ST,PetscInt);

100: SLEPC_EXTERN PetscErrorCode STSetKSP(ST,KSP);
101: SLEPC_EXTERN PetscErrorCode STGetKSP(ST,KSP*);
102: SLEPC_EXTERN PetscErrorCode STSetShift(ST,PetscScalar);
103: SLEPC_EXTERN PetscErrorCode STGetShift(ST,PetscScalar*);
104: SLEPC_EXTERN PetscErrorCode STSetDefaultShift(ST,PetscScalar);
105: SLEPC_EXTERN PetscErrorCode STScaleShift(ST,PetscScalar);
106: SLEPC_EXTERN PetscErrorCode STSetBalanceMatrix(ST,Vec);
107: SLEPC_EXTERN PetscErrorCode STGetBalanceMatrix(ST,Vec*);
108: SLEPC_EXTERN PetscErrorCode STSetTransform(ST,PetscBool);
109: SLEPC_EXTERN PetscErrorCode STGetTransform(ST,PetscBool*);
110: SLEPC_EXTERN PetscErrorCode STSetStructured(ST,PetscBool);
111: SLEPC_EXTERN PetscErrorCode STGetStructured(ST,PetscBool*);

113: SLEPC_EXTERN PetscErrorCode STSetOptionsPrefix(ST,const char*);
114: SLEPC_EXTERN PetscErrorCode STAppendOptionsPrefix(ST,const char*);
115: SLEPC_EXTERN PetscErrorCode STGetOptionsPrefix(ST,const char*[]);

117: SLEPC_EXTERN PetscErrorCode STBackTransform(ST,PetscInt,PetscScalar*,PetscScalar*);
118: SLEPC_EXTERN PetscErrorCode STIsInjective(ST,PetscBool*);

120: SLEPC_EXTERN PetscErrorCode STCheckNullSpace(ST,BV);

122: SLEPC_EXTERN PetscErrorCode STSetPreconditionerMat(ST,Mat);
123: SLEPC_EXTERN PetscErrorCode STGetPreconditionerMat(ST,Mat*);
124: SLEPC_EXTERN PetscErrorCode STSetSplitPreconditioner(ST,PetscInt,Mat[],MatStructure);
125: SLEPC_EXTERN PetscErrorCode STGetSplitPreconditionerTerm(ST,PetscInt,Mat*);
126: SLEPC_EXTERN PetscErrorCode STGetSplitPreconditionerInfo(ST,PetscInt*,MatStructure*);

128: SLEPC_EXTERN PetscErrorCode STMatCreateVecs(ST,Vec*,Vec*);
129: SLEPC_EXTERN PetscErrorCode STMatCreateVecsEmpty(ST,Vec*,Vec*);
130: SLEPC_EXTERN PetscErrorCode STMatGetSize(ST,PetscInt*,PetscInt*);
131: SLEPC_EXTERN PetscErrorCode STMatGetLocalSize(ST,PetscInt*,PetscInt*);

133: /*E
134:     STMatMode - Determines how to handle the coefficient matrix associated
135:     to the spectral transformation

137:     Level: intermediate

139: .seealso: STSetMatMode(), STGetMatMode()
140: E*/
141: typedef enum { ST_MATMODE_COPY,
142:                ST_MATMODE_INPLACE,
143:                ST_MATMODE_SHELL } STMatMode;
144: SLEPC_EXTERN const char *STMatModes[];

146: SLEPC_EXTERN PetscErrorCode STSetMatMode(ST,STMatMode);
147: SLEPC_EXTERN PetscErrorCode STGetMatMode(ST,STMatMode*);
148: SLEPC_EXTERN PetscErrorCode STSetMatStructure(ST,MatStructure);
149: SLEPC_EXTERN PetscErrorCode STGetMatStructure(ST,MatStructure*);

151: SLEPC_EXTERN PetscFunctionList STList;
152: SLEPC_EXTERN PetscErrorCode STRegister(const char[],PetscErrorCode(*)(ST));

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

156: SLEPC_EXTERN PetscErrorCode STShellGetContext(ST,void*);
157: SLEPC_EXTERN PetscErrorCode STShellSetContext(ST,void*);
158: SLEPC_EXTERN PetscErrorCode STShellSetApply(ST,PetscErrorCode (*)(ST,Vec,Vec));
159: SLEPC_EXTERN PetscErrorCode STShellSetApplyTranspose(ST,PetscErrorCode (*)(ST,Vec,Vec));
160: SLEPC_EXTERN PetscErrorCode STShellSetApplyHermitianTranspose(ST,PetscErrorCode (*)(ST,Vec,Vec));
161: SLEPC_EXTERN PetscErrorCode STShellSetBackTransform(ST,PetscErrorCode (*)(ST,PetscInt,PetscScalar*,PetscScalar*));

163: SLEPC_EXTERN PetscErrorCode STCayleyGetAntishift(ST,PetscScalar*);
164: SLEPC_EXTERN PetscErrorCode STCayleySetAntishift(ST,PetscScalar);

166: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "STGetPreconditionerMat()", ) static inline PetscErrorCode STPrecondGetMatForPC(ST st,Mat *A) {return STGetPreconditionerMat(st,A);}
167: PETSC_DEPRECATED_FUNCTION(3, 15, 0, "STSetPreconditionerMat()", ) static inline PetscErrorCode STPrecondSetMatForPC(ST st,Mat A) {return STSetPreconditionerMat(st,A);}
168: SLEPC_EXTERN PetscErrorCode STPrecondGetKSPHasMat(ST,PetscBool*);
169: SLEPC_EXTERN PetscErrorCode STPrecondSetKSPHasMat(ST,PetscBool);

171: SLEPC_EXTERN PetscErrorCode STFilterSetInterval(ST,PetscReal,PetscReal);
172: SLEPC_EXTERN PetscErrorCode STFilterGetInterval(ST,PetscReal*,PetscReal*);
173: SLEPC_EXTERN PetscErrorCode STFilterSetRange(ST,PetscReal,PetscReal);
174: SLEPC_EXTERN PetscErrorCode STFilterGetRange(ST,PetscReal*,PetscReal*);
175: SLEPC_EXTERN PetscErrorCode STFilterSetDegree(ST,PetscInt);
176: SLEPC_EXTERN PetscErrorCode STFilterGetDegree(ST,PetscInt*);
177: SLEPC_EXTERN PetscErrorCode STFilterGetThreshold(ST,PetscReal*);