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