1 /* 2 Defines the interface functions for the Krylov subspace accelerators. 3 */ 4 #ifndef PETSCKSP_H 5 #define PETSCKSP_H 6 7 #include <petscpc.h> 8 9 /* SUBMANSEC = KSP */ 10 11 PETSC_EXTERN PetscErrorCode KSPInitializePackage(void); 12 PETSC_EXTERN PetscErrorCode KSPFinalizePackage(void); 13 14 /*S 15 KSP - Abstract PETSc object that manages all Krylov methods. This is the object that manages the 16 linear solves in PETSc (even those such as direct factorization-based solvers that do no use Krylov accelerators). 17 18 Level: beginner 19 20 Note: 21 When a direct solver is used, but no Krylov solver is used, the `KSP` object is still used but with a 22 `KSPType` of `KSPPREONLY` (or equivalently `KSPNONE`), meaning that only application of the preconditioner is used as the linear solver. 23 24 .seealso: [](doc_linsolve), [](ch_ksp), `KSPCreate()`, `KSPSetType()`, `KSPType`, `SNES`, `TS`, `PC`, `KSP`, `KSPDestroy()`, `KSPCG`, `KSPGMRES` 25 S*/ 26 typedef struct _p_KSP *KSP; 27 28 /*J 29 KSPType - String with the name of a PETSc Krylov method. 30 31 Level: beginner 32 33 .seealso: [](doc_linsolve), [](ch_ksp), `KSPSetType()`, `KSP`, `KSPRegister()`, `KSPCreate()`, `KSPSetFromOptions()` 34 J*/ 35 typedef const char *KSPType; 36 #define KSPRICHARDSON "richardson" 37 #define KSPCHEBYSHEV "chebyshev" 38 #define KSPCG "cg" 39 #define KSPGROPPCG "groppcg" 40 #define KSPPIPECG "pipecg" 41 #define KSPPIPECGRR "pipecgrr" 42 #define KSPPIPELCG "pipelcg" 43 #define KSPPIPEPRCG "pipeprcg" 44 #define KSPPIPECG2 "pipecg2" 45 #define KSPCGNE "cgne" 46 #define KSPNASH "nash" 47 #define KSPSTCG "stcg" 48 #define KSPGLTR "gltr" 49 #define KSPCGNASH PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPNASH", ) "nash" 50 #define KSPCGSTCG PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSTCG", ) "stcg" 51 #define KSPCGGLTR PETSC_DEPRECATED_MACRO(3, 11, 0, "KSPSGLTR", ) "gltr" 52 #define KSPFCG "fcg" 53 #define KSPPIPEFCG "pipefcg" 54 #define KSPGMRES "gmres" 55 #define KSPPIPEFGMRES "pipefgmres" 56 #define KSPFGMRES "fgmres" 57 #define KSPLGMRES "lgmres" 58 #define KSPDGMRES "dgmres" 59 #define KSPPGMRES "pgmres" 60 #define KSPTCQMR "tcqmr" 61 #define KSPBCGS "bcgs" 62 #define KSPIBCGS "ibcgs" 63 #define KSPQMRCGS "qmrcgs" 64 #define KSPFBCGS "fbcgs" 65 #define KSPFBCGSR "fbcgsr" 66 #define KSPBCGSL "bcgsl" 67 #define KSPPIPEBCGS "pipebcgs" 68 #define KSPCGS "cgs" 69 #define KSPTFQMR "tfqmr" 70 #define KSPCR "cr" 71 #define KSPPIPECR "pipecr" 72 #define KSPLSQR "lsqr" 73 #define KSPPREONLY "preonly" 74 #define KSPNONE "none" 75 #define KSPQCG "qcg" 76 #define KSPBICG "bicg" 77 #define KSPMINRES "minres" 78 #define KSPSYMMLQ "symmlq" 79 #define KSPLCD "lcd" 80 #define KSPPYTHON "python" 81 #define KSPGCR "gcr" 82 #define KSPPIPEGCR "pipegcr" 83 #define KSPTSIRM "tsirm" 84 #define KSPCGLS "cgls" 85 #define KSPFETIDP "fetidp" 86 #define KSPHPDDM "hpddm" 87 88 /* Logging support */ 89 PETSC_EXTERN PetscClassId KSP_CLASSID; 90 PETSC_EXTERN PetscClassId KSPGUESS_CLASSID; 91 PETSC_EXTERN PetscClassId DMKSP_CLASSID; 92 93 PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm, KSP *); 94 PETSC_EXTERN PetscErrorCode KSPSetType(KSP, KSPType); 95 PETSC_EXTERN PetscErrorCode KSPGetType(KSP, KSPType *); 96 PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 97 PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 98 PETSC_EXTERN PetscErrorCode KSPSolve(KSP, Vec, Vec); 99 PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP, Vec, Vec); 100 PETSC_EXTERN PetscErrorCode KSPSetUseExplicitTranspose(KSP, PetscBool); 101 PETSC_EXTERN PetscErrorCode KSPMatSolve(KSP, Mat, Mat); 102 PETSC_EXTERN PetscErrorCode KSPMatSolveTranspose(KSP, Mat, Mat); 103 PETSC_EXTERN PetscErrorCode KSPSetMatSolveBatchSize(KSP, PetscInt); 104 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPSetMatSolveBatchSize()", ) static inline PetscErrorCode KSPSetMatSolveBlockSize(KSP ksp, PetscInt n) 105 { 106 return KSPSetMatSolveBatchSize(ksp, n); 107 } 108 PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP, PetscInt *); 109 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPGetMatSolveBatchSize()", ) static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp, PetscInt *n) 110 { 111 return KSPGetMatSolveBatchSize(ksp, n); 112 } 113 PETSC_EXTERN PetscErrorCode KSPReset(KSP); 114 PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP); 115 PETSC_EXTERN PetscErrorCode KSPDestroy(KSP *); 116 PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP, PetscBool); 117 PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP, PetscBool *); 118 PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP, PetscBool); 119 PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP, PC, Vec); 120 121 PETSC_EXTERN PetscFunctionList KSPList; 122 PETSC_EXTERN PetscFunctionList KSPGuessList; 123 PETSC_EXTERN PetscFunctionList KSPMonitorList; 124 PETSC_EXTERN PetscFunctionList KSPMonitorCreateList; 125 PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList; 126 PETSC_EXTERN PetscErrorCode KSPRegister(const char[], PetscErrorCode (*)(KSP)); 127 PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **)); 128 129 PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP, PCSide); 130 PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP, PCSide *); 131 PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP, PetscReal, PetscReal, PetscReal, PetscInt); 132 PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP, PetscReal *, PetscReal *, PetscReal *, PetscInt *); 133 PETSC_EXTERN PetscErrorCode KSPSetMinimumIterations(KSP, PetscInt); 134 PETSC_EXTERN PetscErrorCode KSPGetMinimumIterations(KSP, PetscInt *); 135 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP, PetscBool); 136 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP, PetscBool *); 137 PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP, PetscBool); 138 PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP, PetscBool *); 139 PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP, PetscBool); 140 PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP, PetscBool); 141 PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP, PetscBool *); 142 PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP, PetscBool); 143 PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP, PetscBool *); 144 PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP, Vec *); 145 PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP, Vec *); 146 PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP, PetscReal *); 147 PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP, PetscInt *); 148 PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP, PetscInt *); 149 PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP, PetscInt, Vec **, PetscInt, Vec **); 150 PETSC_DEPRECATED_FUNCTION(3, 6, 0, "KSPCreateVecs()", ) static inline PetscErrorCode KSPGetVecs(KSP ksp, PetscInt n, Vec **x, PetscInt m, Vec **y) 151 { 152 return KSPCreateVecs(ksp, n, x, m, y); 153 } 154 155 PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *); 156 PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP, PetscErrorCode (*)(KSP, Vec, Vec, void *), void *); 157 158 PETSC_EXTERN PetscErrorCode KSPSetPC(KSP, PC); 159 PETSC_EXTERN PetscErrorCode KSPGetPC(KSP, PC *); 160 PETSC_EXTERN PetscErrorCode KSPSetNestLevel(KSP, PetscInt); 161 PETSC_EXTERN PetscErrorCode KSPGetNestLevel(KSP, PetscInt *); 162 163 PETSC_EXTERN PetscErrorCode KSPMonitor(KSP, PetscInt, PetscReal); 164 PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void **)); 165 PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 166 PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP, void *); 167 PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP, const PetscReal *[], PetscInt *); 168 PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP, PetscReal[], PetscInt, PetscBool); 169 PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP, const PetscReal *[], PetscInt *); 170 PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP, PetscReal[], PetscInt, PetscBool); 171 172 PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP, Vec, Vec *); 173 PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP, Vec, Vec, Vec *); 174 PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 175 PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP, PetscInt); 176 177 PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC, KSP *); 178 PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC, KSP); 179 PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 180 PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 181 PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC, PetscInt *, PetscInt *, KSP *[]); 182 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC, PetscInt *, KSP *[]); 183 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC, PetscInt *, KSP *[]); 184 PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC, PetscInt, KSP *); 185 PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC, PetscInt, KSP *); 186 PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC, PetscInt, KSP *); 187 PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC, KSP *); 188 PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC, KSP *); 189 PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC, KSP *); 190 /* 191 PCMGCoarseList contains the list of coarse space constructor currently registered 192 These are added with PCMGRegisterCoarseSpaceConstructor() 193 */ 194 PETSC_EXTERN PetscFunctionList PCMGCoarseList; 195 PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *)); 196 PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat *)); 197 198 PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP, Vec, Vec *); 199 PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP, Vec, Vec, Vec *); 200 201 /*E 202 KSPChebyshevKind - Which kind of Chebyshev polynomial to use 203 204 Values: 205 + `KSP_CHEBYSHEV_FIRST` - "classic" first-kind Chebyshev polynomial 206 . `KSP_CHEBYSHEV_FOURTH` - fourth-kind Chebyshev polynomial 207 - `KSP_CHEBYSHEV_OPT_FOURTH` - optimized fourth-kind Chebyshev polynomial 208 209 Level: intermediate 210 211 .seealso: [](ch_ksp), `KSPCHEBYSHEV`, `KSPChebyshevSetKind()`, `KSPChebyshevGetKind()` 212 E*/ 213 typedef enum { 214 KSP_CHEBYSHEV_FIRST, 215 KSP_CHEBYSHEV_FOURTH, 216 KSP_CHEBYSHEV_OPT_FOURTH 217 } KSPChebyshevKind; 218 219 PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP, PetscReal); 220 PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP, PetscBool); 221 PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP, PetscReal, PetscReal); 222 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP, PetscReal, PetscReal, PetscReal, PetscReal); 223 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP, PetscBool); 224 PETSC_EXTERN PetscErrorCode KSPChebyshevSetKind(KSP, KSPChebyshevKind); 225 PETSC_EXTERN PetscErrorCode KSPChebyshevGetKind(KSP, KSPChebyshevKind *); 226 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP, KSP *); 227 PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP, PetscReal *, PetscReal *); 228 PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP, PetscInt, PetscReal[], PetscReal[], PetscInt *); 229 PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP, PetscInt, PetscReal[], PetscReal[]); 230 PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP, PetscBool, PetscBool, PetscInt *, Vec[], PetscReal[], PetscReal[]); 231 232 /*E 233 234 KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 235 236 Values: 237 + `KSP_FCD_TRUNC_TYPE_STANDARD` - uses all (up to mmax) stored directions 238 - `KSP_FCD_TRUNC_TYPE_NOTAY` - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 239 240 Level: intermediate 241 242 .seealso: [](ch_ksp), `KSP`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()` 243 E*/ 244 typedef enum { 245 KSP_FCD_TRUNC_TYPE_STANDARD, 246 KSP_FCD_TRUNC_TYPE_NOTAY 247 } KSPFCDTruncationType; 248 PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 249 250 PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP, PetscInt); 251 PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP, PetscInt *); 252 PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP, PetscInt); 253 PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP, PetscInt *); 254 PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP, KSPFCDTruncationType); 255 PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP, KSPFCDTruncationType *); 256 257 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP, PetscInt); 258 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP, PetscInt *); 259 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP, PetscInt); 260 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP, PetscInt *); 261 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP, KSPFCDTruncationType); 262 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP, KSPFCDTruncationType *); 263 264 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP, PetscInt); 265 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP, PetscInt *); 266 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP, PetscInt); 267 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP, PetscInt *); 268 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP, KSPFCDTruncationType); 269 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP, KSPFCDTruncationType *); 270 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP, PetscBool); 271 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP, PetscBool *); 272 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 273 274 PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 275 PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *); 276 PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal); 277 PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal); 278 279 PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 280 PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt)); 281 PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt)); 282 PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt); 283 PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt); 284 285 PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt); 286 PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 287 288 PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar); 289 290 PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt); 291 PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *); 292 PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 293 294 PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal); 295 PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *); 296 PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool); 297 298 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *); 299 PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC); 300 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *); 301 PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat); 302 303 PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat); 304 PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *); 305 #if PetscDefined(HAVE_HPDDM) 306 PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMSetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) 307 { 308 return KSPHPDDMSetDeflationMat(ksp, U); 309 } 310 PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMGetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) 311 { 312 return KSPHPDDMGetDeflationMat(ksp, U); 313 } 314 #endif 315 PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPMatSolve()", ) static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) 316 { 317 return KSPMatSolve(ksp, B, X); 318 } 319 /*E 320 KSPHPDDMType - Type of Krylov method used by `KSPHPDDM` 321 322 Values: 323 + `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method 324 . `KSP_HPDDM_TYPE_BGMRES` - block GMRES 325 . `KSP_HPDDM_TYPE_CG` - Conjugate Gradient 326 . `KSP_HPDDM_TYPE_BCG` - block CG 327 . `KSP_HPDDM_TYPE_GCRODR` - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting 328 . `KSP_HPDDM_TYPE_BGCRODR` - block GCRODR 329 . `KSP_HPDDM_TYPE_BFBCG` - breakdown-free BCG 330 - `KSP_HPDDM_TYPE_PREONLY` - apply the preconditioner only 331 332 Level: intermediate 333 334 .seealso: [](ch_ksp), `KSPHPDDM`, `KSPHPDDMSetType()` 335 E*/ 336 typedef enum { 337 KSP_HPDDM_TYPE_GMRES = 0, 338 KSP_HPDDM_TYPE_BGMRES = 1, 339 KSP_HPDDM_TYPE_CG = 2, 340 KSP_HPDDM_TYPE_BCG = 3, 341 KSP_HPDDM_TYPE_GCRODR = 4, 342 KSP_HPDDM_TYPE_BGCRODR = 5, 343 KSP_HPDDM_TYPE_BFBCG = 6, 344 KSP_HPDDM_TYPE_PREONLY = 7 345 } KSPHPDDMType; 346 PETSC_EXTERN const char *const KSPHPDDMTypes[]; 347 348 /*E 349 KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM` 350 351 Values: 352 + `KSP_HPDDM_PRECISION_HALF` - default when PETSc is configured `--with-precision=__fp16` 353 . `KSP_HPDDM_PRECISION_SINGLE` - default when PETSc is configured `--with-precision=single` 354 . `KSP_HPDDM_PRECISION_DOUBLE` - default when PETSc is configured `--with-precision=double` 355 - `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128` 356 357 Level: intermediate 358 359 .seealso: [](ch_ksp), `KSP`, `KSPHPDDM` 360 E*/ 361 typedef enum { 362 KSP_HPDDM_PRECISION_HALF = 0, 363 KSP_HPDDM_PRECISION_SINGLE = 1, 364 KSP_HPDDM_PRECISION_DOUBLE = 2, 365 KSP_HPDDM_PRECISION_QUADRUPLE = 3 366 } KSPHPDDMPrecision; 367 PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType); 368 PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *); 369 370 /*E 371 KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed in the GMRES solvers 372 373 Values: 374 + `KSP_GMRES_CGS_REFINE_NEVER` - one step of classical Gram-Schmidt 375 . `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria 376 - `KSP_GMRES_CGS_REFINE_ALWAYS` - always perform two steps 377 378 Level: advanced 379 380 .seealso: [](ch_ksp), `KSP`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 381 `KSPGMRESGetOrthogonalization()`, 382 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()` 383 E*/ 384 typedef enum { 385 KSP_GMRES_CGS_REFINE_NEVER, 386 KSP_GMRES_CGS_REFINE_IFNEEDED, 387 KSP_GMRES_CGS_REFINE_ALWAYS 388 } KSPGMRESCGSRefinementType; 389 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 390 391 /*MC 392 KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 393 394 Level: advanced 395 396 Note: 397 Possibly unstable, but the fastest to compute 398 399 .seealso: [](ch_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 400 `KSP`, `KSPGMRESGetOrthogonalization()`, 401 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 402 `KSPGMRESModifiedGramSchmidtOrthogonalization()` 403 M*/ 404 405 /*MC 406 KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 407 iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 408 poor orthogonality. 409 410 Level: advanced 411 412 Note: 413 This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to 414 estimate the orthogonality but is more stable. 415 416 .seealso: [](ch_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 417 `KSP`, `KSPGMRESGetOrthogonalization()`, 418 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 419 `KSPGMRESModifiedGramSchmidtOrthogonalization()` 420 M*/ 421 422 /*MC 423 KSP_GMRES_CGS_REFINE_ALWAYS - Do two steps of the classical (unmodified) Gram-Schmidt process. 424 425 Level: advanced 426 427 Notes: 428 This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice 429 but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`. 430 431 You should only use this if you absolutely know that the iterative refinement is needed. 432 433 .seealso: [](ch_ksp), `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 434 `KSP`, `KSPGMRESGetOrthogonalization()`, 435 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 436 `KSPGMRESModifiedGramSchmidtOrthogonalization()` 437 M*/ 438 439 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType); 440 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *); 441 442 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *); 443 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *); 444 PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 445 446 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal); 447 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *); 448 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *); 449 450 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal); 451 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool); 452 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt); 453 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool); 454 455 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 456 PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP); 457 458 PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *); 459 PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm, const char[], const char[], const char[], PetscInt, const char *[], int, int, int, int, PetscDrawLG *); 460 PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 461 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 462 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 463 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 464 PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 465 PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 466 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 467 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 468 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 469 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 470 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 471 PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 472 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 473 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 474 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 475 PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 476 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 477 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 478 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 479 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 480 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 481 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorResidual()", ) static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 482 { 483 return KSPMonitorResidual(ksp, n, rnorm, vf); 484 } 485 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidual()", ) static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 486 { 487 return KSPMonitorTrueResidual(ksp, n, rnorm, vf); 488 } 489 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidualMax()", ) static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 490 { 491 return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf); 492 } 493 494 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *); 495 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *); 496 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **); 497 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *); 498 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal); 499 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *); 500 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **); 501 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **); 502 503 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec); 504 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec); 505 506 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat); 507 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *); 508 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *); 509 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]); 510 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]); 511 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]); 512 513 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool); 514 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *); 515 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool); 516 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *); 517 518 PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer); 519 PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer); 520 PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]); 521 PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer); 522 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **)); 523 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP); 524 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP); 525 PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer); 526 527 PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonView()", ) static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v) 528 { 529 return KSPConvergedReasonView(ksp, v); 530 } 531 PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) 532 { 533 return KSPConvergedReasonViewFromOptions(ksp); 534 } 535 536 #define KSP_FILE_CLASSID 1211223 537 538 PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool); 539 PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool); 540 PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *); 541 PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *); 542 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 543 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 544 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 545 546 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *); 547 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *); 548 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *); 549 550 /*E 551 KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence 552 test routines. 553 554 Values: 555 + `KSP_NORM_DEFAULT` - use the default for the current `KSPType` 556 . `KSP_NORM_NONE` - use no norm calculation 557 . `KSP_NORM_PRECONDITIONED` - use the preconditioned residual norm 558 . `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm 559 - `KSP_NORM_NATURAL` - use the natural norm (the norm induced by the linear operator) 560 561 Level: advanced 562 563 Note: 564 Each solver only supports a subset of these and some may support different ones 565 depending on whether left or right preconditioning is used, see `KSPSetPCSide()` 566 567 Developer Note: 568 This must match the values in petsc/finclude/petscksp.h 569 570 .seealso: [](ch_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`, 571 `KSPSetConvergenceTest()`, `KSPSetPCSide()` 572 E*/ 573 typedef enum { 574 KSP_NORM_DEFAULT = -1, 575 KSP_NORM_NONE = 0, 576 KSP_NORM_PRECONDITIONED = 1, 577 KSP_NORM_UNPRECONDITIONED = 2, 578 KSP_NORM_NATURAL = 3 579 } KSPNormType; 580 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 581 PETSC_EXTERN const char *const *const KSPNormTypes; 582 583 /*MC 584 KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 585 possibly save some computation but means the convergence test cannot 586 be based on a norm of a residual etc. 587 588 Level: advanced 589 590 Note: 591 Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored 592 593 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 594 M*/ 595 596 /*MC 597 KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 598 convergence test routine. 599 600 Level: advanced 601 602 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 603 M*/ 604 605 /*MC 606 KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 607 convergence test routine. 608 609 Level: advanced 610 611 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 612 M*/ 613 614 /*MC 615 KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 616 convergence test routine. This is only supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR` 617 618 Level: advanced 619 620 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()` 621 M*/ 622 623 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType); 624 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *); 625 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt); 626 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt); 627 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool); 628 629 #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_NEG_CURVE", ) 630 #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_STEP_LENGTH", ) 631 #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM(3, 11, 0, "KSP_DIVERGED_PC_FAILED", ) 632 /*E 633 KSPConvergedReason - reason a Krylov method was determined to have converged or diverged 634 635 Values: 636 + `KSP_CONVERGED_RTOL_NORMAL` - requested decrease in the residual for the normal equations 637 . `KSP_CONVERGED_ATOL_NORMAL` - requested absolute value in the residual for the normal equations 638 . `KSP_CONVERGED_RTOL` - requested decrease in the residual 639 . `KSP_CONVERGED_ATOL` - requested absolute value in the residual 640 . `KSP_CONVERGED_ITS` - requested number of iterations 641 . `KSP_CONVERGED_NEG_CURVE` - see note below 642 . `KSP_CONVERGED_STEP_LENGTH` - see note below 643 . `KSP_CONVERGED_HAPPY_BREAKDOWN` - happy breakdown (meaning early convergence of the `KSPType` occurred. 644 . `KSP_DIVERGED_NULL` - breakdown when solving the Hessenberg system within GMRES 645 . `KSP_DIVERGED_ITS` - requested number of iterations 646 . `KSP_DIVERGED_DTOL` - large increase in the residual norm 647 . `KSP_DIVERGED_BREAKDOWN` - breakdown in the Krylov method 648 . `KSP_DIVERGED_BREAKDOWN_BICG` - breakdown in the `KSPBGCS` Krylov method 649 . `KSP_DIVERGED_NONSYMMETRIC` - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry 650 . `KSP_DIVERGED_INDEFINITE_PC` - the preconditioner was indefinite for a `KSPType` that requires it be definite 651 . `KSP_DIVERGED_NANORINF` - a not a number of infinity was detected in a vector during the computation 652 . `KSP_DIVERGED_INDEFINITE_MAT` - the operator was indefinite for a `KSPType` that requires it be definite 653 - `KSP_DIVERGED_PC_FAILED` - the action of the preconditioner failed for some reason 654 655 Level: beginner 656 657 Note: 658 The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by the special `KSPNASH`, 659 `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver. 660 661 Developer Notes: 662 This must match the values in petsc/finclude/petscksp.h 663 664 The string versions of these are `KSPConvergedReasons`; if you change 665 any of the values here also change them that array of names. 666 667 .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()` 668 E*/ 669 typedef enum { /* converged */ 670 KSP_CONVERGED_RTOL_NORMAL = 1, 671 KSP_CONVERGED_ATOL_NORMAL = 9, 672 KSP_CONVERGED_RTOL = 2, 673 KSP_CONVERGED_ATOL = 3, 674 KSP_CONVERGED_ITS = 4, 675 KSP_CONVERGED_NEG_CURVE = 5, 676 KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED = 5, 677 KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6, 678 KSP_CONVERGED_STEP_LENGTH = 6, 679 KSP_CONVERGED_HAPPY_BREAKDOWN = 7, 680 /* diverged */ 681 KSP_DIVERGED_NULL = -2, 682 KSP_DIVERGED_ITS = -3, 683 KSP_DIVERGED_DTOL = -4, 684 KSP_DIVERGED_BREAKDOWN = -5, 685 KSP_DIVERGED_BREAKDOWN_BICG = -6, 686 KSP_DIVERGED_NONSYMMETRIC = -7, 687 KSP_DIVERGED_INDEFINITE_PC = -8, 688 KSP_DIVERGED_NANORINF = -9, 689 KSP_DIVERGED_INDEFINITE_MAT = -10, 690 KSP_DIVERGED_PC_FAILED = -11, 691 KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11, 692 693 KSP_CONVERGED_ITERATING = 0 694 } KSPConvergedReason; 695 PETSC_EXTERN const char *const *KSPConvergedReasons; 696 697 /*MC 698 KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if `KSPConvergedDefaultSetUIRNorm()` was called 699 700 Level: beginner 701 702 See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 703 for left preconditioning it is the 2-norm of the preconditioned residual, and the 704 2-norm of the residual for right preconditioning 705 706 See also `KSP_CONVERGED_ATOL` which may apply before this tolerance. 707 708 .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 709 M*/ 710 711 /*MC 712 KSP_CONVERGED_ATOL - norm(r) <= atol 713 714 Level: beginner 715 716 See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 717 for left preconditioning it is the 2-norm of the preconditioned residual, and the 718 2-norm of the residual for right preconditioning 719 720 See also `KSP_CONVERGED_RTOL` which may apply before this tolerance. 721 722 .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 723 M*/ 724 725 /*MC 726 KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 727 728 Level: beginner 729 730 See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 731 for left preconditioning it is the 2-norm of the preconditioned residual, and the 732 2-norm of the residual for right preconditioning 733 734 Level: beginner 735 736 .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 737 M*/ 738 739 /*MC 740 KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 741 reached 742 743 Level: beginner 744 745 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 746 M*/ 747 748 /*MC 749 KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of 750 the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence 751 test routine is set in KSP. 752 753 Level: beginner 754 755 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 756 M*/ 757 758 /*MC 759 KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 760 method could not continue to enlarge the Krylov space. Could be due to a singular matrix or 761 preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block 762 are colinear. 763 764 Level: beginner 765 766 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 767 M*/ 768 769 /*MC 770 KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the 771 method could not continue to enlarge the Krylov space. 772 773 Level: beginner 774 775 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 776 M*/ 777 778 /*MC 779 KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 780 symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry 781 782 Level: beginner 783 784 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 785 M*/ 786 787 /*MC 788 KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 789 positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to 790 be positive definite 791 792 Level: beginner 793 794 Note: 795 This can happen with the `PCICC` preconditioner, use the options database option `-pc_factor_shift_positive_definite` to force 796 the `PCICC` preconditioner to generate a positive definite preconditioner 797 798 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 799 M*/ 800 801 /*MC 802 KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a 803 zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 804 such as `PCFIELDSPLIT`. 805 806 Level: beginner 807 808 Note: 809 Run with `-ksp_error_if_not_converged` to stop the program when the error is detected and print an error message with details. 810 811 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 812 M*/ 813 814 /*MC 815 KSP_CONVERGED_ITERATING - This flag is returned if `KSPGetConvergedReason()` is called 816 while `KSPSolve()` is still running. 817 818 Level: beginner 819 820 .seealso: [](ch_ksp), `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 821 M*/ 822 823 PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *)); 824 PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *)); 825 PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *)); 826 PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *); 827 PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 828 PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 829 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 830 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 831 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 832 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 833 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool); 834 PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 835 PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *); 836 PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **); 837 PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 838 PETSC_EXTERN PetscErrorCode KSPSetConvergedNegativeCurvature(KSP, PetscBool); 839 PETSC_EXTERN PetscErrorCode KSPGetConvergedNegativeCurvature(KSP, PetscBool *); 840 841 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefault()", ) static inline void KSPDefaultConverged(void) 842 { /* never called */ 843 } 844 #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 845 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultDestroy()", ) static inline void KSPDefaultConvergedDestroy(void) 846 { /* never called */ 847 } 848 #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 849 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultCreate()", ) static inline void KSPDefaultConvergedCreate(void) 850 { /* never called */ 851 } 852 #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 853 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUIRNorm()", ) static inline void KSPDefaultConvergedSetUIRNorm(void) 854 { /* never called */ 855 } 856 #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 857 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedDefaultSetUMIRNorm()", ) static inline void KSPDefaultConvergedSetUMIRNorm(void) 858 { /* never called */ 859 } 860 #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 861 PETSC_DEPRECATED_FUNCTION(3, 5, 0, "KSPConvergedSkip()", ) static inline void KSPSkipConverged(void) 862 { /* never called */ 863 } 864 #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 865 866 PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *); 867 PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPComputeOperator()", ) static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B) 868 { 869 return KSPComputeOperator(A, PETSC_NULLPTR, B); 870 } 871 872 /*E 873 KSPCGType - Determines what type of `KSPCG` to use 874 875 Values: 876 + `KSP_CG_SYMMETRIC` - the matrix is complex symmetric 877 - `KSP_CG_HERMITIAN` - the matrix is complex Hermitian 878 879 Level: beginner 880 881 .seealso: [](ch_ksp), `KSPCG`, `KSP`, `KSPCGSetType()` 882 E*/ 883 typedef enum { 884 KSP_CG_SYMMETRIC = 0, 885 KSP_CG_HERMITIAN = 1 886 } KSPCGType; 887 PETSC_EXTERN const char *const KSPCGTypes[]; 888 889 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType); 890 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool); 891 892 PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal); 893 PETSC_EXTERN PetscErrorCode KSPCGSetObjectiveTarget(KSP, PetscReal); 894 PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *); 895 PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *); 896 897 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *); 898 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *); 899 PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetMinEig()", ) static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x) 900 { 901 return KSPGLTRGetMinEig(ksp, x); 902 } 903 PETSC_DEPRECATED_FUNCTION(3, 12, 0, "KSPGLTRGetLambda()", ) static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x) 904 { 905 return KSPGLTRGetLambda(ksp, x); 906 } 907 908 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]); 909 PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]); 910 911 PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP)); 912 PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP); 913 PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP); 914 915 #include <petscdrawtypes.h> 916 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *); 917 918 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec)); 919 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec)); 920 921 /*S 922 KSPGuess - Abstract PETSc object that manages all initial guess generation methods for Krylov methods. 923 924 Level: intermediate 925 926 .seealso: [](ch_ksp), `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType` 927 S*/ 928 typedef struct _p_KSPGuess *KSPGuess; 929 930 /*J 931 KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods. 932 933 Values: 934 + `KSPGUESSFISCHER` - methodology developed by Paul Fischer 935 - `KSPGUESSPOD` - methodology based on proper orthogonal decomposition 936 937 Level: intermediate 938 939 .seealso: [](ch_ksp), `KSP`, `KSPGuess` 940 J*/ 941 typedef const char *KSPGuessType; 942 #define KSPGUESSFISCHER "fischer" 943 #define KSPGUESSPOD "pod" 944 945 PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess)); 946 PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess); 947 PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *); 948 PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer); 949 PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *); 950 PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *); 951 PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType); 952 PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *); 953 PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal); 954 PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 955 PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec); 956 PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec); 957 PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 958 PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt); 959 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt); 960 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool); 961 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *); 962 963 /*E 964 MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 965 966 Level: intermediate 967 968 .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`, 969 `MatCreateSchurComplementPmat()` 970 E*/ 971 typedef enum { 972 MAT_SCHUR_COMPLEMENT_AINV_DIAG, 973 MAT_SCHUR_COMPLEMENT_AINV_LUMP, 974 MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, 975 MAT_SCHUR_COMPLEMENT_AINV_FULL 976 } MatSchurComplementAinvType; 977 PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 978 979 typedef enum { 980 MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0, 981 MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1, 982 MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2, 983 MAT_LMVM_SYMBROYDEN_SCALE_USER = 3 984 } MatLMVMSymBroydenScaleType; 985 PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[]; 986 987 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *); 988 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *); 989 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP); 990 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 991 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 992 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *); 993 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType); 994 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *); 995 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *); 996 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *); 997 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *); 998 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *); 999 1000 PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *); 1001 PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *); 1002 PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *); 1003 PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1004 PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1005 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1006 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1007 PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 1008 1009 PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 1010 PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *); 1011 PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 1012 PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 1013 PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat); 1014 PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat); 1015 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 1016 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 1017 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 1018 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 1019 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 1020 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec); 1021 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 1022 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *); 1023 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *); 1024 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *); 1025 PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt); 1026 PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *); 1027 PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *); 1028 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar); 1029 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType); 1030 1031 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM); 1032 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool); 1033 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *); 1034 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *); 1035 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *); 1036 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, PetscErrorCode (*func)(KSP, Vec, void *), void *); 1037 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *); 1038 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, PetscErrorCode (*)(KSP, Vec, void *), void *); 1039 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *); 1040 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, PetscErrorCode (**)(KSP, Mat, Mat, void *), void *); 1041 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, PetscErrorCode (*)(KSP, Vec, void *), void *); 1042 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, PetscErrorCode (**)(KSP, Vec, void *), void *); 1043 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, PetscErrorCode (*)(KSP, Vec, void *), void *); 1044 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, PetscErrorCode (**)(KSP, Vec, void *), void *); 1045 1046 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec); 1047 PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, void (**)(PetscInt, PetscInt, PetscInt, const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], PetscReal, const PetscReal[], PetscInt, const PetscScalar[], PetscScalar[]), InsertMode, Vec); 1048 1049 PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *); 1050 PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal); 1051 1052 PETSC_EXTERN PetscErrorCode PCBJKOKKOSSetKSP(PC, KSP); 1053 PETSC_EXTERN PetscErrorCode PCBJKOKKOSGetKSP(PC, KSP *); 1054 #endif 1055