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