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