1f26ada1bSBarry Smith /* 2f26ada1bSBarry Smith Defines the interface functions for the Krylov subspace accelerators. 3f26ada1bSBarry Smith */ 4a4963045SJacob Faibussowitsch #pragma once 5ac09b921SBarry Smith 62c8e378dSBarry Smith #include <petscpc.h> 72eac72dbSBarry Smith 8ac09b921SBarry Smith /* SUBMANSEC = KSP */ 9ac09b921SBarry Smith 10607a6623SBarry Smith PETSC_EXTERN PetscErrorCode KSPInitializePackage(void); 114bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode KSPFinalizePackage(void); 121dbb0a54SBarry Smith 1328ce4d24SBarry Smith /*S 140b4b7b1cSBarry Smith KSP - Abstract PETSc object that manages the linear solves in PETSc (even those such as direct factorization-based solvers that 150b4b7b1cSBarry Smith do not use Krylov accelerators). 1628ce4d24SBarry Smith 1728ce4d24SBarry Smith Level: beginner 1828ce4d24SBarry Smith 190b4b7b1cSBarry Smith Notes: 20a4d1885cSBarry Smith When a direct solver is used, but no Krylov solver is used, the `KSP` object is still used but with a 2116a05f60SBarry Smith `KSPType` of `KSPPREONLY` (or equivalently `KSPNONE`), meaning that only application of the preconditioner is used as the linear solver. 228f6c3df8SBarry Smith 230b4b7b1cSBarry Smith Use `KSPSetType()` or the options database key `-ksp_type` to set the specific Krylov solver algorithm to use with a given `KSP` object 240b4b7b1cSBarry Smith 250b4b7b1cSBarry Smith The `PC` object is used to control preconditioners in PETSc. 260b4b7b1cSBarry Smith 27789736e1SBarry Smith `KSP` can also be used to solve some least squares problems (over or under-determined linear systems), using, for example, `KSPLSQR`, see `PETSCREGRESSORLINEAR` 28789736e1SBarry Smith for additional methods that can be used to solve least squares problems and other linear regressions). 29789736e1SBarry Smith 301cc06b55SBarry Smith .seealso: [](doc_linsolve), [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES` 3128ce4d24SBarry Smith S*/ 32e2a1c21fSSatish Balay typedef struct _p_KSP *KSP; 332eac72dbSBarry Smith 3476bdecfbSBarry Smith /*J 350b4b7b1cSBarry Smith KSPType - String with the name of a PETSc Krylov method. These are all the Krylov solvers that PETSc provides. 3628ce4d24SBarry Smith 3728ce4d24SBarry Smith Level: beginner 3828ce4d24SBarry Smith 391cc06b55SBarry Smith .seealso: [](doc_linsolve), [](ch_ksp), `KSPSetType()`, `KSP`, `KSPRegister()`, `KSPCreate()`, `KSPSetFromOptions()` 4076bdecfbSBarry Smith J*/ 4119fd82e9SBarry Smith typedef const char *KSPType; 4282bf6240SBarry Smith #define KSPRICHARDSON "richardson" 436c9de887SHong Zhang #define KSPCHEBYSHEV "chebyshev" 4482bf6240SBarry Smith #define KSPCG "cg" 452c8fc646SJed Brown #define KSPGROPPCG "groppcg" 462c8fc646SJed Brown #define KSPPIPECG "pipecg" 47901ccb91SSiegfried Cools #define KSPPIPECGRR "pipecgrr" 48265663fdSSiegfried Cools #define KSPPIPELCG "pipelcg" 49b21a8899STyler Chen #define KSPPIPEPRCG "pipeprcg" 50325e8391SManasi #define KSPPIPECG2 "pipecg2" 51df2a969aSvictorle #define KSPCGNE "cgne" 5205de396fSBarry Smith #define KSPNASH "nash" 5305de396fSBarry Smith #define KSPSTCG "stcg" 5405de396fSBarry Smith #define KSPGLTR "gltr" 55edd03b47SJacob Faibussowitsch #define KSPCGNASH PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPNASH", ) "nash" 56edd03b47SJacob Faibussowitsch #define KSPCGSTCG PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSTCG", ) "stcg" 57edd03b47SJacob Faibussowitsch #define KSPCGGLTR PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSGLTR", ) "gltr" 58640f4f53SPatrick Sanan #define KSPFCG "fcg" 59390d8e47SPatrick Sanan #define KSPPIPEFCG "pipefcg" 6082bf6240SBarry Smith #define KSPGMRES "gmres" 61483d6965SPatrick Sanan #define KSPPIPEFGMRES "pipefgmres" 629dcbbd2bSBarry Smith #define KSPFGMRES "fgmres" 639dcbbd2bSBarry Smith #define KSPLGMRES "lgmres" 644f8e6cd9SSatish Balay #define KSPDGMRES "dgmres" 6561b468f9SJed Brown #define KSPPGMRES "pgmres" 6682bf6240SBarry Smith #define KSPTCQMR "tcqmr" 6782bf6240SBarry Smith #define KSPBCGS "bcgs" 68e1c61ce8SBarry Smith #define KSPIBCGS "ibcgs" 69345ecf0bSXiangmin Jiao #define KSPQMRCGS "qmrcgs" 7018ac38e6SHong Zhang #define KSPFBCGS "fbcgs" 71c2b71004SHong Zhang #define KSPFBCGSR "fbcgsr" 72f0037002Svictorle #define KSPBCGSL "bcgsl" 73f154af2dSSiegfried Cools #define KSPPIPEBCGS "pipebcgs" 7482bf6240SBarry Smith #define KSPCGS "cgs" 7582bf6240SBarry Smith #define KSPTFQMR "tfqmr" 7682bf6240SBarry Smith #define KSPCR "cr" 772c8fc646SJed Brown #define KSPPIPECR "pipecr" 7882bf6240SBarry Smith #define KSPLSQR "lsqr" 7982bf6240SBarry Smith #define KSPPREONLY "preonly" 803c2be86cSBarry Smith #define KSPNONE "none" 8182bf6240SBarry Smith #define KSPQCG "qcg" 82c9cf9b11SBarry Smith #define KSPBICG "bicg" 83b4ac9ba4SBarry Smith #define KSPMINRES "minres" 8401c5daebSSatish Balay #define KSPSYMMLQ "symmlq" 8521356919SSatish Balay #define KSPLCD "lcd" 861d6018f0SLisandro Dalcin #define KSPPYTHON "python" 8758ad63f4SBarry Smith #define KSPGCR "gcr" 88fad47a0aSPatrick Sanan #define KSPPIPEGCR "pipegcr" 89e4d80e07Szianekhodja #define KSPTSIRM "tsirm" 90e4d80e07Szianekhodja #define KSPCGLS "cgls" 91329cd39dSStefano Zampini #define KSPFETIDP "fetidp" 9238cfc46eSPierre Jolivet #define KSPHPDDM "hpddm" 932eac72dbSBarry Smith 948ba1e511SMatthew Knepley /* Logging support */ 95014dd563SJed Brown PETSC_EXTERN PetscClassId KSP_CLASSID; 96ba36735cSStefano Zampini PETSC_EXTERN PetscClassId KSPGUESS_CLASSID; 975358d0d4SBarry Smith PETSC_EXTERN PetscClassId DMKSP_CLASSID; 988ba1e511SMatthew Knepley 99014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm, KSP *); 10019fd82e9SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetType(KSP, KSPType); 101c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetType(KSP, KSPType *); 102014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 103014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 104014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolve(KSP, Vec, Vec); 105014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP, Vec, Vec); 10675f8e9bdSHong Zhang PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP, PetscBool); 1075ea7661aSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP, Mat, Mat); 108bbbebc2cSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPMatSolveTranspose(KSP, Mat, Mat); 1093e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt); 110edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPSetMatSolveBatchSize()", ) static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n) 111d71ae5a4SJacob Faibussowitsch { 1129371c9d4SSatish Balay return KSPSetMatSolveBatchSize(ksp, n); 1139371c9d4SSatish Balay } 1143e00ce70SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *); 115edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPGetMatSolveBatchSize()", ) static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n) 116d71ae5a4SJacob Faibussowitsch { 1179371c9d4SSatish Balay return KSPGetMatSolveBatchSize(ksp, n); 1189371c9d4SSatish Balay } 119014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPReset(KSP); 120ef17e8b6SStefano Zampini PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP); 121014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *); 12223ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool); 1237d85ae06SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *); 12425c41539SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool); 125c0decd05SBarry Smith PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec); 1262eac72dbSBarry Smith 127140e18c1SBarry Smith PETSC_EXTERN PetscFunctionList KSPList; 128ba36735cSStefano Zampini PETSC_EXTERN PetscFunctionList KSPGuessList; 129798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorList; 130798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorCreateList; 131798534f6SMatthew G. Knepley PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList; 132bdf89e91SBarry Smith PETSC_EXTERN PetscErrorCode KSPRegister(const char[], PetscErrorCode (*)(KSP)); 1334d4d2bdcSBarry Smith 1344d4d2bdcSBarry Smith /*S 1354d4d2bdcSBarry Smith KSPMonitorRegisterFn - A function prototype for functions provided to `KSPMonitorRegister()` 1364d4d2bdcSBarry Smith 1374d4d2bdcSBarry Smith Calling Sequence: 1384d4d2bdcSBarry Smith + ksp - iterative solver obtained from `KSPCreate()` 1394d4d2bdcSBarry Smith . it - iteration number 1404d4d2bdcSBarry Smith . rnorm - (estimated) 2-norm of (preconditioned) residual 1414d4d2bdcSBarry Smith - ctx - `PetscViewerAndFormat` object 1424d4d2bdcSBarry Smith 1434d4d2bdcSBarry Smith Level: beginner 1444d4d2bdcSBarry Smith 1454d4d2bdcSBarry Smith Note: 1464d4d2bdcSBarry Smith This is a `KSPMonitorFn` specialized for a context of `PetscViewerAndFormat` 1474d4d2bdcSBarry Smith 1484d4d2bdcSBarry Smith .seealso: [](ch_snes), `KSP`, `KSPMonitorSet()`, `KSPMonitorRegister()`, `KSPMonitorFn`, `KSPMonitorRegisterCreateFn`, `KSPMonitorRegisterDestroyFn` 1494d4d2bdcSBarry Smith S*/ 1504d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPMonitorRegisterFn(KSP ksp, PetscInt it, PetscReal rnorm, PetscViewerAndFormat *ctx); 1514d4d2bdcSBarry Smith 1524d4d2bdcSBarry Smith /*S 1534d4d2bdcSBarry Smith KSPMonitorRegisterCreateFn - A function prototype for functions that do the creation when provided to `KSPMonitorRegister()` 1544d4d2bdcSBarry Smith 1554d4d2bdcSBarry Smith Calling Sequence: 1564d4d2bdcSBarry Smith + viewer - the viewer to be used with the `KSPMonitorRegisterFn` 1574d4d2bdcSBarry Smith . format - the format of the viewer 1584d4d2bdcSBarry Smith . ctx - a context for the monitor 1594d4d2bdcSBarry Smith - result - a `PetscViewerAndFormat` object 1604d4d2bdcSBarry Smith 1614d4d2bdcSBarry Smith Level: beginner 1624d4d2bdcSBarry Smith 1634d4d2bdcSBarry Smith .seealso: [](ch_snes), `KSPMonitorRegisterFn`, `KSP`, `KSPMonitorSet()`, `KSPMonitorRegister()`, `KSPMonitorFn`, `KSPMonitorRegisterDestroyFn` 1644d4d2bdcSBarry Smith S*/ 1652a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPMonitorRegisterCreateFn(PetscViewer viewer, PetscViewerFormat format, PetscCtx ctx, PetscViewerAndFormat **result); 1664d4d2bdcSBarry Smith 1674d4d2bdcSBarry Smith /*S 1684d4d2bdcSBarry Smith KSPMonitorRegisterDestroyFn - A function prototype for functions that do the after use destruction when provided to `KSPMonitorRegister()` 1694d4d2bdcSBarry Smith 1704d4d2bdcSBarry Smith Calling Sequence: 1714d4d2bdcSBarry Smith . vf - a `PetscViewerAndFormat` object to be destroyed, including any context 1724d4d2bdcSBarry Smith 1734d4d2bdcSBarry Smith Level: beginner 1744d4d2bdcSBarry Smith 1754d4d2bdcSBarry Smith .seealso: [](ch_snes), `KSPMonitorRegisterFn`, `KSP`, `KSPMonitorSet()`, `KSPMonitorRegister()`, `KSPMonitorFn`, `KSPMonitorRegisterCreateFn` 1764d4d2bdcSBarry Smith S*/ 1774d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPMonitorRegisterDestroyFn(PetscViewerAndFormat **result); 1784d4d2bdcSBarry Smith 1794d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, KSPMonitorRegisterFn *, KSPMonitorRegisterCreateFn *, KSPMonitorRegisterDestroyFn *); 18030de9b25SBarry Smith 181014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide); 182014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *); 183014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt); 184c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *); 18525c92fe2SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetMinimumIterations(KSP, PetscInt); 18625c92fe2SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetMinimumIterations(KSP, PetscInt *); 187014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool); 188014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *); 189014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool); 190014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *); 191014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool); 1927d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool); 193c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *); 194014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool); 195c60c7ad4SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *); 196014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *); 197014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *); 198014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *); 199014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *); 200734794cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *); 2012a7a6963SBarry Smith PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **); 202edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 6, 0, "KSPCreateVecs()", ) static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y) 203d71ae5a4SJacob Faibussowitsch { 2049371c9d4SSatish Balay return KSPCreateVecs(ksp, n, x, m, y); 2059371c9d4SSatish Balay } 2062eac72dbSBarry Smith 2074d4d2bdcSBarry Smith /*S 2084d4d2bdcSBarry Smith KSPPSolveFn - A function prototype for functions provided to `KSPSetPreSolve()` and `KSPSetPostSolve()` 2094d4d2bdcSBarry Smith 2104d4d2bdcSBarry Smith Calling Sequence: 2114d4d2bdcSBarry Smith + ksp - the `KSP` context 2124d4d2bdcSBarry Smith . rhs - the right-hand side vector 2134d4d2bdcSBarry Smith . x - the solution vector 2144d4d2bdcSBarry Smith - ctx - optional context that was provided with `KSPSetPreSolve()` or `KSPSetPostSolve()` 2154d4d2bdcSBarry Smith 2164d4d2bdcSBarry Smith Level: intermediate 2174d4d2bdcSBarry Smith 2184d4d2bdcSBarry Smith .seealso: [](ch_snes), `KSP`, `KSPSetPreSolve()`, `KSPSetPostSolve()`, `PCShellPSolveFn` 2194d4d2bdcSBarry Smith S*/ 2202a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPPSolveFn(KSP ksp, Vec rhs, Vec x, PetscCtx ctx); 2214d4d2bdcSBarry Smith 2222a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, KSPPSolveFn *, PetscCtx); 2232a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, KSPPSolveFn *, PetscCtx); 224d4211eb9SBarry Smith 225014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC); 226014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *); 2273821be0aSBarry Smith PETSC_EXTERN PetscErrorCode KSPSetNestLevel(KSP, PetscInt); 2283821be0aSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetNestLevel(KSP, PetscInt *); 229aabeff55SBarry Smith 2304d4d2bdcSBarry Smith /*S 2314d4d2bdcSBarry Smith KSPMonitorFn - A function prototype for functions provided to `KSPMonitorSet()` 2324d4d2bdcSBarry Smith 2334d4d2bdcSBarry Smith Calling Sequence: 2344d4d2bdcSBarry Smith + ksp - iterative solver obtained from `KSPCreate()` 2354d4d2bdcSBarry Smith . it - iteration number 2364d4d2bdcSBarry Smith . rnorm - (estimated) 2-norm of (preconditioned) residual 2374d4d2bdcSBarry Smith - ctx - optional monitoring context, as provided with `KSPMonitorSet()` 2384d4d2bdcSBarry Smith 2394d4d2bdcSBarry Smith Level: beginner 2404d4d2bdcSBarry Smith 2414d4d2bdcSBarry Smith .seealso: [](ch_snes), `KSP`, `KSPMonitorSet()` 2424d4d2bdcSBarry Smith S*/ 2432a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPMonitorFn(KSP ksp, PetscInt it, PetscReal rnorm, PetscCtx ctx); 2444d4d2bdcSBarry Smith 245014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal); 2462a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, KSPMonitorFn *, PetscCtx, PetscCtxDestroyFn *); 247014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 2482a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, PetscCtxRt); 24994a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *); 2506497c311SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscCount, PetscBool); 25194a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *); 2526497c311SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscCount, PetscBool); 2534b0e389bSBarry Smith 254fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *); 255fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *); 256fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 257fa0ddf94SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt); 258fa0ddf94SBarry Smith 259014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *); 260cfebe74eSStefano Zampini PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP); 261014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 262014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 263014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 264fb80e629SPablo Brubeck PETSC_EXTERN PetscErrorCode PCPatchGetSubKSP(PC, PetscInt *, KSP *[]); 265014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]); 266285fb4e2SStefano Zampini PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]); 267b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *); 268b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *); 269b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *); 270b4876bcbSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *); 271014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *); 2724a99276eSJakub Kruzik PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *); 2734d4d2bdcSBarry Smith 2744d4d2bdcSBarry Smith /*S 2754d4d2bdcSBarry Smith PCMGCoarseSpaceConstructorFn - A function prototype for functions registered with `PCMGRegisterCoarseSpaceConstructor()` 2764d4d2bdcSBarry Smith 2774d4d2bdcSBarry Smith Calling Sequence: 2784d4d2bdcSBarry Smith + pc - The `PC` object 2794d4d2bdcSBarry Smith . l - The multigrid level, 0 is the coarse level 2804d4d2bdcSBarry Smith . dm - The `DM` for this level 2814d4d2bdcSBarry Smith . smooth - The level smoother 2824d4d2bdcSBarry Smith . Nc - The size of the coarse space 2834d4d2bdcSBarry Smith . initGuess - Basis for an initial guess for the space 2844d4d2bdcSBarry Smith - coarseSp - A basis for the computed coarse space 2854d4d2bdcSBarry Smith 2864d4d2bdcSBarry Smith Level: beginner 2874d4d2bdcSBarry Smith 2884d4d2bdcSBarry Smith .seealso: [](ch_ksp), `PCMGRegisterCoarseSpaceConstructor()`, `PCMGGetCoarseSpaceConstructor()` 2894d4d2bdcSBarry Smith S*/ 2904d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PCMGCoarseSpaceConstructorFn(PC pc, PetscInt l, DM dm, KSP smooth, PetscInt Nc, Mat initGuess, Mat *coarseSp); 2914d4d2bdcSBarry Smith 292f3b08a26SMatthew G. Knepley PETSC_EXTERN PetscFunctionList PCMGCoarseList; 2934d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PCMGCoarseSpaceConstructorFn *); 2944d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PCMGCoarseSpaceConstructorFn **); 295f3b08a26SMatthew G. Knepley 296014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *); 297014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *); 2982eac72dbSBarry Smith 299f2edd1f0SMalachi Phillips /*E 30095bd0b28SBarry Smith KSPChebyshevKind - Which kind of Chebyshev polynomial to use with `KSPCHEBYSHEV` 301f2edd1f0SMalachi Phillips 302f2edd1f0SMalachi Phillips Values: 303f2edd1f0SMalachi Phillips + `KSP_CHEBYSHEV_FIRST` - "classic" first-kind Chebyshev polynomial 304f2edd1f0SMalachi Phillips . `KSP_CHEBYSHEV_FOURTH` - fourth-kind Chebyshev polynomial 305f2edd1f0SMalachi Phillips - `KSP_CHEBYSHEV_OPT_FOURTH` - optimized fourth-kind Chebyshev polynomial 306f2edd1f0SMalachi Phillips 307f2edd1f0SMalachi Phillips Level: intermediate 308f2edd1f0SMalachi Phillips 3091cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPCHEBYSHEV`, `KSPChebyshevSetKind()`, `KSPChebyshevGetKind()` 310f2edd1f0SMalachi Phillips E*/ 311f2edd1f0SMalachi Phillips typedef enum { 312f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_FIRST, 313f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_FOURTH, 314f2edd1f0SMalachi Phillips KSP_CHEBYSHEV_OPT_FOURTH 315f2edd1f0SMalachi Phillips } KSPChebyshevKind; 316f2edd1f0SMalachi Phillips 317014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal); 318014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool); 319014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal); 32058450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal); 321b65fd476SMark Adams PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool); 322f2edd1f0SMalachi Phillips PETSC_EXTERN PetscErrorCode KSPChebyshevSetKind(KSP, KSPChebyshevKind); 32316a05f60SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevGetKind(KSP, KSPChebyshevKind *); 32458450487SBarry Smith PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *); 325014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *); 326d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *); 327d38fc84dSBarry Smith PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]); 3287d7fa35eSSylvain Mercier PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]); 3294b0e389bSBarry Smith 330640f4f53SPatrick Sanan /*E 331640f4f53SPatrick Sanan 3326364ff28SBarry Smith KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate gradient/directions methods 333640f4f53SPatrick Sanan 33467b8a455SSatish Balay Values: 3356364ff28SBarry Smith + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to `mmax`) stored directions 3366364ff28SBarry Smith - `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last `max(1,mod(i,mmax))` stored directions at iteration i = 0, 1, ... 337640f4f53SPatrick Sanan 3382b26979fSBarry Smith Level: intermediate 339640f4f53SPatrick Sanan 3406364ff28SBarry Smith Note: 3416364ff28SBarry Smith Function such as `KSPFCGSetMmax()`, `KSPPIPEGCRSetNMax()`, `KSPPIPEGCRSetNMax()`, and `KSPPIPEFCGSetNMax()` may be 3426364ff28SBarry Smith used to provide `nmax` or they may be provided with the option database. 3436364ff28SBarry Smith 3446364ff28SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()`, 3456364ff28SBarry Smith `KSPPIPEGCRSetTruncationType()`, `KSPPIPEGCRGetTruncationType()`, 3466364ff28SBarry Smith `KSPFCGSetMmax()`, `KSPPIPEGCRSetNMax()`, `KSPPIPEGCRGetNMax()`, `KSPPIPEFCGGetNMax()` 347640f4f53SPatrick Sanan E*/ 3489371c9d4SSatish Balay typedef enum { 3499371c9d4SSatish Balay KSP_FCD_TRUNC_TYPE_STANDARD, 3509371c9d4SSatish Balay KSP_FCD_TRUNC_TYPE_NOTAY 3519371c9d4SSatish Balay } KSPFCDTruncationType; 35206137d0aSPatrick Sanan PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 353640f4f53SPatrick Sanan 354640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt); 355640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *); 356640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt); 357640f4f53SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *); 35806137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType); 35906137d0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *); 360640f4f53SPatrick Sanan 361390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt); 362390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *); 363390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt); 364390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *); 365390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType); 366390d8e47SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *); 367390d8e47SPatrick Sanan 368fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt); 369fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *); 370fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt); 371fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *); 372fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType); 373fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *); 374fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool); 375fad47a0aSPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *); 3764d4d2bdcSBarry Smith 3774d4d2bdcSBarry Smith /*S 3784d4d2bdcSBarry Smith KSPFlexibleModifyPCFn - A prototype of a function used to modify the preconditioner during the use of flexible `KSP` methods, such as `KSPFGMRES` 3794d4d2bdcSBarry Smith 3804d4d2bdcSBarry Smith Calling Sequence: 3814d4d2bdcSBarry Smith + ksp - the `KSP` context being used. 3824d4d2bdcSBarry Smith . total_its - the total number of iterations that have occurred. 3834d4d2bdcSBarry Smith . local_its - the number of iterations since last restart if applicable 3844d4d2bdcSBarry Smith . res_norm - the current residual norm 3856364ff28SBarry Smith - ctx - optional context variable set with `KSPFlexibleSetModifyPC()` 3864d4d2bdcSBarry Smith 3874d4d2bdcSBarry Smith Level: beginner 3884d4d2bdcSBarry Smith 3896364ff28SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPFlexibleSetModifyPC()` 3904d4d2bdcSBarry Smith S*/ 3912a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPFlexibleModifyPCFn(KSP ksp, PetscInt total_its, PetscInt local_its, PetscReal res_norm, PetscCtx ctx); 3924d4d2bdcSBarry Smith 3932a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPFlexibleSetModifyPC(KSP, KSPFlexibleModifyPCFn *, PetscCtx, PetscCtxDestroyFn *); 3946364ff28SBarry Smith 3956364ff28SBarry Smith PETSC_DEPRECATED_FUNCTION(3, 25, 0, "KSPFlexibleSetModifyPC()", ) 3966364ff28SBarry Smith static inline PetscErrorCode KSPPIPEGCRSetModifyPC(KSP ksp, KSPFlexibleModifyPCFn *fun, PetscCtx ctx, PetscCtxDestroyFn *dfun) 3976364ff28SBarry Smith { 3986364ff28SBarry Smith return KSPFlexibleSetModifyPC(ksp, fun, ctx, dfun); 3996364ff28SBarry Smith } 40083f0b094SBarry Smith 401014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 402014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *); 403014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal); 404e3729115SFande Kong PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal); 4059f236934SBarry Smith 406014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 407014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt)); 408014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt)); 409014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt); 410014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt); 4111d73ed98SBarry Smith 412014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt); 413014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 4141d73ed98SBarry Smith 415483d6965SPatrick Sanan PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar); 416483d6965SPatrick Sanan 417014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt); 418014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *); 4196364ff28SBarry Smith 4206364ff28SBarry Smith PETSC_DEPRECATED_FUNCTION(3, 25, 0, "KSPFlexibleSetModifyPC()", ) 4216364ff28SBarry Smith static inline PetscErrorCode KSPGCRSetModifyPC(KSP ksp, KSPFlexibleModifyPCFn *fun, PetscCtx ctx, PetscCtxDestroyFn *dfun) 4226364ff28SBarry Smith { 4236364ff28SBarry Smith return KSPFlexibleSetModifyPC(ksp, fun, ctx, dfun); 4246364ff28SBarry Smith } 42558ad63f4SBarry Smith 426f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal); 427f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *); 428f87a0b54SStefano Zampini PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool); 429f87a0b54SStefano Zampini 43004d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *); 43104d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC); 43204d8bde8SStefano Zampini PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *); 433568c15deSstefano_zampini PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat); 434e82af88dSprj- 4356cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat); 4366cac28cbSPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *); 4376cac28cbSPierre Jolivet #if PetscDefined(HAVE_HPDDM) 438edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMSetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) 439d71ae5a4SJacob Faibussowitsch { 4409371c9d4SSatish Balay return KSPHPDDMSetDeflationMat(ksp, U); 4419371c9d4SSatish Balay } 442edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMGetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) 443d71ae5a4SJacob Faibussowitsch { 4449371c9d4SSatish Balay return KSPHPDDMGetDeflationMat(ksp, U); 4459371c9d4SSatish Balay } 4466cac28cbSPierre Jolivet #endif 447edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPMatSolve()", ) static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) 448d71ae5a4SJacob Faibussowitsch { 4499371c9d4SSatish Balay return KSPMatSolve(ksp, B, X); 4509371c9d4SSatish Balay } 451d55ab951SPierre Jolivet /*E 45287497f52SBarry Smith KSPHPDDMType - Type of Krylov method used by `KSPHPDDM` 453d55ab951SPierre Jolivet 454d55ab951SPierre Jolivet Values: 455a4d1885cSBarry Smith + `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method 456a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BGMRES` - block GMRES 457a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_CG` - Conjugate Gradient 458a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BCG` - block CG 459a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_GCRODR` - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting 460a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BGCRODR` - block GCRODR 461a4d1885cSBarry Smith . `KSP_HPDDM_TYPE_BFBCG` - breakdown-free BCG 462a4d1885cSBarry Smith - `KSP_HPDDM_TYPE_PREONLY` - apply the preconditioner only 463d55ab951SPierre Jolivet 46416a05f60SBarry Smith Level: intermediate 46516a05f60SBarry Smith 4661cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPHPDDM`, `KSPHPDDMSetType()` 467d55ab951SPierre Jolivet E*/ 468d55ab951SPierre Jolivet typedef enum { 469d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GMRES = 0, 470d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGMRES = 1, 471d55ab951SPierre Jolivet KSP_HPDDM_TYPE_CG = 2, 472d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BCG = 3, 473d55ab951SPierre Jolivet KSP_HPDDM_TYPE_GCRODR = 4, 474d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BGCRODR = 5, 475d55ab951SPierre Jolivet KSP_HPDDM_TYPE_BFBCG = 6, 476d55ab951SPierre Jolivet KSP_HPDDM_TYPE_PREONLY = 7 477d55ab951SPierre Jolivet } KSPHPDDMType; 478d55ab951SPierre Jolivet PETSC_EXTERN const char *const KSPHPDDMTypes[]; 479a4d1885cSBarry Smith 480d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType); 481d55ab951SPierre Jolivet PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *); 482e82af88dSprj- 483b49fd9e1SBarry Smith /*E 48416a05f60SBarry Smith KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed in the GMRES solvers 485b49fd9e1SBarry Smith 486a4d1885cSBarry Smith Values: 487a4d1885cSBarry Smith + `KSP_GMRES_CGS_REFINE_NEVER` - one step of classical Gram-Schmidt 488a4d1885cSBarry Smith . `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria 489a4d1885cSBarry Smith - `KSP_GMRES_CGS_REFINE_ALWAYS` - always perform two steps 490b49fd9e1SBarry Smith 49116a05f60SBarry Smith Level: advanced 49216a05f60SBarry Smith 49395bd0b28SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPGMRES`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 494a4d1885cSBarry Smith `KSPGMRESGetOrthogonalization()`, 495a4d1885cSBarry Smith `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()` 496b49fd9e1SBarry Smith E*/ 4979371c9d4SSatish Balay typedef enum { 4989371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_NEVER, 4999371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_IFNEEDED, 5009371c9d4SSatish Balay KSP_GMRES_CGS_REFINE_ALWAYS 5019371c9d4SSatish Balay } KSPGMRESCGSRefinementType; 5026a6fc655SJed Brown PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 50316a05f60SBarry Smith 5041f7e983dSSatish Balay /*MC 5051957e957SBarry Smith KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 5068c5b8ba0SBarry Smith 5078c5b8ba0SBarry Smith Level: advanced 5088c5b8ba0SBarry Smith 50987497f52SBarry Smith Note: 51087497f52SBarry Smith Possibly unstable, but the fastest to compute 5118c5b8ba0SBarry Smith 51295bd0b28SBarry Smith .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 513a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 514db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 515db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 5168c5b8ba0SBarry Smith M*/ 5178c5b8ba0SBarry Smith 5181f7e983dSSatish Balay /*MC 5198c5b8ba0SBarry Smith KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 5208c5b8ba0SBarry Smith iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 5218c5b8ba0SBarry Smith poor orthogonality. 5228c5b8ba0SBarry Smith 5238c5b8ba0SBarry Smith Level: advanced 5248c5b8ba0SBarry Smith 525a4d1885cSBarry Smith Note: 526a4d1885cSBarry Smith This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to 5278c5b8ba0SBarry Smith estimate the orthogonality but is more stable. 5288c5b8ba0SBarry Smith 52995bd0b28SBarry Smith .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 530a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 531db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 532db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 5338c5b8ba0SBarry Smith M*/ 5348c5b8ba0SBarry Smith 5351f7e983dSSatish Balay /*MC 53616a05f60SBarry Smith KSP_GMRES_CGS_REFINE_ALWAYS - Do two steps of the classical (unmodified) Gram-Schmidt process. 5378c5b8ba0SBarry Smith 5388c5b8ba0SBarry Smith Level: advanced 5398c5b8ba0SBarry Smith 54087497f52SBarry Smith Notes: 54187497f52SBarry Smith This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice 54287497f52SBarry Smith but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`. 5438c5b8ba0SBarry Smith 5448c5b8ba0SBarry Smith You should only use this if you absolutely know that the iterative refinement is needed. 5458c5b8ba0SBarry Smith 54695bd0b28SBarry Smith .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 547a4d1885cSBarry Smith `KSP`, `KSPGMRESGetOrthogonalization()`, 548db781477SPatrick Sanan `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 549db781477SPatrick Sanan `KSPGMRESModifiedGramSchmidtOrthogonalization()` 5508c5b8ba0SBarry Smith M*/ 5518c5b8ba0SBarry Smith 552014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType); 553014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *); 55408480c60SBarry Smith 5556364ff28SBarry Smith PETSC_EXTERN KSPFlexibleModifyPCFn KSPFlexibleModifyPCNoChange; 5566364ff28SBarry Smith PETSC_EXTERN KSPFlexibleModifyPCFn KSPFlexibleModifyPCKSP; 5576364ff28SBarry Smith 5586364ff28SBarry Smith PETSC_DEPRECATED_FUNCTION(3, 25, 0, "KSPFlexibleModifyPCNoChange()", ) 5596364ff28SBarry Smith static inline PetscErrorCode KSPFGMRESModifyPCNoChange(KSP ksp, PetscInt total_its, PetscInt loc_its, PetscReal res_norm, PetscCtx ctx) 5606364ff28SBarry Smith { 5616364ff28SBarry Smith return KSPFlexibleModifyPCNoChange(ksp, total_its, loc_its, res_norm, ctx); 5626364ff28SBarry Smith } 5636364ff28SBarry Smith 5646364ff28SBarry Smith PETSC_DEPRECATED_FUNCTION(3, 25, 0, "KSPFlexibleModifyPCKSP()", ) 5656364ff28SBarry Smith static inline PetscErrorCode KSPFGMRESModifyPCKSP(KSP ksp, PetscInt total_its, PetscInt loc_its, PetscReal res_norm, PetscCtx ctx) 5666364ff28SBarry Smith { 5676364ff28SBarry Smith return KSPFlexibleModifyPCKSP(ksp, total_its, loc_its, res_norm, ctx); 5686364ff28SBarry Smith } 569c38d4ed2SBarry Smith 570014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal); 571014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *); 572014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *); 573121fd945SKris Buschelman 574014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal); 575014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool); 576014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt); 577e4b508b2SJed Brown PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool); 578d9492815SBarry Smith 579014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 58087d6299eSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP); 5812eac72dbSBarry Smith 5822a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], PetscCtx); 5834d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidual; 58452e0830fSMatthew Knepley PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualView; 58552e0830fSMatthew Knepley PETSC_DEPRECATED_FUNCTION(3, 23, 0, "KSPMonitorResidualDraw()", ) static inline PetscErrorCode KSPMonitorResidualDraw(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 58652e0830fSMatthew Knepley { 58752e0830fSMatthew Knepley return KSPMonitorResidualView(ksp, n, rnorm, vf); 58852e0830fSMatthew Knepley } 5894d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualDrawLG; 5902a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **); 5914d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualShort; 5924d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorResidualRange; 5934d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidual; 59452e0830fSMatthew Knepley PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualView; 59552e0830fSMatthew Knepley PETSC_DEPRECATED_FUNCTION(3, 23, 0, "KSPMonitorTrueResidualDraw()", ) static inline PetscErrorCode KSPMonitorTrueResidualDraw(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 59652e0830fSMatthew Knepley { 59752e0830fSMatthew Knepley return KSPMonitorTrueResidualView(ksp, n, rnorm, vf); 59852e0830fSMatthew Knepley } 5994d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualDrawLG; 6002a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **); 6014d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorTrueResidualMax; 6024d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorError; 6034d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorErrorDraw; 6044d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorErrorDrawLG; 6052a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **); 6064d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolution; 6074d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolutionDraw; 6084d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSolutionDrawLG; 6092a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **); 6104d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPMonitorSingularValue; 6112a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, PetscCtx, PetscViewerAndFormat **); 612edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorResidual()", ) static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 613d71ae5a4SJacob Faibussowitsch { 6149371c9d4SSatish Balay return KSPMonitorResidual(ksp, n, rnorm, vf); 6159371c9d4SSatish Balay } 616edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidual()", ) static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 617d71ae5a4SJacob Faibussowitsch { 6189371c9d4SSatish Balay return KSPMonitorTrueResidual(ksp, n, rnorm, vf); 6199371c9d4SSatish Balay } 620edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidualMax()", ) static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 621d71ae5a4SJacob Faibussowitsch { 6229371c9d4SSatish Balay return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf); 6239371c9d4SSatish Balay } 624798534f6SMatthew G. Knepley 625798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *); 6262a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, PetscCtx); 6272a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(PetscCtxRt); 628341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *); 629341b8579SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal); 630e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *); 631e04113cfSBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **); 6322a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(PetscCtxRt); 63384cb2905SBarry Smith 634014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec); 635014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec); 636c01c455dSBarry Smith 63723ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat); 63823ee1639SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *); 639014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *); 640014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]); 641014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]); 642014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]); 6431eb62cbbSBarry Smith 644014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool); 645014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *); 646014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool); 647014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *); 6481f7f0c4fSBarry Smith 6494d4d2bdcSBarry Smith /*S 6504d4d2bdcSBarry Smith KSPConvergedReasonViewFn - A prototype of a function used with `KSPConvergedReasonViewSet()` 6514d4d2bdcSBarry Smith 6524d4d2bdcSBarry Smith Calling Sequence: 6534d4d2bdcSBarry Smith + ksp - the `KSP` object whose `KSPConvergedReason` is to be viewed 6544d4d2bdcSBarry Smith - ctx - context used by the function, set with `KSPConvergedReasonViewSet()` 6554d4d2bdcSBarry Smith 6564d4d2bdcSBarry Smith Level: beginner 6574d4d2bdcSBarry Smith 6584d4d2bdcSBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPConvergedReasonView()`, `KSPConvergedReasonViewSet()`, `KSPConvergedReasonViewFromOptions()`, `KSPView()` 6594d4d2bdcSBarry Smith S*/ 6602a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPConvergedReasonViewFn(KSP ksp, PetscCtx ctx); 6614d4d2bdcSBarry Smith 662014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer); 66355849f57SBarry Smith PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer); 664fe2efc57SMark PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]); 66519a666eeSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer); 6664d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, KSPConvergedReasonViewFn *, void *, PetscCtxDestroyFn *); 6671b2b9847SBarry Smith PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP); 668c4421ceaSFande Kong PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP); 66994a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer); 6701b2b9847SBarry Smith 671edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonView()", ) static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v) 672d71ae5a4SJacob Faibussowitsch { 6739371c9d4SSatish Balay return KSPConvergedReasonView(ksp, v); 6749371c9d4SSatish Balay } 675edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) 676d71ae5a4SJacob Faibussowitsch { 6779371c9d4SSatish Balay return KSPConvergedReasonViewFromOptions(ksp); 6789371c9d4SSatish Balay } 67955849f57SBarry Smith 68055849f57SBarry Smith #define KSP_FILE_CLASSID 1211223 6811eb62cbbSBarry Smith 682dcc87044SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool); 6830e827bf4SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool); 684014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *); 685884807c7SVaclav Hapla PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *); 6864d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPLSQRMonitorResidual; 6874d4d2bdcSBarry Smith PETSC_EXTERN KSPMonitorRegisterFn KSPLSQRMonitorResidualDrawLG; 688798534f6SMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 689db9b2ab1SHong Zhang 690014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *); 691014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *); 69268ddcbeaSBarry Smith PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *); 6939f0612e4SBarry Smith PETSC_EXTERN PetscErrorCode PCMPIGetKSP(PC, KSP *); 69483ab6a24SBarry Smith 69528ce4d24SBarry Smith /*E 696a4d1885cSBarry Smith KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence 6978a4b9c5cSBarry Smith test routines. 6988a4b9c5cSBarry Smith 699a4d1885cSBarry Smith Values: 700a4d1885cSBarry Smith + `KSP_NORM_DEFAULT` - use the default for the current `KSPType` 701a4d1885cSBarry Smith . `KSP_NORM_NONE` - use no norm calculation 702a4d1885cSBarry Smith . `KSP_NORM_PRECONDITIONED` - use the preconditioned residual norm 703a4d1885cSBarry Smith . `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm 704a4d1885cSBarry Smith - `KSP_NORM_NATURAL` - use the natural norm (the norm induced by the linear operator) 705a4d1885cSBarry Smith 7068a4b9c5cSBarry Smith Level: advanced 7078a4b9c5cSBarry Smith 708a4d1885cSBarry Smith Note: 709a3f661c8SBarry Smith Each solver only supports a subset of these and some may support different ones 71016a05f60SBarry Smith depending on whether left or right preconditioning is used, see `KSPSetPCSide()` 711a3f661c8SBarry Smith 7121cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`, 71395bd0b28SBarry Smith `KSPSetConvergenceTest()`, `KSPSetPCSide()`, `KSP_NORM_DEFAULT`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 7148a4b9c5cSBarry Smith E*/ 7159371c9d4SSatish Balay typedef enum { 7169371c9d4SSatish Balay KSP_NORM_DEFAULT = -1, 7179371c9d4SSatish Balay KSP_NORM_NONE = 0, 7189371c9d4SSatish Balay KSP_NORM_PRECONDITIONED = 1, 7199371c9d4SSatish Balay KSP_NORM_UNPRECONDITIONED = 2, 7209371c9d4SSatish Balay KSP_NORM_NATURAL = 3 7219371c9d4SSatish Balay } KSPNormType; 7229e568555SJed Brown #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 723014dd563SJed Brown PETSC_EXTERN const char *const *const KSPNormTypes; 724a21b2a99SBarry Smith 7251f7e983dSSatish Balay /*MC 7269793b452SBarry Smith KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 7278c5b8ba0SBarry Smith possibly save some computation but means the convergence test cannot 7288c5b8ba0SBarry Smith be based on a norm of a residual etc. 7298c5b8ba0SBarry Smith 7308c5b8ba0SBarry Smith Level: advanced 7318c5b8ba0SBarry Smith 732a4d1885cSBarry Smith Note: 733a4d1885cSBarry Smith Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored 7348c5b8ba0SBarry Smith 7351cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 7368c5b8ba0SBarry Smith M*/ 7378c5b8ba0SBarry Smith 7381f7e983dSSatish Balay /*MC 7391957e957SBarry Smith KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 7408c5b8ba0SBarry Smith convergence test routine. 7418c5b8ba0SBarry Smith 7428c5b8ba0SBarry Smith Level: advanced 7438c5b8ba0SBarry Smith 7441cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 7458c5b8ba0SBarry Smith M*/ 7468c5b8ba0SBarry Smith 7471f7e983dSSatish Balay /*MC 748ce9499c7SBarry Smith KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 7498c5b8ba0SBarry Smith convergence test routine. 7508c5b8ba0SBarry Smith 7518c5b8ba0SBarry Smith Level: advanced 7528c5b8ba0SBarry Smith 7531cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 7548c5b8ba0SBarry Smith M*/ 7558c5b8ba0SBarry Smith 7561f7e983dSSatish Balay /*MC 757ce9499c7SBarry Smith KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 75887497f52SBarry Smith convergence test routine. This is only supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR` 7598c5b8ba0SBarry Smith 7608c5b8ba0SBarry Smith Level: advanced 7618c5b8ba0SBarry Smith 7621cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()` 7638c5b8ba0SBarry Smith M*/ 7648c5b8ba0SBarry Smith 765014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType); 766014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *); 767ce78bad3SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP, KSPNormType, PCSide, PetscInt); 768014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt); 769014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool); 7708a4b9c5cSBarry Smith 771edd03b47SJacob Faibussowitsch #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_NEG_CURVE", ) 772edd03b47SJacob Faibussowitsch #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_STEP_LENGTH", ) 773b5c1ed96SAlexander #define KSP_CONVERGED_RTOL_NORMAL_DEPRECATED KSP_CONVERGED_RTOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_RTOL_NORMAL_EQUATIONS", ) 774b5c1ed96SAlexander #define KSP_CONVERGED_ATOL_NORMAL_DEPRECATED KSP_CONVERGED_ATOL_NORMAL PETSC_DEPRECATED_ENUM(3, 24, 0, "KSP_CONVERGED_ATOL_NORMAL_EQUATIONS", ) 7758a4b9c5cSBarry Smith /*E 776a4d1885cSBarry Smith KSPConvergedReason - reason a Krylov method was determined to have converged or diverged 77728ce4d24SBarry Smith 778a4d1885cSBarry Smith Values: 77978daaec8SBarry Smith + `KSP_CONVERGED_RTOL_NORMAL_EQUATIONS` - requested decrease in the residual of the normal equations, for `KSPLSQR` 78078daaec8SBarry Smith . `KSP_CONVERGED_ATOL_NORMAL_EQUATIONS` - requested absolute value in the residual of the normal equations, for `KSPLSQR` 781a4d1885cSBarry Smith . `KSP_CONVERGED_RTOL` - requested decrease in the residual 782a4d1885cSBarry Smith . `KSP_CONVERGED_ATOL` - requested absolute value in the residual 783a4d1885cSBarry Smith . `KSP_CONVERGED_ITS` - requested number of iterations 7844a221d59SStefano Zampini . `KSP_CONVERGED_NEG_CURVE` - see note below 785a4d1885cSBarry Smith . `KSP_CONVERGED_STEP_LENGTH` - see note below 7866e7835fbSStefano Zampini . `KSP_CONVERGED_HAPPY_BREAKDOWN` - happy breakdown (meaning early convergence of the `KSPType` occurred). 787cc8c0b42SAlex Lindsay . `KSP_CONVERGED_USER` - the user has indicated convergence for an arbitrary reason 78878daaec8SBarry Smith . `KSP_DIVERGED_NULL` - breakdown when solving the Hessenberg system within `KSPGMRES` 789a4d1885cSBarry Smith . `KSP_DIVERGED_ITS` - requested number of iterations 79078daaec8SBarry Smith . `KSP_DIVERGED_DTOL` - large increase in the residual norm indicating the solution is diverging 791a4d1885cSBarry Smith . `KSP_DIVERGED_BREAKDOWN` - breakdown in the Krylov method 7928f47816eSPierre Jolivet . `KSP_DIVERGED_BREAKDOWN_BICG` - breakdown in the `KSPBCGS` Krylov method 793a4d1885cSBarry Smith . `KSP_DIVERGED_NONSYMMETRIC` - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry 79478daaec8SBarry Smith . `KSP_DIVERGED_INDEFINITE_PC` - the preconditioner was indefinite for a `KSPType` that requires it be definite, such as `KSPCG` 795a4d1885cSBarry Smith . `KSP_DIVERGED_NANORINF` - a not a number of infinity was detected in a vector during the computation 79678daaec8SBarry Smith . `KSP_DIVERGED_INDEFINITE_MAT` - the operator was indefinite for a `KSPType` that requires it be definite, such as `KSPCG` 797cc8c0b42SAlex Lindsay . `KSP_DIVERGED_PC_FAILED` - the action of the preconditioner failed for some reason 798cc8c0b42SAlex Lindsay - `KSP_DIVERGED_USER` - the user has indicated divergence for an arbitrary reason 799a4d1885cSBarry Smith 80016a05f60SBarry Smith Level: beginner 80116a05f60SBarry Smith 802a4d1885cSBarry Smith Note: 8036e7835fbSStefano Zampini The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by `KSPCG`, `KSPMINRES` and by 8046e7835fbSStefano Zampini the special `KSPNASH`, `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver. 80528ce4d24SBarry Smith 80695bd0b28SBarry Smith Developer Note: 80787497f52SBarry Smith The string versions of these are `KSPConvergedReasons`; if you change 8084d0a8057SBarry Smith any of the values here also change them that array of names. 80986c02ca4SBarry Smith 8101cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()` 81128ce4d24SBarry Smith E*/ 812d15094e1SBarry Smith typedef enum { /* converged */ 81378daaec8SBarry Smith KSP_CONVERGED_RTOL_NORMAL_DEPRECATED = 1, 81478daaec8SBarry Smith KSP_CONVERGED_RTOL_NORMAL_EQUATIONS = 1, 81578daaec8SBarry Smith KSP_CONVERGED_ATOL_NORMAL_DEPRECATED = 9, 81678daaec8SBarry Smith KSP_CONVERGED_ATOL_NORMAL_EQUATIONS = 9, 817d15094e1SBarry Smith KSP_CONVERGED_RTOL = 2, 818d15094e1SBarry Smith KSP_CONVERGED_ATOL = 3, 819b335793eSSatish Balay KSP_CONVERGED_ITS = 4, 8204a221d59SStefano Zampini KSP_CONVERGED_NEG_CURVE = 5, 8214a221d59SStefano Zampini KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED = 5, 8224a221d59SStefano Zampini KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6, 8234a221d59SStefano Zampini KSP_CONVERGED_STEP_LENGTH = 6, 8244a221d59SStefano Zampini KSP_CONVERGED_HAPPY_BREAKDOWN = 7, 825cc8c0b42SAlex Lindsay KSP_CONVERGED_USER = 8, 826d15094e1SBarry Smith /* diverged */ 827b3cc6726SBarry Smith KSP_DIVERGED_NULL = -2, 828d15094e1SBarry Smith KSP_DIVERGED_ITS = -3, 829d15094e1SBarry Smith KSP_DIVERGED_DTOL = -4, 830d15094e1SBarry Smith KSP_DIVERGED_BREAKDOWN = -5, 831b4ac9ba4SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG = -6, 832b4ac9ba4SBarry Smith KSP_DIVERGED_NONSYMMETRIC = -7, 833b4ac9ba4SBarry Smith KSP_DIVERGED_INDEFINITE_PC = -8, 8344d51c080SBarry Smith KSP_DIVERGED_NANORINF = -9, 8356aee1118SBarry Smith KSP_DIVERGED_INDEFINITE_MAT = -10, 836c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED = -11, 837aa4d2078SSatish Balay KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11, 838cc8c0b42SAlex Lindsay KSP_DIVERGED_USER = -12, 839d15094e1SBarry Smith 8409371c9d4SSatish Balay KSP_CONVERGED_ITERATING = 0 8419371c9d4SSatish Balay } KSPConvergedReason; 842014dd563SJed Brown PETSC_EXTERN const char *const *KSPConvergedReasons; 843d15094e1SBarry Smith 844c838673bSBarry Smith /*MC 845af27ebaaSBarry Smith KSP_CONVERGED_RTOL - $||r|| \le rtol*||b||$ or $rtol*||b - A*x_0||$ if `KSPConvergedDefaultSetUIRNorm()` was called 846c838673bSBarry Smith 847c838673bSBarry Smith Level: beginner 848c838673bSBarry Smith 84995bd0b28SBarry Smith Notes: 85087497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 851c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 852c838673bSBarry Smith 2-norm of the residual for right preconditioning 853c838673bSBarry Smith 85487497f52SBarry Smith See also `KSP_CONVERGED_ATOL` which may apply before this tolerance. 855f9fed41fSBarry Smith 8561cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 857c838673bSBarry Smith M*/ 858c838673bSBarry Smith 859c838673bSBarry Smith /*MC 860af27ebaaSBarry Smith KSP_CONVERGED_ATOL - $||r|| \le atol$ 861c838673bSBarry Smith 862c838673bSBarry Smith Level: beginner 863c838673bSBarry Smith 86495bd0b28SBarry Smith Notes: 86587497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 866c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 867c838673bSBarry Smith 2-norm of the residual for right preconditioning 868c838673bSBarry Smith 86987497f52SBarry Smith See also `KSP_CONVERGED_RTOL` which may apply before this tolerance. 870c838673bSBarry Smith 8711cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 872c838673bSBarry Smith M*/ 873c838673bSBarry Smith 874c838673bSBarry Smith /*MC 875af27ebaaSBarry Smith KSP_DIVERGED_DTOL - $||r|| \ge dtol*||b||$ 876c838673bSBarry Smith 877c838673bSBarry Smith Level: beginner 878c838673bSBarry Smith 87995bd0b28SBarry Smith Note: 88087497f52SBarry Smith See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 881c838673bSBarry Smith for left preconditioning it is the 2-norm of the preconditioned residual, and the 882c838673bSBarry Smith 2-norm of the residual for right preconditioning 883c838673bSBarry Smith 8840241b274SPierre Jolivet .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_CONVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 885c838673bSBarry Smith M*/ 886c838673bSBarry Smith 887c838673bSBarry Smith /*MC 888c838673bSBarry Smith KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 889c838673bSBarry Smith reached 890c838673bSBarry Smith 891c838673bSBarry Smith Level: beginner 892c838673bSBarry Smith 8931cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 894c838673bSBarry Smith M*/ 895c838673bSBarry Smith 896c838673bSBarry Smith /*MC 89787497f52SBarry Smith KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of 89887497f52SBarry Smith the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence 899af27ebaaSBarry Smith test routine is set in `KSP`. 900c838673bSBarry Smith 901c838673bSBarry Smith Level: beginner 902c838673bSBarry Smith 9031cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 904c838673bSBarry Smith M*/ 905c838673bSBarry Smith 906c838673bSBarry Smith /*MC 907c838673bSBarry Smith KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 9081de96524SPierre Jolivet method could not continue to enlarge the Krylov space. Could be due to a singular matrix or 90987497f52SBarry Smith preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block 91046091a0eSPierre Jolivet are collinear. 911c838673bSBarry Smith 912c838673bSBarry Smith Level: beginner 913c838673bSBarry Smith 9141cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 915c838673bSBarry Smith M*/ 916c838673bSBarry Smith 917c838673bSBarry Smith /*MC 91887497f52SBarry Smith KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the 919c838673bSBarry Smith method could not continue to enlarge the Krylov space. 920c838673bSBarry Smith 921c838673bSBarry Smith Level: beginner 922c838673bSBarry Smith 9231cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 924c838673bSBarry Smith M*/ 925c838673bSBarry Smith 926c838673bSBarry Smith /*MC 927c838673bSBarry Smith KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 92887497f52SBarry Smith symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry 929c838673bSBarry Smith 930c838673bSBarry Smith Level: beginner 931c838673bSBarry Smith 9321cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 933c838673bSBarry Smith M*/ 934c838673bSBarry Smith 935c838673bSBarry Smith /*MC 936c838673bSBarry Smith KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 93787497f52SBarry Smith positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to 9380b4b7b1cSBarry Smith be symmetric positive definite (SPD). 939c838673bSBarry Smith 940c838673bSBarry Smith Level: beginner 941c838673bSBarry Smith 94287497f52SBarry Smith Note: 943a4d1885cSBarry Smith This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force 94487497f52SBarry Smith the `PCICC` preconditioner to generate a positive definite preconditioner 945c838673bSBarry Smith 9461cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 947c838673bSBarry Smith M*/ 948c838673bSBarry Smith 949c838673bSBarry Smith /*MC 950c0decd05SBarry Smith KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a 9519fc87aa7SBarry Smith zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 95287497f52SBarry Smith such as `PCFIELDSPLIT`. 9539fc87aa7SBarry Smith 9549fc87aa7SBarry Smith Level: beginner 9559fc87aa7SBarry Smith 956a4d1885cSBarry Smith Note: 957a4d1885cSBarry Smith Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details. 9589fc87aa7SBarry Smith 9591cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 9609fc87aa7SBarry Smith M*/ 9619fc87aa7SBarry Smith 9629fc87aa7SBarry Smith /*MC 963a4d1885cSBarry Smith KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called 964a4d1885cSBarry Smith while `KSPSolve()` is still running. 965c838673bSBarry Smith 966c838673bSBarry Smith Level: beginner 967c838673bSBarry Smith 9681cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 969c838673bSBarry Smith M*/ 970c838673bSBarry Smith 9714d4d2bdcSBarry Smith /*S 9724d4d2bdcSBarry Smith KSPConvergenceTestFn - A prototype of a function used with `KSPSetConvergenceTest()` 9734d4d2bdcSBarry Smith 9744d4d2bdcSBarry Smith Calling Sequence: 9754d4d2bdcSBarry Smith + ksp - iterative solver obtained from `KSPCreate()` 9764d4d2bdcSBarry Smith . it - iteration number 9774d4d2bdcSBarry Smith . rnorm - (estimated) 2-norm of (preconditioned) residual 9784d4d2bdcSBarry Smith . reason - the reason why it has converged or diverged 9794d4d2bdcSBarry Smith - ctx - optional convergence context, as set by `KSPSetConvergenceTest()` 9804d4d2bdcSBarry Smith 9814d4d2bdcSBarry Smith Level: beginner 9824d4d2bdcSBarry Smith 9834d4d2bdcSBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPSetConvergenceTest()`, `KSPGetConvergenceTest()` 9844d4d2bdcSBarry Smith S*/ 9852a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPConvergenceTestFn(KSP ksp, PetscInt it, PetscReal rnorm, KSPConvergedReason *reason, PetscCtx ctx); 9864d4d2bdcSBarry Smith 9874d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, KSPConvergenceTestFn *, void *, PetscCtxDestroyFn *); 9882a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, KSPConvergenceTestFn **, PetscCtxRt, PetscCtxDestroyFn **); 9892a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, KSPConvergenceTestFn **, PetscCtxRt, PetscCtxDestroyFn **); 9902a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, PetscCtxRt); 9914d4d2bdcSBarry Smith PETSC_EXTERN KSPConvergenceTestFn KSPConvergedDefault; 9924d4d2bdcSBarry Smith PETSC_EXTERN KSPConvergenceTestFn KSPLSQRConvergedDefault; 9934d4d2bdcSBarry Smith PETSC_EXTERN PetscCtxDestroyFn KSPConvergedDefaultDestroy; 9948de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 9958de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 9968de4e6dcSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 99754b05d9cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool); 9980059c7bdSJed Brown PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 999014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *); 10004d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char *[]); 100194a4d4ebSMatthew G. Knepley PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 10022a28d964SStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetConvergedNegativeCurvature(KSP, PetscBool); 10032a28d964SStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetConvergedNegativeCurvature(KSP, PetscBool *); 1004abef13c0SSatish Balay 1005edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefault()", ) static inline void KSPDefaultConverged(void) 1006d71ae5a4SJacob Faibussowitsch { /* never called */ 10079371c9d4SSatish Balay } 10088ea1b3e6SJed Brown #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 1009edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultDestroy()", ) static inline void KSPDefaultConvergedDestroy(void) 1010d71ae5a4SJacob Faibussowitsch { /* never called */ 10119371c9d4SSatish Balay } 10128ea1b3e6SJed Brown #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 1013edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultCreate()", ) static inline void KSPDefaultConvergedCreate(void) 1014d71ae5a4SJacob Faibussowitsch { /* never called */ 10159371c9d4SSatish Balay } 10168ea1b3e6SJed Brown #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 1017edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUIRNorm()", ) static inline void KSPDefaultConvergedSetUIRNorm(void) 1018d71ae5a4SJacob Faibussowitsch { /* never called */ 10199371c9d4SSatish Balay } 10208ea1b3e6SJed Brown #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 1021edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUMIRNorm()", ) static inline void KSPDefaultConvergedSetUMIRNorm(void) 1022d71ae5a4SJacob Faibussowitsch { /* never called */ 10239371c9d4SSatish Balay } 10248ea1b3e6SJed Brown #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 1025edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedSkip()", ) static inline void KSPSkipConverged(void) 1026d71ae5a4SJacob Faibussowitsch { /* never called */ 10279371c9d4SSatish Balay } 10288ea1b3e6SJed Brown #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 10298ea1b3e6SJed Brown 10300bacdadaSStefano Zampini PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *); 1031edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPComputeOperator()", ) static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B) 1032d71ae5a4SJacob Faibussowitsch { 1033f22e26b7SPierre Jolivet return KSPComputeOperator(A, PETSC_NULLPTR, B); 10349371c9d4SSatish Balay } 1035d4fbbf0eSBarry Smith 103628ce4d24SBarry Smith /*E 1037a4d1885cSBarry Smith KSPCGType - Determines what type of `KSPCG` to use 103828ce4d24SBarry Smith 1039a4d1885cSBarry Smith Values: 1040a4d1885cSBarry Smith + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric 1041a4d1885cSBarry Smith - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian 1042a4d1885cSBarry Smith 104316a05f60SBarry Smith Level: beginner 104416a05f60SBarry Smith 10451cc06b55SBarry Smith .seealso: [](ch_ksp), `KSPCG`, `KSP`, `KSPCGSetType()` 104628ce4d24SBarry Smith E*/ 10479371c9d4SSatish Balay typedef enum { 10489371c9d4SSatish Balay KSP_CG_SYMMETRIC = 0, 10499371c9d4SSatish Balay KSP_CG_HERMITIAN = 1 10509371c9d4SSatish Balay } KSPCGType; 10516a6fc655SJed Brown PETSC_EXTERN const char *const KSPCGTypes[]; 105228ce4d24SBarry Smith 1053014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType); 1054014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool); 10558031f4adStmunson 1056dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal); 1057fb01098fSStefano Zampini PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal); 1058dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *); 1059dad1adc6STodd Munson PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *); 1060fcae7a14Stmunson 106105de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *); 106205de396fSBarry Smith PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *); 1063edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetMinEig()", ) static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x) 1064d71ae5a4SJacob Faibussowitsch { 10659371c9d4SSatish Balay return KSPGLTRGetMinEig(ksp, x); 10669371c9d4SSatish Balay } 1067edd03b47SJacob Faibussowitsch PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetLambda()", ) static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x) 1068d71ae5a4SJacob Faibussowitsch { 10699371c9d4SSatish Balay return KSPGLTRGetLambda(ksp, x); 10709371c9d4SSatish Balay } 10718031f4adStmunson 1072014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]); 1073ebead697SStefano Zampini PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]); 10741d6018f0SLisandro Dalcin 1075014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP); 1076014dd563SJed Brown PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP); 10773369ce9aSBarry Smith 1078014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *); 10792f2e5d10SKris Buschelman 10804d4d2bdcSBarry Smith /*S 10814d4d2bdcSBarry Smith PCShellPSolveFn - A function prototype for functions provided to `PCShellSetPreSolve()` and `PCShellSetPostSolve()` 10824d4d2bdcSBarry Smith 10834d4d2bdcSBarry Smith Calling Sequence: 10844d4d2bdcSBarry Smith + pc - the preconditioner `PC` context 10854d4d2bdcSBarry Smith . ksp - the `KSP` context 10864d4d2bdcSBarry Smith . xin - input vector 10874d4d2bdcSBarry Smith - xout - output vector 10884d4d2bdcSBarry Smith 10894d4d2bdcSBarry Smith Level: intermediate 10904d4d2bdcSBarry Smith 10914d4d2bdcSBarry Smith .seealso: [](ch_snes), `KSPPSolveFn`, `KSP`, `PCShellSetPreSolve()`, `PCShellSetPostSolve()` 10924d4d2bdcSBarry Smith S*/ 10934d4d2bdcSBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode PCShellPSolveFn(PC pc, KSP ksp, Vec xim, Vec xout); 10944d4d2bdcSBarry Smith 10954d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PCShellPSolveFn *); 10964d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PCShellPSolveFn *); 109703919abeSBarry Smith 1098ba36735cSStefano Zampini /*S 1099a4d1885cSBarry Smith KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods. 1100f8a50e2bSBarry Smith 1101a4d1885cSBarry Smith Level: intermediate 1102f8a50e2bSBarry Smith 11030b4b7b1cSBarry Smith Note: 11040b4b7b1cSBarry Smith These methods generate initial guesses based on a series of previous, related, linear solves. For example, 11050b4b7b1cSBarry Smith in implicit time-stepping with `TS`. 11060b4b7b1cSBarry Smith 11079168452cSPierre Jolivet .seealso: [](ch_ksp), `KSPCreate()`, `KSPGuessSetType()`, `KSPGuessType` 1108ba36735cSStefano Zampini S*/ 1109ba36735cSStefano Zampini typedef struct _p_KSPGuess *KSPGuess; 111016a05f60SBarry Smith 1111ba36735cSStefano Zampini /*J 111287497f52SBarry Smith KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods. 1113ba36735cSStefano Zampini 1114a4d1885cSBarry Smith Values: 1115a4d1885cSBarry Smith + `KSPGUESSFISCHER` - methodology developed by Paul Fischer 11160b4b7b1cSBarry Smith - `KSPGUESSPOD` - methodology based on proper orthogonal decomposition (POD) 1117a4d1885cSBarry Smith 111816a05f60SBarry Smith Level: intermediate 111916a05f60SBarry Smith 11201cc06b55SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPGuess` 1121ba36735cSStefano Zampini J*/ 1122ba36735cSStefano Zampini typedef const char *KSPGuessType; 1123ba36735cSStefano Zampini #define KSPGUESSFISCHER "fischer" 1124ba36735cSStefano Zampini #define KSPGUESSPOD "pod" 1125a4d1885cSBarry Smith 11261d36bdfdSBarry Smith PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess)); 1127ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess); 1128ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *); 1129ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer); 1130ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *); 1131ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *); 1132ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType); 1133ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *); 11348410009bSDavid Wells PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal); 1135ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 1136ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec); 1137ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec); 1138ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 1139ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt); 1140014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt); 1141ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool); 1142ba36735cSStefano Zampini PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *); 1143f8a50e2bSBarry Smith 1144470b340bSDmitry Karpeev /*E 11457addb90fSBarry Smith MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement matrix assembly routines 1146470b340bSDmitry Karpeev 1147470b340bSDmitry Karpeev Level: intermediate 1148470b340bSDmitry Karpeev 1149a4d1885cSBarry Smith .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`, 11500b4b7b1cSBarry Smith `MatCreateSchurComplementPmat()`, `MatCreateSchurComplement()` 1151470b340bSDmitry Karpeev E*/ 11529371c9d4SSatish Balay typedef enum { 11539371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_DIAG, 11549371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_LUMP, 11559371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, 11569371c9d4SSatish Balay MAT_SCHUR_COMPLEMENT_AINV_FULL 11579371c9d4SSatish Balay } MatSchurComplementAinvType; 1158470b340bSDmitry Karpeev PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 1159470b340bSDmitry Karpeev 1160014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *); 1161014dd563SJed Brown PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *); 1162d2768c1eSMatthew G Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP); 1163bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 1164aa6c7ce3SBarry Smith PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 1165bee83525SDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *); 1166470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType); 1167470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *); 11685bf40c1dSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *); 11695a1c817aSMatthew G. Knepley PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *); 1170470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *); 1171470b340bSDmitry Karpeev PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *); 11723f22127dSBarry Smith 117378e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *); 117478e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *); 1175bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatCreateLMVMDBFGS(MPI_Comm, PetscInt, PetscInt, Mat *); 1176bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatCreateLMVMDDFP(MPI_Comm, PetscInt, PetscInt, Mat *); 1177bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatCreateLMVMDQN(MPI_Comm, PetscInt, PetscInt, Mat *); 117878e4361aSAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *); 1179864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1180864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1181864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1182864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1183864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1184cd929ea3SAlp Dener 1185cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 1186b2d8c577SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *); 1187cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 1188cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 1189e0ed867bSAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat); 11900ad3a497SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat); 1191cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 1192cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 1193cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 1194cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 1195cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 11962d5e3849SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec); 1197cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 11981ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMGetLastUpdate(Mat, Vec *, Vec *); 1199cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *); 1200cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *); 1201cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *); 120292f76d53SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt); 1203bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatLMVMGetHistorySize(Mat, PetscInt *); 1204cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *); 1205cd929ea3SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *); 1206864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar); 1207e96a0efeSStefano Zampini 1208bbb72809SHansol Suh /*E 12091ca5963aSToby Isaac MatLMVMMultAlgorithm - The type of algorithm used for matrix-vector products and solves used internally by a `MatLMVM` matrix 12101ca5963aSToby Isaac 12111ca5963aSToby Isaac Values: 12121ca5963aSToby Isaac + `MAT_LMVM_MULT_RECURSIVE` - Use recursive formulas for products and solves 12131ca5963aSToby Isaac . `MAT_LMVM_MULT_DENSE` - Use dense formulas for products and solves when possible 12141ca5963aSToby Isaac - `MAT_LMVM_MULT_COMPACT_DENSE` - The same as `MATLMVM_MULT_DENSE`, but go further and ensure products and solves are computed in compact low-rank update form 12151ca5963aSToby Isaac 12161ca5963aSToby Isaac Level: advanced 12171ca5963aSToby Isaac 12181ca5963aSToby Isaac Options Database Keys: 12191ca5963aSToby Isaac . -mat_lmvm_mult_algorithm - the algorithm to use for multiplication (recursive, dense, compact_dense) 12201ca5963aSToby Isaac 12211ca5963aSToby Isaac .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSetMultAlgorithm()`, `MatLMVMGetMultAlgorithm()` 12221ca5963aSToby Isaac E*/ 12231ca5963aSToby Isaac typedef enum { 12241ca5963aSToby Isaac MAT_LMVM_MULT_RECURSIVE, 12251ca5963aSToby Isaac MAT_LMVM_MULT_DENSE, 12261ca5963aSToby Isaac MAT_LMVM_MULT_COMPACT_DENSE, 12271ca5963aSToby Isaac } MatLMVMMultAlgorithm; 12281ca5963aSToby Isaac 12291ca5963aSToby Isaac PETSC_EXTERN const char *const MatLMVMMultAlgorithms[]; 12301ca5963aSToby Isaac 12311ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSetMultAlgorithm(Mat, MatLMVMMultAlgorithm); 12321ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMGetMultAlgorithm(Mat, MatLMVMMultAlgorithm *); 12331ca5963aSToby Isaac 12341ca5963aSToby Isaac /*E 12351ca5963aSToby Isaac MatLMVMSymBroydenScaleType - Rescaling type for the initial Hessian of a symmetric Broyden matrix. 1236bbb72809SHansol Suh 1237bbb72809SHansol Suh Values: 1238b6982945SToby Isaac + `MAT_LMVM_SYMBROYDEN_SCALE_NONE` - no rescaling 1239b6982945SToby Isaac . `MAT_LMVM_SYMBROYDEN_SCALE_SCALAR` - scalar rescaling 1240b6982945SToby Isaac . `MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL` - diagonal rescaling 1241b6982945SToby Isaac . `MAT_LMVM_SYMBROYDEN_SCALE_USER` - same as `MAT_LMVM_SYMBROYDN_SCALE_NONE` 1242b6982945SToby Isaac - `MAT_LMVM_SYMBROYDEN_SCALE_DECIDE` - let PETSc decide rescaling 1243bbb72809SHansol Suh 1244bbb72809SHansol Suh Level: intermediate 1245bbb72809SHansol Suh 1246bbb72809SHansol Suh .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSymBroydenSetScaleType()` 1247bbb72809SHansol Suh E*/ 1248e96a0efeSStefano Zampini typedef enum { 1249e96a0efeSStefano Zampini MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0, 1250e96a0efeSStefano Zampini MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1, 1251e96a0efeSStefano Zampini MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2, 1252b6982945SToby Isaac MAT_LMVM_SYMBROYDEN_SCALE_USER = 3, 1253b6982945SToby Isaac MAT_LMVM_SYMBROYDEN_SCALE_DECIDE = 4 1254e96a0efeSStefano Zampini } MatLMVMSymBroydenScaleType; 1255e96a0efeSStefano Zampini PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[]; 1256e96a0efeSStefano Zampini 1257864588a7SAlp Dener PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType); 12581ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenGetPhi(Mat, PetscReal *); 12591ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetPhi(Mat, PetscReal); 12601ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenGetPsi(Mat, PetscReal *); 12611ca5963aSToby Isaac PETSC_EXTERN PetscErrorCode MatLMVMSymBadBroydenSetPsi(Mat, PetscReal); 1262cd929ea3SAlp Dener 1263bbb72809SHansol Suh /*E 12640b4b7b1cSBarry Smith MatLMVMDenseType - Memory storage strategy for dense variants of `MATLMVM`. 1265bbb72809SHansol Suh 1266bbb72809SHansol Suh Values: 1267bbb72809SHansol Suh + `MAT_LMVM_DENSE_REORDER` - reorders memory to minimize kernel launch 1268bbb72809SHansol Suh - `MAT_LMVM_DENSE_INPLACE` - computes inplace to minimize memory movement 1269bbb72809SHansol Suh 1270bbb72809SHansol Suh Level: intermediate 1271bbb72809SHansol Suh 1272bbb72809SHansol Suh .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMDenseSetType()` 1273bbb72809SHansol Suh E*/ 1274bbb72809SHansol Suh typedef enum { 1275bbb72809SHansol Suh MAT_LMVM_DENSE_REORDER, 1276bbb72809SHansol Suh MAT_LMVM_DENSE_INPLACE 1277bbb72809SHansol Suh } MatLMVMDenseType; 1278bbb72809SHansol Suh PETSC_EXTERN const char *const MatLMVMDenseTypes[]; 1279bbb72809SHansol Suh 1280bbb72809SHansol Suh PETSC_EXTERN PetscErrorCode MatLMVMDenseSetType(Mat, MatLMVMDenseType); 1281bbb72809SHansol Suh 1282014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM); 1283*bf0c7fc2SBarry Smith 1284*bf0c7fc2SBarry Smith /*E 1285*bf0c7fc2SBarry Smith KSPDMActive - Indicates if the `DM` attached to the `KSP` should be used to compute the operator, the right-hand side, or the initial guess 1286*bf0c7fc2SBarry Smith 1287*bf0c7fc2SBarry Smith Values: 1288*bf0c7fc2SBarry Smith + `KSP_DMACTIVE_OPERATOR` - compute the operator 1289*bf0c7fc2SBarry Smith . `KSP_DMACTIVE_RHS` - compute the right-hand side 1290*bf0c7fc2SBarry Smith . `KSP_DMACTIVE_INITIAL_GUESS` - compute the initial guess 1291*bf0c7fc2SBarry Smith - `KSP_DMACTIVE_ALL` - compute all of them 1292*bf0c7fc2SBarry Smith 1293*bf0c7fc2SBarry Smith Level: intermediate 1294*bf0c7fc2SBarry Smith 1295*bf0c7fc2SBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPSetDMActive()`, `KSPSetDM()` 1296*bf0c7fc2SBarry Smith E*/ 1297*bf0c7fc2SBarry Smith typedef enum { 1298*bf0c7fc2SBarry Smith KSP_DMACTIVE_OPERATOR = 1, 1299*bf0c7fc2SBarry Smith KSP_DMACTIVE_RHS = 2, 1300*bf0c7fc2SBarry Smith KSP_DMACTIVE_INITIAL_GUESS = 4, 1301*bf0c7fc2SBarry Smith KSP_DMACTIVE_ALL = 1 + 2 + 4 1302*bf0c7fc2SBarry Smith } KSPDMActive; 1303*bf0c7fc2SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, KSPDMActive, PetscBool); 1304*bf0c7fc2SBarry Smith 1305014dd563SJed Brown PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *); 13062a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, PetscCtx); 13072a8381b2SBarry Smith PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, PetscCtxRt); 13087b578ef6SBarry Smith 13097b578ef6SBarry Smith /*S 13108434afd1SBarry Smith KSPComputeRHSFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeRHS()` 13117b578ef6SBarry Smith 13127b578ef6SBarry Smith Calling Sequence: 13137b578ef6SBarry Smith + ksp - `ksp` context 13147b578ef6SBarry Smith . b - output vector 13157b578ef6SBarry Smith - ctx - [optional] user-defined function context 13167b578ef6SBarry Smith 13177b578ef6SBarry Smith Level: beginner 13187b578ef6SBarry Smith 13194d4d2bdcSBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeInitialGuessFn`, `KSPComputeOperatorsFn` 13207b578ef6SBarry Smith S*/ 13212a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeRHSFn(KSP ksp, Vec b, PetscCtx ctx); 13227b578ef6SBarry Smith 13238434afd1SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, KSPComputeRHSFn *, void *); 13247b578ef6SBarry Smith 13257b578ef6SBarry Smith /*S 13268434afd1SBarry Smith KSPComputeOperatorsFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeOperators()` 13277b578ef6SBarry Smith 13287b578ef6SBarry Smith Calling Sequence: 13297b578ef6SBarry Smith + ksp - `KSP` context 13307b578ef6SBarry Smith . A - the operator that defines the linear system 13317b578ef6SBarry Smith . P - an operator from which to build the preconditioner (often the same as `A`) 13327b578ef6SBarry Smith - ctx - [optional] user-defined function context 13337b578ef6SBarry Smith 13347b578ef6SBarry Smith Level: beginner 13357b578ef6SBarry Smith 13364d4d2bdcSBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeInitialGuessFn` 13377b578ef6SBarry Smith S*/ 13382a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeOperatorsFn(KSP ksp, Mat A, Mat P, PetscCtx ctx); 13397b578ef6SBarry Smith 13408434afd1SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, KSPComputeOperatorsFn, void *); 13417b578ef6SBarry Smith 13427b578ef6SBarry Smith /*S 13438434afd1SBarry Smith KSPComputeInitialGuessFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeInitialGuess()` 13447b578ef6SBarry Smith 13457b578ef6SBarry Smith Calling Sequence: 13467b578ef6SBarry Smith + ksp - `ksp` context 13477b578ef6SBarry Smith . x - output vector 13487b578ef6SBarry Smith - ctx - [optional] user-defined function context 13497b578ef6SBarry Smith 13507b578ef6SBarry Smith Level: beginner 13517b578ef6SBarry Smith 13524d4d2bdcSBarry Smith .seealso: [](ch_ksp), `KSP`, `KSPSetComputeInitialGuess()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeOperatorsFn` 13537b578ef6SBarry Smith S*/ 13542a8381b2SBarry Smith PETSC_EXTERN_TYPEDEF typedef PetscErrorCode KSPComputeInitialGuessFn(KSP ksp, Vec x, PetscCtx ctx); 13557b578ef6SBarry Smith 13568434afd1SBarry Smith PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, KSPComputeInitialGuessFn *, void *); 13578434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, KSPComputeOperatorsFn *, void *); 13588434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, KSPComputeOperatorsFn **, void *); 13598434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, KSPComputeRHSFn *, void *); 13608434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, KSPComputeRHSFn **, void *); 13618434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, KSPComputeInitialGuessFn *, void *); 13628434afd1SBarry Smith PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, KSPComputeInitialGuessFn **, void *); 13636c699258SBarry Smith 136402b02e71SToby Isaac PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec); 13651898fd5cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM, DM, PetscInt, const char *[], Vec[], ScatterMode); 13661898fd5cSMatthew G. Knepley PETSC_EXTERN PetscErrorCode DMSwarmProjectGradientFields(DM, DM, PetscInt, const char *[], Vec[], ScatterMode); 1367557cf195SMatthew G. Knepley 13682b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *); 13692b3cbbdaSStefano Zampini PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal); 13704bf303faSJacob Faibussowitsch 13714bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCBJKOKKOSSetKSP(PC, KSP); 13724bf303faSJacob Faibussowitsch PETSC_EXTERN PetscErrorCode PCBJKOKKOSGetKSP(PC, KSP *); 13735d83a8b1SBarry Smith 13745d83a8b1SBarry Smith PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM); 13754d4d2bdcSBarry Smith 13764d4d2bdcSBarry Smith #include <petscdstypes.h> 13774d4d2bdcSBarry Smith PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, PetscPointFn **, InsertMode, Vec); 1378