1 /* 2 Defines the interface functions for the Krylov subspace accelerators. 3 */ 4 #pragma once 5 6 #include <petscpc.h> 7 8 /* SUBMANSEC = KSP */ 9 10 PETSC_EXTERN PetscErrorCode KSPInitializePackage(void); 11 PETSC_EXTERN PetscErrorCode KSPFinalizePackage(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 with `KSPCHEBYSHEV` 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 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 272 273 PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 274 PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt *); 275 PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP, PetscReal); 276 PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP, PetscReal); 277 278 PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 279 PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP, PetscErrorCode (*)(KSP, PetscInt)); 280 PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP, PetscErrorCode (**)(KSP, PetscInt)); 281 PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP, PetscInt); 282 PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP, PetscInt); 283 284 PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP, PetscInt); 285 PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 286 287 PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP, PetscScalar); 288 289 PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP, PetscInt); 290 PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP, PetscInt *); 291 PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 292 293 PETSC_EXTERN PetscErrorCode KSPMINRESSetRadius(KSP, PetscReal); 294 PETSC_EXTERN PetscErrorCode KSPMINRESGetUseQLP(KSP, PetscBool *); 295 PETSC_EXTERN PetscErrorCode KSPMINRESSetUseQLP(KSP, PetscBool); 296 297 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP, PC *); 298 PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC); 299 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *); 300 PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat); 301 302 PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat); 303 PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *); 304 #if PetscDefined(HAVE_HPDDM) 305 PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMSetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) 306 { 307 return KSPHPDDMSetDeflationMat(ksp, U); 308 } 309 PETSC_DEPRECATED_FUNCTION(3, 18, 0, "KSPHPDDMGetDeflationMat()", ) static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) 310 { 311 return KSPHPDDMGetDeflationMat(ksp, U); 312 } 313 #endif 314 PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPMatSolve()", ) static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) 315 { 316 return KSPMatSolve(ksp, B, X); 317 } 318 /*E 319 KSPHPDDMType - Type of Krylov method used by `KSPHPDDM` 320 321 Values: 322 + `KSP_HPDDM_TYPE_GMRES` (default) - Generalized Minimal Residual method 323 . `KSP_HPDDM_TYPE_BGMRES` - block GMRES 324 . `KSP_HPDDM_TYPE_CG` - Conjugate Gradient 325 . `KSP_HPDDM_TYPE_BCG` - block CG 326 . `KSP_HPDDM_TYPE_GCRODR` - Generalized Conjugate Residual method with inner Orthogonalization and Deflated Restarting 327 . `KSP_HPDDM_TYPE_BGCRODR` - block GCRODR 328 . `KSP_HPDDM_TYPE_BFBCG` - breakdown-free BCG 329 - `KSP_HPDDM_TYPE_PREONLY` - apply the preconditioner only 330 331 Level: intermediate 332 333 .seealso: [](ch_ksp), `KSPHPDDM`, `KSPHPDDMSetType()` 334 E*/ 335 typedef enum { 336 KSP_HPDDM_TYPE_GMRES = 0, 337 KSP_HPDDM_TYPE_BGMRES = 1, 338 KSP_HPDDM_TYPE_CG = 2, 339 KSP_HPDDM_TYPE_BCG = 3, 340 KSP_HPDDM_TYPE_GCRODR = 4, 341 KSP_HPDDM_TYPE_BGCRODR = 5, 342 KSP_HPDDM_TYPE_BFBCG = 6, 343 KSP_HPDDM_TYPE_PREONLY = 7 344 } KSPHPDDMType; 345 PETSC_EXTERN const char *const KSPHPDDMTypes[]; 346 347 /*E 348 KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM` 349 350 Values: 351 + `KSP_HPDDM_PRECISION_HALF` - default when PETSc is configured `--with-precision=__fp16` 352 . `KSP_HPDDM_PRECISION_SINGLE` - default when PETSc is configured `--with-precision=single` 353 . `KSP_HPDDM_PRECISION_DOUBLE` - default when PETSc is configured `--with-precision=double` 354 - `KSP_HPDDM_PRECISION_QUADRUPLE` - default when PETSc is configured `--with-precision=__float128` 355 356 Level: intermediate 357 358 .seealso: [](ch_ksp), `KSP`, `KSPHPDDM` 359 E*/ 360 typedef enum { 361 KSP_HPDDM_PRECISION_HALF = 0, 362 KSP_HPDDM_PRECISION_SINGLE = 1, 363 KSP_HPDDM_PRECISION_DOUBLE = 2, 364 KSP_HPDDM_PRECISION_QUADRUPLE = 3 365 } KSPHPDDMPrecision; 366 PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType); 367 PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *); 368 369 /*E 370 KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed in the GMRES solvers 371 372 Values: 373 + `KSP_GMRES_CGS_REFINE_NEVER` - one step of classical Gram-Schmidt 374 . `KSP_GMRES_CGS_REFINE_IFNEEDED` - a second step is performed if the first step does not satisfy some criteria 375 - `KSP_GMRES_CGS_REFINE_ALWAYS` - always perform two steps 376 377 Level: advanced 378 379 .seealso: [](ch_ksp), `KSP`, `KSPGMRES`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 380 `KSPGMRESGetOrthogonalization()`, 381 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()` 382 E*/ 383 typedef enum { 384 KSP_GMRES_CGS_REFINE_NEVER, 385 KSP_GMRES_CGS_REFINE_IFNEEDED, 386 KSP_GMRES_CGS_REFINE_ALWAYS 387 } KSPGMRESCGSRefinementType; 388 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 389 390 /*MC 391 KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 392 393 Level: advanced 394 395 Note: 396 Possibly unstable, but the fastest to compute 397 398 .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 399 `KSP`, `KSPGMRESGetOrthogonalization()`, 400 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 401 `KSPGMRESModifiedGramSchmidtOrthogonalization()` 402 M*/ 403 404 /*MC 405 KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 406 iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 407 poor orthogonality. 408 409 Level: advanced 410 411 Note: 412 This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to 413 estimate the orthogonality but is more stable. 414 415 .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 416 `KSP`, `KSPGMRESGetOrthogonalization()`, 417 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 418 `KSPGMRESModifiedGramSchmidtOrthogonalization()` 419 M*/ 420 421 /*MC 422 KSP_GMRES_CGS_REFINE_ALWAYS - Do two steps of the classical (unmodified) Gram-Schmidt process. 423 424 Level: advanced 425 426 Notes: 427 This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice 428 but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`. 429 430 You should only use this if you absolutely know that the iterative refinement is needed. 431 432 .seealso: [](ch_ksp), `KSPGMRES`, `KSPGMRESCGSRefinementType`, `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, 433 `KSP`, `KSPGMRESGetOrthogonalization()`, 434 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 435 `KSPGMRESModifiedGramSchmidtOrthogonalization()` 436 M*/ 437 438 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType); 439 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *); 440 441 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *); 442 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *); 443 PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 444 445 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal); 446 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *); 447 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *); 448 449 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal); 450 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool); 451 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt); 452 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool); 453 454 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 455 PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP); 456 457 PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *); 458 PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm, const char[], const char[], const char[], PetscInt, const char *[], int, int, int, int, PetscDrawLG *); 459 PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 460 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 461 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 462 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 463 PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 464 PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 465 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 466 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 467 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 468 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 469 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 470 PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 471 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 472 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 473 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 474 PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 475 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 476 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 477 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 478 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 479 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 480 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorResidual()", ) static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 481 { 482 return KSPMonitorResidual(ksp, n, rnorm, vf); 483 } 484 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidual()", ) static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 485 { 486 return KSPMonitorTrueResidual(ksp, n, rnorm, vf); 487 } 488 PETSC_DEPRECATED_FUNCTION(3, 15, 0, "KSPMonitorTrueResidualMax()", ) static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 489 { 490 return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf); 491 } 492 493 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *); 494 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *); 495 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **); 496 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *); 497 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal); 498 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *); 499 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **); 500 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **); 501 502 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec); 503 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec); 504 505 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat); 506 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *); 507 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *); 508 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]); 509 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]); 510 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]); 511 512 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool); 513 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *); 514 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool); 515 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *); 516 517 PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer); 518 PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer); 519 PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]); 520 PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer); 521 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **)); 522 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP); 523 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP); 524 PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer); 525 526 PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonView()", ) static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v) 527 { 528 return KSPConvergedReasonView(ksp, v); 529 } 530 PETSC_DEPRECATED_FUNCTION(3, 14, 0, "KSPConvergedReasonViewFromOptions()", ) static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) 531 { 532 return KSPConvergedReasonViewFromOptions(ksp); 533 } 534 535 #define KSP_FILE_CLASSID 1211223 536 537 PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool); 538 PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool); 539 PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *); 540 PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *); 541 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 542 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 543 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 544 545 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *); 546 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *); 547 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *); 548 PETSC_EXTERN PetscErrorCode PCMPIGetKSP(PC, KSP *); 549 550 /*E 551 KSPNormType - Norm calculated by the `KSP` and passed in the Krylov convergence 552 test routines. 553 554 Values: 555 + `KSP_NORM_DEFAULT` - use the default for the current `KSPType` 556 . `KSP_NORM_NONE` - use no norm calculation 557 . `KSP_NORM_PRECONDITIONED` - use the preconditioned residual norm 558 . `KSP_NORM_UNPRECONDITIONED` - use the unpreconditioned residual norm 559 - `KSP_NORM_NATURAL` - use the natural norm (the norm induced by the linear operator) 560 561 Level: advanced 562 563 Note: 564 Each solver only supports a subset of these and some may support different ones 565 depending on whether left or right preconditioning is used, see `KSPSetPCSide()` 566 567 .seealso: [](ch_ksp), `KSP`, `PCSide`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`, 568 `KSPSetConvergenceTest()`, `KSPSetPCSide()`, `KSP_NORM_DEFAULT`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 569 E*/ 570 typedef enum { 571 KSP_NORM_DEFAULT = -1, 572 KSP_NORM_NONE = 0, 573 KSP_NORM_PRECONDITIONED = 1, 574 KSP_NORM_UNPRECONDITIONED = 2, 575 KSP_NORM_NATURAL = 3 576 } KSPNormType; 577 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 578 PETSC_EXTERN const char *const *const KSPNormTypes; 579 580 /*MC 581 KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 582 possibly save some computation but means the convergence test cannot 583 be based on a norm of a residual etc. 584 585 Level: advanced 586 587 Note: 588 Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored 589 590 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 591 M*/ 592 593 /*MC 594 KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 595 convergence test routine. 596 597 Level: advanced 598 599 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 600 M*/ 601 602 /*MC 603 KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 604 convergence test routine. 605 606 Level: advanced 607 608 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 609 M*/ 610 611 /*MC 612 KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 613 convergence test routine. This is only supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR` 614 615 Level: advanced 616 617 .seealso: [](ch_ksp), `KSPNormType`, `KSP`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()` 618 M*/ 619 620 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType); 621 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *); 622 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt); 623 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt); 624 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool); 625 626 #define KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED KSP_CONVERGED_CG_NEG_CURVE PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_NEG_CURVE", ) 627 #define KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED KSP_CONVERGED_CG_CONSTRAINED PETSC_DEPRECATED_ENUM(3, 19, 0, "KSP_CONVERGED_STEP_LENGTH", ) 628 #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM(3, 11, 0, "KSP_DIVERGED_PC_FAILED", ) 629 /*E 630 KSPConvergedReason - reason a Krylov method was determined to have converged or diverged 631 632 Values: 633 + `KSP_CONVERGED_RTOL_NORMAL` - requested decrease in the residual for the normal equations 634 . `KSP_CONVERGED_ATOL_NORMAL` - requested absolute value in the residual for the normal equations 635 . `KSP_CONVERGED_RTOL` - requested decrease in the residual 636 . `KSP_CONVERGED_ATOL` - requested absolute value in the residual 637 . `KSP_CONVERGED_ITS` - requested number of iterations 638 . `KSP_CONVERGED_NEG_CURVE` - see note below 639 . `KSP_CONVERGED_STEP_LENGTH` - see note below 640 . `KSP_CONVERGED_HAPPY_BREAKDOWN` - happy breakdown (meaning early convergence of the `KSPType` occurred). 641 . `KSP_DIVERGED_NULL` - breakdown when solving the Hessenberg system within GMRES 642 . `KSP_DIVERGED_ITS` - requested number of iterations 643 . `KSP_DIVERGED_DTOL` - large increase in the residual norm 644 . `KSP_DIVERGED_BREAKDOWN` - breakdown in the Krylov method 645 . `KSP_DIVERGED_BREAKDOWN_BICG` - breakdown in the `KSPBGCS` Krylov method 646 . `KSP_DIVERGED_NONSYMMETRIC` - the operator or preonditioner was not symmetric for a `KSPType` that requires symmetry 647 . `KSP_DIVERGED_INDEFINITE_PC` - the preconditioner was indefinite for a `KSPType` that requires it be definite 648 . `KSP_DIVERGED_NANORINF` - a not a number of infinity was detected in a vector during the computation 649 . `KSP_DIVERGED_INDEFINITE_MAT` - the operator was indefinite for a `KSPType` that requires it be definite 650 - `KSP_DIVERGED_PC_FAILED` - the action of the preconditioner failed for some reason 651 652 Level: beginner 653 654 Note: 655 The values `KSP_CONVERGED_NEG_CURVE`, and `KSP_CONVERGED_STEP_LENGTH` are returned only by `KSPCG`, `KSPMINRES` and by 656 the special `KSPNASH`, `KSPSTCG`, and `KSPGLTR` solvers which are used by the `SNESNEWTONTR` (trust region) solver. 657 658 Developer Note: 659 The string versions of these are `KSPConvergedReasons`; if you change 660 any of the values here also change them that array of names. 661 662 .seealso: [](ch_ksp), `KSP`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()` 663 E*/ 664 typedef enum { /* converged */ 665 KSP_CONVERGED_RTOL_NORMAL = 1, 666 KSP_CONVERGED_ATOL_NORMAL = 9, 667 KSP_CONVERGED_RTOL = 2, 668 KSP_CONVERGED_ATOL = 3, 669 KSP_CONVERGED_ITS = 4, 670 KSP_CONVERGED_NEG_CURVE = 5, 671 KSP_CONVERGED_CG_NEG_CURVE_DEPRECATED = 5, 672 KSP_CONVERGED_CG_CONSTRAINED_DEPRECATED = 6, 673 KSP_CONVERGED_STEP_LENGTH = 6, 674 KSP_CONVERGED_HAPPY_BREAKDOWN = 7, 675 /* diverged */ 676 KSP_DIVERGED_NULL = -2, 677 KSP_DIVERGED_ITS = -3, 678 KSP_DIVERGED_DTOL = -4, 679 KSP_DIVERGED_BREAKDOWN = -5, 680 KSP_DIVERGED_BREAKDOWN_BICG = -6, 681 KSP_DIVERGED_NONSYMMETRIC = -7, 682 KSP_DIVERGED_INDEFINITE_PC = -8, 683 KSP_DIVERGED_NANORINF = -9, 684 KSP_DIVERGED_INDEFINITE_MAT = -10, 685 KSP_DIVERGED_PC_FAILED = -11, 686 KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11, 687 688 KSP_CONVERGED_ITERATING = 0 689 } KSPConvergedReason; 690 PETSC_EXTERN const char *const *KSPConvergedReasons; 691 692 /*MC 693 KSP_CONVERGED_RTOL - $||r|| \le rtol*||b||$ or $rtol*||b - A*x_0||$ if `KSPConvergedDefaultSetUIRNorm()` was called 694 695 Level: beginner 696 697 Notes: 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 - $||r|| \le atol$ 709 710 Level: beginner 711 712 Notes: 713 See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 714 for left preconditioning it is the 2-norm of the preconditioned residual, and the 715 2-norm of the residual for right preconditioning 716 717 See also `KSP_CONVERGED_RTOL` which may apply before this tolerance. 718 719 .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 720 M*/ 721 722 /*MC 723 KSP_DIVERGED_DTOL - $||r|| \ge dtol*||b||$ 724 725 Level: beginner 726 727 Note: 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 .seealso: [](ch_ksp), `KSPNormType`, `KSP_CONVERGED_ATOL`, `KSP_CONVERGED_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()`, `KSPGuessSetType()`, `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 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *); 976 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *); 977 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP); 978 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 979 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 980 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *); 981 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType); 982 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *); 983 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *); 984 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *); 985 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *); 986 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *); 987 988 PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *); 989 PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *); 990 PETSC_EXTERN PetscErrorCode MatCreateLMVMDBFGS(MPI_Comm, PetscInt, PetscInt, Mat *); 991 PETSC_EXTERN PetscErrorCode MatCreateLMVMDDFP(MPI_Comm, PetscInt, PetscInt, Mat *); 992 PETSC_EXTERN PetscErrorCode MatCreateLMVMDQN(MPI_Comm, PetscInt, PetscInt, Mat *); 993 PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *); 994 PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 995 PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 996 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 997 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 998 PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 999 1000 PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 1001 PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *); 1002 PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 1003 PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 1004 PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat); 1005 PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat); 1006 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 1007 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 1008 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 1009 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 1010 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 1011 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec); 1012 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 1013 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *); 1014 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *); 1015 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *); 1016 PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt); 1017 PETSC_EXTERN PetscErrorCode MatLMVMGetHistorySize(Mat, PetscInt *); 1018 PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *); 1019 PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *); 1020 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar); 1021 1022 /*E 1023 MatLMVMSymBroydenScaleType - Scaling type for symmetric Broyden. 1024 1025 Values: 1026 + `MAT_LMVM_SYMBROYDEN_SCALE_NONE` - No scaling 1027 . `MAT_LMVM_SYMBROYDEN_SCALE_SCALAR` - scalar scaling 1028 . `MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL` - diagonal scaling 1029 - `MAT_LMVM_SYMBROYDEN_SCALE_USER` - user-provided scale option 1030 1031 Level: intermediate 1032 1033 .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMSymBroydenSetScaleType()` 1034 E*/ 1035 typedef enum { 1036 MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0, 1037 MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1, 1038 MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2, 1039 MAT_LMVM_SYMBROYDEN_SCALE_USER = 3 1040 } MatLMVMSymBroydenScaleType; 1041 PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[]; 1042 1043 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType); 1044 1045 /*E 1046 MatLMVMDenseType - Memory storage strategy for dense variants `MATLMVM`. 1047 1048 Values: 1049 + `MAT_LMVM_DENSE_REORDER` - reorders memory to minimize kernel launch 1050 - `MAT_LMVM_DENSE_INPLACE` - computes inplace to minimize memory movement 1051 1052 Level: intermediate 1053 1054 .seealso: [](ch_matrices), `MatLMVM`, `MatLMVMDenseSetType()` 1055 E*/ 1056 typedef enum { 1057 MAT_LMVM_DENSE_REORDER, 1058 MAT_LMVM_DENSE_INPLACE 1059 } MatLMVMDenseType; 1060 PETSC_EXTERN const char *const MatLMVMDenseTypes[]; 1061 1062 PETSC_EXTERN PetscErrorCode MatLMVMDenseSetType(Mat, MatLMVMDenseType); 1063 1064 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM); 1065 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool); 1066 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *); 1067 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *); 1068 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *); 1069 1070 /*S 1071 KSPComputeRHSFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeRHS()` 1072 1073 Calling Sequence: 1074 + ksp - `ksp` context 1075 . b - output vector 1076 - ctx - [optional] user-defined function context 1077 1078 Level: beginner 1079 1080 .seealso: [](ch_snes), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeInitialGuessFn`, `KSPComputeOperatorsFn` 1081 S*/ 1082 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeRHSFn)(KSP ksp, Vec b, void *ctx); 1083 1084 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, KSPComputeRHSFn *, void *); 1085 1086 /*S 1087 KSPComputeOperatorsFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeOperators()` 1088 1089 Calling Sequence: 1090 + ksp - `KSP` context 1091 . A - the operator that defines the linear system 1092 . P - an operator from which to build the preconditioner (often the same as `A`) 1093 - ctx - [optional] user-defined function context 1094 1095 Level: beginner 1096 1097 .seealso: [](ch_snes), `KSP`, `KSPSetComputeRHS()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeInitialGuessFn` 1098 S*/ 1099 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeOperatorsFn)(KSP ksp, Mat A, Mat P, void *ctx); 1100 1101 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, KSPComputeOperatorsFn, void *); 1102 1103 /*S 1104 KSPComputeInitialGuessFn - A prototype of a `KSP` evaluation function that would be passed to `KSPSetComputeInitialGuess()` 1105 1106 Calling Sequence: 1107 + ksp - `ksp` context 1108 . x - output vector 1109 - ctx - [optional] user-defined function context 1110 1111 Level: beginner 1112 1113 .seealso: [](ch_snes), `KSP`, `KSPSetComputeInitialGuess()`, `SNESGetFunction()`, `KSPComputeRHSFn`, `KSPComputeOperatorsFn` 1114 S*/ 1115 PETSC_EXTERN_TYPEDEF typedef PetscErrorCode(KSPComputeInitialGuessFn)(KSP ksp, Vec x, void *ctx); 1116 1117 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, KSPComputeInitialGuessFn *, void *); 1118 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, KSPComputeOperatorsFn *, void *); 1119 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, KSPComputeOperatorsFn **, void *); 1120 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, KSPComputeRHSFn *, void *); 1121 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, KSPComputeRHSFn **, void *); 1122 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, KSPComputeInitialGuessFn *, void *); 1123 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, KSPComputeInitialGuessFn **, void *); 1124 1125 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec); 1126 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); 1127 PETSC_EXTERN PetscErrorCode DMSwarmProjectFields(DM, DM, PetscInt, const char **, Vec[], ScatterMode mode); 1128 1129 PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *); 1130 PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal); 1131 1132 PETSC_EXTERN PetscErrorCode PCBJKOKKOSSetKSP(PC, KSP); 1133 PETSC_EXTERN PetscErrorCode PCBJKOKKOSGetKSP(PC, KSP *); 1134 1135 PETSC_EXTERN PetscErrorCode DMCopyDMKSP(DM, DM); 1136