1f26ada1bSBarry Smith /* 2f26ada1bSBarry Smith Defines the interface functions for the Krylov subspace accelerators. 3f26ada1bSBarry Smith */ 426bd1501SBarry Smith #ifndef PETSCKSP_H 526bd1501SBarry Smith #define PETSCKSP_H 6ac09b921SBarry Smith 72c8e378dSBarry Smith #include <petscpc.h> 82eac72dbSBarry Smith 9ac09b921SBarry Smith /* SUBMANSEC = KSP */ 10ac09b921SBarry Smith 11607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void); 121dbb0a54SBarry Smith 1328ce4d24SBarry Smith /*S 148f6c3df8SBarry Smith KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the 158f6c3df8SBarry Smith linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators). 1628ce4d24SBarry Smith 1728ce4d24SBarry Smith Level: beginner 1828ce4d24SBarry Smith 1987497f52SBarry Smith Note: 20a4d1885cSBarry Smith When a direct solver is used, but no Krylov solver is used, the `KSP` object is still used but with a 2187497f52SBarry Smith `KSPType` of `KSPPREONLY`, meaning that only application of the preconditioner is used as the linear solver. 228f6c3df8SBarry Smith 23a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES` 2428ce4d24SBarry Smith S*/ 25e2a1c21fSSatish Balay typedef struct _p_KSP *KSP; 262eac72dbSBarry Smith 2776bdecfbSBarry Smith /*J 288f6c3df8SBarry Smith KSPType - String with the name of a PETSc Krylov method. 2928ce4d24SBarry Smith 3028ce4d24SBarry Smith Level: beginner 3128ce4d24SBarry Smith 32a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSetType()`, `KSP`, `KSPRegister()`, `KSPCreate()`, `KSPSetFromOptions()` 3376bdecfbSBarry Smith J*/ 3419fd82e9SBarry Smith typedef const char *KSPType; 3582bf6240SBarry Smith #define KSPRICHARDSON "richardson" 366c9de887SHong Zhang #define KSPCHEBYSHEV "chebyshev" 3782bf6240SBarry Smith #define KSPCG "cg" 382c8fc646SJed Brown #define KSPGROPPCG "groppcg" 392c8fc646SJed Brown #define KSPPIPECG "pipecg" 40901ccb91SSiegfried Cools #define KSPPIPECGRR "pipecgrr" 41265663fdSSiegfried Cools #define KSPPIPELCG "pipelcg" 42b21a8899STyler Chen #define KSPPIPEPRCG "pipeprcg" 43325e8391SManasi #define KSPPIPECG2 "pipecg2" 44df2a969aSvictorle #define KSPCGNE "cgne" 4505de396fSBarry Smith #define KSPNASH "nash" 4605de396fSBarry Smith #define KSPSTCG "stcg" 4705de396fSBarry Smith #define KSPGLTR "gltr" 4805de396fSBarry Smith #define KSPCGNASH PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"") "nash" 4905de396fSBarry Smith #define KSPCGSTCG PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"") "stcg" 5005de396fSBarry Smith #define KSPCGGLTR PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr" 51640f4f53SPatrick Sanan #define KSPFCG "fcg" 52390d8e47SPatrick Sanan #define KSPPIPEFCG "pipefcg" 5382bf6240SBarry Smith #define KSPGMRES "gmres" 54483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres" 559dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 569dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 574f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 5861b468f9SJed Brown #define KSPPGMRES "pgmres" 5982bf6240SBarry Smith #define KSPTCQMR "tcqmr" 6082bf6240SBarry Smith #define KSPBCGS "bcgs" 61e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 62345ecf0bSXiangmin Jiao #define KSPQMRCGS "qmrcgs" 6318ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 64c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 65f0037002Svictorle #define KSPBCGSL "bcgsl" 66f154af2dSSiegfried Cools #define KSPPIPEBCGS "pipebcgs" 6782bf6240SBarry Smith #define KSPCGS "cgs" 6882bf6240SBarry Smith #define KSPTFQMR "tfqmr" 6982bf6240SBarry Smith #define KSPCR "cr" 702c8fc646SJed Brown #define KSPPIPECR "pipecr" 7182bf6240SBarry Smith #define KSPLSQR "lsqr" 7282bf6240SBarry Smith #define KSPPREONLY "preonly" 733c2be86cSBarry Smith #define KSPNONE "none" 7482bf6240SBarry Smith #define KSPQCG "qcg" 75c9cf9b11SBarry Smith #define KSPBICG "bicg" 76b4ac9ba4SBarry Smith #define KSPMINRES "minres" 7701c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 7821356919SSatish Balay #define KSPLCD "lcd" 791d6018f0SLisandro Dalcin #define KSPPYTHON "python" 8058ad63f4SBarry Smith #define KSPGCR "gcr" 81fad47a0aSPatrick Sanan #define KSPPIPEGCR "pipegcr" 82e4d80e07Szianekhodja #define KSPTSIRM "tsirm" 83e4d80e07Szianekhodja #define KSPCGLS "cgls" 84329cd39dSStefano Zampini #define KSPFETIDP "fetidp" 8538cfc46eSPierre Jolivet #define KSPHPDDM "hpddm" 862eac72dbSBarry Smith 878ba1e511SMatthew Knepley /* Logging support */ 88014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 89ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID; 905358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 918ba1e511SMatthew Knepley 92014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm, KSP *); 9319fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP, KSPType); 94c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP, KSPType *); 95014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 97014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP, Vec, Vec); 98014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP, Vec, Vec); 9975f8e9bdSHong Zhang PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP, PetscBool); 1005ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP, Mat, Mat); 101*bbbebc2cSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolveTranspose(KSP, Mat, Mat); 1023e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt); 103d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPSetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n) 104d71ae5a4SJacob Faibussowitsch { 1059371c9d4SSatish Balay return KSPSetMatSolveBatchSize(ksp, n); 1069371c9d4SSatish Balay } 1073e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *); 108d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n) 109d71ae5a4SJacob Faibussowitsch { 1109371c9d4SSatish Balay return KSPGetMatSolveBatchSize(ksp, n); 1119371c9d4SSatish Balay } 112014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 113ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP); 114014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *); 11523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool); 1167d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *); 11725c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool); 118c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec); 1192eac72dbSBarry Smith 120140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 121ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList; 122798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorList; 123798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorCreateList; 124798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList; 125bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[], PetscErrorCode (*)(KSP)); 126798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **)); 12730de9b25SBarry Smith 128014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide); 129014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *); 130014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt); 131c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *); 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool); 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *); 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool); 135014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *); 136014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool); 1377d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool); 138c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *); 139014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool); 140c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *); 141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *); 142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *); 143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *); 144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *); 145734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *); 1462a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **); 147d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y) 148d71ae5a4SJacob Faibussowitsch { 1499371c9d4SSatish Balay return KSPCreateVecs(ksp, n, x, m, y); 1509371c9d4SSatish Balay } 1512eac72dbSBarry Smith 152d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *); 153d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *); 154d4211eb9SBarry Smith 155014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC); 156014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *); 157aabeff55SBarry Smith 158014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal); 159014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void **)); 160014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 1613ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *); 16294a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *); 163014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscInt, PetscBool); 16494a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *); 16594a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscInt, PetscBool); 1664b0e389bSBarry Smith 167fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *); 168fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *); 169fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 170fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt); 171fa0ddf94SBarry Smith 172014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *); 173cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP); 174014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 175014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 177014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]); 178285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]); 179b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *); 180b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *); 181b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *); 182b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *); 1844a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *); 185f3b08a26SMatthew G. Knepley /* 186f3b08a26SMatthew G. Knepley PCMGCoarseList contains the list of coarse space constructor currently registered 187f3b08a26SMatthew G. Knepley These are added with PCMGRegisterCoarseSpaceConstructor() 188f3b08a26SMatthew G. Knepley */ 189f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList; 1902b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *)); 1912b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *)); 192f3b08a26SMatthew G. Knepley 193014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *); 194014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *); 1952eac72dbSBarry Smith 196f2edd1f0SMalachi Phillips /*E 197f2edd1f0SMalachi Phillips KSPChebyshevKind - Which Chebyshev polynomial to use 198f2edd1f0SMalachi Phillips 199f2edd1f0SMalachi Phillips Values: 200f2edd1f0SMalachi Phillips + `KSP_CHEBYSHEV_FIRST` - "classic" first-kind Chebyshev polynomial 201f2edd1f0SMalachi Phillips . `KSP_CHEBYSHEV_FOURTH` - fourth-kind Chebyshev polynomial 202f2edd1f0SMalachi Phillips - `KSP_CHEBYSHEV_OPT_FOURTH` - optimized fourth-kind Chebyshev polynomial 203f2edd1f0SMalachi Phillips 204f2edd1f0SMalachi Phillips Level: intermediate 205f2edd1f0SMalachi Phillips 206f2edd1f0SMalachi Phillips .seealso: [](chapter_ksp), `KSPCHEBYSHEV`, `KSPChebyshevSetKind` 207f2edd1f0SMalachi Phillips E*/ 208f2edd1f0SMalachi Phillips typedef enum { 209f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_FIRST, 210f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_FOURTH, 211f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_OPT_FOURTH 212f2edd1f0SMalachi Phillips } KSPChebyshevKind; 213f2edd1f0SMalachi Phillips 214014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal); 215014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool); 216014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal); 21758450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal); 218b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool); 219f2edd1f0SMalachi Phillips PETSC_EXTERN PetscErrorCode KSPChebyshevSetKind(KSP, KSPChebyshevKind); 22058450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *); 221014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *); 222d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *); 223d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]); 2247d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]); 2254b0e389bSBarry Smith 226640f4f53SPatrick Sanan /*E 227640f4f53SPatrick Sanan 22806137d0aSPatrick Sanan KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 229640f4f53SPatrick Sanan 23067b8a455SSatish Balay Values: 231a1cb98faSBarry Smith + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions 232a1cb98faSBarry Smith - `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 233640f4f53SPatrick Sanan 2342b26979fSBarry Smith Level: intermediate 235640f4f53SPatrick Sanan 236a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()` 237640f4f53SPatrick Sanan E*/ 2389371c9d4SSatish Balay typedef enum { 2399371c9d4SSatish Balay KSP_FCD_TRUNC_TYPE_STANDARD, 2409371c9d4SSatish Balay KSP_FCD_TRUNC_TYPE_NOTAY 2419371c9d4SSatish Balay } KSPFCDTruncationType; 24206137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 243640f4f53SPatrick Sanan 244640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt); 245640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *); 246640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt); 247640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *); 24806137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType); 24906137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *); 250640f4f53SPatrick Sanan 251390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt); 252390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *); 253390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt); 254390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *); 255390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType); 256390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *); 257390d8e47SPatrick Sanan 258fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt); 259fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *); 260fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt); 261fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *); 262fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType); 263fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *); 264fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool); 265fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *); 26683f0b094SBarry Smith 267014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 268014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *); 269014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal); 270e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal); 2719f236934SBarry Smith 272014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 273014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt)); 274014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt)); 275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt); 276014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt); 2771d73ed98SBarry Smith 278014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt); 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2801d73ed98SBarry Smith 281483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar); 282483d6965SPatrick Sanan 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt); 284014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *); 285014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 28658ad63f4SBarry Smith 287f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal); 288f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *); 289f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool); 290f87a0b54SStefano Zampini 29104d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *); 29204d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC); 29304d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *); 294568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat); 295e82af88dSprj- 2966cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat); 2976cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *); 2986cac28cbSPierre Jolivet #if PetscDefined(HAVE_HPDDM) 299d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMSetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) 300d71ae5a4SJacob Faibussowitsch { 3019371c9d4SSatish Balay return KSPHPDDMSetDeflationMat(ksp, U); 3029371c9d4SSatish Balay } 303d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMGetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) 304d71ae5a4SJacob Faibussowitsch { 3059371c9d4SSatish Balay return KSPHPDDMGetDeflationMat(ksp, U); 3069371c9d4SSatish Balay } 3076cac28cbSPierre Jolivet #endif 308d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) 309d71ae5a4SJacob Faibussowitsch { 3109371c9d4SSatish Balay return KSPMatSolve(ksp, B, X); 3119371c9d4SSatish Balay } 312d55ab951SPierre Jolivet /*E 31387497f52SBarry Smith KSPHPDDMType - Type of Krylov method used by `KSPHPDDM` 314d55ab951SPierre Jolivet 315d55ab951SPierre Jolivet Level: intermediate 316d55ab951SPierre Jolivet 317d55ab951SPierre Jolivet Values: 318a4d1885cSBarry Smith + `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method 319a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BGMRES` - block GMRES 320a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_CG` - Conjugate Gradient 321a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BCG` - block CG 322a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_GCRODR` - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting 323a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BGCRODR` - block GCRODR 324a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BFBCG` - breakdown-free BCG 325a4d1885cSBarry Smith - `KSP_HPDDM_TYPE_PREONLY` - apply the preconditioner only 326d55ab951SPierre Jolivet 327a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPHPDDM`, `KSPHPDDMSetType()` 328d55ab951SPierre Jolivet E*/ 329d55ab951SPierre Jolivet typedef enum { 330d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GMRES = 0, 331d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGMRES = 1, 332d55ab951SPierre Jolivet KSP_HPDDM_TYPE_CG = 2, 333d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BCG = 3, 334d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GCRODR = 4, 335d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGCRODR = 5, 336d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BFBCG = 6, 337d55ab951SPierre Jolivet KSP_HPDDM_TYPE_PREONLY = 7 338d55ab951SPierre Jolivet } KSPHPDDMType; 339d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[]; 340a4d1885cSBarry Smith 3412dd49c90SPierre Jolivet /*E 34287497f52SBarry Smith KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM` 3432dd49c90SPierre Jolivet 3442dd49c90SPierre Jolivet Level: intermediate 3452dd49c90SPierre Jolivet 3462dd49c90SPierre Jolivet Values: 347a4d1885cSBarry Smith + `KSP_HPDDM_PRECISION_HALF` - default when PETSc is configured `--with-precision=__fp16` 348a4d1885cSBarry Smith . `KSP_HPDDM_PRECISION_SINGLE` - default when PETSc is configured `--with-precision=single` 349a4d1885cSBarry Smith . `KSP_HPDDM_PRECISION_DOUBLE` - default when PETSc is configured `--with-precision=double` 350a4d1885cSBarry Smith - `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128` 3512dd49c90SPierre Jolivet 352a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSP`, `KSPHPDDM` 3532dd49c90SPierre Jolivet E*/ 3542dd49c90SPierre Jolivet typedef enum { 3552dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_HALF = 0, 3562dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_SINGLE = 1, 3572dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_DOUBLE = 2, 3582dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_QUADRUPLE = 3 3592dd49c90SPierre Jolivet } KSPHPDDMPrecision; 360d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType); 361d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *); 362e82af88dSprj- 363b49fd9e1SBarry Smith /*E 364b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 365b49fd9e1SBarry Smith 366b49fd9e1SBarry Smith Level: advanced 367b49fd9e1SBarry Smith 368a4d1885cSBarry Smith Values: 369a4d1885cSBarry Smith + `KSP_GMRES_CGS_REFINE_NEVER` - one step of classical Gram-Schmidt 370a4d1885cSBarry Smith . `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria 371a4d1885cSBarry Smith - `KSP_GMRES_CGS_REFINE_ALWAYS` - always perform two steps 372b49fd9e1SBarry Smith 373a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSP`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 374a4d1885cSBarry Smith `KSPGMRESGetOrthogonalization()`, 375a4d1885cSBarry Smith `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()` 376b49fd9e1SBarry Smith E*/ 3779371c9d4SSatish Balay typedef enum { 3789371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_NEVER, 3799371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_IFNEEDED, 3809371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_ALWAYS 3819371c9d4SSatish Balay } KSPGMRESCGSRefinementType; 3826a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 3831f7e983dSSatish Balay /*MC 3841957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 3858c5b8ba0SBarry Smith 3868c5b8ba0SBarry Smith Level: advanced 3878c5b8ba0SBarry Smith 38887497f52SBarry Smith Note: 38987497f52SBarry Smith Possibly unstable, but the fastest to compute 3908c5b8ba0SBarry Smith 391a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 392a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 393db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 394db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 3958c5b8ba0SBarry Smith M*/ 3968c5b8ba0SBarry Smith 3971f7e983dSSatish Balay /*MC 3988c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 3998c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 4008c5b8ba0SBarry Smith poor orthogonality. 4018c5b8ba0SBarry Smith 4028c5b8ba0SBarry Smith Level: advanced 4038c5b8ba0SBarry Smith 404a4d1885cSBarry Smith Note: 405a4d1885cSBarry Smith This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to 4068c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 4078c5b8ba0SBarry Smith 408a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 409a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 410db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 411db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 4128c5b8ba0SBarry Smith M*/ 4138c5b8ba0SBarry Smith 4141f7e983dSSatish Balay /*MC 4158c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 4168c5b8ba0SBarry Smith 4178c5b8ba0SBarry Smith Level: advanced 4188c5b8ba0SBarry Smith 41987497f52SBarry Smith Notes: 42087497f52SBarry Smith This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice 42187497f52SBarry Smith but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`. 4228c5b8ba0SBarry Smith 4238c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 4248c5b8ba0SBarry Smith 425a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 426a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 427db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 428db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 4298c5b8ba0SBarry Smith M*/ 4308c5b8ba0SBarry Smith 431014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType); 432014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *); 43308480c60SBarry Smith 434014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *); 435014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *); 436014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 437c38d4ed2SBarry Smith 438014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal); 439014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *); 440014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *); 441121fd945SKris Buschelman 442014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal); 443014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool); 444014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt); 445e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool); 446d9492815SBarry Smith 447014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 44887d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP); 449014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 4502eac72dbSBarry Smith 451798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *); 452798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm, const char[], const char[], const char[], PetscInt, const char *[], int, int, int, int, PetscDrawLG *); 453798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 454798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 455798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 456798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 457798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 458798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 459798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 460798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 461798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 462798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 463798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 464798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 465798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 466798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 467798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 468798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 469798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 470798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 471798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 472fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 473798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 474d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 475d71ae5a4SJacob Faibussowitsch { 4769371c9d4SSatish Balay return KSPMonitorResidual(ksp, n, rnorm, vf); 4779371c9d4SSatish Balay } 478d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 479d71ae5a4SJacob Faibussowitsch { 4809371c9d4SSatish Balay return KSPMonitorTrueResidual(ksp, n, rnorm, vf); 4819371c9d4SSatish Balay } 482d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 483d71ae5a4SJacob Faibussowitsch { 4849371c9d4SSatish Balay return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf); 4859371c9d4SSatish Balay } 486798534f6SMatthew G. Knepley 487798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *); 488341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *); 489341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **); 490341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *); 491341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal); 492e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *); 493e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **); 494e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **); 49584cb2905SBarry Smith 496014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec); 497014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec); 498c01c455dSBarry Smith 49923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat); 50023ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *); 501014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *); 502014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]); 503014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]); 504014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]); 5051eb62cbbSBarry Smith 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool); 507014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *); 508014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool); 509014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *); 5101f7f0c4fSBarry Smith 511014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer); 51255849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer); 513fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]); 51419a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer); 515c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **)); 5161b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP); 517c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP); 51894a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer); 5191b2b9847SBarry Smith 520d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v) 521d71ae5a4SJacob Faibussowitsch { 5229371c9d4SSatish Balay return KSPConvergedReasonView(ksp, v); 5239371c9d4SSatish Balay } 524d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) 525d71ae5a4SJacob Faibussowitsch { 5269371c9d4SSatish Balay return KSPConvergedReasonViewFromOptions(ksp); 5279371c9d4SSatish Balay } 52855849f57SBarry Smith 52955849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 5301eb62cbbSBarry Smith 531dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool); 5320e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool); 533014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *); 534884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *); 535798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 536798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 537798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 538db9b2ab1SHong Zhang 539014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *); 540014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *); 54168ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *); 54283ab6a24SBarry Smith 54328ce4d24SBarry Smith /*E 544a4d1885cSBarry Smith KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence 5458a4b9c5cSBarry Smith test routines. 5468a4b9c5cSBarry Smith 547a4d1885cSBarry Smith Values: 548a4d1885cSBarry Smith + `KSP_NORM_DEFAULT` - use the default for the current `KSPType` 549a4d1885cSBarry Smith . `KSP_NORM_NONE` - use no norm calculation 550a4d1885cSBarry Smith . `KSP_NORM_PRECONDITIONED` - use the preconditioned residual norm 551a4d1885cSBarry Smith . `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm 552a4d1885cSBarry Smith - `KSP_NORM_NATURAL` - use the natural norm (the norm induced by the linear operator) 553a4d1885cSBarry Smith 5548a4b9c5cSBarry Smith Level: advanced 5558a4b9c5cSBarry Smith 556a4d1885cSBarry Smith Note: 557a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 55887497f52SBarry Smith depending on left or right preconditioning, see `KSPSetPCSide()` 559a3f661c8SBarry Smith 56087497f52SBarry Smith Developer Note: 561a4d1885cSBarry Smith This must match the values in petsc/finclude/petscksp.h 5628a4b9c5cSBarry Smith 563a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`, 564db781477SPatrick Sanan `KSPSetConvergenceTest()`, `KSPSetPCSide()` 5658a4b9c5cSBarry Smith E*/ 5669371c9d4SSatish Balay typedef enum { 5679371c9d4SSatish Balay KSP_NORM_DEFAULT = -1, 5689371c9d4SSatish Balay KSP_NORM_NONE = 0, 5699371c9d4SSatish Balay KSP_NORM_PRECONDITIONED = 1, 5709371c9d4SSatish Balay KSP_NORM_UNPRECONDITIONED = 2, 5719371c9d4SSatish Balay KSP_NORM_NATURAL = 3 5729371c9d4SSatish Balay } KSPNormType; 5739e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 574014dd563SJed Brown PETSC_EXTERN const char *const *const KSPNormTypes; 575a21b2a99SBarry Smith 5761f7e983dSSatish Balay /*MC 5779793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 5788c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 5798c5b8ba0SBarry Smith be based on a norm of a residual etc. 5808c5b8ba0SBarry Smith 5818c5b8ba0SBarry Smith Level: advanced 5828c5b8ba0SBarry Smith 583a4d1885cSBarry Smith Note: 584a4d1885cSBarry Smith Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored 5858c5b8ba0SBarry Smith 586a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 5878c5b8ba0SBarry Smith M*/ 5888c5b8ba0SBarry Smith 5891f7e983dSSatish Balay /*MC 5901957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 5918c5b8ba0SBarry Smith convergence test routine. 5928c5b8ba0SBarry Smith 5938c5b8ba0SBarry Smith Level: advanced 5948c5b8ba0SBarry Smith 595a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 5968c5b8ba0SBarry Smith M*/ 5978c5b8ba0SBarry Smith 5981f7e983dSSatish Balay /*MC 599ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 6008c5b8ba0SBarry Smith convergence test routine. 6018c5b8ba0SBarry Smith 6028c5b8ba0SBarry Smith Level: advanced 6038c5b8ba0SBarry Smith 604a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 6058c5b8ba0SBarry Smith M*/ 6068c5b8ba0SBarry Smith 6071f7e983dSSatish Balay /*MC 608ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 60987497f52SBarry Smith convergence test routine. This is only supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR` 6108c5b8ba0SBarry Smith 6118c5b8ba0SBarry Smith Level: advanced 6128c5b8ba0SBarry Smith 613a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()` 6148c5b8ba0SBarry Smith M*/ 6158c5b8ba0SBarry Smith 616014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType); 617014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *); 618014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt); 619014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt); 620014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool); 6218a4b9c5cSBarry Smith 6224a221d59SStefano Zampini #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM("Use KSP_CONVERGED_NEG_CURVE (since version 3.19)") 6234a221d59SStefano Zampini #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM("Use KSP_CONVERGED_STEP_LENGTH (since version 3.19)") 624f0fc11ceSJed Brown #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)") 6258a4b9c5cSBarry Smith /*E 626a4d1885cSBarry Smith KSPConvergedReason - reason a Krylov method was determined to have converged or diverged 62728ce4d24SBarry Smith 62828ce4d24SBarry Smith Level: beginner 62928ce4d24SBarry Smith 630a4d1885cSBarry Smith Values: 631a4d1885cSBarry Smith + `KSP_CONVERGED_RTOL_NORMAL` - requested decrease in the residual for the normal equations 632a4d1885cSBarry Smith . `KSP_CONVERGED_ATOL_NORMAL` - requested absolute value in the residual for the normal equations 633a4d1885cSBarry Smith . `KSP_CONVERGED_RTOL` - requested decrease in the residual 634a4d1885cSBarry Smith . `KSP_CONVERGED_ATOL` - requested absolute value in the residual 635a4d1885cSBarry Smith . `KSP_CONVERGED_ITS` - requested number of iterations 6364a221d59SStefano Zampini . `KSP_CONVERGED_NEG_CURVE` - see note below 637a4d1885cSBarry Smith . `KSP_CONVERGED_STEP_LENGTH` - see note below 638a4d1885cSBarry Smith . `KSP_CONVERGED_HAPPY_BREAKDOWN` - happy breakdown (meaning early convergence of the `KSPType` occurred. 639a4d1885cSBarry Smith . `KSP_DIVERGED_NULL` - breakdown when solving the Hessenberg system within GMRES 640a4d1885cSBarry Smith . `KSP_DIVERGED_ITS` - requested number of iterations 641a4d1885cSBarry Smith . `KSP_DIVERGED_DTOL` - large increase in the residual norm 642a4d1885cSBarry Smith . `KSP_DIVERGED_BREAKDOWN` - breakdown in the Krylov method 643a4d1885cSBarry Smith . `KSP_DIVERGED_BREAKDOWN_BICG` - breakdown in the `KSPBGCS` Krylov method 644a4d1885cSBarry Smith . `KSP_DIVERGED_NONSYMMETRIC` - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry 645a4d1885cSBarry Smith . `KSP_DIVERGED_INDEFINITE_PC` - the preconditioner was indefinite for a `KSPType` that requires it be definite 646a4d1885cSBarry Smith . `KSP_DIVERGED_NANORINF` - a not a number of infinity was detected in a vector during the computation 647a4d1885cSBarry Smith . `KSP_DIVERGED_INDEFINITE_MAT` - the operator was indefinite for a `KSPType` that requires it be definite 648a4d1885cSBarry Smith - `KSP_DIVERGED_PC_FAILED` - the action of the preconditioner failed for some reason 649a4d1885cSBarry Smith 650a4d1885cSBarry Smith Note: 6514a221d59SStefano Zampini The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by the special `KSPNASH`, 652a4d1885cSBarry Smith `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver. 65328ce4d24SBarry Smith 65495452b02SPatrick Sanan Developer Notes: 655a4d1885cSBarry Smith This must match the values in petsc/finclude/petscksp.h 6564d0a8057SBarry Smith 65787497f52SBarry Smith The string versions of these are `KSPConvergedReasons`; if you change 6584d0a8057SBarry Smith any of the values here also change them that array of names. 65986c02ca4SBarry Smith 660a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()` 66128ce4d24SBarry Smith E*/ 662d15094e1SBarry Smith typedef enum { /* converged */ 6639ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 6649ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 665d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 666d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 667b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 6684a221d59SStefano Zampini KSP_CONVERGED_NEG_CURVE = 5, 6694a221d59SStefano Zampini KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED = 5, 6704a221d59SStefano Zampini KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6, 6714a221d59SStefano Zampini KSP_CONVERGED_STEP_LENGTH = 6, 6724a221d59SStefano Zampini KSP_CONVERGED_HAPPY_BREAKDOWN = 7, 673d15094e1SBarry Smith /* diverged */ 674b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 675d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 676d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 677d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 678b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 679b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 680b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 6814d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 6826aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 683c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED = -11, 684aa4d2078SSatish Balay KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11, 685d15094e1SBarry Smith 6869371c9d4SSatish Balay KSP_CONVERGED_ITERATING = 0 6879371c9d4SSatish Balay } KSPConvergedReason; 688014dd563SJed Brown PETSC_EXTERN const char *const *KSPConvergedReasons; 689d15094e1SBarry Smith 690c838673bSBarry Smith /*MC 69187497f52SBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if `KSPConvergedDefaultSetUIRNorm()` was called 692c838673bSBarry Smith 693c838673bSBarry Smith Level: beginner 694c838673bSBarry Smith 69587497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 696c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 697c838673bSBarry Smith 2-norm of the residual for right preconditioning 698c838673bSBarry Smith 69987497f52SBarry Smith See also `KSP_CONVERGED_ATOL` which may apply before this tolerance. 700f9fed41fSBarry Smith 701a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 702c838673bSBarry Smith 703c838673bSBarry Smith M*/ 704c838673bSBarry Smith 705c838673bSBarry Smith /*MC 706c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 707c838673bSBarry Smith 708c838673bSBarry Smith Level: beginner 709c838673bSBarry Smith 71087497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 711c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 712c838673bSBarry Smith 2-norm of the residual for right preconditioning 713c838673bSBarry Smith 71487497f52SBarry Smith See also `KSP_CONVERGED_RTOL` which may apply before this tolerance. 715c838673bSBarry Smith 716a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 717c838673bSBarry Smith 718c838673bSBarry Smith M*/ 719c838673bSBarry Smith 720c838673bSBarry Smith /*MC 721c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 722c838673bSBarry Smith 723c838673bSBarry Smith Level: beginner 724c838673bSBarry Smith 72587497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 726c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 727c838673bSBarry Smith 2-norm of the residual for right preconditioning 728c838673bSBarry Smith 729c838673bSBarry Smith Level: beginner 730c838673bSBarry Smith 731a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 732c838673bSBarry Smith 733c838673bSBarry Smith M*/ 734c838673bSBarry Smith 735c838673bSBarry Smith /*MC 736c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 737c838673bSBarry Smith reached 738c838673bSBarry Smith 739c838673bSBarry Smith Level: beginner 740c838673bSBarry Smith 741a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 742c838673bSBarry Smith 743c838673bSBarry Smith M*/ 744c838673bSBarry Smith 745c838673bSBarry Smith /*MC 74687497f52SBarry Smith KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of 74787497f52SBarry Smith the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence 7484d0a8057SBarry Smith test routine is set in KSP. 749c838673bSBarry Smith 750c838673bSBarry Smith Level: beginner 751c838673bSBarry Smith 752a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 753c838673bSBarry Smith 754c838673bSBarry Smith M*/ 755c838673bSBarry Smith 756c838673bSBarry Smith /*MC 757c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 7581de96524SPierre Jolivet method could not continue to enlarge the Krylov space. Could be due to a singular matrix or 75987497f52SBarry Smith preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block 7601de96524SPierre Jolivet are colinear. 761c838673bSBarry Smith 762c838673bSBarry Smith Level: beginner 763c838673bSBarry Smith 764a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 765c838673bSBarry Smith 766c838673bSBarry Smith M*/ 767c838673bSBarry Smith 768c838673bSBarry Smith /*MC 76987497f52SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the 770c838673bSBarry Smith method could not continue to enlarge the Krylov space. 771c838673bSBarry Smith 772c838673bSBarry Smith Level: beginner 773c838673bSBarry Smith 774a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 775c838673bSBarry Smith 776c838673bSBarry Smith M*/ 777c838673bSBarry Smith 778c838673bSBarry Smith /*MC 779c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 78087497f52SBarry Smith symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry 781c838673bSBarry Smith 782c838673bSBarry Smith Level: beginner 783c838673bSBarry Smith 784a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 785c838673bSBarry Smith 786c838673bSBarry Smith M*/ 787c838673bSBarry Smith 788c838673bSBarry Smith /*MC 789c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 79087497f52SBarry Smith positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to 791c838673bSBarry Smith be positive definite 792c838673bSBarry Smith 793c838673bSBarry Smith Level: beginner 794c838673bSBarry Smith 79587497f52SBarry Smith Note: 796a4d1885cSBarry Smith This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force 79787497f52SBarry Smith the `PCICC` preconditioner to generate a positive definite preconditioner 798c838673bSBarry Smith 799a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 800c838673bSBarry Smith 801c838673bSBarry Smith M*/ 802c838673bSBarry Smith 803c838673bSBarry Smith /*MC 804c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a 8059fc87aa7SBarry Smith zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 80687497f52SBarry Smith such as `PCFIELDSPLIT`. 8079fc87aa7SBarry Smith 8089fc87aa7SBarry Smith Level: beginner 8099fc87aa7SBarry Smith 810a4d1885cSBarry Smith Note: 811a4d1885cSBarry Smith Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details. 8129fc87aa7SBarry Smith 813a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 8149fc87aa7SBarry Smith 8159fc87aa7SBarry Smith M*/ 8169fc87aa7SBarry Smith 8179fc87aa7SBarry Smith /*MC 818a4d1885cSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called 819a4d1885cSBarry Smith while `KSPSolve()` is still running. 820c838673bSBarry Smith 821c838673bSBarry Smith Level: beginner 822c838673bSBarry Smith 823a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 824c838673bSBarry Smith 825c838673bSBarry Smith M*/ 826c838673bSBarry Smith 827014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *)); 828df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *)); 829df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *)); 8303ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *); 8318de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 8323eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 8338de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 8348de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 8358de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 8368de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 83754b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool); 8380059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 839014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *); 840c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **); 84194a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 8422a28d964SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetConvergedNegativeCurvature(KSP, PetscBool); 8432a28d964SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergedNegativeCurvature(KSP, PetscBool *); 844abef13c0SSatish Balay 845d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") static inline void KSPDefaultConverged(void) 846d71ae5a4SJacob Faibussowitsch { /* never called */ 8479371c9d4SSatish Balay } 8488ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 849d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") static inline void KSPDefaultConvergedDestroy(void) 850d71ae5a4SJacob Faibussowitsch { /* never called */ 8519371c9d4SSatish Balay } 8528ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 853d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") static inline void KSPDefaultConvergedCreate(void) 854d71ae5a4SJacob Faibussowitsch { /* never called */ 8559371c9d4SSatish Balay } 8568ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 857d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUIRNorm(void) 858d71ae5a4SJacob Faibussowitsch { /* never called */ 8599371c9d4SSatish Balay } 8608ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 861d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUMIRNorm(void) 862d71ae5a4SJacob Faibussowitsch { /* never called */ 8639371c9d4SSatish Balay } 8648ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 865d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") static inline void KSPSkipConverged(void) 866d71ae5a4SJacob Faibussowitsch { /* never called */ 8679371c9d4SSatish Balay } 8688ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 8698ea1b3e6SJed Brown 8700bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *); 871d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B) 872d71ae5a4SJacob Faibussowitsch { 8739371c9d4SSatish Balay return KSPComputeOperator(A, NULL, B); 8749371c9d4SSatish Balay } 875d4fbbf0eSBarry Smith 87628ce4d24SBarry Smith /*E 877a4d1885cSBarry Smith KSPCGType - Determines what type of `KSPCG` to use 87828ce4d24SBarry Smith 87928ce4d24SBarry Smith Level: beginner 88028ce4d24SBarry Smith 881a4d1885cSBarry Smith Values: 882a4d1885cSBarry Smith + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric 883a4d1885cSBarry Smith - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian 884a4d1885cSBarry Smith 885a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSP`, `KSPCGSetType()` 88628ce4d24SBarry Smith E*/ 8879371c9d4SSatish Balay typedef enum { 8889371c9d4SSatish Balay KSP_CG_SYMMETRIC = 0, 8899371c9d4SSatish Balay KSP_CG_HERMITIAN = 1 8909371c9d4SSatish Balay } KSPCGType; 8916a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 89228ce4d24SBarry Smith 893014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType); 894014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool); 8958031f4adStmunson 896dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal); 897dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *); 898dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *); 899fcae7a14Stmunson 90005de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *); 90105de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *); 902d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x) 903d71ae5a4SJacob Faibussowitsch { 9049371c9d4SSatish Balay return KSPGLTRGetMinEig(ksp, x); 9059371c9d4SSatish Balay } 906d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x) 907d71ae5a4SJacob Faibussowitsch { 9089371c9d4SSatish Balay return KSPGLTRGetLambda(ksp, x); 9099371c9d4SSatish Balay } 9108031f4adStmunson 911014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]); 912ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]); 9131d6018f0SLisandro Dalcin 914f560b561SHong Zhang PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP)); 915014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP); 916014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP); 9173369ce9aSBarry Smith 9189804daf3SBarry Smith #include <petscdrawtypes.h> 919014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *); 9202f2e5d10SKris Buschelman 921014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec)); 922014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec)); 92303919abeSBarry Smith 924ba36735cSStefano Zampini /*S 925a4d1885cSBarry Smith KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods. 926f8a50e2bSBarry Smith 927a4d1885cSBarry Smith Level: intermediate 928f8a50e2bSBarry Smith 929a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType` 930ba36735cSStefano Zampini S*/ 931ba36735cSStefano Zampini typedef struct _p_KSPGuess *KSPGuess; 932ba36735cSStefano Zampini /*J 93387497f52SBarry Smith KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods. 934ba36735cSStefano Zampini 935a4d1885cSBarry Smith Level: intermediate 936ba36735cSStefano Zampini 937a4d1885cSBarry Smith Values: 938a4d1885cSBarry Smith + `KSPGUESSFISCHER` - methodology developed by Paul Fischer 939a4d1885cSBarry Smith - `KSPGUESSPOD` - methodology based on proper orthogonal decomposition 940a4d1885cSBarry Smith 941a4d1885cSBarry Smith .seealso: [](chapter_ksp), `KSP`, `KSPGuess` 942ba36735cSStefano Zampini J*/ 943ba36735cSStefano Zampini typedef const char *KSPGuessType; 944ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer" 945ba36735cSStefano Zampini #define KSPGUESSPOD "pod" 946a4d1885cSBarry Smith 9471d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess)); 948ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess); 949ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *); 950ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer); 951ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *); 952ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *); 953ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType); 954ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *); 9558410009bSDavid Wells PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal); 956ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 957ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec); 958ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec); 959ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 960ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt); 961014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt); 962ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool); 963ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *); 964f8a50e2bSBarry Smith 965470b340bSDmitry Karpeev /*E 966470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 967470b340bSDmitry Karpeev 968470b340bSDmitry Karpeev Level: intermediate 969470b340bSDmitry Karpeev 970a4d1885cSBarry Smith .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`, 971a4d1885cSBarry Smith `MatCreateSchurComplementPmat()` 972470b340bSDmitry Karpeev E*/ 9739371c9d4SSatish Balay typedef enum { 9749371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_DIAG, 9759371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_LUMP, 9769371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, 9779371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_FULL 9789371c9d4SSatish Balay } MatSchurComplementAinvType; 979470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 980470b340bSDmitry Karpeev 9819371c9d4SSatish Balay typedef enum { 9829371c9d4SSatish Balay MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0, 983864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1, 984864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2, 9859371c9d4SSatish Balay MAT_LMVM_SYMBROYDEN_SCALE_USER = 3 9869371c9d4SSatish Balay } MatLMVMSymBroydenScaleType; 987864588a7SAlp Dener PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[]; 98892f76d53SAlp Dener 989014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *); 990014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *); 991d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP); 992bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 993aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 994bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *); 995470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType); 996470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *); 9975bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *); 9985a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *); 999470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *); 1000470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *); 10013f22127dSBarry Smith 100278e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *); 100378e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *); 100478e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *); 1005864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1006864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1007864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1008864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1009864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1010cd929ea3SAlp Dener 1011cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 1012b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *); 1013cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 1014cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 1015e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat); 10160ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat); 1017cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 1018cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 1019cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 1020cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 1021cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 10222d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec); 1023cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 1024cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *); 1025cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *); 1026cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *); 102792f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt); 1028cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *); 1029cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *); 1030864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar); 1031864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType); 1032cd929ea3SAlp Dener 1033014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM); 1034014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool); 1035014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *); 1036014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *); 1037014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *); 1038fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, PetscErrorCode (*func)(KSP, Vec, void *), void *); 103923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *); 1040fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, PetscErrorCode (*)(KSP, Vec, void *), void *); 104123ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *); 104223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, PetscErrorCode (**)(KSP, Mat, Mat, void *), void *); 1043014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, PetscErrorCode (*)(KSP, Vec, void *), void *); 1044014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, PetscErrorCode (**)(KSP, Vec, void *), void *); 1045fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, PetscErrorCode (*)(KSP, Vec, void *), void *); 1046fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, PetscErrorCode (**)(KSP, Vec, void *), void *); 10476c699258SBarry Smith 104802b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec); 10499371c9d4SSatish Balay PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode, Vec); 1050557cf195SMatthew G. Knepley 10512b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *); 10522b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal); 10532eac72dbSBarry Smith #endif 1054