1f26ada1bSBarry Smith /* 2f26ada1bSBarry Smith Defines the interface functions for the Krylov subspace accelerators. 3f26ada1bSBarry Smith */ 426bd1501SBarry Smith #ifndef PETSCKSP_H 526bd1501SBarry Smith #define PETSCKSP_H 6ac09b921SBarry Smith 72c8e378dSBarry Smith #include <petscpc.h> 82eac72dbSBarry Smith 9ac09b921SBarry Smith /* SUBMANSEC = KSP */ 10ac09b921SBarry Smith 11607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void); 121dbb0a54SBarry Smith 1328ce4d24SBarry Smith /*S 148f6c3df8SBarry Smith KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the 158f6c3df8SBarry Smith linear solves in PETSc (even those such as direct solvers that do no use Krylov accelerators). 1628ce4d24SBarry Smith 1728ce4d24SBarry Smith Level: beginner 1828ce4d24SBarry Smith 1987497f52SBarry Smith Note: 20db3b2ab5SMatthew Knepley When a direct solver is used, but no Krylov solver is used, the KSP object is still used but with a 2187497f52SBarry Smith `KSPType` of `KSPPREONLY`, meaning that only application of the preconditioner is used as the linear solver. 228f6c3df8SBarry Smith 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" 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); 1013e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt); 102d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPSetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n) 103d71ae5a4SJacob Faibussowitsch { 1049371c9d4SSatish Balay return KSPSetMatSolveBatchSize(ksp, n); 1059371c9d4SSatish Balay } 1063e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *); 107d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n) 108d71ae5a4SJacob Faibussowitsch { 1099371c9d4SSatish Balay return KSPGetMatSolveBatchSize(ksp, n); 1109371c9d4SSatish Balay } 111014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 112ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP); 113014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *); 11423ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool); 1157d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *); 11625c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool); 117c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec); 1182eac72dbSBarry Smith 119140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 120ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList; 121798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorList; 122798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorCreateList; 123798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList; 124bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[], PetscErrorCode (*)(KSP)); 125798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **)); 12630de9b25SBarry Smith 127014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide); 128014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *); 129014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt); 130c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *); 131014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool); 132014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *); 133014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool); 134014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *); 135014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool); 1367d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool); 137c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *); 138014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool); 139c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *); 140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *); 141014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *); 142014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *); 143014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *); 144734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *); 1452a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **); 146d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y) 147d71ae5a4SJacob Faibussowitsch { 1489371c9d4SSatish Balay return KSPCreateVecs(ksp, n, x, m, y); 1499371c9d4SSatish Balay } 1502eac72dbSBarry Smith 151d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *); 152d4211eb9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *); 153d4211eb9SBarry Smith 154014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC); 155014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *); 156aabeff55SBarry Smith 157014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal); 158014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void **)); 159014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 1603ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *); 16194a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *); 162014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscInt, PetscBool); 16394a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *); 16494a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscInt, PetscBool); 1654b0e389bSBarry Smith 166fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *); 167fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *); 168fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 169fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt); 170fa0ddf94SBarry Smith 171014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *); 172cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP); 173014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 174014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 175014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 176014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]); 177285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]); 178b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *); 179b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *); 180b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *); 181b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *); 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *); 1834a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *); 184f3b08a26SMatthew G. Knepley /* 185f3b08a26SMatthew G. Knepley PCMGCoarseList contains the list of coarse space constructor currently registered 186f3b08a26SMatthew G. Knepley These are added with PCMGRegisterCoarseSpaceConstructor() 187f3b08a26SMatthew G. Knepley */ 188f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList; 1892b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *)); 1902b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *)); 191f3b08a26SMatthew G. Knepley 192014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *); 193014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *); 1942eac72dbSBarry Smith 195014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal); 196014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool); 197014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal); 19858450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal); 199b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool); 20058450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *); 201014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *); 202d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *); 203d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]); 2047d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]); 2054b0e389bSBarry Smith 206640f4f53SPatrick Sanan /*E 207640f4f53SPatrick Sanan 20806137d0aSPatrick Sanan KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 209640f4f53SPatrick Sanan 21067b8a455SSatish Balay Values: 21187497f52SBarry Smith $ `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions 21287497f52SBarry Smith $ `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 213640f4f53SPatrick Sanan 2142b26979fSBarry Smith Level: intermediate 215640f4f53SPatrick Sanan 21687497f52SBarry Smith .seealso `:` `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()` 217640f4f53SPatrick Sanan E*/ 2189371c9d4SSatish Balay typedef enum { 2199371c9d4SSatish Balay KSP_FCD_TRUNC_TYPE_STANDARD, 2209371c9d4SSatish Balay KSP_FCD_TRUNC_TYPE_NOTAY 2219371c9d4SSatish Balay } KSPFCDTruncationType; 22206137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 223640f4f53SPatrick Sanan 224640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt); 225640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *); 226640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt); 227640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *); 22806137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType); 22906137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *); 230640f4f53SPatrick Sanan 231390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt); 232390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *); 233390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt); 234390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *); 235390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType); 236390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *); 237390d8e47SPatrick Sanan 238fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt); 239fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *); 240fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt); 241fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *); 242fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType); 243fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *); 244fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool); 245fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *); 24683f0b094SBarry Smith 247014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 248014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *); 249014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal); 250e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal); 2519f236934SBarry Smith 252014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 253014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt)); 254014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt)); 255014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt); 256014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt); 2571d73ed98SBarry Smith 258014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt); 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 2601d73ed98SBarry Smith 261483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar); 262483d6965SPatrick Sanan 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt); 264014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 26658ad63f4SBarry Smith 26704d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *); 26804d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC); 26904d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *); 270568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat); 271e82af88dSprj- 2726cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat); 2736cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *); 2746cac28cbSPierre Jolivet #if PetscDefined(HAVE_HPDDM) 275d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMSetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) 276d71ae5a4SJacob Faibussowitsch { 2779371c9d4SSatish Balay return KSPHPDDMSetDeflationMat(ksp, U); 2789371c9d4SSatish Balay } 279d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMGetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) 280d71ae5a4SJacob Faibussowitsch { 2819371c9d4SSatish Balay return KSPHPDDMGetDeflationMat(ksp, U); 2829371c9d4SSatish Balay } 2836cac28cbSPierre Jolivet #endif 284d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) 285d71ae5a4SJacob Faibussowitsch { 2869371c9d4SSatish Balay return KSPMatSolve(ksp, B, X); 2879371c9d4SSatish Balay } 288d55ab951SPierre Jolivet /*E 28987497f52SBarry Smith KSPHPDDMType - Type of Krylov method used by `KSPHPDDM` 290d55ab951SPierre Jolivet 291d55ab951SPierre Jolivet Level: intermediate 292d55ab951SPierre Jolivet 293d55ab951SPierre Jolivet Values: 29487497f52SBarry Smith $ `KSP_HPDDM_TYPE_GMRES` (default) 29587497f52SBarry Smith $ `KSP_HPDDM_TYPE_BGMRES` 29687497f52SBarry Smith $ `KSP_HPDDM_TYPE_CG` 29787497f52SBarry Smith $ `KSP_HPDDM_TYPE_BCG` 29887497f52SBarry Smith $ `KSP_HPDDM_TYPE_GCRODR` 29987497f52SBarry Smith $ `KSP_HPDDM_TYPE_BGCRODR` 30087497f52SBarry Smith $ `KSP_HPDDM_TYPE_BFBCG` 30187497f52SBarry Smith $ `KSP_HPDDM_TYPE_PREONLY` 302d55ab951SPierre Jolivet 303db781477SPatrick Sanan .seealso: `KSPHPDDM`, `KSPHPDDMSetType()` 304d55ab951SPierre Jolivet E*/ 305d55ab951SPierre Jolivet typedef enum { 306d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GMRES = 0, 307d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGMRES = 1, 308d55ab951SPierre Jolivet KSP_HPDDM_TYPE_CG = 2, 309d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BCG = 3, 310d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GCRODR = 4, 311d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGCRODR = 5, 312d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BFBCG = 6, 313d55ab951SPierre Jolivet KSP_HPDDM_TYPE_PREONLY = 7 314d55ab951SPierre Jolivet } KSPHPDDMType; 315d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[]; 3162dd49c90SPierre Jolivet /*E 31787497f52SBarry Smith KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM` 3182dd49c90SPierre Jolivet 3192dd49c90SPierre Jolivet Level: intermediate 3202dd49c90SPierre Jolivet 3212dd49c90SPierre Jolivet Values: 3224c92ab2fSPierre Jolivet $ `KSP_HPDDM_PRECISION_HALF` 32387497f52SBarry Smith $ `KSP_HPDDM_PRECISION_SINGLE` (default when PETSc is configured --with-precision=single) 32487497f52SBarry Smith $ `KSP_HPDDM_PRECISION_DOUBLE` (default when PETSc is configured --with-precision=double) 32587497f52SBarry Smith $ `KSP_HPDDM_PRECISION_QUADRUPLE` (default when PETSc is configured --with-precision=__float128) 3262dd49c90SPierre Jolivet 327db781477SPatrick Sanan .seealso: `KSPHPDDM` 3282dd49c90SPierre Jolivet E*/ 3292dd49c90SPierre Jolivet typedef enum { 3302dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_HALF = 0, 3312dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_SINGLE = 1, 3322dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_DOUBLE = 2, 3332dd49c90SPierre Jolivet KSP_HPDDM_PRECISION_QUADRUPLE = 3 3342dd49c90SPierre Jolivet } KSPHPDDMPrecision; 335d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType); 336d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *); 337e82af88dSprj- 338b49fd9e1SBarry Smith /*E 339b49fd9e1SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 340b49fd9e1SBarry Smith 341b49fd9e1SBarry Smith Level: advanced 342b49fd9e1SBarry Smith 343db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 344db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()` 345b49fd9e1SBarry Smith 346b49fd9e1SBarry Smith E*/ 3479371c9d4SSatish Balay typedef enum { 3489371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_NEVER, 3499371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_IFNEEDED, 3509371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_ALWAYS 3519371c9d4SSatish Balay } KSPGMRESCGSRefinementType; 3526a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 3531f7e983dSSatish Balay /*MC 3541957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 3558c5b8ba0SBarry Smith 3568c5b8ba0SBarry Smith Level: advanced 3578c5b8ba0SBarry Smith 35887497f52SBarry Smith Note: 35987497f52SBarry Smith Possibly unstable, but the fastest to compute 3608c5b8ba0SBarry Smith 361db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 362db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 363db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 3648c5b8ba0SBarry Smith M*/ 3658c5b8ba0SBarry Smith 3661f7e983dSSatish Balay /*MC 3678c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 3688c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 3698c5b8ba0SBarry Smith poor orthogonality. 3708c5b8ba0SBarry Smith 3718c5b8ba0SBarry Smith Level: advanced 3728c5b8ba0SBarry Smith 37387497f52SBarry Smith Note: This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to 3748c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 3758c5b8ba0SBarry Smith 376db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 377db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 378db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 3798c5b8ba0SBarry Smith M*/ 3808c5b8ba0SBarry Smith 3811f7e983dSSatish Balay /*MC 3828c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 3838c5b8ba0SBarry Smith 3848c5b8ba0SBarry Smith Level: advanced 3858c5b8ba0SBarry Smith 38687497f52SBarry Smith Notes: 38787497f52SBarry Smith This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice 38887497f52SBarry Smith but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`. 3898c5b8ba0SBarry Smith 3908c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 3918c5b8ba0SBarry Smith 392db781477SPatrick Sanan .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 393db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 394db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 3958c5b8ba0SBarry Smith M*/ 3968c5b8ba0SBarry Smith 397014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType); 398014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *); 39908480c60SBarry Smith 400014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *); 401014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *); 402014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 403c38d4ed2SBarry Smith 404014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal); 405014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *); 406014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *); 407121fd945SKris Buschelman 408014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal); 409014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool); 410014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt); 411e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool); 412d9492815SBarry Smith 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 41487d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP); 415014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 4162eac72dbSBarry Smith 417798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *); 418798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm, const char[], const char[], const char[], PetscInt, const char *[], int, int, int, int, PetscDrawLG *); 419798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 420798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 421798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 422798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 423798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 424798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 425798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 426798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 427798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 428798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 429798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 430798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 431798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 432798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 433798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 434798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 435798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 436798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 437798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 438fe01d993SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 439798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 440d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 441d71ae5a4SJacob Faibussowitsch { 4429371c9d4SSatish Balay return KSPMonitorResidual(ksp, n, rnorm, vf); 4439371c9d4SSatish Balay } 444d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 445d71ae5a4SJacob Faibussowitsch { 4469371c9d4SSatish Balay return KSPMonitorTrueResidual(ksp, n, rnorm, vf); 4479371c9d4SSatish Balay } 448d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 449d71ae5a4SJacob Faibussowitsch { 4509371c9d4SSatish Balay return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf); 4519371c9d4SSatish Balay } 452798534f6SMatthew G. Knepley 453798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *); 454*341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *); 455*341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **); 456*341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *); 457*341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal); 458e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *); 459e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **); 460e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **); 46184cb2905SBarry Smith 462014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec); 463014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec); 464c01c455dSBarry Smith 46523ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat); 46623ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *); 467014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *); 468014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]); 469014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]); 470014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]); 4711eb62cbbSBarry Smith 472014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool); 473014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *); 474014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool); 475014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *); 4761f7f0c4fSBarry Smith 477014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer); 47855849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer); 479fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]); 48019a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer); 481c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **)); 4821b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP); 483c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP); 48494a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer); 4851b2b9847SBarry Smith 486d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v) 487d71ae5a4SJacob Faibussowitsch { 4889371c9d4SSatish Balay return KSPConvergedReasonView(ksp, v); 4899371c9d4SSatish Balay } 490d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) 491d71ae5a4SJacob Faibussowitsch { 4929371c9d4SSatish Balay return KSPConvergedReasonViewFromOptions(ksp); 4939371c9d4SSatish Balay } 49455849f57SBarry Smith 49555849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 4961eb62cbbSBarry Smith 497dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool); 4980e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool); 499014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *); 500884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *); 501798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 502798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 503798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 504db9b2ab1SHong Zhang 505014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *); 506014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *); 50768ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *); 50883ab6a24SBarry Smith 50928ce4d24SBarry Smith /*E 5108a4b9c5cSBarry Smith KSPNormType - Norm that is passed in the Krylov convergence 5118a4b9c5cSBarry Smith test routines. 5128a4b9c5cSBarry Smith 5138a4b9c5cSBarry Smith Level: advanced 5148a4b9c5cSBarry Smith 515a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 51687497f52SBarry Smith depending on left or right preconditioning, see `KSPSetPCSide()` 517a3f661c8SBarry Smith 51887497f52SBarry Smith Developer Note: 51995452b02SPatrick Sanan this must match petsc/finclude/petscksp.h 5208a4b9c5cSBarry Smith 521db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`, 522db781477SPatrick Sanan `KSPSetConvergenceTest()`, `KSPSetPCSide()` 5238a4b9c5cSBarry Smith E*/ 5249371c9d4SSatish Balay typedef enum { 5259371c9d4SSatish Balay KSP_NORM_DEFAULT = -1, 5269371c9d4SSatish Balay KSP_NORM_NONE = 0, 5279371c9d4SSatish Balay KSP_NORM_PRECONDITIONED = 1, 5289371c9d4SSatish Balay KSP_NORM_UNPRECONDITIONED = 2, 5299371c9d4SSatish Balay KSP_NORM_NATURAL = 3 5309371c9d4SSatish Balay } KSPNormType; 5319e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 532014dd563SJed Brown PETSC_EXTERN const char *const *const KSPNormTypes; 533a21b2a99SBarry Smith 5341f7e983dSSatish Balay /*MC 5359793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 5368c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 5378c5b8ba0SBarry Smith be based on a norm of a residual etc. 5388c5b8ba0SBarry Smith 5398c5b8ba0SBarry Smith Level: advanced 5408c5b8ba0SBarry Smith 54187497f52SBarry Smith Note: Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored 5428c5b8ba0SBarry Smith 543db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 5448c5b8ba0SBarry Smith M*/ 5458c5b8ba0SBarry Smith 5461f7e983dSSatish Balay /*MC 5471957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 5488c5b8ba0SBarry Smith convergence test routine. 5498c5b8ba0SBarry Smith 5508c5b8ba0SBarry Smith Level: advanced 5518c5b8ba0SBarry Smith 552db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 5538c5b8ba0SBarry Smith M*/ 5548c5b8ba0SBarry Smith 5551f7e983dSSatish Balay /*MC 556ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 5578c5b8ba0SBarry Smith convergence test routine. 5588c5b8ba0SBarry Smith 5598c5b8ba0SBarry Smith Level: advanced 5608c5b8ba0SBarry Smith 561db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 5628c5b8ba0SBarry Smith M*/ 5638c5b8ba0SBarry Smith 5641f7e983dSSatish Balay /*MC 565ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 56687497f52SBarry Smith convergence test routine. This is only supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR` 5678c5b8ba0SBarry Smith 5688c5b8ba0SBarry Smith Level: advanced 5698c5b8ba0SBarry Smith 570db781477SPatrick Sanan .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()` 5718c5b8ba0SBarry Smith M*/ 5728c5b8ba0SBarry Smith 573014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType); 574014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *); 575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt); 576014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt); 577014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool); 5788a4b9c5cSBarry Smith 579f0fc11ceSJed Brown #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)") 5808a4b9c5cSBarry Smith /*E 5811957e957SBarry Smith KSPConvergedReason - reason a Krylov method was said to have converged or diverged 58228ce4d24SBarry Smith 58328ce4d24SBarry Smith Level: beginner 58428ce4d24SBarry Smith 58595452b02SPatrick Sanan Notes: 58687497f52SBarry Smith See `KSPGetConvergedReason()` for explanation of each value 58728ce4d24SBarry Smith 58895452b02SPatrick Sanan Developer Notes: 58987497f52SBarry Smith This must match petsc/finclude/petscksp.h 5904d0a8057SBarry Smith 59187497f52SBarry Smith The string versions of these are `KSPConvergedReasons`; if you change 5924d0a8057SBarry Smith any of the values here also change them that array of names. 59386c02ca4SBarry Smith 594db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()` 59528ce4d24SBarry Smith E*/ 596d15094e1SBarry Smith typedef enum { /* converged */ 5979ecb1286SBarry Smith KSP_CONVERGED_RTOL_NORMAL = 1, 5989ecb1286SBarry Smith KSP_CONVERGED_ATOL_NORMAL = 9, 599d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 600d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 601b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 6028031f4adStmunson KSP_CONVERGED_CG_NEG_CURVE = 5, 6038031f4adStmunson KSP_CONVERGED_CG_CONSTRAINED = 6, 604329f5518SBarry Smith KSP_CONVERGED_STEP_LENGTH = 7, 605af3d8002SBarry Smith KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 606d15094e1SBarry Smith /* diverged */ 607b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 608d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 609d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 610d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 611b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 612b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 613b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 6144d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 6156aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 616c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED = -11, 617aa4d2078SSatish Balay KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11, 618d15094e1SBarry Smith 6199371c9d4SSatish Balay KSP_CONVERGED_ITERATING = 0 6209371c9d4SSatish Balay } KSPConvergedReason; 621014dd563SJed Brown PETSC_EXTERN const char *const *KSPConvergedReasons; 622d15094e1SBarry Smith 623c838673bSBarry Smith /*MC 62487497f52SBarry Smith KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if `KSPConvergedDefaultSetUIRNorm()` was called 625c838673bSBarry Smith 626c838673bSBarry Smith Level: beginner 627c838673bSBarry Smith 62887497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 629c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 630c838673bSBarry Smith 2-norm of the residual for right preconditioning 631c838673bSBarry Smith 63287497f52SBarry Smith See also `KSP_CONVERGED_ATOL` which may apply before this tolerance. 633f9fed41fSBarry Smith 634db781477SPatrick Sanan .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 635c838673bSBarry Smith 636c838673bSBarry Smith M*/ 637c838673bSBarry Smith 638c838673bSBarry Smith /*MC 639c838673bSBarry Smith KSP_CONVERGED_ATOL - norm(r) <= atol 640c838673bSBarry Smith 641c838673bSBarry Smith Level: beginner 642c838673bSBarry Smith 64387497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 644c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 645c838673bSBarry Smith 2-norm of the residual for right preconditioning 646c838673bSBarry Smith 64787497f52SBarry Smith See also `KSP_CONVERGED_RTOL` which may apply before this tolerance. 648c838673bSBarry Smith 649db781477SPatrick Sanan .seealso: `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 650c838673bSBarry Smith 651c838673bSBarry Smith M*/ 652c838673bSBarry Smith 653c838673bSBarry Smith /*MC 654c838673bSBarry Smith KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 655c838673bSBarry Smith 656c838673bSBarry Smith Level: beginner 657c838673bSBarry Smith 65887497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 659c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 660c838673bSBarry Smith 2-norm of the residual for right preconditioning 661c838673bSBarry Smith 662c838673bSBarry Smith Level: beginner 663c838673bSBarry Smith 664db781477SPatrick Sanan .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 665c838673bSBarry Smith 666c838673bSBarry Smith M*/ 667c838673bSBarry Smith 668c838673bSBarry Smith /*MC 669c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 670c838673bSBarry Smith reached 671c838673bSBarry Smith 672c838673bSBarry Smith Level: beginner 673c838673bSBarry Smith 674db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 675c838673bSBarry Smith 676c838673bSBarry Smith M*/ 677c838673bSBarry Smith 678c838673bSBarry Smith /*MC 67987497f52SBarry Smith KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of 68087497f52SBarry Smith the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence 6814d0a8057SBarry Smith test routine is set in KSP. 682c838673bSBarry Smith 683c838673bSBarry Smith Level: beginner 684c838673bSBarry Smith 685db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 686c838673bSBarry Smith 687c838673bSBarry Smith M*/ 688c838673bSBarry Smith 689c838673bSBarry Smith /*MC 690c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 6911de96524SPierre Jolivet method could not continue to enlarge the Krylov space. Could be due to a singular matrix or 69287497f52SBarry Smith preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block 6931de96524SPierre Jolivet are colinear. 694c838673bSBarry Smith 695c838673bSBarry Smith Level: beginner 696c838673bSBarry Smith 697db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 698c838673bSBarry Smith 699c838673bSBarry Smith M*/ 700c838673bSBarry Smith 701c838673bSBarry Smith /*MC 70287497f52SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the 703c838673bSBarry Smith method could not continue to enlarge the Krylov space. 704c838673bSBarry Smith 705c838673bSBarry Smith Level: beginner 706c838673bSBarry Smith 707db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 708c838673bSBarry Smith 709c838673bSBarry Smith M*/ 710c838673bSBarry Smith 711c838673bSBarry Smith /*MC 712c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 71387497f52SBarry Smith symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry 714c838673bSBarry Smith 715c838673bSBarry Smith Level: beginner 716c838673bSBarry Smith 717db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 718c838673bSBarry Smith 719c838673bSBarry Smith M*/ 720c838673bSBarry Smith 721c838673bSBarry Smith /*MC 722c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 72387497f52SBarry Smith positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to 724c838673bSBarry Smith be positive definite 725c838673bSBarry Smith 726c838673bSBarry Smith Level: beginner 727c838673bSBarry Smith 72887497f52SBarry Smith Note: 72987497f52SBarry Smith This can happen with the `PCICC` preconditioner, use -pc_factor_shift_positive_definite to force 73087497f52SBarry Smith the `PCICC` preconditioner to generate a positive definite preconditioner 731c838673bSBarry Smith 732db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 733c838673bSBarry Smith 734c838673bSBarry Smith M*/ 735c838673bSBarry Smith 736c838673bSBarry Smith /*MC 737c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a 7389fc87aa7SBarry Smith zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 73987497f52SBarry Smith such as `PCFIELDSPLIT`. 7409fc87aa7SBarry Smith 7419fc87aa7SBarry Smith Level: beginner 7429fc87aa7SBarry Smith 74395452b02SPatrick Sanan Notes: 74495452b02SPatrick Sanan Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details. 7459fc87aa7SBarry Smith 746db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 7479fc87aa7SBarry Smith 7489fc87aa7SBarry Smith M*/ 7499fc87aa7SBarry Smith 7509fc87aa7SBarry Smith /*MC 75187497f52SBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if you call `KSPGetConvergedReason()` 75287497f52SBarry Smith while the `KSPSolve()` is still running. 753c838673bSBarry Smith 754c838673bSBarry Smith Level: beginner 755c838673bSBarry Smith 756db781477SPatrick Sanan .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 757c838673bSBarry Smith 758c838673bSBarry Smith M*/ 759c838673bSBarry Smith 760014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *)); 761df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *)); 762df8705c3SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *)); 7633ec1f749SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *); 7648de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 7653eeb4429SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 7668de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 7678de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 7688de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 7698de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 77054b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool); 7710059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 772014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *); 773c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **); 77494a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 775abef13c0SSatish Balay 776d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") static inline void KSPDefaultConverged(void) 777d71ae5a4SJacob Faibussowitsch { /* never called */ 7789371c9d4SSatish Balay } 7798ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 780d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") static inline void KSPDefaultConvergedDestroy(void) 781d71ae5a4SJacob Faibussowitsch { /* never called */ 7829371c9d4SSatish Balay } 7838ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 784d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") static inline void KSPDefaultConvergedCreate(void) 785d71ae5a4SJacob Faibussowitsch { /* never called */ 7869371c9d4SSatish Balay } 7878ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 788d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUIRNorm(void) 789d71ae5a4SJacob Faibussowitsch { /* never called */ 7909371c9d4SSatish Balay } 7918ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 792d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUMIRNorm(void) 793d71ae5a4SJacob Faibussowitsch { /* never called */ 7949371c9d4SSatish Balay } 7958ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 796d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") static inline void KSPSkipConverged(void) 797d71ae5a4SJacob Faibussowitsch { /* never called */ 7989371c9d4SSatish Balay } 7998ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 8008ea1b3e6SJed Brown 8010bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *); 802d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B) 803d71ae5a4SJacob Faibussowitsch { 8049371c9d4SSatish Balay return KSPComputeOperator(A, NULL, B); 8059371c9d4SSatish Balay } 806d4fbbf0eSBarry Smith 80728ce4d24SBarry Smith /*E 80828ce4d24SBarry Smith KSPCGType - Determines what type of CG to use 80928ce4d24SBarry Smith 81028ce4d24SBarry Smith Level: beginner 81128ce4d24SBarry Smith 812db781477SPatrick Sanan .seealso: `KSPCGSetType()` 81328ce4d24SBarry Smith E*/ 8149371c9d4SSatish Balay typedef enum { 8159371c9d4SSatish Balay KSP_CG_SYMMETRIC = 0, 8169371c9d4SSatish Balay KSP_CG_HERMITIAN = 1 8179371c9d4SSatish Balay } KSPCGType; 8186a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 81928ce4d24SBarry Smith 820014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType); 821014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool); 8228031f4adStmunson 823dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal); 824dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *); 825dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *); 826fcae7a14Stmunson 82705de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *); 82805de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *); 829d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x) 830d71ae5a4SJacob Faibussowitsch { 8319371c9d4SSatish Balay return KSPGLTRGetMinEig(ksp, x); 8329371c9d4SSatish Balay } 833d71ae5a4SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x) 834d71ae5a4SJacob Faibussowitsch { 8359371c9d4SSatish Balay return KSPGLTRGetLambda(ksp, x); 8369371c9d4SSatish Balay } 8378031f4adStmunson 838014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]); 839ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]); 8401d6018f0SLisandro Dalcin 841f560b561SHong Zhang PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP)); 842014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP); 843014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP); 8443369ce9aSBarry Smith 8459804daf3SBarry Smith #include <petscdrawtypes.h> 846014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *); 8472f2e5d10SKris Buschelman 848014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec)); 849014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec)); 85003919abeSBarry Smith 851ba36735cSStefano Zampini /*S 852ba36735cSStefano Zampini KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods. 853f8a50e2bSBarry Smith 854ba36735cSStefano Zampini Level: beginner 855f8a50e2bSBarry Smith 856db781477SPatrick Sanan .seealso: `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType` 857ba36735cSStefano Zampini S*/ 858ba36735cSStefano Zampini typedef struct _p_KSPGuess *KSPGuess; 859ba36735cSStefano Zampini /*J 86087497f52SBarry Smith KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods. 861ba36735cSStefano Zampini 862ba36735cSStefano Zampini Level: beginner 863ba36735cSStefano Zampini 864db781477SPatrick Sanan .seealso: `KSPGuess` 865ba36735cSStefano Zampini J*/ 866ba36735cSStefano Zampini typedef const char *KSPGuessType; 867ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer" 868ba36735cSStefano Zampini #define KSPGUESSPOD "pod" 8691d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess)); 870ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess); 871ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *); 872ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer); 873ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *); 874ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *); 875ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType); 876ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *); 8778410009bSDavid Wells PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal); 878ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 879ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec); 880ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec); 881ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 882ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt); 883014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt); 884ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool); 885ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *); 886f8a50e2bSBarry Smith 887470b340bSDmitry Karpeev /*E 888470b340bSDmitry Karpeev MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 889470b340bSDmitry Karpeev 890470b340bSDmitry Karpeev Level: intermediate 891470b340bSDmitry Karpeev 892db781477SPatrick Sanan .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`, `MatCreateSchurComplementPmat()` 893470b340bSDmitry Karpeev E*/ 8949371c9d4SSatish Balay typedef enum { 8959371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_DIAG, 8969371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_LUMP, 8979371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, 8989371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_FULL 8999371c9d4SSatish Balay } MatSchurComplementAinvType; 900470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 901470b340bSDmitry Karpeev 9029371c9d4SSatish Balay typedef enum { 9039371c9d4SSatish Balay MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0, 904864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1, 905864588a7SAlp Dener MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2, 9069371c9d4SSatish Balay MAT_LMVM_SYMBROYDEN_SCALE_USER = 3 9079371c9d4SSatish Balay } MatLMVMSymBroydenScaleType; 908864588a7SAlp Dener PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[]; 90992f76d53SAlp Dener 910014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *); 911014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *); 912d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP); 913bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 914aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 915bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *); 916470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType); 917470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *); 9185bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *); 9195a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *); 920470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *); 921470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *); 9223f22127dSBarry Smith 92378e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *); 92478e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *); 92578e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *); 926864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 927864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 928864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 929864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 930864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 931cd929ea3SAlp Dener 932cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 933b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *); 934cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 935cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 936e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat); 9370ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat); 938cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 939cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 940cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 941cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 942cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 9432d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec); 944cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 945cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *); 946cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *); 947cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *); 94892f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt); 949cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *); 950cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *); 951864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar); 952864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType); 953cd929ea3SAlp Dener 954014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM); 955014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool); 956014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *); 957014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *); 958014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *); 959fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, PetscErrorCode (*func)(KSP, Vec, void *), void *); 96023ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *); 961fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, PetscErrorCode (*)(KSP, Vec, void *), void *); 96223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *); 96323ee1639SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, PetscErrorCode (**)(KSP, Mat, Mat, void *), void *); 964014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, PetscErrorCode (*)(KSP, Vec, void *), void *); 965014dd563SJed Brown PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, PetscErrorCode (**)(KSP, Vec, void *), void *); 966fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, PetscErrorCode (*)(KSP, Vec, void *), void *); 967fbd07b44SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, PetscErrorCode (**)(KSP, Vec, void *), void *); 9686c699258SBarry Smith 96902b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec); 9709371c9d4SSatish 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); 971557cf195SMatthew G. Knepley 9722b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *); 9732b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal); 9742eac72dbSBarry Smith #endif 975