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 1516a05f60SBarry Smith linear solves in PETSc (even those such as direct factorization-based 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 2116a05f60SBarry Smith `KSPType` of `KSPPREONLY` (or equivalently `KSPNONE`), meaning that only application of the preconditioner is used as the linear solver. 228f6c3df8SBarry Smith 231cc06b55SBarry Smith .seealso: [](doc_linsolve), [](ch_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 321cc06b55SBarry Smith .seealso: [](doc_linsolve), [](ch_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" 48*edd03b47SJacob Faibussowitsch #define KSPCGNASH PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPNASH", ) "nash" 49*edd03b47SJacob Faibussowitsch #define KSPCGSTCG PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSTCG", ) "stcg" 50*edd03b47SJacob Faibussowitsch #define KSPCGGLTR PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSGLTR", ) "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); 101bbbebc2cSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolveTranspose(KSP, Mat, Mat); 1023e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt); 103*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPSetMatSolveBatchSize()", ) 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 *); 108*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPGetMatSolveBatchSize()", ) 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 *); 13225c92fe2SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetMinimumIterations(KSP, PetscInt); 13325c92fe2SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetMinimumIterations(KSP, PetscInt *); 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool); 135014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *); 136014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool); 137014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *); 138014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool); 1397d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool); 140c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *); 141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool); 142c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *); 143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *); 144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *); 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *); 146014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *); 147734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *); 1482a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **); 149*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 6, 0, "KSPCreateVecs()", ) static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y) 150d71ae5a4SJacob Faibussowitsch { 1519371c9d4SSatish Balay return KSPCreateVecs(ksp, n, x, m, y); 1529371c9d4SSatish Balay } 1532eac72dbSBarry Smith 154d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *); 155d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *); 156d4211eb9SBarry Smith 157014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC); 158014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *); 159aabeff55SBarry Smith 160014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal); 161014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void **)); 162014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 1633ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *); 16494a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *); 165014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscInt, PetscBool); 16694a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *); 16794a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscInt, PetscBool); 1684b0e389bSBarry Smith 169fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *); 170fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *); 171fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 172fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt); 173fa0ddf94SBarry Smith 174014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *); 175cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP); 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 177014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 178014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 179014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]); 180285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]); 181b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *); 182b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *); 183b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *); 184b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *); 185014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *); 1864a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *); 187f3b08a26SMatthew G. Knepley /* 188f3b08a26SMatthew G. Knepley PCMGCoarseList contains the list of coarse space constructor currently registered 189f3b08a26SMatthew G. Knepley These are added with PCMGRegisterCoarseSpaceConstructor() 190f3b08a26SMatthew G. Knepley */ 191f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList; 1922b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *)); 1932b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *)); 194f3b08a26SMatthew G. Knepley 195014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *); 196014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *); 1972eac72dbSBarry Smith 198f2edd1f0SMalachi Phillips /*E 19916a05f60SBarry Smith KSPChebyshevKind - Which kind of Chebyshev polynomial to use 200f2edd1f0SMalachi Phillips 201f2edd1f0SMalachi Phillips Values: 202f2edd1f0SMalachi Phillips + `KSP_CHEBYSHEV_FIRST` - "classic" first-kind Chebyshev polynomial 203f2edd1f0SMalachi Phillips . `KSP_CHEBYSHEV_FOURTH` - fourth-kind Chebyshev polynomial 204f2edd1f0SMalachi Phillips - `KSP_CHEBYSHEV_OPT_FOURTH` - optimized fourth-kind Chebyshev polynomial 205f2edd1f0SMalachi Phillips 206f2edd1f0SMalachi Phillips Level: intermediate 207f2edd1f0SMalachi Phillips 2081cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPCHEBYSHEV`, `KSPChebyshevSetKind()`, `KSPChebyshevGetKind()` 209f2edd1f0SMalachi Phillips E*/ 210f2edd1f0SMalachi Phillips typedef enum { 211f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_FIRST, 212f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_FOURTH, 213f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_OPT_FOURTH 214f2edd1f0SMalachi Phillips } KSPChebyshevKind; 215f2edd1f0SMalachi Phillips 216014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal); 217014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool); 218014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal); 21958450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal); 220b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool); 221f2edd1f0SMalachi Phillips PETSC_EXTERN PetscErrorCode KSPChebyshevSetKind(KSP, KSPChebyshevKind); 22216a05f60SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevGetKind(KSP, KSPChebyshevKind *); 22358450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *); 224014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *); 225d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *); 226d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]); 2277d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]); 2284b0e389bSBarry Smith 229640f4f53SPatrick Sanan /*E 230640f4f53SPatrick Sanan 23106137d0aSPatrick Sanan KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 232640f4f53SPatrick Sanan 23367b8a455SSatish Balay Values: 234a1cb98faSBarry Smith + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions 235a1cb98faSBarry Smith - `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 236640f4f53SPatrick Sanan 2372b26979fSBarry Smith Level: intermediate 238640f4f53SPatrick Sanan 2391cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()` 240640f4f53SPatrick Sanan E*/ 2419371c9d4SSatish Balay typedef enum { 2429371c9d4SSatish Balay KSP_FCD_TRUNC_TYPE_STANDARD, 2439371c9d4SSatish Balay KSP_FCD_TRUNC_TYPE_NOTAY 2449371c9d4SSatish Balay } KSPFCDTruncationType; 24506137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 246640f4f53SPatrick Sanan 247640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt); 248640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *); 249640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt); 250640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *); 25106137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType); 25206137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *); 253640f4f53SPatrick Sanan 254390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt); 255390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *); 256390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt); 257390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *); 258390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType); 259390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *); 260390d8e47SPatrick Sanan 261fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt); 262fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *); 263fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt); 264fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *); 265fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType); 266fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *); 267fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool); 268fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *); 26983f0b094SBarry Smith 270014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *); 272014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal); 273e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal); 2749f236934SBarry Smith 275014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 276014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt)); 277014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt)); 278014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt); 279014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt); 2801d73ed98SBarry Smith 281014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt); 282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2831d73ed98SBarry Smith 284483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar); 285483d6965SPatrick Sanan 286014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt); 287014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *); 288014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 28958ad63f4SBarry Smith 290f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal); 291f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *); 292f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool); 293f87a0b54SStefano Zampini 29404d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *); 29504d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC); 29604d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *); 297568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat); 298e82af88dSprj- 2996cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat); 3006cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *); 3016cac28cbSPierre Jolivet #if PetscDefined(HAVE_HPDDM) 302*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMSetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) 303d71ae5a4SJacob Faibussowitsch { 3049371c9d4SSatish Balay return KSPHPDDMSetDeflationMat(ksp, U); 3059371c9d4SSatish Balay } 306*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMGetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) 307d71ae5a4SJacob Faibussowitsch { 3089371c9d4SSatish Balay return KSPHPDDMGetDeflationMat(ksp, U); 3099371c9d4SSatish Balay } 3106cac28cbSPierre Jolivet #endif 311*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPMatSolve()", ) static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) 312d71ae5a4SJacob Faibussowitsch { 3139371c9d4SSatish Balay return KSPMatSolve(ksp, B, X); 3149371c9d4SSatish Balay } 315d55ab951SPierre Jolivet /*E 31687497f52SBarry Smith KSPHPDDMType - Type of Krylov method used by `KSPHPDDM` 317d55ab951SPierre Jolivet 318d55ab951SPierre Jolivet Values: 319a4d1885cSBarry Smith + `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method 320a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BGMRES` - block GMRES 321a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_CG` - Conjugate Gradient 322a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BCG` - block CG 323a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_GCRODR` - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting 324a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BGCRODR` - block GCRODR 325a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BFBCG` - breakdown-free BCG 326a4d1885cSBarry Smith - `KSP_HPDDM_TYPE_PREONLY` - apply the preconditioner only 327d55ab951SPierre Jolivet 32816a05f60SBarry Smith Level: intermediate 32916a05f60SBarry Smith 3301cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPHPDDM`, `KSPHPDDMSetType()` 331d55ab951SPierre Jolivet E*/ 332d55ab951SPierre Jolivet typedef enum { 333d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GMRES = 0, 334d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGMRES = 1, 335d55ab951SPierre Jolivet KSP_HPDDM_TYPE_CG = 2, 336d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BCG = 3, 337d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GCRODR = 4, 338d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGCRODR = 5, 339d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BFBCG = 6, 340d55ab951SPierre Jolivet KSP_HPDDM_TYPE_PREONLY = 7 341d55ab951SPierre Jolivet } KSPHPDDMType; 342d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[]; 343a4d1885cSBarry Smith 3442dd49c90SPierre Jolivet /*E 34587497f52SBarry Smith KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM` 3462dd49c90SPierre Jolivet 3472dd49c90SPierre Jolivet Values: 348a4d1885cSBarry Smith + `KSP_HPDDM_PRECISION_HALF` - default when PETSc is configured `--with-precision=__fp16` 349a4d1885cSBarry Smith . `KSP_HPDDM_PRECISION_SINGLE` - default when PETSc is configured `--with-precision=single` 350a4d1885cSBarry Smith . `KSP_HPDDM_PRECISION_DOUBLE` - default when PETSc is configured `--with-precision=double` 351a4d1885cSBarry Smith - `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128` 3522dd49c90SPierre Jolivet 35316a05f60SBarry Smith Level: intermediate 35416a05f60SBarry Smith 3551cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPHPDDM` 3562dd49c90SPierre Jolivet E*/ 3572dd49c90SPierre Jolivet typedef enum { 3582dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_HALF = 0, 3592dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_SINGLE = 1, 3602dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_DOUBLE = 2, 3612dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_QUADRUPLE = 3 3622dd49c90SPierre Jolivet } KSPHPDDMPrecision; 363d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType); 364d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *); 365e82af88dSprj- 366b49fd9e1SBarry Smith /*E 36716a05f60SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed in the GMRES solvers 368b49fd9e1SBarry Smith 369a4d1885cSBarry Smith Values: 370a4d1885cSBarry Smith + `KSP_GMRES_CGS_REFINE_NEVER` - one step of classical Gram-Schmidt 371a4d1885cSBarry Smith . `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria 372a4d1885cSBarry Smith - `KSP_GMRES_CGS_REFINE_ALWAYS` - always perform two steps 373b49fd9e1SBarry Smith 37416a05f60SBarry Smith Level: advanced 37516a05f60SBarry Smith 3761cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 377a4d1885cSBarry Smith `KSPGMRESGetOrthogonalization()`, 378a4d1885cSBarry Smith `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()` 379b49fd9e1SBarry Smith E*/ 3809371c9d4SSatish Balay typedef enum { 3819371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_NEVER, 3829371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_IFNEEDED, 3839371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_ALWAYS 3849371c9d4SSatish Balay } KSPGMRESCGSRefinementType; 3856a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 38616a05f60SBarry Smith 3871f7e983dSSatish Balay /*MC 3881957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 3898c5b8ba0SBarry Smith 3908c5b8ba0SBarry Smith Level: advanced 3918c5b8ba0SBarry Smith 39287497f52SBarry Smith Note: 39387497f52SBarry Smith Possibly unstable, but the fastest to compute 3948c5b8ba0SBarry Smith 3951cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 396a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 397db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 398db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 3998c5b8ba0SBarry Smith M*/ 4008c5b8ba0SBarry Smith 4011f7e983dSSatish Balay /*MC 4028c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 4038c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 4048c5b8ba0SBarry Smith poor orthogonality. 4058c5b8ba0SBarry Smith 4068c5b8ba0SBarry Smith Level: advanced 4078c5b8ba0SBarry Smith 408a4d1885cSBarry Smith Note: 409a4d1885cSBarry Smith This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to 4108c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 4118c5b8ba0SBarry Smith 4121cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 413a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 414db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 415db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 4168c5b8ba0SBarry Smith M*/ 4178c5b8ba0SBarry Smith 4181f7e983dSSatish Balay /*MC 41916a05f60SBarry Smith KSP_GMRES_CGS_REFINE_ALWAYS - Do two steps of the classical (unmodified) Gram-Schmidt process. 4208c5b8ba0SBarry Smith 4218c5b8ba0SBarry Smith Level: advanced 4228c5b8ba0SBarry Smith 42387497f52SBarry Smith Notes: 42487497f52SBarry Smith This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice 42587497f52SBarry Smith but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`. 4268c5b8ba0SBarry Smith 4278c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 4288c5b8ba0SBarry Smith 4291cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 430a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 431db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 432db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 4338c5b8ba0SBarry Smith M*/ 4348c5b8ba0SBarry Smith 435014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType); 436014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *); 43708480c60SBarry Smith 438014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *); 439014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *); 440014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 441c38d4ed2SBarry Smith 442014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal); 443014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *); 444014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *); 445121fd945SKris Buschelman 446014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal); 447014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool); 448014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt); 449e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool); 450d9492815SBarry Smith 451014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 45287d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP); 453014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 4542eac72dbSBarry Smith 455798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *); 456798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm, const char[], const char[], const char[], PetscInt, const char *[], int, int, int, int, PetscDrawLG *); 457798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 458798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 459798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 460798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 461798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 462798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 463798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 464798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 465798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 466798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 467798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 468798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 469798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 470798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 471798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 472798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 473798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 474798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 475798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 476fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 477798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 478*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorResidual()", ) static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 479d71ae5a4SJacob Faibussowitsch { 4809371c9d4SSatish Balay return KSPMonitorResidual(ksp, n, rnorm, vf); 4819371c9d4SSatish Balay } 482*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidual()", ) static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 483d71ae5a4SJacob Faibussowitsch { 4849371c9d4SSatish Balay return KSPMonitorTrueResidual(ksp, n, rnorm, vf); 4859371c9d4SSatish Balay } 486*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidualMax()", ) static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 487d71ae5a4SJacob Faibussowitsch { 4889371c9d4SSatish Balay return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf); 4899371c9d4SSatish Balay } 490798534f6SMatthew G. Knepley 491798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *); 492341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *); 493341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **); 494341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *); 495341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal); 496e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *); 497e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **); 498e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **); 49984cb2905SBarry Smith 500014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec); 501014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec); 502c01c455dSBarry Smith 50323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat); 50423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *); 505014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *); 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]); 507014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]); 508014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]); 5091eb62cbbSBarry Smith 510014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool); 511014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *); 512014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool); 513014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *); 5141f7f0c4fSBarry Smith 515014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer); 51655849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer); 517fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]); 51819a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer); 519c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **)); 5201b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP); 521c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP); 52294a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer); 5231b2b9847SBarry Smith 524*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonView()", ) static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v) 525d71ae5a4SJacob Faibussowitsch { 5269371c9d4SSatish Balay return KSPConvergedReasonView(ksp, v); 5279371c9d4SSatish Balay } 528*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) 529d71ae5a4SJacob Faibussowitsch { 5309371c9d4SSatish Balay return KSPConvergedReasonViewFromOptions(ksp); 5319371c9d4SSatish Balay } 53255849f57SBarry Smith 53355849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 5341eb62cbbSBarry Smith 535dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool); 5360e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool); 537014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *); 538884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *); 539798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 540798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 541798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 542db9b2ab1SHong Zhang 543014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *); 544014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *); 54568ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *); 54683ab6a24SBarry Smith 54728ce4d24SBarry Smith /*E 548a4d1885cSBarry Smith KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence 5498a4b9c5cSBarry Smith test routines. 5508a4b9c5cSBarry Smith 551a4d1885cSBarry Smith Values: 552a4d1885cSBarry Smith + `KSP_NORM_DEFAULT` - use the default for the current `KSPType` 553a4d1885cSBarry Smith . `KSP_NORM_NONE` - use no norm calculation 554a4d1885cSBarry Smith . `KSP_NORM_PRECONDITIONED` - use the preconditioned residual norm 555a4d1885cSBarry Smith . `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm 556a4d1885cSBarry Smith - `KSP_NORM_NATURAL` - use the natural norm (the norm induced by the linear operator) 557a4d1885cSBarry Smith 5588a4b9c5cSBarry Smith Level: advanced 5598a4b9c5cSBarry Smith 560a4d1885cSBarry Smith Note: 561a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 56216a05f60SBarry Smith depending on whether left or right preconditioning is used, see `KSPSetPCSide()` 563a3f661c8SBarry Smith 56487497f52SBarry Smith Developer Note: 565a4d1885cSBarry Smith This must match the values in petsc/finclude/petscksp.h 5668a4b9c5cSBarry Smith 5671cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`, 568db781477SPatrick Sanan `KSPSetConvergenceTest()`, `KSPSetPCSide()` 5698a4b9c5cSBarry Smith E*/ 5709371c9d4SSatish Balay typedef enum { 5719371c9d4SSatish Balay KSP_NORM_DEFAULT = -1, 5729371c9d4SSatish Balay KSP_NORM_NONE = 0, 5739371c9d4SSatish Balay KSP_NORM_PRECONDITIONED = 1, 5749371c9d4SSatish Balay KSP_NORM_UNPRECONDITIONED = 2, 5759371c9d4SSatish Balay KSP_NORM_NATURAL = 3 5769371c9d4SSatish Balay } KSPNormType; 5779e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 578014dd563SJed Brown PETSC_EXTERN const char *const *const KSPNormTypes; 579a21b2a99SBarry Smith 5801f7e983dSSatish Balay /*MC 5819793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 5828c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 5838c5b8ba0SBarry Smith be based on a norm of a residual etc. 5848c5b8ba0SBarry Smith 5858c5b8ba0SBarry Smith Level: advanced 5868c5b8ba0SBarry Smith 587a4d1885cSBarry Smith Note: 588a4d1885cSBarry Smith Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored 5898c5b8ba0SBarry Smith 5901cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 5918c5b8ba0SBarry Smith M*/ 5928c5b8ba0SBarry Smith 5931f7e983dSSatish Balay /*MC 5941957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 5958c5b8ba0SBarry Smith convergence test routine. 5968c5b8ba0SBarry Smith 5978c5b8ba0SBarry Smith Level: advanced 5988c5b8ba0SBarry Smith 5991cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 6008c5b8ba0SBarry Smith M*/ 6018c5b8ba0SBarry Smith 6021f7e983dSSatish Balay /*MC 603ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 6048c5b8ba0SBarry Smith convergence test routine. 6058c5b8ba0SBarry Smith 6068c5b8ba0SBarry Smith Level: advanced 6078c5b8ba0SBarry Smith 6081cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 6098c5b8ba0SBarry Smith M*/ 6108c5b8ba0SBarry Smith 6111f7e983dSSatish Balay /*MC 612ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 61387497f52SBarry Smith convergence test routine. This is only supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR` 6148c5b8ba0SBarry Smith 6158c5b8ba0SBarry Smith Level: advanced 6168c5b8ba0SBarry Smith 6171cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()` 6188c5b8ba0SBarry Smith M*/ 6198c5b8ba0SBarry Smith 620014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType); 621014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *); 622014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt); 623014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt); 624014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool); 6258a4b9c5cSBarry Smith 626*edd03b47SJacob Faibussowitsch #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_NEG_CURVE", ) 627*edd03b47SJacob Faibussowitsch #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_STEP_LENGTH", ) 628*edd03b47SJacob Faibussowitsch #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM(3, 11, 0, "KSP_DIVERGED_PC_FAILED", ) 6298a4b9c5cSBarry Smith /*E 630a4d1885cSBarry Smith KSPConvergedReason - reason a Krylov method was determined to have converged or diverged 63128ce4d24SBarry Smith 632a4d1885cSBarry Smith Values: 633a4d1885cSBarry Smith + `KSP_CONVERGED_RTOL_NORMAL` - requested decrease in the residual for the normal equations 634a4d1885cSBarry Smith . `KSP_CONVERGED_ATOL_NORMAL` - requested absolute value in the residual for the normal equations 635a4d1885cSBarry Smith . `KSP_CONVERGED_RTOL` - requested decrease in the residual 636a4d1885cSBarry Smith . `KSP_CONVERGED_ATOL` - requested absolute value in the residual 637a4d1885cSBarry Smith . `KSP_CONVERGED_ITS` - requested number of iterations 6384a221d59SStefano Zampini . `KSP_CONVERGED_NEG_CURVE` - see note below 639a4d1885cSBarry Smith . `KSP_CONVERGED_STEP_LENGTH` - see note below 640a4d1885cSBarry Smith . `KSP_CONVERGED_HAPPY_BREAKDOWN` - happy breakdown (meaning early convergence of the `KSPType` occurred. 641a4d1885cSBarry Smith . `KSP_DIVERGED_NULL` - breakdown when solving the Hessenberg system within GMRES 642a4d1885cSBarry Smith . `KSP_DIVERGED_ITS` - requested number of iterations 643a4d1885cSBarry Smith . `KSP_DIVERGED_DTOL` - large increase in the residual norm 644a4d1885cSBarry Smith . `KSP_DIVERGED_BREAKDOWN` - breakdown in the Krylov method 645a4d1885cSBarry Smith . `KSP_DIVERGED_BREAKDOWN_BICG` - breakdown in the `KSPBGCS` Krylov method 646a4d1885cSBarry Smith . `KSP_DIVERGED_NONSYMMETRIC` - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry 647a4d1885cSBarry Smith . `KSP_DIVERGED_INDEFINITE_PC` - the preconditioner was indefinite for a `KSPType` that requires it be definite 648a4d1885cSBarry Smith . `KSP_DIVERGED_NANORINF` - a not a number of infinity was detected in a vector during the computation 649a4d1885cSBarry Smith . `KSP_DIVERGED_INDEFINITE_MAT` - the operator was indefinite for a `KSPType` that requires it be definite 650a4d1885cSBarry Smith - `KSP_DIVERGED_PC_FAILED` - the action of the preconditioner failed for some reason 651a4d1885cSBarry Smith 65216a05f60SBarry Smith Level: beginner 65316a05f60SBarry Smith 654a4d1885cSBarry Smith Note: 6554a221d59SStefano Zampini The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by the special `KSPNASH`, 656a4d1885cSBarry Smith `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver. 65728ce4d24SBarry Smith 65895452b02SPatrick Sanan Developer Notes: 659a4d1885cSBarry Smith This must match the values in petsc/finclude/petscksp.h 6604d0a8057SBarry Smith 66187497f52SBarry Smith The string versions of these are `KSPConvergedReasons`; if you change 6624d0a8057SBarry Smith any of the values here also change them that array of names. 66386c02ca4SBarry Smith 6641cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()` 66528ce4d24SBarry Smith E*/ 666d15094e1SBarry Smith typedef enum { /* converged */ 6679ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 6689ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 669d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 670d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 671b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 6724a221d59SStefano Zampini KSP_CONVERGED_NEG_CURVE = 5, 6734a221d59SStefano Zampini KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED = 5, 6744a221d59SStefano Zampini KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6, 6754a221d59SStefano Zampini KSP_CONVERGED_STEP_LENGTH = 6, 6764a221d59SStefano Zampini KSP_CONVERGED_HAPPY_BREAKDOWN = 7, 677d15094e1SBarry Smith /* diverged */ 678b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 679d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 680d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 681d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 682b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 683b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 684b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 6854d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 6866aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 687c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED = -11, 688aa4d2078SSatish Balay KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11, 689d15094e1SBarry Smith 6909371c9d4SSatish Balay KSP_CONVERGED_ITERATING = 0 6919371c9d4SSatish Balay } KSPConvergedReason; 692014dd563SJed Brown PETSC_EXTERN const char *const *KSPConvergedReasons; 693d15094e1SBarry Smith 694c838673bSBarry Smith /*MC 69587497f52SBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if `KSPConvergedDefaultSetUIRNorm()` was called 696c838673bSBarry Smith 697c838673bSBarry Smith Level: beginner 698c838673bSBarry Smith 69987497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 700c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 701c838673bSBarry Smith 2-norm of the residual for right preconditioning 702c838673bSBarry Smith 70387497f52SBarry Smith See also `KSP_CONVERGED_ATOL` which may apply before this tolerance. 704f9fed41fSBarry Smith 7051cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 706c838673bSBarry Smith M*/ 707c838673bSBarry Smith 708c838673bSBarry Smith /*MC 709c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 710c838673bSBarry Smith 711c838673bSBarry Smith Level: beginner 712c838673bSBarry Smith 71387497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 714c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 715c838673bSBarry Smith 2-norm of the residual for right preconditioning 716c838673bSBarry Smith 71787497f52SBarry Smith See also `KSP_CONVERGED_RTOL` which may apply before this tolerance. 718c838673bSBarry Smith 7191cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 720c838673bSBarry Smith M*/ 721c838673bSBarry Smith 722c838673bSBarry Smith /*MC 723c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 724c838673bSBarry Smith 725c838673bSBarry Smith Level: beginner 726c838673bSBarry Smith 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 731c838673bSBarry Smith Level: beginner 732c838673bSBarry Smith 7331cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 734c838673bSBarry Smith M*/ 735c838673bSBarry Smith 736c838673bSBarry Smith /*MC 737c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 738c838673bSBarry Smith reached 739c838673bSBarry Smith 740c838673bSBarry Smith Level: beginner 741c838673bSBarry Smith 7421cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 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 7521cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 753c838673bSBarry Smith M*/ 754c838673bSBarry Smith 755c838673bSBarry Smith /*MC 756c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 7571de96524SPierre Jolivet method could not continue to enlarge the Krylov space. Could be due to a singular matrix or 75887497f52SBarry Smith preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block 7591de96524SPierre Jolivet are colinear. 760c838673bSBarry Smith 761c838673bSBarry Smith Level: beginner 762c838673bSBarry Smith 7631cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 764c838673bSBarry Smith M*/ 765c838673bSBarry Smith 766c838673bSBarry Smith /*MC 76787497f52SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the 768c838673bSBarry Smith method could not continue to enlarge the Krylov space. 769c838673bSBarry Smith 770c838673bSBarry Smith Level: beginner 771c838673bSBarry Smith 7721cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 773c838673bSBarry Smith M*/ 774c838673bSBarry Smith 775c838673bSBarry Smith /*MC 776c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 77787497f52SBarry Smith symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry 778c838673bSBarry Smith 779c838673bSBarry Smith Level: beginner 780c838673bSBarry Smith 7811cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 782c838673bSBarry Smith M*/ 783c838673bSBarry Smith 784c838673bSBarry Smith /*MC 785c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 78687497f52SBarry Smith positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to 787c838673bSBarry Smith be positive definite 788c838673bSBarry Smith 789c838673bSBarry Smith Level: beginner 790c838673bSBarry Smith 79187497f52SBarry Smith Note: 792a4d1885cSBarry Smith This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force 79387497f52SBarry Smith the `PCICC` preconditioner to generate a positive definite preconditioner 794c838673bSBarry Smith 7951cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 796c838673bSBarry Smith M*/ 797c838673bSBarry Smith 798c838673bSBarry Smith /*MC 799c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a 8009fc87aa7SBarry Smith zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 80187497f52SBarry Smith such as `PCFIELDSPLIT`. 8029fc87aa7SBarry Smith 8039fc87aa7SBarry Smith Level: beginner 8049fc87aa7SBarry Smith 805a4d1885cSBarry Smith Note: 806a4d1885cSBarry Smith Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details. 8079fc87aa7SBarry Smith 8081cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 8099fc87aa7SBarry Smith M*/ 8109fc87aa7SBarry Smith 8119fc87aa7SBarry Smith /*MC 812a4d1885cSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called 813a4d1885cSBarry Smith while `KSPSolve()` is still running. 814c838673bSBarry Smith 815c838673bSBarry Smith Level: beginner 816c838673bSBarry Smith 8171cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 818c838673bSBarry Smith M*/ 819c838673bSBarry Smith 820014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *)); 821df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *)); 822df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *)); 8233ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *); 8248de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 8253eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 8268de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 8278de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 8288de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 8298de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 83054b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool); 8310059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 832014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *); 833c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **); 83494a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 8352a28d964SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetConvergedNegativeCurvature(KSP, PetscBool); 8362a28d964SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergedNegativeCurvature(KSP, PetscBool *); 837abef13c0SSatish Balay 838*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefault()", ) static inline void KSPDefaultConverged(void) 839d71ae5a4SJacob Faibussowitsch { /* never called */ 8409371c9d4SSatish Balay } 8418ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 842*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultDestroy()", ) static inline void KSPDefaultConvergedDestroy(void) 843d71ae5a4SJacob Faibussowitsch { /* never called */ 8449371c9d4SSatish Balay } 8458ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 846*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultCreate()", ) static inline void KSPDefaultConvergedCreate(void) 847d71ae5a4SJacob Faibussowitsch { /* never called */ 8489371c9d4SSatish Balay } 8498ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 850*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUIRNorm()", ) static inline void KSPDefaultConvergedSetUIRNorm(void) 851d71ae5a4SJacob Faibussowitsch { /* never called */ 8529371c9d4SSatish Balay } 8538ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 854*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUMIRNorm()", ) static inline void KSPDefaultConvergedSetUMIRNorm(void) 855d71ae5a4SJacob Faibussowitsch { /* never called */ 8569371c9d4SSatish Balay } 8578ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 858*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedSkip()", ) static inline void KSPSkipConverged(void) 859d71ae5a4SJacob Faibussowitsch { /* never called */ 8609371c9d4SSatish Balay } 8618ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 8628ea1b3e6SJed Brown 8630bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *); 864*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPComputeOperator()", ) static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B) 865d71ae5a4SJacob Faibussowitsch { 866f22e26b7SPierre Jolivet return KSPComputeOperator(A, PETSC_NULLPTR, B); 8679371c9d4SSatish Balay } 868d4fbbf0eSBarry Smith 86928ce4d24SBarry Smith /*E 870a4d1885cSBarry Smith KSPCGType - Determines what type of `KSPCG` to use 87128ce4d24SBarry Smith 872a4d1885cSBarry Smith Values: 873a4d1885cSBarry Smith + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric 874a4d1885cSBarry Smith - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian 875a4d1885cSBarry Smith 87616a05f60SBarry Smith Level: beginner 87716a05f60SBarry Smith 8781cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPCG`, `KSP`, `KSPCGSetType()` 87928ce4d24SBarry Smith E*/ 8809371c9d4SSatish Balay typedef enum { 8819371c9d4SSatish Balay KSP_CG_SYMMETRIC = 0, 8829371c9d4SSatish Balay KSP_CG_HERMITIAN = 1 8839371c9d4SSatish Balay } KSPCGType; 8846a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 88528ce4d24SBarry Smith 886014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType); 887014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool); 8888031f4adStmunson 889dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal); 890fb01098fSStefano Zampini PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal); 891dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *); 892dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *); 893fcae7a14Stmunson 89405de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *); 89505de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *); 896*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetMinEig()", ) static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x) 897d71ae5a4SJacob Faibussowitsch { 8989371c9d4SSatish Balay return KSPGLTRGetMinEig(ksp, x); 8999371c9d4SSatish Balay } 900*edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetLambda()", ) static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x) 901d71ae5a4SJacob Faibussowitsch { 9029371c9d4SSatish Balay return KSPGLTRGetLambda(ksp, x); 9039371c9d4SSatish Balay } 9048031f4adStmunson 905014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]); 906ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]); 9071d6018f0SLisandro Dalcin 908f560b561SHong Zhang PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP)); 909014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP); 910014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP); 9113369ce9aSBarry Smith 9129804daf3SBarry Smith #include <petscdrawtypes.h> 913014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *); 9142f2e5d10SKris Buschelman 915014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec)); 916014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec)); 91703919abeSBarry Smith 918ba36735cSStefano Zampini /*S 919a4d1885cSBarry Smith KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods. 920f8a50e2bSBarry Smith 921a4d1885cSBarry Smith Level: intermediate 922f8a50e2bSBarry Smith 9231cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType` 924ba36735cSStefano Zampini S*/ 925ba36735cSStefano Zampini typedef struct _p_KSPGuess *KSPGuess; 92616a05f60SBarry Smith 927ba36735cSStefano Zampini /*J 92887497f52SBarry Smith KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods. 929ba36735cSStefano Zampini 930a4d1885cSBarry Smith Values: 931a4d1885cSBarry Smith + `KSPGUESSFISCHER` - methodology developed by Paul Fischer 932a4d1885cSBarry Smith - `KSPGUESSPOD` - methodology based on proper orthogonal decomposition 933a4d1885cSBarry Smith 93416a05f60SBarry Smith Level: intermediate 93516a05f60SBarry Smith 9361cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPGuess` 937ba36735cSStefano Zampini J*/ 938ba36735cSStefano Zampini typedef const char *KSPGuessType; 939ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer" 940ba36735cSStefano Zampini #define KSPGUESSPOD "pod" 941a4d1885cSBarry Smith 9421d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess)); 943ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess); 944ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *); 945ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer); 946ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *); 947ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *); 948ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType); 949ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *); 9508410009bSDavid Wells PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal); 951ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 952ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec); 953ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec); 954ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 955ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt); 956014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt); 957ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool); 958ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *); 959f8a50e2bSBarry Smith 960470b340bSDmitry Karpeev /*E 961470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 962470b340bSDmitry Karpeev 963470b340bSDmitry Karpeev Level: intermediate 964470b340bSDmitry Karpeev 965a4d1885cSBarry Smith .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`, 966a4d1885cSBarry Smith `MatCreateSchurComplementPmat()` 967470b340bSDmitry Karpeev E*/ 9689371c9d4SSatish Balay typedef enum { 9699371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_DIAG, 9709371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_LUMP, 9719371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, 9729371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_FULL 9739371c9d4SSatish Balay } MatSchurComplementAinvType; 974470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 975470b340bSDmitry Karpeev 9769371c9d4SSatish Balay typedef enum { 9779371c9d4SSatish Balay MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0, 978864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1, 979864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2, 9809371c9d4SSatish Balay MAT_LMVM_SYMBROYDEN_SCALE_USER = 3 9819371c9d4SSatish Balay } MatLMVMSymBroydenScaleType; 982864588a7SAlp Dener PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[]; 98392f76d53SAlp Dener 984014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *); 985014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *); 986d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP); 987bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 988aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 989bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *); 990470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType); 991470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *); 9925bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *); 9935a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *); 994470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *); 995470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *); 9963f22127dSBarry Smith 99778e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *); 99878e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *); 99978e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *); 1000864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1001864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1002864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1003864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1004864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1005cd929ea3SAlp Dener 1006cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 1007b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *); 1008cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 1009cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 1010e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat); 10110ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat); 1012cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 1013cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 1014cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 1015cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 1016cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 10172d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec); 1018cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 1019cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *); 1020cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *); 1021cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *); 102292f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt); 1023cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *); 1024cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *); 1025864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar); 1026864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType); 1027cd929ea3SAlp Dener 1028014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM); 1029014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool); 1030014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *); 1031014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *); 1032014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *); 1033fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, PetscErrorCode (*func)(KSP, Vec, void *), void *); 103423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *); 1035fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, PetscErrorCode (*)(KSP, Vec, void *), void *); 103623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *); 103723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, PetscErrorCode (**)(KSP, Mat, Mat, void *), void *); 1038014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, PetscErrorCode (*)(KSP, Vec, void *), void *); 1039014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, PetscErrorCode (**)(KSP, Vec, void *), void *); 1040fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, PetscErrorCode (*)(KSP, Vec, void *), void *); 1041fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, PetscErrorCode (**)(KSP, Vec, void *), void *); 10426c699258SBarry Smith 104302b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec); 10449371c9d4SSatish 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); 1045557cf195SMatthew G. Knepley 10462b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *); 10472b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal); 10482eac72dbSBarry Smith #endif 1049