1f26ada1bSBarry Smith /* 2f26ada1bSBarry Smith Defines the interface functions for the Krylov subspace accelerators. 3f26ada1bSBarry Smith */ 4a4963045SJacob Faibussowitsch #pragma once 5ac09b921SBarry Smith 62c8e378dSBarry Smith #include <petscpc.h> 72eac72dbSBarry Smith 8ac09b921SBarry Smith /* SUBMANSEC = KSP */ 9ac09b921SBarry Smith 10607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void); 114bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode KSPFinalizePackage(void); 121dbb0a54SBarry Smith 1328ce4d24SBarry Smith /*S 140b4b7b1cSBarry Smith KSP - Abstract PETSc object that manages the linear solves in PETSc (even those such as direct factorization-based solvers that 150b4b7b1cSBarry Smith do not use Krylov accelerators). 1628ce4d24SBarry Smith 1728ce4d24SBarry Smith Level: beginner 1828ce4d24SBarry Smith 190b4b7b1cSBarry Smith Notes: 20a4d1885cSBarry Smith When a direct solver is used, but no Krylov solver is used, the `KSP` object is still used but with a 2116a05f60SBarry Smith `KSPType` of `KSPPREONLY` (or equivalently `KSPNONE`), meaning that only application of the preconditioner is used as the linear solver. 228f6c3df8SBarry Smith 230b4b7b1cSBarry Smith Use `KSPSetType()` or the options database key `-ksp_type` to set the specific Krylov solver algorithm to use with a given `KSP` object 240b4b7b1cSBarry Smith 250b4b7b1cSBarry Smith The `PC` object is used to control preconditioners in PETSc. 260b4b7b1cSBarry Smith 27789736e1SBarry Smith `KSP` can also be used to solve some least squares problems (over or under-determined linear systems), using, for example, `KSPLSQR`, see `PETSCREGRESSORLINEAR` 28789736e1SBarry Smith for additional methods that can be used to solve least squares problems and other linear regressions). 29789736e1SBarry Smith 301cc06b55SBarry Smith .seealso: [](doc_linsolve), [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES` 3128ce4d24SBarry Smith S*/ 32e2a1c21fSSatish Balay typedef struct _p_KSP *KSP; 332eac72dbSBarry Smith 3476bdecfbSBarry Smith /*J 350b4b7b1cSBarry Smith KSPType - String with the name of a PETSc Krylov method. These are all the Krylov solvers that PETSc provides. 3628ce4d24SBarry Smith 3728ce4d24SBarry Smith Level: beginner 3828ce4d24SBarry Smith 391cc06b55SBarry Smith .seealso: [](doc_linsolve), [](ch_ksp), `KSPSetType()`, `KSP`, `KSPRegister()`, `KSPCreate()`, `KSPSetFromOptions()` 4076bdecfbSBarry Smith J*/ 4119fd82e9SBarry Smith typedef const char *KSPType; 4282bf6240SBarry Smith #define KSPRICHARDSON "richardson" 436c9de887SHong Zhang #define KSPCHEBYSHEV "chebyshev" 4482bf6240SBarry Smith #define KSPCG "cg" 452c8fc646SJed Brown #define KSPGROPPCG "groppcg" 462c8fc646SJed Brown #define KSPPIPECG "pipecg" 47901ccb91SSiegfried Cools #define KSPPIPECGRR "pipecgrr" 48265663fdSSiegfried Cools #define KSPPIPELCG "pipelcg" 49b21a8899STyler Chen #define KSPPIPEPRCG "pipeprcg" 50325e8391SManasi #define KSPPIPECG2 "pipecg2" 51df2a969aSvictorle #define KSPCGNE "cgne" 5205de396fSBarry Smith #define KSPNASH "nash" 5305de396fSBarry Smith #define KSPSTCG "stcg" 5405de396fSBarry Smith #define KSPGLTR "gltr" 55edd03b47SJacob Faibussowitsch #define KSPCGNASH PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPNASH", ) "nash" 56edd03b47SJacob Faibussowitsch #define KSPCGSTCG PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSTCG", ) "stcg" 57edd03b47SJacob Faibussowitsch #define KSPCGGLTR PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSGLTR", ) "gltr" 58640f4f53SPatrick Sanan #define KSPFCG "fcg" 59390d8e47SPatrick Sanan #define KSPPIPEFCG "pipefcg" 6082bf6240SBarry Smith #define KSPGMRES "gmres" 61483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres" 629dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 639dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 644f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 6561b468f9SJed Brown #define KSPPGMRES "pgmres" 6682bf6240SBarry Smith #define KSPTCQMR "tcqmr" 6782bf6240SBarry Smith #define KSPBCGS "bcgs" 68e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 69345ecf0bSXiangmin Jiao #define KSPQMRCGS "qmrcgs" 7018ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 71c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 72f0037002Svictorle #define KSPBCGSL "bcgsl" 73f154af2dSSiegfried Cools #define KSPPIPEBCGS "pipebcgs" 7482bf6240SBarry Smith #define KSPCGS "cgs" 7582bf6240SBarry Smith #define KSPTFQMR "tfqmr" 7682bf6240SBarry Smith #define KSPCR "cr" 772c8fc646SJed Brown #define KSPPIPECR "pipecr" 7882bf6240SBarry Smith #define KSPLSQR "lsqr" 7982bf6240SBarry Smith #define KSPPREONLY "preonly" 803c2be86cSBarry Smith #define KSPNONE "none" 8182bf6240SBarry Smith #define KSPQCG "qcg" 82c9cf9b11SBarry Smith #define KSPBICG "bicg" 83b4ac9ba4SBarry Smith #define KSPMINRES "minres" 8401c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 8521356919SSatish Balay #define KSPLCD "lcd" 861d6018f0SLisandro Dalcin #define KSPPYTHON "python" 8758ad63f4SBarry Smith #define KSPGCR "gcr" 88fad47a0aSPatrick Sanan #define KSPPIPEGCR "pipegcr" 89e4d80e07Szianekhodja #define KSPTSIRM "tsirm" 90e4d80e07Szianekhodja #define KSPCGLS "cgls" 91329cd39dSStefano Zampini #define KSPFETIDP "fetidp" 9238cfc46eSPierre Jolivet #define KSPHPDDM "hpddm" 932eac72dbSBarry Smith 948ba1e511SMatthew Knepley /* Logging support */ 95014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 96ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID; 975358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 988ba1e511SMatthew Knepley 99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm, KSP *); 10019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP, KSPType); 101c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP, KSPType *); 102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 103014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP, Vec, Vec); 105014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP, Vec, Vec); 10675f8e9bdSHong Zhang PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP, PetscBool); 1075ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP, Mat, Mat); 108bbbebc2cSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolveTranspose(KSP, Mat, Mat); 1093e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt); 110edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPSetMatSolveBatchSize()", ) static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n) 111d71ae5a4SJacob Faibussowitsch { 1129371c9d4SSatish Balay return KSPSetMatSolveBatchSize(ksp, n); 1139371c9d4SSatish Balay } 1143e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *); 115edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPGetMatSolveBatchSize()", ) static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n) 116d71ae5a4SJacob Faibussowitsch { 1179371c9d4SSatish Balay return KSPGetMatSolveBatchSize(ksp, n); 1189371c9d4SSatish Balay } 119014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 120ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP); 121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *); 12223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool); 1237d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *); 12425c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool); 125c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec); 1262eac72dbSBarry Smith 127140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 128ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList; 129798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorList; 130798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorCreateList; 131798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList; 132bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[], PetscErrorCode (*)(KSP)); 133798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **)); 13430de9b25SBarry Smith 135014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide); 136014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *); 137014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt); 138c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *); 13925c92fe2SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetMinimumIterations(KSP, PetscInt); 14025c92fe2SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetMinimumIterations(KSP, PetscInt *); 141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool); 142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *); 143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool); 144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *); 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool); 1467d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool); 147c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *); 148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool); 149c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *); 150014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *); 151014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *); 152014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *); 153014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *); 154734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *); 1552a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **); 156edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 6, 0, "KSPCreateVecs()", ) static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y) 157d71ae5a4SJacob Faibussowitsch { 1589371c9d4SSatish Balay return KSPCreateVecs(ksp, n, x, m, y); 1599371c9d4SSatish Balay } 1602eac72dbSBarry Smith 161d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *); 162d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *); 163d4211eb9SBarry Smith 164014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC); 165014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *); 1663821be0aSBarry Smith PETSC_EXTERN PetscErrorCode KSPSetNestLevel(KSP, PetscInt); 1673821be0aSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetNestLevel(KSP, PetscInt *); 168aabeff55SBarry Smith 169014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal); 17049abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscCtxDestroyFn *); 171014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 1723ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *); 17394a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *); 1746497c311SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscCount, PetscBool); 17594a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *); 1766497c311SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscCount, PetscBool); 1774b0e389bSBarry Smith 178fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *); 179fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *); 180fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 181fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt); 182fa0ddf94SBarry Smith 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *); 184cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP); 185014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 188fb80e629SPablo Brubeck PETSC_EXTERN PetscErrorCode PCPatchGetSubKSP(PC, PetscInt *, KSP *[]); 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]); 190285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]); 191b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *); 192b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *); 193b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *); 194b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *); 195014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *); 1964a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *); 197f3b08a26SMatthew G. Knepley /* 198f3b08a26SMatthew G. Knepley PCMGCoarseList contains the list of coarse space constructor currently registered 199f3b08a26SMatthew G. Knepley These are added with PCMGRegisterCoarseSpaceConstructor() 200f3b08a26SMatthew G. Knepley */ 201f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList; 2022b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *)); 2032b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *)); 204f3b08a26SMatthew G. Knepley 205014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *); 206014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *); 2072eac72dbSBarry Smith 208f2edd1f0SMalachi Phillips /*E 20995bd0b28SBarry Smith KSPChebyshevKind - Which kind of Chebyshev polynomial to use with `KSPCHEBYSHEV` 210f2edd1f0SMalachi Phillips 211f2edd1f0SMalachi Phillips Values: 212f2edd1f0SMalachi Phillips + `KSP_CHEBYSHEV_FIRST` - "classic" first-kind Chebyshev polynomial 213f2edd1f0SMalachi Phillips . `KSP_CHEBYSHEV_FOURTH` - fourth-kind Chebyshev polynomial 214f2edd1f0SMalachi Phillips - `KSP_CHEBYSHEV_OPT_FOURTH` - optimized fourth-kind Chebyshev polynomial 215f2edd1f0SMalachi Phillips 216f2edd1f0SMalachi Phillips Level: intermediate 217f2edd1f0SMalachi Phillips 2181cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPCHEBYSHEV`, `KSPChebyshevSetKind()`, `KSPChebyshevGetKind()` 219f2edd1f0SMalachi Phillips E*/ 220f2edd1f0SMalachi Phillips typedef enum { 221f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_FIRST, 222f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_FOURTH, 223f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_OPT_FOURTH 224f2edd1f0SMalachi Phillips } KSPChebyshevKind; 225f2edd1f0SMalachi Phillips 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal); 227014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool); 228014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal); 22958450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal); 230b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool); 231f2edd1f0SMalachi Phillips PETSC_EXTERN PetscErrorCode KSPChebyshevSetKind(KSP, KSPChebyshevKind); 23216a05f60SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevGetKind(KSP, KSPChebyshevKind *); 23358450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *); 234014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *); 235d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *); 236d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]); 2377d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]); 2384b0e389bSBarry Smith 239640f4f53SPatrick Sanan /*E 240640f4f53SPatrick Sanan 24106137d0aSPatrick Sanan KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 242640f4f53SPatrick Sanan 24367b8a455SSatish Balay Values: 244a1cb98faSBarry Smith + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions 245a1cb98faSBarry Smith - `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 246640f4f53SPatrick Sanan 2472b26979fSBarry Smith Level: intermediate 248640f4f53SPatrick Sanan 2491cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()` 250640f4f53SPatrick Sanan E*/ 2519371c9d4SSatish Balay typedef enum { 2529371c9d4SSatish Balay KSP_FCD_TRUNC_TYPE_STANDARD, 2539371c9d4SSatish Balay KSP_FCD_TRUNC_TYPE_NOTAY 2549371c9d4SSatish Balay } KSPFCDTruncationType; 25506137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 256640f4f53SPatrick Sanan 257640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt); 258640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *); 259640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt); 260640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *); 26106137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType); 26206137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *); 263640f4f53SPatrick Sanan 264390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt); 265390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *); 266390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt); 267390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *); 268390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType); 269390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *); 270390d8e47SPatrick Sanan 271fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt); 272fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *); 273fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt); 274fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *); 275fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType); 276fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *); 277fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool); 278fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *); 2794bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 28083f0b094SBarry Smith 281014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *); 283014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal); 284e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal); 2859f236934SBarry Smith 286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt)); 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt)); 289014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt); 290014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt); 2911d73ed98SBarry Smith 292014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt); 293014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2941d73ed98SBarry Smith 295483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar); 296483d6965SPatrick Sanan 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt); 298014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *); 299014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 30058ad63f4SBarry Smith 301f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal); 302f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *); 303f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool); 304f87a0b54SStefano Zampini 30504d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *); 30604d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC); 30704d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *); 308568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat); 309e82af88dSprj- 3106cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat); 3116cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *); 3126cac28cbSPierre Jolivet #if PetscDefined(HAVE_HPDDM) 313edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMSetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) 314d71ae5a4SJacob Faibussowitsch { 3159371c9d4SSatish Balay return KSPHPDDMSetDeflationMat(ksp, U); 3169371c9d4SSatish Balay } 317edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMGetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) 318d71ae5a4SJacob Faibussowitsch { 3199371c9d4SSatish Balay return KSPHPDDMGetDeflationMat(ksp, U); 3209371c9d4SSatish Balay } 3216cac28cbSPierre Jolivet #endif 322edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPMatSolve()", ) static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) 323d71ae5a4SJacob Faibussowitsch { 3249371c9d4SSatish Balay return KSPMatSolve(ksp, B, X); 3259371c9d4SSatish Balay } 326d55ab951SPierre Jolivet /*E 32787497f52SBarry Smith KSPHPDDMType - Type of Krylov method used by `KSPHPDDM` 328d55ab951SPierre Jolivet 329d55ab951SPierre Jolivet Values: 330a4d1885cSBarry Smith + `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method 331a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BGMRES` - block GMRES 332a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_CG` - Conjugate Gradient 333a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BCG` - block CG 334a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_GCRODR` - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting 335a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BGCRODR` - block GCRODR 336a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BFBCG` - breakdown-free BCG 337a4d1885cSBarry Smith - `KSP_HPDDM_TYPE_PREONLY` - apply the preconditioner only 338d55ab951SPierre Jolivet 33916a05f60SBarry Smith Level: intermediate 34016a05f60SBarry Smith 3411cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPHPDDM`, `KSPHPDDMSetType()` 342d55ab951SPierre Jolivet E*/ 343d55ab951SPierre Jolivet typedef enum { 344d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GMRES = 0, 345d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGMRES = 1, 346d55ab951SPierre Jolivet KSP_HPDDM_TYPE_CG = 2, 347d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BCG = 3, 348d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GCRODR = 4, 349d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGCRODR = 5, 350d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BFBCG = 6, 351d55ab951SPierre Jolivet KSP_HPDDM_TYPE_PREONLY = 7 352d55ab951SPierre Jolivet } KSPHPDDMType; 353d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[]; 354a4d1885cSBarry Smith 3552dd49c90SPierre Jolivet /*E 35687497f52SBarry Smith KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM` 3572dd49c90SPierre Jolivet 3582dd49c90SPierre Jolivet Values: 359a4d1885cSBarry Smith + `KSP_HPDDM_PRECISION_HALF` - default when PETSc is configured `--with-precision=__fp16` 360a4d1885cSBarry Smith . `KSP_HPDDM_PRECISION_SINGLE` - default when PETSc is configured `--with-precision=single` 361a4d1885cSBarry Smith . `KSP_HPDDM_PRECISION_DOUBLE` - default when PETSc is configured `--with-precision=double` 362a4d1885cSBarry Smith - `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128` 3632dd49c90SPierre Jolivet 36416a05f60SBarry Smith Level: intermediate 36516a05f60SBarry Smith 3661cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPHPDDM` 3672dd49c90SPierre Jolivet E*/ 3682dd49c90SPierre Jolivet typedef enum { 3692dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_HALF = 0, 3702dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_SINGLE = 1, 3712dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_DOUBLE = 2, 3722dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_QUADRUPLE = 3 3732dd49c90SPierre Jolivet } KSPHPDDMPrecision; 374d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType); 375d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *); 376e82af88dSprj- 377b49fd9e1SBarry Smith /*E 37816a05f60SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed in the GMRES solvers 379b49fd9e1SBarry Smith 380a4d1885cSBarry Smith Values: 381a4d1885cSBarry Smith + `KSP_GMRES_CGS_REFINE_NEVER` - one step of classical Gram-Schmidt 382a4d1885cSBarry Smith . `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria 383a4d1885cSBarry Smith - `KSP_GMRES_CGS_REFINE_ALWAYS` - always perform two steps 384b49fd9e1SBarry Smith 38516a05f60SBarry Smith Level: advanced 38616a05f60SBarry Smith 38795bd0b28SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPGMRES`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 388a4d1885cSBarry Smith `KSPGMRESGetOrthogonalization()`, 389a4d1885cSBarry Smith `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()` 390b49fd9e1SBarry Smith E*/ 3919371c9d4SSatish Balay typedef enum { 3929371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_NEVER, 3939371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_IFNEEDED, 3949371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_ALWAYS 3959371c9d4SSatish Balay } KSPGMRESCGSRefinementType; 3966a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 39716a05f60SBarry Smith 3981f7e983dSSatish Balay /*MC 3991957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 4008c5b8ba0SBarry Smith 4018c5b8ba0SBarry Smith Level: advanced 4028c5b8ba0SBarry Smith 40387497f52SBarry Smith Note: 40487497f52SBarry Smith Possibly unstable, but the fastest to compute 4058c5b8ba0SBarry Smith 40695bd0b28SBarry Smith .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 407a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 408db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 409db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 4108c5b8ba0SBarry Smith M*/ 4118c5b8ba0SBarry Smith 4121f7e983dSSatish Balay /*MC 4138c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 4148c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 4158c5b8ba0SBarry Smith poor orthogonality. 4168c5b8ba0SBarry Smith 4178c5b8ba0SBarry Smith Level: advanced 4188c5b8ba0SBarry Smith 419a4d1885cSBarry Smith Note: 420a4d1885cSBarry Smith This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to 4218c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 4228c5b8ba0SBarry Smith 42395bd0b28SBarry Smith .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 424a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 425db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 426db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 4278c5b8ba0SBarry Smith M*/ 4288c5b8ba0SBarry Smith 4291f7e983dSSatish Balay /*MC 43016a05f60SBarry Smith KSP_GMRES_CGS_REFINE_ALWAYS - Do two steps of the classical (unmodified) Gram-Schmidt process. 4318c5b8ba0SBarry Smith 4328c5b8ba0SBarry Smith Level: advanced 4338c5b8ba0SBarry Smith 43487497f52SBarry Smith Notes: 43587497f52SBarry Smith This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice 43687497f52SBarry Smith but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`. 4378c5b8ba0SBarry Smith 4388c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 4398c5b8ba0SBarry Smith 44095bd0b28SBarry Smith .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 441a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 442db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 443db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 4448c5b8ba0SBarry Smith M*/ 4458c5b8ba0SBarry Smith 446014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType); 447014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *); 44808480c60SBarry Smith 449014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *); 450014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *); 451014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 452c38d4ed2SBarry Smith 453014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal); 454014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *); 455014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *); 456121fd945SKris Buschelman 457014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal); 458014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool); 459014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt); 460e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool); 461d9492815SBarry Smith 462014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 46387d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP); 4642eac72dbSBarry Smith 465798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *); 466798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 467798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 468798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 469798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 470798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 471798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 472798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 473798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 474798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 475798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 476798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 477798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 478798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 479798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 480798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 481798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 482798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 483798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 484798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 485fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 486798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 487edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorResidual()", ) static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 488d71ae5a4SJacob Faibussowitsch { 4899371c9d4SSatish Balay return KSPMonitorResidual(ksp, n, rnorm, vf); 4909371c9d4SSatish Balay } 491edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidual()", ) static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 492d71ae5a4SJacob Faibussowitsch { 4939371c9d4SSatish Balay return KSPMonitorTrueResidual(ksp, n, rnorm, vf); 4949371c9d4SSatish Balay } 495edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidualMax()", ) static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 496d71ae5a4SJacob Faibussowitsch { 4979371c9d4SSatish Balay return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf); 4989371c9d4SSatish Balay } 499798534f6SMatthew G. Knepley 500798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *); 501341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *); 502341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **); 503341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *); 504341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal); 505e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *); 506e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **); 507e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **); 50884cb2905SBarry Smith 509014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec); 510014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec); 511c01c455dSBarry Smith 51223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat); 51323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *); 514014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *); 515014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]); 516014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]); 517014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]); 5181eb62cbbSBarry Smith 519014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool); 520014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *); 521014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool); 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *); 5231f7f0c4fSBarry Smith 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer); 52555849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer); 526fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]); 52719a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer); 52849abdd8aSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *, PetscCtxDestroyFn *); 5291b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP); 530c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP); 53194a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer); 5321b2b9847SBarry Smith 533edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonView()", ) static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v) 534d71ae5a4SJacob Faibussowitsch { 5359371c9d4SSatish Balay return KSPConvergedReasonView(ksp, v); 5369371c9d4SSatish Balay } 537edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) 538d71ae5a4SJacob Faibussowitsch { 5399371c9d4SSatish Balay return KSPConvergedReasonViewFromOptions(ksp); 5409371c9d4SSatish Balay } 54155849f57SBarry Smith 54255849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 5431eb62cbbSBarry Smith 544dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool); 5450e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool); 546014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *); 547884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *); 548798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 549798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 550798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 551db9b2ab1SHong Zhang 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *); 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *); 55468ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *); 5559f0612e4SBarry Smith PETSC_EXTERN PetscErrorCode PCMPIGetKSP(PC, KSP *); 55683ab6a24SBarry Smith 55728ce4d24SBarry Smith /*E 558a4d1885cSBarry Smith KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence 5598a4b9c5cSBarry Smith test routines. 5608a4b9c5cSBarry Smith 561a4d1885cSBarry Smith Values: 562a4d1885cSBarry Smith + `KSP_NORM_DEFAULT` - use the default for the current `KSPType` 563a4d1885cSBarry Smith . `KSP_NORM_NONE` - use no norm calculation 564a4d1885cSBarry Smith . `KSP_NORM_PRECONDITIONED` - use the preconditioned residual norm 565a4d1885cSBarry Smith . `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm 566a4d1885cSBarry Smith - `KSP_NORM_NATURAL` - use the natural norm (the norm induced by the linear operator) 567a4d1885cSBarry Smith 5688a4b9c5cSBarry Smith Level: advanced 5698a4b9c5cSBarry Smith 570a4d1885cSBarry Smith Note: 571a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 57216a05f60SBarry Smith depending on whether left or right preconditioning is used, see `KSPSetPCSide()` 573a3f661c8SBarry Smith 5741cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`, 57595bd0b28SBarry Smith `KSPSetConvergenceTest()`, `KSPSetPCSide()`, `KSP_NORM_DEFAULT`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 5768a4b9c5cSBarry Smith E*/ 5779371c9d4SSatish Balay typedef enum { 5789371c9d4SSatish Balay KSP_NORM_DEFAULT = -1, 5799371c9d4SSatish Balay KSP_NORM_NONE = 0, 5809371c9d4SSatish Balay KSP_NORM_PRECONDITIONED = 1, 5819371c9d4SSatish Balay KSP_NORM_UNPRECONDITIONED = 2, 5829371c9d4SSatish Balay KSP_NORM_NATURAL = 3 5839371c9d4SSatish Balay } KSPNormType; 5849e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 585014dd563SJed Brown PETSC_EXTERN const char *const *const KSPNormTypes; 586a21b2a99SBarry Smith 5871f7e983dSSatish Balay /*MC 5889793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 5898c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 5908c5b8ba0SBarry Smith be based on a norm of a residual etc. 5918c5b8ba0SBarry Smith 5928c5b8ba0SBarry Smith Level: advanced 5938c5b8ba0SBarry Smith 594a4d1885cSBarry Smith Note: 595a4d1885cSBarry Smith Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored 5968c5b8ba0SBarry Smith 5971cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 5988c5b8ba0SBarry Smith M*/ 5998c5b8ba0SBarry Smith 6001f7e983dSSatish Balay /*MC 6011957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 6028c5b8ba0SBarry Smith convergence test routine. 6038c5b8ba0SBarry Smith 6048c5b8ba0SBarry Smith Level: advanced 6058c5b8ba0SBarry Smith 6061cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 6078c5b8ba0SBarry Smith M*/ 6088c5b8ba0SBarry Smith 6091f7e983dSSatish Balay /*MC 610ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 6118c5b8ba0SBarry Smith convergence test routine. 6128c5b8ba0SBarry Smith 6138c5b8ba0SBarry Smith Level: advanced 6148c5b8ba0SBarry Smith 6151cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 6168c5b8ba0SBarry Smith M*/ 6178c5b8ba0SBarry Smith 6181f7e983dSSatish Balay /*MC 619ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 62087497f52SBarry Smith convergence test routine. This is only supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR` 6218c5b8ba0SBarry Smith 6228c5b8ba0SBarry Smith Level: advanced 6238c5b8ba0SBarry Smith 6241cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()` 6258c5b8ba0SBarry Smith M*/ 6268c5b8ba0SBarry Smith 627014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType); 628014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *); 629ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP, KSPNormType, PCSide, PetscInt); 630014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt); 631014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool); 6328a4b9c5cSBarry Smith 633edd03b47SJacob Faibussowitsch #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_NEG_CURVE", ) 634edd03b47SJacob Faibussowitsch #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_STEP_LENGTH", ) 635b5c1ed96SAlexander #define KSP_CONVERGED_RTOL_NORMAL_DEPRECATED KSP_CONVERGED_RTOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_RTOL_NORMAL_EQUATIONS", ) 636b5c1ed96SAlexander #define KSP_CONVERGED_ATOL_NORMAL_DEPRECATED KSP_CONVERGED_ATOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_ATOL_NORMAL_EQUATIONS", ) 6378a4b9c5cSBarry Smith /*E 638a4d1885cSBarry Smith KSPConvergedReason - reason a Krylov method was determined to have converged or diverged 63928ce4d24SBarry Smith 640a4d1885cSBarry Smith Values: 64178daaec8SBarry Smith + `KSP_CONVERGED_RTOL_NORMAL_EQUATIONS` - requested decrease in the residual of the normal equations, for `KSPLSQR` 64278daaec8SBarry Smith . `KSP_CONVERGED_ATOL_NORMAL_EQUATIONS` - requested absolute value in the residual of the normal equations, for `KSPLSQR` 643a4d1885cSBarry Smith . `KSP_CONVERGED_RTOL` - requested decrease in the residual 644a4d1885cSBarry Smith . `KSP_CONVERGED_ATOL` - requested absolute value in the residual 645a4d1885cSBarry Smith . `KSP_CONVERGED_ITS` - requested number of iterations 6464a221d59SStefano Zampini . `KSP_CONVERGED_NEG_CURVE` - see note below 647a4d1885cSBarry Smith . `KSP_CONVERGED_STEP_LENGTH` - see note below 6486e7835fbSStefano Zampini . `KSP_CONVERGED_HAPPY_BREAKDOWN` - happy breakdown (meaning early convergence of the `KSPType` occurred). 649*cc8c0b42SAlex Lindsay . `KSP_CONVERGED_USER` - the user has indicated convergence for an arbitrary reason 65078daaec8SBarry Smith . `KSP_DIVERGED_NULL` - breakdown when solving the Hessenberg system within `KSPGMRES` 651a4d1885cSBarry Smith . `KSP_DIVERGED_ITS` - requested number of iterations 65278daaec8SBarry Smith . `KSP_DIVERGED_DTOL` - large increase in the residual norm indicating the solution is diverging 653a4d1885cSBarry Smith . `KSP_DIVERGED_BREAKDOWN` - breakdown in the Krylov method 6548f47816eSPierre Jolivet . `KSP_DIVERGED_BREAKDOWN_BICG` - breakdown in the `KSPBCGS` Krylov method 655a4d1885cSBarry Smith . `KSP_DIVERGED_NONSYMMETRIC` - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry 65678daaec8SBarry Smith . `KSP_DIVERGED_INDEFINITE_PC` - the preconditioner was indefinite for a `KSPType` that requires it be definite, such as `KSPCG` 657a4d1885cSBarry Smith . `KSP_DIVERGED_NANORINF` - a not a number of infinity was detected in a vector during the computation 65878daaec8SBarry Smith . `KSP_DIVERGED_INDEFINITE_MAT` - the operator was indefinite for a `KSPType` that requires it be definite, such as `KSPCG` 659*cc8c0b42SAlex Lindsay . `KSP_DIVERGED_PC_FAILED` - the action of the preconditioner failed for some reason 660*cc8c0b42SAlex Lindsay - `KSP_DIVERGED_USER` - the user has indicated divergence for an arbitrary reason 661a4d1885cSBarry Smith 66216a05f60SBarry Smith Level: beginner 66316a05f60SBarry Smith 664a4d1885cSBarry Smith Note: 6656e7835fbSStefano Zampini The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by `KSPCG`, `KSPMINRES` and by 6666e7835fbSStefano Zampini the special `KSPNASH`, `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver. 66728ce4d24SBarry Smith 66895bd0b28SBarry Smith Developer Note: 66987497f52SBarry Smith The string versions of these are `KSPConvergedReasons`; if you change 6704d0a8057SBarry Smith any of the values here also change them that array of names. 67186c02ca4SBarry Smith 6721cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()` 67328ce4d24SBarry Smith E*/ 674d15094e1SBarry Smith typedef enum { /* converged */ 67578daaec8SBarry Smith KSP_CONVERGED_RTOL_NORMAL_DEPRECATED = 1, 67678daaec8SBarry Smith KSP_CONVERGED_RTOL_NORMAL_EQUATIONS = 1, 67778daaec8SBarry Smith KSP_CONVERGED_ATOL_NORMAL_DEPRECATED = 9, 67878daaec8SBarry Smith KSP_CONVERGED_ATOL_NORMAL_EQUATIONS = 9, 679d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 680d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 681b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 6824a221d59SStefano Zampini KSP_CONVERGED_NEG_CURVE = 5, 6834a221d59SStefano Zampini KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED = 5, 6844a221d59SStefano Zampini KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6, 6854a221d59SStefano Zampini KSP_CONVERGED_STEP_LENGTH = 6, 6864a221d59SStefano Zampini KSP_CONVERGED_HAPPY_BREAKDOWN = 7, 687*cc8c0b42SAlex Lindsay KSP_CONVERGED_USER = 8, 688d15094e1SBarry Smith /* diverged */ 689b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 690d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 691d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 692d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 693b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 694b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 695b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 6964d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 6976aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 698c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED = -11, 699aa4d2078SSatish Balay KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11, 700*cc8c0b42SAlex Lindsay KSP_DIVERGED_USER = -12, 701d15094e1SBarry Smith 7029371c9d4SSatish Balay KSP_CONVERGED_ITERATING = 0 7039371c9d4SSatish Balay } KSPConvergedReason; 704014dd563SJed Brown PETSC_EXTERN const char *const *KSPConvergedReasons; 705d15094e1SBarry Smith 706c838673bSBarry Smith /*MC 707af27ebaaSBarry Smith KSP_CONVERGED_RTOL - $||r|| \le rtol*||b||$ or $rtol*||b - A*x_0||$ if `KSPConvergedDefaultSetUIRNorm()` was called 708c838673bSBarry Smith 709c838673bSBarry Smith Level: beginner 710c838673bSBarry Smith 71195bd0b28SBarry Smith Notes: 71287497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 713c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 714c838673bSBarry Smith 2-norm of the residual for right preconditioning 715c838673bSBarry Smith 71687497f52SBarry Smith See also `KSP_CONVERGED_ATOL` which may apply before this tolerance. 717f9fed41fSBarry Smith 7181cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 719c838673bSBarry Smith M*/ 720c838673bSBarry Smith 721c838673bSBarry Smith /*MC 722af27ebaaSBarry Smith KSP_CONVERGED_ATOL - $||r|| \le atol$ 723c838673bSBarry Smith 724c838673bSBarry Smith Level: beginner 725c838673bSBarry Smith 72695bd0b28SBarry Smith Notes: 72787497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 728c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 729c838673bSBarry Smith 2-norm of the residual for right preconditioning 730c838673bSBarry Smith 73187497f52SBarry Smith See also `KSP_CONVERGED_RTOL` which may apply before this tolerance. 732c838673bSBarry Smith 7331cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 734c838673bSBarry Smith M*/ 735c838673bSBarry Smith 736c838673bSBarry Smith /*MC 737af27ebaaSBarry Smith KSP_DIVERGED_DTOL - $||r|| \ge dtol*||b||$ 738c838673bSBarry Smith 739c838673bSBarry Smith Level: beginner 740c838673bSBarry Smith 74195bd0b28SBarry Smith Note: 74287497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 743c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 744c838673bSBarry Smith 2-norm of the residual for right preconditioning 745c838673bSBarry Smith 7460241b274SPierre Jolivet .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_CONVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 747c838673bSBarry Smith M*/ 748c838673bSBarry Smith 749c838673bSBarry Smith /*MC 750c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 751c838673bSBarry Smith reached 752c838673bSBarry Smith 753c838673bSBarry Smith Level: beginner 754c838673bSBarry Smith 7551cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 756c838673bSBarry Smith M*/ 757c838673bSBarry Smith 758c838673bSBarry Smith /*MC 75987497f52SBarry Smith KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of 76087497f52SBarry Smith the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence 761af27ebaaSBarry Smith test routine is set in `KSP`. 762c838673bSBarry Smith 763c838673bSBarry Smith Level: beginner 764c838673bSBarry Smith 7651cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 766c838673bSBarry Smith M*/ 767c838673bSBarry Smith 768c838673bSBarry Smith /*MC 769c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 7701de96524SPierre Jolivet method could not continue to enlarge the Krylov space. Could be due to a singular matrix or 77187497f52SBarry Smith preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block 77246091a0eSPierre Jolivet are collinear. 773c838673bSBarry Smith 774c838673bSBarry Smith Level: beginner 775c838673bSBarry Smith 7761cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 777c838673bSBarry Smith M*/ 778c838673bSBarry Smith 779c838673bSBarry Smith /*MC 78087497f52SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the 781c838673bSBarry Smith method could not continue to enlarge the Krylov space. 782c838673bSBarry Smith 783c838673bSBarry Smith Level: beginner 784c838673bSBarry Smith 7851cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 786c838673bSBarry Smith M*/ 787c838673bSBarry Smith 788c838673bSBarry Smith /*MC 789c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 79087497f52SBarry Smith symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry 791c838673bSBarry Smith 792c838673bSBarry Smith Level: beginner 793c838673bSBarry Smith 7941cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 795c838673bSBarry Smith M*/ 796c838673bSBarry Smith 797c838673bSBarry Smith /*MC 798c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 79987497f52SBarry Smith positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to 8000b4b7b1cSBarry Smith be symmetric positive definite (SPD). 801c838673bSBarry Smith 802c838673bSBarry Smith Level: beginner 803c838673bSBarry Smith 80487497f52SBarry Smith Note: 805a4d1885cSBarry Smith This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force 80687497f52SBarry Smith the `PCICC` preconditioner to generate a positive definite preconditioner 807c838673bSBarry Smith 8081cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 809c838673bSBarry Smith M*/ 810c838673bSBarry Smith 811c838673bSBarry Smith /*MC 812c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a 8139fc87aa7SBarry Smith zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 81487497f52SBarry Smith such as `PCFIELDSPLIT`. 8159fc87aa7SBarry Smith 8169fc87aa7SBarry Smith Level: beginner 8179fc87aa7SBarry Smith 818a4d1885cSBarry Smith Note: 819a4d1885cSBarry Smith Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details. 8209fc87aa7SBarry Smith 8211cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 8229fc87aa7SBarry Smith M*/ 8239fc87aa7SBarry Smith 8249fc87aa7SBarry Smith /*MC 825a4d1885cSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called 826a4d1885cSBarry Smith while `KSPSolve()` is still running. 827c838673bSBarry Smith 828c838673bSBarry Smith Level: beginner 829c838673bSBarry Smith 8301cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 831c838673bSBarry Smith M*/ 832c838673bSBarry Smith 833014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *)); 834df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *)); 835df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *)); 8363ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *); 8378de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 8383eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 8398de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 8408de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 8418de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 8428de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 84354b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool); 8440059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 845014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *); 846c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **); 84794a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 8482a28d964SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetConvergedNegativeCurvature(KSP, PetscBool); 8492a28d964SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergedNegativeCurvature(KSP, PetscBool *); 850abef13c0SSatish Balay 851edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefault()", ) static inline void KSPDefaultConverged(void) 852d71ae5a4SJacob Faibussowitsch { /* never called */ 8539371c9d4SSatish Balay } 8548ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 855edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultDestroy()", ) static inline void KSPDefaultConvergedDestroy(void) 856d71ae5a4SJacob Faibussowitsch { /* never called */ 8579371c9d4SSatish Balay } 8588ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 859edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultCreate()", ) static inline void KSPDefaultConvergedCreate(void) 860d71ae5a4SJacob Faibussowitsch { /* never called */ 8619371c9d4SSatish Balay } 8628ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 863edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUIRNorm()", ) static inline void KSPDefaultConvergedSetUIRNorm(void) 864d71ae5a4SJacob Faibussowitsch { /* never called */ 8659371c9d4SSatish Balay } 8668ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 867edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUMIRNorm()", ) static inline void KSPDefaultConvergedSetUMIRNorm(void) 868d71ae5a4SJacob Faibussowitsch { /* never called */ 8699371c9d4SSatish Balay } 8708ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 871edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedSkip()", ) static inline void KSPSkipConverged(void) 872d71ae5a4SJacob Faibussowitsch { /* never called */ 8739371c9d4SSatish Balay } 8748ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 8758ea1b3e6SJed Brown 8760bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *); 877edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPComputeOperator()", ) static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B) 878d71ae5a4SJacob Faibussowitsch { 879f22e26b7SPierre Jolivet return KSPComputeOperator(A, PETSC_NULLPTR, B); 8809371c9d4SSatish Balay } 881d4fbbf0eSBarry Smith 88228ce4d24SBarry Smith /*E 883a4d1885cSBarry Smith KSPCGType - Determines what type of `KSPCG` to use 88428ce4d24SBarry Smith 885a4d1885cSBarry Smith Values: 886a4d1885cSBarry Smith + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric 887a4d1885cSBarry Smith - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian 888a4d1885cSBarry Smith 88916a05f60SBarry Smith Level: beginner 89016a05f60SBarry Smith 8911cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPCG`, `KSP`, `KSPCGSetType()` 89228ce4d24SBarry Smith E*/ 8939371c9d4SSatish Balay typedef enum { 8949371c9d4SSatish Balay KSP_CG_SYMMETRIC = 0, 8959371c9d4SSatish Balay KSP_CG_HERMITIAN = 1 8969371c9d4SSatish Balay } KSPCGType; 8976a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 89828ce4d24SBarry Smith 899014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType); 900014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool); 9018031f4adStmunson 902dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal); 903fb01098fSStefano Zampini PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal); 904dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *); 905dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *); 906fcae7a14Stmunson 90705de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *); 90805de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *); 909edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetMinEig()", ) static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x) 910d71ae5a4SJacob Faibussowitsch { 9119371c9d4SSatish Balay return KSPGLTRGetMinEig(ksp, x); 9129371c9d4SSatish Balay } 913edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetLambda()", ) static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x) 914d71ae5a4SJacob Faibussowitsch { 9159371c9d4SSatish Balay return KSPGLTRGetLambda(ksp, x); 9169371c9d4SSatish Balay } 9178031f4adStmunson 918014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]); 919ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]); 9201d6018f0SLisandro Dalcin 921f560b561SHong Zhang PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP)); 922014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP); 923014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP); 9243369ce9aSBarry Smith 925014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *); 9262f2e5d10SKris Buschelman 927014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec)); 928014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec)); 92903919abeSBarry Smith 930ba36735cSStefano Zampini /*S 931a4d1885cSBarry Smith KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods. 932f8a50e2bSBarry Smith 933a4d1885cSBarry Smith Level: intermediate 934f8a50e2bSBarry Smith 9350b4b7b1cSBarry Smith Note: 9360b4b7b1cSBarry Smith These methods generate initial guesses based on a series of previous, related, linear solves. For example, 9370b4b7b1cSBarry Smith in implicit time-stepping with `TS`. 9380b4b7b1cSBarry Smith 9399168452cSPierre Jolivet .seealso: [](ch_ksp), `KSPCreate()`, `KSPGuessSetType()`, `KSPGuessType` 940ba36735cSStefano Zampini S*/ 941ba36735cSStefano Zampini typedef struct _p_KSPGuess *KSPGuess; 94216a05f60SBarry Smith 943ba36735cSStefano Zampini /*J 94487497f52SBarry Smith KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods. 945ba36735cSStefano Zampini 946a4d1885cSBarry Smith Values: 947a4d1885cSBarry Smith + `KSPGUESSFISCHER` - methodology developed by Paul Fischer 9480b4b7b1cSBarry Smith - `KSPGUESSPOD` - methodology based on proper orthogonal decomposition (POD) 949a4d1885cSBarry Smith 95016a05f60SBarry Smith Level: intermediate 95116a05f60SBarry Smith 9521cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPGuess` 953ba36735cSStefano Zampini J*/ 954ba36735cSStefano Zampini typedef const char *KSPGuessType; 955ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer" 956ba36735cSStefano Zampini #define KSPGUESSPOD "pod" 957a4d1885cSBarry Smith 9581d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess)); 959ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess); 960ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *); 961ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer); 962ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *); 963ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *); 964ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType); 965ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *); 9668410009bSDavid Wells PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal); 967ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 968ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec); 969ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec); 970ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 971ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt); 972014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt); 973ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool); 974ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *); 975f8a50e2bSBarry Smith 976470b340bSDmitry Karpeev /*E 9777addb90fSBarry Smith MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement matrix assembly routines 978470b340bSDmitry Karpeev 979470b340bSDmitry Karpeev Level: intermediate 980470b340bSDmitry Karpeev 981a4d1885cSBarry Smith .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`, 9820b4b7b1cSBarry Smith `MatCreateSchurComplementPmat()`, `MatCreateSchurComplement()` 983470b340bSDmitry Karpeev E*/ 9849371c9d4SSatish Balay typedef enum { 9859371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_DIAG, 9869371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_LUMP, 9879371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, 9889371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_FULL 9899371c9d4SSatish Balay } MatSchurComplementAinvType; 990470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 991470b340bSDmitry Karpeev 992014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *); 993014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *); 994d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP); 995bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 996aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 997bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *); 998470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType); 999470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *); 10005bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *); 10015a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *); 1002470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *); 1003470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *); 10043f22127dSBarry Smith 100578e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *); 100678e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *); 1007bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatCreateLMVMDBFGS(MPI_Comm, PetscInt, PetscInt, Mat *); 1008bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatCreateLMVMDDFP(MPI_Comm, PetscInt, PetscInt, Mat *); 1009bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatCreateLMVMDQN(MPI_Comm, PetscInt, PetscInt, Mat *); 101078e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *); 1011864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1012864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1013864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1014864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1015864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1016cd929ea3SAlp Dener 1017cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 1018b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *); 1019cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 1020cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 1021e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat); 10220ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat); 1023cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 1024cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 1025cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 1026cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 1027cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 10282d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec); 1029cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 10301ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMGetLastUpdate(Mat, Vec *, Vec *); 1031cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *); 1032cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *); 1033cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *); 103492f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt); 1035bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatLMVMGetHistorySize(Mat, PetscInt *); 1036cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *); 1037cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *); 1038864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar); 1039e96a0efeSStefano Zampini 1040bbb72809SHansol Suh /*E 10411ca5963aSToby Isaac MatLMVMMultAlgorithm - The type of algorithm used for matrix-vector products and solves used internally by a `MatLMVM` matrix 10421ca5963aSToby Isaac 10431ca5963aSToby Isaac Values: 10441ca5963aSToby Isaac + `MAT_LMVM_MULT_RECURSIVE` - Use recursive formulas for products and solves 10451ca5963aSToby Isaac . `MAT_LMVM_MULT_DENSE` - Use dense formulas for products and solves when possible 10461ca5963aSToby Isaac - `MAT_LMVM_MULT_COMPACT_DENSE` - The same as `MATLMVM_MULT_DENSE`, but go further and ensure products and solves are computed in compact low-rank update form 10471ca5963aSToby Isaac 10481ca5963aSToby Isaac Level: advanced 10491ca5963aSToby Isaac 10501ca5963aSToby Isaac Options Database Keys: 10511ca5963aSToby Isaac . -mat_lmvm_mult_algorithm - the algorithm to use for multiplication (recursive, dense, compact_dense) 10521ca5963aSToby Isaac 10531ca5963aSToby Isaac .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSetMultAlgorithm()`, `MatLMVMGetMultAlgorithm()` 10541ca5963aSToby Isaac E*/ 10551ca5963aSToby Isaac typedef enum { 10561ca5963aSToby Isaac MAT_LMVM_MULT_RECURSIVE, 10571ca5963aSToby Isaac MAT_LMVM_MULT_DENSE, 10581ca5963aSToby Isaac MAT_LMVM_MULT_COMPACT_DENSE, 10591ca5963aSToby Isaac } MatLMVMMultAlgorithm; 10601ca5963aSToby Isaac 10611ca5963aSToby Isaac PETSC_EXTERN const char *const MatLMVMMultAlgorithms[]; 10621ca5963aSToby Isaac 10631ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSetMultAlgorithm(Mat, MatLMVMMultAlgorithm); 10641ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMGetMultAlgorithm(Mat, MatLMVMMultAlgorithm *); 10651ca5963aSToby Isaac 10661ca5963aSToby Isaac /*E 10671ca5963aSToby Isaac MatLMVMSymBroydenScaleType - Rescaling type for the initial Hessian of a symmetric Broyden matrix. 1068bbb72809SHansol Suh 1069bbb72809SHansol Suh Values: 1070b6982945SToby Isaac + `MAT_LMVM_SYMBROYDEN_SCALE_NONE` - no rescaling 1071b6982945SToby Isaac . `MAT_LMVM_SYMBROYDEN_SCALE_SCALAR` - scalar rescaling 1072b6982945SToby Isaac . `MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL` - diagonal rescaling 1073b6982945SToby Isaac . `MAT_LMVM_SYMBROYDEN_SCALE_USER` - same as `MAT_LMVM_SYMBROYDN_SCALE_NONE` 1074b6982945SToby Isaac - `MAT_LMVM_SYMBROYDEN_SCALE_DECIDE` - let PETSc decide rescaling 1075bbb72809SHansol Suh 1076bbb72809SHansol Suh Level: intermediate 1077bbb72809SHansol Suh 1078bbb72809SHansol Suh .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSymBroydenSetScaleType()` 1079bbb72809SHansol Suh E*/ 1080e96a0efeSStefano Zampini typedef enum { 1081e96a0efeSStefano Zampini MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0, 1082e96a0efeSStefano Zampini MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1, 1083e96a0efeSStefano Zampini MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2, 1084b6982945SToby Isaac MAT_LMVM_SYMBROYDEN_SCALE_USER = 3, 1085b6982945SToby Isaac MAT_LMVM_SYMBROYDEN_SCALE_DECIDE = 4 1086e96a0efeSStefano Zampini } MatLMVMSymBroydenScaleType; 1087e96a0efeSStefano Zampini PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[]; 1088e96a0efeSStefano Zampini 1089864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType); 10901ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenGetPhi(Mat, PetscReal *); 10911ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetPhi(Mat, PetscReal); 10921ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenGetPsi(Mat, PetscReal *); 10931ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenSetPsi(Mat, PetscReal); 1094cd929ea3SAlp Dener 1095bbb72809SHansol Suh /*E 10960b4b7b1cSBarry Smith MatLMVMDenseType - Memory storage strategy for dense variants of `MATLMVM`. 1097bbb72809SHansol Suh 1098bbb72809SHansol Suh Values: 1099bbb72809SHansol Suh + `MAT_LMVM_DENSE_REORDER` - reorders memory to minimize kernel launch 1100bbb72809SHansol Suh - `MAT_LMVM_DENSE_INPLACE` - computes inplace to minimize memory movement 1101bbb72809SHansol Suh 1102bbb72809SHansol Suh Level: intermediate 1103bbb72809SHansol Suh 1104bbb72809SHansol Suh .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMDenseSetType()` 1105bbb72809SHansol Suh E*/ 1106bbb72809SHansol Suh typedef enum { 1107bbb72809SHansol Suh MAT_LMVM_DENSE_REORDER, 1108bbb72809SHansol Suh MAT_LMVM_DENSE_INPLACE 1109bbb72809SHansol Suh } MatLMVMDenseType; 1110bbb72809SHansol Suh PETSC_EXTERN const char *const MatLMVMDenseTypes[]; 1111bbb72809SHansol Suh 1112bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatLMVMDenseSetType(Mat, MatLMVMDenseType); 1113bbb72809SHansol Suh 1114014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM); 1115014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool); 1116014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *); 1117014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *); 1118014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *); 11197b578ef6SBarry Smith 11207b578ef6SBarry Smith /*S 11218434afd1SBarry Smith KSPComputeRHSFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeRHS()` 11227b578ef6SBarry Smith 11237b578ef6SBarry Smith Calling Sequence: 11247b578ef6SBarry Smith + ksp - `ksp` context 11257b578ef6SBarry Smith . b - output vector 11267b578ef6SBarry Smith - ctx - [optional] user-defined function context 11277b578ef6SBarry Smith 11287b578ef6SBarry Smith Level: beginner 11297b578ef6SBarry Smith 11308434afd1SBarry Smith .seealso: [](ch_snes), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeInitialGuessFn`, `KSPComputeOperatorsFn` 11317b578ef6SBarry Smith S*/ 11328434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeRHSFn)(KSP ksp, Vec b, void *ctx); 11337b578ef6SBarry Smith 11348434afd1SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, KSPComputeRHSFn *, void *); 11357b578ef6SBarry Smith 11367b578ef6SBarry Smith /*S 11378434afd1SBarry Smith KSPComputeOperatorsFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeOperators()` 11387b578ef6SBarry Smith 11397b578ef6SBarry Smith Calling Sequence: 11407b578ef6SBarry Smith + ksp - `KSP` context 11417b578ef6SBarry Smith . A - the operator that defines the linear system 11427b578ef6SBarry Smith . P - an operator from which to build the preconditioner (often the same as `A`) 11437b578ef6SBarry Smith - ctx - [optional] user-defined function context 11447b578ef6SBarry Smith 11457b578ef6SBarry Smith Level: beginner 11467b578ef6SBarry Smith 11478434afd1SBarry Smith .seealso: [](ch_snes), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeInitialGuessFn` 11487b578ef6SBarry Smith S*/ 11498434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeOperatorsFn)(KSP ksp, Mat A, Mat P, void *ctx); 11507b578ef6SBarry Smith 11518434afd1SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, KSPComputeOperatorsFn, void *); 11527b578ef6SBarry Smith 11537b578ef6SBarry Smith /*S 11548434afd1SBarry Smith KSPComputeInitialGuessFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeInitialGuess()` 11557b578ef6SBarry Smith 11567b578ef6SBarry Smith Calling Sequence: 11577b578ef6SBarry Smith + ksp - `ksp` context 11587b578ef6SBarry Smith . x - output vector 11597b578ef6SBarry Smith - ctx - [optional] user-defined function context 11607b578ef6SBarry Smith 11617b578ef6SBarry Smith Level: beginner 11627b578ef6SBarry Smith 11638434afd1SBarry Smith .seealso: [](ch_snes), `KSP`, `KSPSetComputeInitialGuess()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeOperatorsFn` 11647b578ef6SBarry Smith S*/ 11658434afd1SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeInitialGuessFn)(KSP ksp, Vec x, void *ctx); 11667b578ef6SBarry Smith 11678434afd1SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, KSPComputeInitialGuessFn *, void *); 11688434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, KSPComputeOperatorsFn *, void *); 11698434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, KSPComputeOperatorsFn **, void *); 11708434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, KSPComputeRHSFn *, void *); 11718434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, KSPComputeRHSFn **, void *); 11728434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, KSPComputeInitialGuessFn *, void *); 11738434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, KSPComputeInitialGuessFn **, void *); 11746c699258SBarry Smith 117502b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec); 11769371c9d4SSatish 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); 1177953494dbSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM, DM, PetscInt, const char **, Vec[], ScatterMode mode); 1178557cf195SMatthew G. Knepley 11792b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *); 11802b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal); 11814bf303faSJacob Faibussowitsch 11824bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCBJKOKKOSSetKSP(PC, KSP); 11834bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCBJKOKKOSGetKSP(PC, KSP *); 11845d83a8b1SBarry Smith 11855d83a8b1SBarry Smith PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM); 1186