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