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 6*ac09b921SBarry Smith 72c8e378dSBarry Smith #include <petscpc.h> 82eac72dbSBarry Smith 9*ac09b921SBarry Smith /* SUBMANSEC = KSP */ 10*ac09b921SBarry Smith 11607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void); 121dbb0a54SBarry Smith 1328ce4d24SBarry Smith /*S 148f6c3df8SBarry Smith KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the 158f6c3df8SBarry Smith linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators). 1628ce4d24SBarry Smith 1728ce4d24SBarry Smith Level: beginner 1828ce4d24SBarry Smith 1995452b02SPatrick Sanan Notes: 20db3b2ab5SMatthew Knepley When a direct solver is used, but no Krylov solver is used, the KSP object is still used but with a 21db3b2ab5SMatthew Knepley KSPType of KSPPREONLY, meaning that only application of the preconditioner is used as the linear solver. 228f6c3df8SBarry Smith 23db781477SPatrick Sanan .seealso: `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 32db781477SPatrick Sanan .seealso: `KSPSetType()`, `KSP`, `KSPRegister()`, `KSPCreate()`, `KSPSetFromOptions()` 3376bdecfbSBarry Smith J*/ 3419fd82e9SBarry Smith typedef const char* KSPType; 3582bf6240SBarry Smith #define KSPRICHARDSON "richardson" 366c9de887SHong Zhang #define KSPCHEBYSHEV "chebyshev" 3782bf6240SBarry Smith #define KSPCG "cg" 382c8fc646SJed Brown #define KSPGROPPCG "groppcg" 392c8fc646SJed Brown #define KSPPIPECG "pipecg" 40901ccb91SSiegfried Cools #define KSPPIPECGRR "pipecgrr" 41265663fdSSiegfried Cools #define KSPPIPELCG "pipelcg" 42b21a8899STyler Chen #define KSPPIPEPRCG "pipeprcg" 43325e8391SManasi #define KSPPIPECG2 "pipecg2" 44df2a969aSvictorle #define KSPCGNE "cgne" 4505de396fSBarry Smith #define KSPNASH "nash" 4605de396fSBarry Smith #define KSPSTCG "stcg" 4705de396fSBarry Smith #define KSPGLTR "gltr" 4805de396fSBarry Smith #define KSPCGNASH PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGNASH macro is deprecated use KSPNASH (since version 3.11)\"") "nash" 4905de396fSBarry Smith #define KSPCGSTCG PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGSTCG macro is deprecated use KSPSTCG (since version 3.11)\"") "stcg" 5005de396fSBarry Smith #define KSPCGGLTR PETSC_DEPRECATED_MACRO("GCC warning \"KSPCGGLTR macro is deprecated use KSPSGLTR (since version 3.11)\"") "gltr" 51640f4f53SPatrick Sanan #define KSPFCG "fcg" 52390d8e47SPatrick Sanan #define KSPPIPEFCG "pipefcg" 5382bf6240SBarry Smith #define KSPGMRES "gmres" 54483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres" 559dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 569dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 574f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 5861b468f9SJed Brown #define KSPPGMRES "pgmres" 5982bf6240SBarry Smith #define KSPTCQMR "tcqmr" 6082bf6240SBarry Smith #define KSPBCGS "bcgs" 61e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 62345ecf0bSXiangmin Jiao #define KSPQMRCGS "qmrcgs" 6318ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 64c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 65f0037002Svictorle #define KSPBCGSL "bcgsl" 66f154af2dSSiegfried Cools #define KSPPIPEBCGS "pipebcgs" 6782bf6240SBarry Smith #define KSPCGS "cgs" 6882bf6240SBarry Smith #define KSPTFQMR "tfqmr" 6982bf6240SBarry Smith #define KSPCR "cr" 702c8fc646SJed Brown #define KSPPIPECR "pipecr" 7182bf6240SBarry Smith #define KSPLSQR "lsqr" 7282bf6240SBarry Smith #define KSPPREONLY "preonly" 7382bf6240SBarry Smith #define KSPQCG "qcg" 74c9cf9b11SBarry Smith #define KSPBICG "bicg" 75b4ac9ba4SBarry Smith #define KSPMINRES "minres" 7601c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 7721356919SSatish Balay #define KSPLCD "lcd" 781d6018f0SLisandro Dalcin #define KSPPYTHON "python" 7958ad63f4SBarry Smith #define KSPGCR "gcr" 80fad47a0aSPatrick Sanan #define KSPPIPEGCR "pipegcr" 81e4d80e07Szianekhodja #define KSPTSIRM "tsirm" 82e4d80e07Szianekhodja #define KSPCGLS "cgls" 83329cd39dSStefano Zampini #define KSPFETIDP "fetidp" 8438cfc46eSPierre Jolivet #define KSPHPDDM "hpddm" 852eac72dbSBarry Smith 868ba1e511SMatthew Knepley /* Logging support */ 87014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 88ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID; 895358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 908ba1e511SMatthew Knepley 91014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP*); 9219fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 93c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType*); 94014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 95014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 96014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 97014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 9875f8e9bdSHong Zhang PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP,PetscBool); 995ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP,Mat,Mat); 1003e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP,PetscInt); 1019fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPSetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp,PetscInt n) {return KSPSetMatSolveBatchSize(ksp,n);} 1023e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP,PetscInt*); 1039fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp,PetscInt *n) {return KSPGetMatSolveBatchSize(ksp,n);} 104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 105ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP); 106014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 10723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 1087d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP,PetscBool*); 10925c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool); 110c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec); 1112eac72dbSBarry Smith 112140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 113ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList; 114798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorList; 115798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorCreateList; 116798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList; 117bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 118798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **)); 11930de9b25SBarry Smith 120014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 122014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 123c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 124014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool); 125014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*); 126014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool); 127014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*); 128014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool); 1297d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool); 130c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*); 131014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool); 132c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*); 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*); 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*); 135014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 136014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 137734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*); 1382a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 1399fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") static inline PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);} 1402eac72dbSBarry Smith 141d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 142d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 143d4211eb9SBarry Smith 144014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 145014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 146aabeff55SBarry Smith 147014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 148014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**)); 149014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 1503ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void*); 15194a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,const PetscReal*[],PetscInt*); 152014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool); 15394a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP,const PetscReal*[],PetscInt*); 15494a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP,PetscReal[],PetscInt,PetscBool); 1554b0e389bSBarry Smith 156fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 157fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*); 158fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 159fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 160fa0ddf94SBarry Smith 161014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 162cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP); 163014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 164014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 165014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 166014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 167285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]); 168b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 169b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 170b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 171b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 172014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*); 1734a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC,KSP*); 174f3b08a26SMatthew G. Knepley /* 175f3b08a26SMatthew G. Knepley PCMGCoarseList contains the list of coarse space constructor currently registered 176f3b08a26SMatthew G. Knepley These are added with PCMGRegisterCoarseSpaceConstructor() 177f3b08a26SMatthew G. Knepley */ 178f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList; 179f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **)); 180f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, const Vec[], Vec **)); 181f3b08a26SMatthew G. Knepley 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*); 1842eac72dbSBarry Smith 185014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 186014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 18858450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 189b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool); 19058450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*); 191014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 192d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*); 193d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]); 1947d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]); 1954b0e389bSBarry Smith 196640f4f53SPatrick Sanan /*E 197640f4f53SPatrick Sanan 19806137d0aSPatrick Sanan KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 199640f4f53SPatrick Sanan 20067b8a455SSatish Balay Values: 20167b8a455SSatish Balay $ KSP_FCD_TRUNC_TYPE_STANDARD - uses all (up to mmax) stored directions 20267b8a455SSatish Balay $ KSP_FCD_TRUNC_TYPE_NOTAY - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 203640f4f53SPatrick Sanan 2042b26979fSBarry Smith Level: intermediate 205db781477SPatrick Sanan .seealso `:` `KSPFCG,KSPPIPEFCG,KSPPIPEGCR,KSPFCGSetTruncationType(),KSPFCGGetTruncationType()` 206640f4f53SPatrick Sanan 207640f4f53SPatrick Sanan E*/ 20893f1e87bSPatrick Sanan typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType; 20906137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 210640f4f53SPatrick Sanan 211640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 212640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 213640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 214640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 21506137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType); 21606137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*); 217640f4f53SPatrick Sanan 218390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt); 219390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*); 220390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt); 221390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*); 222390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType); 223390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*); 224390d8e47SPatrick Sanan 225fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt); 226fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*); 227fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt); 228fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*); 229fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType); 230fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*); 231fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool); 232fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*); 23383f0b094SBarry Smith 234014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 235014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 236014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 237e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP,PetscReal); 2389f236934SBarry Smith 239014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 240014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 241014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 242014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 243014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 2441d73ed98SBarry Smith 245014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 246014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2471d73ed98SBarry Smith 248483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar); 249483d6965SPatrick Sanan 250014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 251014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 25358ad63f4SBarry Smith 25404d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*); 25504d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC); 25604d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*); 257568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat); 258e82af88dSprj- 2596cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP,Mat); 2606cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP,Mat*); 2616cac28cbSPierre Jolivet #if PetscDefined(HAVE_HPDDM) 2626cac28cbSPierre Jolivet PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMSetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) { return KSPHPDDMSetDeflationMat(ksp, U); } 2636cac28cbSPierre Jolivet PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMGetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) { return KSPHPDDMGetDeflationMat(ksp, U); } 2646cac28cbSPierre Jolivet #endif 2659fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) { return KSPMatSolve(ksp, B, X); } 266d55ab951SPierre Jolivet /*E 267d55ab951SPierre Jolivet KSPHPDDMType - Type of Krylov method used by KSPHPDDM 268d55ab951SPierre Jolivet 269d55ab951SPierre Jolivet Level: intermediate 270d55ab951SPierre Jolivet 271d55ab951SPierre Jolivet Values: 27267b8a455SSatish Balay $ KSP_HPDDM_TYPE_GMRES (default) 27367b8a455SSatish Balay $ KSP_HPDDM_TYPE_BGMRES 27467b8a455SSatish Balay $ KSP_HPDDM_TYPE_CG 27567b8a455SSatish Balay $ KSP_HPDDM_TYPE_BCG 27667b8a455SSatish Balay $ KSP_HPDDM_TYPE_GCRODR 27767b8a455SSatish Balay $ KSP_HPDDM_TYPE_BGCRODR 27867b8a455SSatish Balay $ KSP_HPDDM_TYPE_BFBCG 27967b8a455SSatish Balay $ KSP_HPDDM_TYPE_PREONLY 280d55ab951SPierre Jolivet 281db781477SPatrick Sanan .seealso: `KSPHPDDM`, `KSPHPDDMSetType()` 282d55ab951SPierre Jolivet E*/ 283d55ab951SPierre Jolivet typedef enum { 284d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GMRES = 0, 285d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGMRES = 1, 286d55ab951SPierre Jolivet KSP_HPDDM_TYPE_CG = 2, 287d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BCG = 3, 288d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GCRODR = 4, 289d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGCRODR = 5, 290d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BFBCG = 6, 291d55ab951SPierre Jolivet KSP_HPDDM_TYPE_PREONLY = 7 292d55ab951SPierre Jolivet } KSPHPDDMType; 293d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[]; 2942dd49c90SPierre Jolivet /*E 2952dd49c90SPierre Jolivet KSPHPDDMPrecision - Precision of Krylov bases used by KSPHPDDM 2962dd49c90SPierre Jolivet 2972dd49c90SPierre Jolivet Level: intermediate 2982dd49c90SPierre Jolivet 2992dd49c90SPierre Jolivet Values: 30067b8a455SSatish Balay $ KSP_HPDDM_PRECISION_HALF (currently unsupported) 30167b8a455SSatish Balay $ KSP_HPDDM_PRECISION_SINGLE (default when PETSc is configured --with-precision=single) 30267b8a455SSatish Balay $ KSP_HPDDM_PRECISION_DOUBLE (default when PETSc is configured --with-precision=double) 30367b8a455SSatish Balay $ KSP_HPDDM_PRECISION_QUADRUPLE (currently unsupported) 3042dd49c90SPierre Jolivet 305db781477SPatrick Sanan .seealso: `KSPHPDDM` 3062dd49c90SPierre Jolivet E*/ 3072dd49c90SPierre Jolivet typedef enum { 3082dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_HALF = 0, 3092dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_SINGLE = 1, 3102dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_DOUBLE = 2, 3112dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_QUADRUPLE = 3 3122dd49c90SPierre Jolivet } KSPHPDDMPrecision; 313d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP,KSPHPDDMType); 314d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP,KSPHPDDMType*); 315e82af88dSprj- 316b49fd9e1SBarry Smith /*E 317b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 318b49fd9e1SBarry Smith 319b49fd9e1SBarry Smith Level: advanced 320b49fd9e1SBarry Smith 321db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 322db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()` 323b49fd9e1SBarry Smith 324b49fd9e1SBarry Smith E*/ 32578d1dd23SKris Buschelman typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 3266a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 3271f7e983dSSatish Balay /*MC 3281957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 3298c5b8ba0SBarry Smith 3308c5b8ba0SBarry Smith Level: advanced 3318c5b8ba0SBarry Smith 3328c5b8ba0SBarry Smith Note: Possible unstable, but the fastest to compute 3338c5b8ba0SBarry Smith 334db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 335db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 336db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 3378c5b8ba0SBarry Smith M*/ 3388c5b8ba0SBarry Smith 3391f7e983dSSatish Balay /*MC 3408c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 3418c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 3428c5b8ba0SBarry Smith poor orthogonality. 3438c5b8ba0SBarry Smith 3448c5b8ba0SBarry Smith Level: advanced 3458c5b8ba0SBarry Smith 3468c5b8ba0SBarry Smith Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 3478c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 3488c5b8ba0SBarry Smith 349db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 350db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 351db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 3528c5b8ba0SBarry Smith M*/ 3538c5b8ba0SBarry Smith 3541f7e983dSSatish Balay /*MC 3558c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 3568c5b8ba0SBarry Smith 3578c5b8ba0SBarry Smith Level: advanced 3588c5b8ba0SBarry Smith 3598c5b8ba0SBarry Smith Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 3608c5b8ba0SBarry Smith but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 3618c5b8ba0SBarry Smith 3628c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 3638c5b8ba0SBarry Smith 364db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 365db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 366db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 3678c5b8ba0SBarry Smith M*/ 3688c5b8ba0SBarry Smith 369014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 370014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 37108480c60SBarry Smith 372014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 373014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 374014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 375c38d4ed2SBarry Smith 376014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 377014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 378014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 379121fd945SKris Buschelman 380014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 381014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool); 382014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 383e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 384d9492815SBarry Smith 385014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 38687d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP); 387014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 3882eac72dbSBarry Smith 389798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],void*); 390798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm,const char[],const char[],const char[],PetscInt,const char *[],int,int,int,int,PetscDrawLG *); 391798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 392798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 393798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 394798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 395798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 396798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 397798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 398798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 399798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 400798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 401798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 402798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 403798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 404798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 405798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 406798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 407798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 408798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 409798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 410fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 411798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 4129fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorResidual(ksp,n,rnorm,vf);} 4139fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorTrueResidual(ksp,n,rnorm,vf);} 4149fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorTrueResidualMax(ksp,n,rnorm,vf);} 415798534f6SMatthew G. Knepley 416798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*); 417390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 418390160ecSHong Zhang PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 419e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 420e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 421e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 42284cb2905SBarry Smith 423014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 424014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 425c01c455dSBarry Smith 42623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 42723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 428014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*); 429014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 430014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 431014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 4321eb62cbbSBarry Smith 433014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool); 434014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*); 435014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool); 436014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*); 4371f7f0c4fSBarry Smith 438014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 43955849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 440fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]); 44119a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP,PetscViewer); 442c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP,PetscErrorCode (*)(KSP,void*),void *vctx,PetscErrorCode (*)(void**)); 4431b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP); 444c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP); 44594a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP,PetscViewer); 4461b2b9847SBarry Smith 4479fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") static inline PetscErrorCode KSPReasonView(KSP ksp,PetscViewer v) {return KSPConvergedReasonView(ksp,v);} 4489fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) {return KSPConvergedReasonViewFromOptions(ksp);} 44955849f57SBarry Smith 45055849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 4511eb62cbbSBarry Smith 452dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool); 4530e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool); 454014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 455884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*); 456798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 457798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 458798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 459db9b2ab1SHong Zhang 460014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 461014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 46268ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*); 46383ab6a24SBarry Smith 46428ce4d24SBarry Smith /*E 4658a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 4668a4b9c5cSBarry Smith test routines. 4678a4b9c5cSBarry Smith 4688a4b9c5cSBarry Smith Level: advanced 4698a4b9c5cSBarry Smith 470a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 4711718446fSBarry Smith depending on left or right preconditioning, see KSPSetPCSide() 472a3f661c8SBarry Smith 47395452b02SPatrick Sanan Notes: 47495452b02SPatrick Sanan this must match petsc/finclude/petscksp.h 4758a4b9c5cSBarry Smith 476db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`, 477db781477SPatrick Sanan `KSPSetConvergenceTest()`, `KSPSetPCSide()` 4788a4b9c5cSBarry Smith E*/ 4799e568555SJed Brown typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 4809e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 481014dd563SJed Brown PETSC_EXTERN const char *const*const KSPNormTypes; 482a21b2a99SBarry Smith 4831f7e983dSSatish Balay /*MC 4849793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 4858c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 4868c5b8ba0SBarry Smith be based on a norm of a residual etc. 4878c5b8ba0SBarry Smith 4888c5b8ba0SBarry Smith Level: advanced 4898c5b8ba0SBarry Smith 4901957e957SBarry Smith Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 4918c5b8ba0SBarry Smith 492db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 4938c5b8ba0SBarry Smith M*/ 4948c5b8ba0SBarry Smith 4951f7e983dSSatish Balay /*MC 4961957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 4978c5b8ba0SBarry Smith convergence test routine. 4988c5b8ba0SBarry Smith 4998c5b8ba0SBarry Smith Level: advanced 5008c5b8ba0SBarry Smith 501db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 5028c5b8ba0SBarry Smith M*/ 5038c5b8ba0SBarry Smith 5041f7e983dSSatish Balay /*MC 505ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 5068c5b8ba0SBarry Smith convergence test routine. 5078c5b8ba0SBarry Smith 5088c5b8ba0SBarry Smith Level: advanced 5098c5b8ba0SBarry Smith 510db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 5118c5b8ba0SBarry Smith M*/ 5128c5b8ba0SBarry Smith 5131f7e983dSSatish Balay /*MC 514ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 515390d8e47SPatrick Sanan convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR 5168c5b8ba0SBarry Smith 5178c5b8ba0SBarry Smith Level: advanced 5188c5b8ba0SBarry Smith 519db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()` 5208c5b8ba0SBarry Smith M*/ 5218c5b8ba0SBarry Smith 522014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 523014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 524014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 525014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 526014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 5278a4b9c5cSBarry Smith 528f0fc11ceSJed Brown #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)") 5298a4b9c5cSBarry Smith /*E 5301957e957SBarry Smith KSPConvergedReason - reason a Krylov method was said to have converged or diverged 53128ce4d24SBarry Smith 53228ce4d24SBarry Smith Level: beginner 53328ce4d24SBarry Smith 53495452b02SPatrick Sanan Notes: 53595452b02SPatrick Sanan See KSPGetConvergedReason() for explanation of each value 53628ce4d24SBarry Smith 53795452b02SPatrick Sanan Developer Notes: 53895452b02SPatrick Sanan this must match petsc/finclude/petscksp.h 5394d0a8057SBarry Smith 5404d0a8057SBarry Smith The string versions of these are KSPConvergedReasons; if you change 5414d0a8057SBarry Smith any of the values here also change them that array of names. 54286c02ca4SBarry Smith 543db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()` 54428ce4d24SBarry Smith E*/ 545d15094e1SBarry Smith typedef enum {/* converged */ 5469ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 5479ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 548d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 549d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 550b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 5518031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 5528031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 553329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 554af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 555d15094e1SBarry Smith /* diverged */ 556b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 557d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 558d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 559d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 560b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 561b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 562b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 5634d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 5646aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 565c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED = -11, 566aa4d2078SSatish Balay KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11, 567d15094e1SBarry Smith 568d15094e1SBarry Smith KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 569014dd563SJed Brown PETSC_EXTERN const char *const*KSPConvergedReasons; 570d15094e1SBarry Smith 571c838673bSBarry Smith /*MC 572f9fed41fSBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called 573c838673bSBarry Smith 574c838673bSBarry Smith Level: beginner 575c838673bSBarry Smith 576c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 577c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 578c838673bSBarry Smith 2-norm of the residual for right preconditioning 579c838673bSBarry Smith 580f9fed41fSBarry Smith See also KSP_CONVERGED_ATOL which may apply before this tolerance. 581f9fed41fSBarry Smith 582db781477SPatrick Sanan .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 583c838673bSBarry Smith 584c838673bSBarry Smith M*/ 585c838673bSBarry Smith 586c838673bSBarry Smith /*MC 587c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 588c838673bSBarry Smith 589c838673bSBarry Smith Level: beginner 590c838673bSBarry Smith 591c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 592c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 593c838673bSBarry Smith 2-norm of the residual for right preconditioning 594c838673bSBarry Smith 595f9fed41fSBarry Smith See also KSP_CONVERGED_RTOL which may apply before this tolerance. 596c838673bSBarry Smith 597db781477SPatrick Sanan .seealso: `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 598c838673bSBarry Smith 599c838673bSBarry Smith M*/ 600c838673bSBarry Smith 601c838673bSBarry Smith /*MC 602c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 603c838673bSBarry Smith 604c838673bSBarry Smith Level: beginner 605c838673bSBarry Smith 606c838673bSBarry Smith See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 607c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 608c838673bSBarry Smith 2-norm of the residual for right preconditioning 609c838673bSBarry Smith 610c838673bSBarry Smith Level: beginner 611c838673bSBarry Smith 612db781477SPatrick Sanan .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 613c838673bSBarry Smith 614c838673bSBarry Smith M*/ 615c838673bSBarry Smith 616c838673bSBarry Smith /*MC 617c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 618c838673bSBarry Smith reached 619c838673bSBarry Smith 620c838673bSBarry Smith Level: beginner 621c838673bSBarry Smith 622db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 623c838673bSBarry Smith 624c838673bSBarry Smith M*/ 625c838673bSBarry Smith 626c838673bSBarry Smith /*MC 6278c9406f6SLisandro Dalcin KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 6280059c7bdSJed Brown the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 6294d0a8057SBarry Smith test routine is set in KSP. 630c838673bSBarry Smith 631c838673bSBarry Smith Level: beginner 632c838673bSBarry Smith 633db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 634c838673bSBarry Smith 635c838673bSBarry Smith M*/ 636c838673bSBarry Smith 637c838673bSBarry Smith /*MC 638c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 6391de96524SPierre Jolivet method could not continue to enlarge the Krylov space. Could be due to a singular matrix or 6401de96524SPierre Jolivet preconditioner. In KSPHPDDM, this is also returned when some search directions within a block 6411de96524SPierre Jolivet are colinear. 642c838673bSBarry Smith 643c838673bSBarry Smith Level: beginner 644c838673bSBarry Smith 645db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 646c838673bSBarry Smith 647c838673bSBarry Smith M*/ 648c838673bSBarry Smith 649c838673bSBarry Smith /*MC 650c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 651c838673bSBarry Smith method could not continue to enlarge the Krylov space. 652c838673bSBarry Smith 653c838673bSBarry Smith Level: beginner 654c838673bSBarry Smith 655db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 656c838673bSBarry Smith 657c838673bSBarry Smith M*/ 658c838673bSBarry Smith 659c838673bSBarry Smith /*MC 660c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 661c838673bSBarry Smith symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 662c838673bSBarry Smith 663c838673bSBarry Smith Level: beginner 664c838673bSBarry Smith 665db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 666c838673bSBarry Smith 667c838673bSBarry Smith M*/ 668c838673bSBarry Smith 669c838673bSBarry Smith /*MC 670c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 671c838673bSBarry Smith positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 672c838673bSBarry Smith be positive definite 673c838673bSBarry Smith 674c838673bSBarry Smith Level: beginner 675c838673bSBarry Smith 67695452b02SPatrick Sanan Notes: 67795452b02SPatrick Sanan This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 678c838673bSBarry Smith the PCICC preconditioner to generate a positive definite preconditioner 679c838673bSBarry Smith 680db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 681c838673bSBarry Smith 682c838673bSBarry Smith M*/ 683c838673bSBarry Smith 684c838673bSBarry Smith /*MC 685c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a 6869fc87aa7SBarry Smith zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 6879fc87aa7SBarry Smith such as PCFIELDSPLIT. 6889fc87aa7SBarry Smith 6899fc87aa7SBarry Smith Level: beginner 6909fc87aa7SBarry Smith 69195452b02SPatrick Sanan Notes: 69295452b02SPatrick Sanan Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details. 6939fc87aa7SBarry Smith 694db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 6959fc87aa7SBarry Smith 6969fc87aa7SBarry Smith M*/ 6979fc87aa7SBarry Smith 6989fc87aa7SBarry Smith /*MC 699c838673bSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 700c838673bSBarry Smith while the KSPSolve() is still running. 701c838673bSBarry Smith 702c838673bSBarry Smith Level: beginner 703c838673bSBarry Smith 704db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 705c838673bSBarry Smith 706c838673bSBarry Smith M*/ 707c838673bSBarry Smith 708014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*)); 709df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*)); 710df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*)); 7113ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void*); 7128de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 7133eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 7148de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*); 7158de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**); 7168de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 7178de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 71854b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP,PetscBool); 7190059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 720014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*); 721c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP,const char**); 72294a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP,PetscReal*,PetscReal*,PetscReal*,PetscReal*); 723abef13c0SSatish Balay 7249fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") static inline void KSPDefaultConverged(void) { /* never called */ } 7258ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 7269fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") static inline void KSPDefaultConvergedDestroy(void) { /* never called */ } 7278ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 7289fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") static inline void KSPDefaultConvergedCreate(void) { /* never called */ } 7298ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 7309fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 7318ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 7329fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 7338ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 7349fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") static inline void KSPSkipConverged(void) { /* never called */ } 7358ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 7368ea1b3e6SJed Brown 7370bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*); 7389fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") static inline PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); } 739d4fbbf0eSBarry Smith 74028ce4d24SBarry Smith /*E 74128ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 74228ce4d24SBarry Smith 74328ce4d24SBarry Smith Level: beginner 74428ce4d24SBarry Smith 745db781477SPatrick Sanan .seealso: `KSPCGSetType()` 74628ce4d24SBarry Smith E*/ 7479dcbbd2bSBarry Smith typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 7486a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 74928ce4d24SBarry Smith 750014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 751014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool); 7528031f4adStmunson 753dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal); 754dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*); 755dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*); 756fcae7a14Stmunson 75705de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*); 75805de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*); 7599fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);} 7609fbee547SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);} 7618031f4adStmunson 762014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 7631d6018f0SLisandro Dalcin 764f560b561SHong Zhang PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC,PetscErrorCode (*)(PC,KSP)); 765014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 766014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 7673369ce9aSBarry Smith 7689804daf3SBarry Smith #include <petscdrawtypes.h> 769014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 7702f2e5d10SKris Buschelman 771014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 772014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 77303919abeSBarry Smith 774ba36735cSStefano Zampini /*S 775ba36735cSStefano Zampini KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods. 776f8a50e2bSBarry Smith 777ba36735cSStefano Zampini Level: beginner 778f8a50e2bSBarry Smith 779db781477SPatrick Sanan .seealso: `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType` 780ba36735cSStefano Zampini S*/ 781ba36735cSStefano Zampini typedef struct _p_KSPGuess* KSPGuess; 782ba36735cSStefano Zampini /*J 783ba36735cSStefano Zampini KSPGuessType - String with the name of a PETSc initial guess for Krylov methods. 784ba36735cSStefano Zampini 785ba36735cSStefano Zampini Level: beginner 786ba36735cSStefano Zampini 787db781477SPatrick Sanan .seealso: `KSPGuess` 788ba36735cSStefano Zampini J*/ 789ba36735cSStefano Zampini typedef const char* KSPGuessType; 790ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer" 791ba36735cSStefano Zampini #define KSPGUESSPOD "pod" 7921d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess)); 793ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess); 794ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*); 795ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer); 796ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*); 797ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*); 798ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType); 799ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*); 8008410009bSDavid Wells PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess,PetscReal); 801ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 802ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec); 803ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec); 804ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 805ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt); 806014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 807ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool); 808ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*); 809f8a50e2bSBarry Smith 810470b340bSDmitry Karpeev /*E 811470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 812470b340bSDmitry Karpeev 813470b340bSDmitry Karpeev Level: intermediate 814470b340bSDmitry Karpeev 815db781477SPatrick Sanan .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`, `MatCreateSchurComplementPmat()` 816470b340bSDmitry Karpeev E*/ 8170c42975cSTristan Konolige typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType; 818470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 819470b340bSDmitry Karpeev 820864588a7SAlp Dener typedef enum {MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0, 821864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1, 822864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2, 823864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_USER = 3} MatLMVMSymBroydenScaleType; 824864588a7SAlp Dener PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[]; 82592f76d53SAlp Dener 826014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 827014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 828d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 829bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 830aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 831bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 832470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 833470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 8345bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 8355a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 836470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*); 837470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 8383f22127dSBarry Smith 83978e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*); 84078e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*); 84178e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*); 842864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 843864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 844864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 845864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 846864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 847cd929ea3SAlp Dener 848cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 849b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*); 850cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 851cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 852e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat); 8530ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat); 854cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 855cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 856cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 857cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 858cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 8592d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec); 860cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 861cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*); 862cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*); 863cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*); 86492f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt); 865cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*); 866cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*); 867864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar); 868864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType); 869cd929ea3SAlp Dener 870014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 871014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool); 872014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 873014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 874014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 875fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*); 87623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 877fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 87823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 87923ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 880014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 881014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 882fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 883fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 8846c699258SBarry Smith 88502b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 8868c6c5593SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, 887e29c0834SMatthew G. Knepley void (**)(PetscInt, PetscInt, PetscInt, 888e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 889e29c0834SMatthew G. Knepley const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 890191494d9SMatthew G. Knepley PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec); 891557cf195SMatthew G. Knepley 892557cf195SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, PetscInt, Vec[], Vec[], Mat *, void *); 893557cf195SMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, PetscInt, Vec[], Vec[], PetscReal); 8942eac72dbSBarry Smith #endif 895