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