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: `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: `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: `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 KSPFETIDPGetInnerBDDC(KSP, PC *); 268 PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP, PC); 269 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP, KSP *); 270 PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP, Mat); 271 272 PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP, Mat); 273 PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP, Mat *); 274 #if PetscDefined(HAVE_HPDDM) 275 PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMSetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) 276 { 277 return KSPHPDDMSetDeflationMat(ksp, U); 278 } 279 PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMGetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) 280 { 281 return KSPHPDDMGetDeflationMat(ksp, U); 282 } 283 #endif 284 PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) 285 { 286 return KSPMatSolve(ksp, B, X); 287 } 288 /*E 289 KSPHPDDMType - Type of Krylov method used by `KSPHPDDM` 290 291 Level: intermediate 292 293 Values: 294 $ `KSP_HPDDM_TYPE_GMRES` (default) 295 $ `KSP_HPDDM_TYPE_BGMRES` 296 $ `KSP_HPDDM_TYPE_CG` 297 $ `KSP_HPDDM_TYPE_BCG` 298 $ `KSP_HPDDM_TYPE_GCRODR` 299 $ `KSP_HPDDM_TYPE_BGCRODR` 300 $ `KSP_HPDDM_TYPE_BFBCG` 301 $ `KSP_HPDDM_TYPE_PREONLY` 302 303 .seealso: `KSPHPDDM`, `KSPHPDDMSetType()` 304 E*/ 305 typedef enum { 306 KSP_HPDDM_TYPE_GMRES = 0, 307 KSP_HPDDM_TYPE_BGMRES = 1, 308 KSP_HPDDM_TYPE_CG = 2, 309 KSP_HPDDM_TYPE_BCG = 3, 310 KSP_HPDDM_TYPE_GCRODR = 4, 311 KSP_HPDDM_TYPE_BGCRODR = 5, 312 KSP_HPDDM_TYPE_BFBCG = 6, 313 KSP_HPDDM_TYPE_PREONLY = 7 314 } KSPHPDDMType; 315 PETSC_EXTERN const char *const KSPHPDDMTypes[]; 316 /*E 317 KSPHPDDMPrecision - Precision of Krylov bases used by `KSPHPDDM` 318 319 Level: intermediate 320 321 Values: 322 $ `KSP_HPDDM_PRECISION_HALF` 323 $ `KSP_HPDDM_PRECISION_SINGLE` (default when PETSc is configured --with-precision=single) 324 $ `KSP_HPDDM_PRECISION_DOUBLE` (default when PETSc is configured --with-precision=double) 325 $ `KSP_HPDDM_PRECISION_QUADRUPLE` (default when PETSc is configured --with-precision=__float128) 326 327 .seealso: `KSPHPDDM` 328 E*/ 329 typedef enum { 330 KSP_HPDDM_PRECISION_HALF = 0, 331 KSP_HPDDM_PRECISION_SINGLE = 1, 332 KSP_HPDDM_PRECISION_DOUBLE = 2, 333 KSP_HPDDM_PRECISION_QUADRUPLE = 3 334 } KSPHPDDMPrecision; 335 PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP, KSPHPDDMType); 336 PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP, KSPHPDDMType *); 337 338 /*E 339 KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 340 341 Level: advanced 342 343 .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 344 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()` 345 346 E*/ 347 typedef enum { 348 KSP_GMRES_CGS_REFINE_NEVER, 349 KSP_GMRES_CGS_REFINE_IFNEEDED, 350 KSP_GMRES_CGS_REFINE_ALWAYS 351 } KSPGMRESCGSRefinementType; 352 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 353 /*MC 354 KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 355 356 Level: advanced 357 358 Note: 359 Possibly unstable, but the fastest to compute 360 361 .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 362 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 363 `KSPGMRESModifiedGramSchmidtOrthogonalization()` 364 M*/ 365 366 /*MC 367 KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 368 iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 369 poor orthogonality. 370 371 Level: advanced 372 373 Note: This is slower than `KSP_GMRES_CGS_REFINE_NEVER` because it requires an extra norm computation to 374 estimate the orthogonality but is more stable. 375 376 .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 377 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 378 `KSPGMRESModifiedGramSchmidtOrthogonalization()` 379 M*/ 380 381 /*MC 382 KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 383 384 Level: advanced 385 386 Notes: 387 This is roughly twice the cost of `KSP_GMRES_CGS_REFINE_NEVER` because it performs the process twice 388 but it saves the extra norm calculation needed by `KSP_GMRES_CGS_REFINE_IFNEEDED`. 389 390 You should only use this if you absolutely know that the iterative refinement is needed. 391 392 .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 393 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 394 `KSPGMRESModifiedGramSchmidtOrthogonalization()` 395 M*/ 396 397 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP, KSPGMRESCGSRefinementType); 398 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP, KSPGMRESCGSRefinementType *); 399 400 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP, PetscInt, PetscInt, PetscReal, void *); 401 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP, PetscInt, PetscInt, PetscReal, void *); 402 PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscInt, PetscReal, void *), void *, PetscErrorCode (*)(void *)); 403 404 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP, PetscReal); 405 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP, PetscReal *); 406 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP, PetscReal *); 407 408 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP, PetscReal); 409 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP, PetscBool); 410 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP, PetscInt); 411 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP, PetscBool); 412 413 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 414 PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP); 415 PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 416 417 PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP, const char[], const char[], void *); 418 PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm, const char[], const char[], const char[], PetscInt, const char *[], int, int, int, int, PetscDrawLG *); 419 PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 420 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 421 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 422 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 423 PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 424 PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 425 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 426 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 427 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 428 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 429 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 430 PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 431 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 432 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 433 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 434 PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 435 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 436 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 437 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 438 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 439 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 440 PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 441 { 442 return KSPMonitorResidual(ksp, n, rnorm, vf); 443 } 444 PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 445 { 446 return KSPMonitorTrueResidual(ksp, n, rnorm, vf); 447 } 448 PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) 449 { 450 return KSPMonitorTrueResidualMax(ksp, n, rnorm, vf); 451 } 452 453 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP, PetscInt, PetscReal, void *); 454 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP, PetscInt, PetscReal, void *); 455 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **); 456 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceCreate(void *); 457 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceSetCoefficient(void *, PetscReal); 458 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP, PetscInt, PetscReal, void *); 459 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP, void **); 460 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void **); 461 462 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP, Vec, Vec); 463 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP, Vec, Vec, Vec, Vec, Vec); 464 465 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP, Mat, Mat); 466 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP, Mat *, Mat *); 467 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP, PetscBool *, PetscBool *); 468 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP, const char[]); 469 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP, const char[]); 470 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP, const char *[]); 471 472 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP, PetscBool); 473 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP, PetscBool *); 474 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP, PetscBool); 475 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP, PetscBool *); 476 477 PETSC_EXTERN PetscErrorCode KSPView(KSP, PetscViewer); 478 PETSC_EXTERN PetscErrorCode KSPLoad(KSP, PetscViewer); 479 PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP, PetscObject, const char[]); 480 PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP, PetscViewer); 481 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP, PetscErrorCode (*)(KSP, void *), void *vctx, PetscErrorCode (*)(void **)); 482 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP); 483 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP); 484 PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP, PetscViewer); 485 486 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") static inline PetscErrorCode KSPReasonView(KSP ksp, PetscViewer v) 487 { 488 return KSPConvergedReasonView(ksp, v); 489 } 490 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) 491 { 492 return KSPConvergedReasonViewFromOptions(ksp); 493 } 494 495 #define KSP_FILE_CLASSID 1211223 496 497 PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP, PetscBool); 498 PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP, PetscBool); 499 PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP, Vec *); 500 PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP, PetscReal *, PetscReal *); 501 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 502 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP, PetscInt, PetscReal, PetscViewerAndFormat *); 503 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **); 504 505 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC, KSP *); 506 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC, KSP *); 507 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC, KSP *); 508 509 /*E 510 KSPNormType - Norm that is passed in the Krylov convergence 511 test routines. 512 513 Level: advanced 514 515 Each solver only supports a subset of these and some may support different ones 516 depending on left or right preconditioning, see `KSPSetPCSide()` 517 518 Developer Note: 519 this must match petsc/finclude/petscksp.h 520 521 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`, 522 `KSPSetConvergenceTest()`, `KSPSetPCSide()` 523 E*/ 524 typedef enum { 525 KSP_NORM_DEFAULT = -1, 526 KSP_NORM_NONE = 0, 527 KSP_NORM_PRECONDITIONED = 1, 528 KSP_NORM_UNPRECONDITIONED = 2, 529 KSP_NORM_NATURAL = 3 530 } KSPNormType; 531 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 532 PETSC_EXTERN const char *const *const KSPNormTypes; 533 534 /*MC 535 KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 536 possibly save some computation but means the convergence test cannot 537 be based on a norm of a residual etc. 538 539 Level: advanced 540 541 Note: Some Krylov methods need to compute a residual norm (such as `KPSGMRES`) and then this option is ignored 542 543 .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 544 M*/ 545 546 /*MC 547 KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 548 convergence test routine. 549 550 Level: advanced 551 552 .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 553 M*/ 554 555 /*MC 556 KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 557 convergence test routine. 558 559 Level: advanced 560 561 .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 562 M*/ 563 564 /*MC 565 KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 566 convergence test routine. This is only supported by `KSPCG`, `KSPCR`, `KSPCGNE`, `KSPCGS`, `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR` 567 568 Level: advanced 569 570 .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()` 571 M*/ 572 573 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP, KSPNormType); 574 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP, KSPNormType *); 575 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp, KSPNormType, PCSide, PetscInt); 576 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP, PetscInt); 577 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP, PetscBool); 578 579 #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)") 580 /*E 581 KSPConvergedReason - reason a Krylov method was said to have converged or diverged 582 583 Level: beginner 584 585 Notes: 586 See `KSPGetConvergedReason()` for explanation of each value 587 588 Developer Notes: 589 This must match petsc/finclude/petscksp.h 590 591 The string versions of these are `KSPConvergedReasons`; if you change 592 any of the values here also change them that array of names. 593 594 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()` 595 E*/ 596 typedef enum { /* converged */ 597 KSP_CONVERGED_RTOL_NORMAL = 1, 598 KSP_CONVERGED_ATOL_NORMAL = 9, 599 KSP_CONVERGED_RTOL = 2, 600 KSP_CONVERGED_ATOL = 3, 601 KSP_CONVERGED_ITS = 4, 602 KSP_CONVERGED_CG_NEG_CURVE = 5, 603 KSP_CONVERGED_CG_CONSTRAINED = 6, 604 KSP_CONVERGED_STEP_LENGTH = 7, 605 KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 606 /* diverged */ 607 KSP_DIVERGED_NULL = -2, 608 KSP_DIVERGED_ITS = -3, 609 KSP_DIVERGED_DTOL = -4, 610 KSP_DIVERGED_BREAKDOWN = -5, 611 KSP_DIVERGED_BREAKDOWN_BICG = -6, 612 KSP_DIVERGED_NONSYMMETRIC = -7, 613 KSP_DIVERGED_INDEFINITE_PC = -8, 614 KSP_DIVERGED_NANORINF = -9, 615 KSP_DIVERGED_INDEFINITE_MAT = -10, 616 KSP_DIVERGED_PC_FAILED = -11, 617 KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11, 618 619 KSP_CONVERGED_ITERATING = 0 620 } KSPConvergedReason; 621 PETSC_EXTERN const char *const *KSPConvergedReasons; 622 623 /*MC 624 KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if `KSPConvergedDefaultSetUIRNorm()` was called 625 626 Level: beginner 627 628 See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 629 for left preconditioning it is the 2-norm of the preconditioned residual, and the 630 2-norm of the residual for right preconditioning 631 632 See also `KSP_CONVERGED_ATOL` which may apply before this tolerance. 633 634 .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 635 636 M*/ 637 638 /*MC 639 KSP_CONVERGED_ATOL - norm(r) <= atol 640 641 Level: beginner 642 643 See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 644 for left preconditioning it is the 2-norm of the preconditioned residual, and the 645 2-norm of the residual for right preconditioning 646 647 See also `KSP_CONVERGED_RTOL` which may apply before this tolerance. 648 649 .seealso: `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 650 651 M*/ 652 653 /*MC 654 KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 655 656 Level: beginner 657 658 See `KSPNormType` and `KSPSetNormType()` for possible norms that may be used. By default 659 for left preconditioning it is the 2-norm of the preconditioned residual, and the 660 2-norm of the residual for right preconditioning 661 662 Level: beginner 663 664 .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 665 666 M*/ 667 668 /*MC 669 KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 670 reached 671 672 Level: beginner 673 674 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 675 676 M*/ 677 678 /*MC 679 KSP_CONVERGED_ITS - Used by the `KSPPREONLY` solver after the single iteration of 680 the preconditioner is applied. Also used when the `KSPConvergedSkip()` convergence 681 test routine is set in KSP. 682 683 Level: beginner 684 685 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 686 687 M*/ 688 689 /*MC 690 KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 691 method could not continue to enlarge the Krylov space. Could be due to a singular matrix or 692 preconditioner. In `KSPHPDDM`, this is also returned when some search directions within a block 693 are colinear. 694 695 Level: beginner 696 697 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 698 699 M*/ 700 701 /*MC 702 KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the `KSPBICG` method was detected so the 703 method could not continue to enlarge the Krylov space. 704 705 Level: beginner 706 707 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 708 709 M*/ 710 711 /*MC 712 KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 713 symmetric and this Krylov method (`KSPCG`, `KSPMINRES`, `KSPCR`) requires symmetry 714 715 Level: beginner 716 717 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 718 719 M*/ 720 721 /*MC 722 KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 723 positive and negative eigenvalues) and this Krylov method (`KSPCG`) requires it to 724 be positive definite 725 726 Level: beginner 727 728 Note: 729 This can happen with the `PCICC` preconditioner, use -pc_factor_shift_positive_definite to force 730 the `PCICC` preconditioner to generate a positive definite preconditioner 731 732 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 733 734 M*/ 735 736 /*MC 737 KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a 738 zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 739 such as `PCFIELDSPLIT`. 740 741 Level: beginner 742 743 Notes: 744 Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details. 745 746 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 747 748 M*/ 749 750 /*MC 751 KSP_CONVERGED_ITERATING - This flag is returned if you call `KSPGetConvergedReason()` 752 while the `KSPSolve()` is still running. 753 754 Level: beginner 755 756 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 757 758 M*/ 759 760 PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP, PetscErrorCode (*)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void *, PetscErrorCode (*)(void *)); 761 PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *)); 762 PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP, PetscErrorCode (**)(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *), void **, PetscErrorCode (**)(void *)); 763 PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP, void *); 764 PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 765 PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 766 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 767 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 768 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 769 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 770 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP, PetscBool); 771 PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP, PetscInt, PetscReal, KSPConvergedReason *, void *); 772 PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP, KSPConvergedReason *); 773 PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP, const char **); 774 PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP, PetscReal *, PetscReal *, PetscReal *, PetscReal *); 775 776 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") static inline void KSPDefaultConverged(void) 777 { /* never called */ 778 } 779 #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 780 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") static inline void KSPDefaultConvergedDestroy(void) 781 { /* never called */ 782 } 783 #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 784 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") static inline void KSPDefaultConvergedCreate(void) 785 { /* never called */ 786 } 787 #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 788 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUIRNorm(void) 789 { /* never called */ 790 } 791 #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 792 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUMIRNorm(void) 793 { /* never called */ 794 } 795 #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 796 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") static inline void KSPSkipConverged(void) 797 { /* never called */ 798 } 799 #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 800 801 PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP, MatType, Mat *); 802 PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") static inline PetscErrorCode KSPComputeExplicitOperator(KSP A, Mat *B) 803 { 804 return KSPComputeOperator(A, NULL, B); 805 } 806 807 /*E 808 KSPCGType - Determines what type of CG to use 809 810 Level: beginner 811 812 .seealso: `KSPCGSetType()` 813 E*/ 814 typedef enum { 815 KSP_CG_SYMMETRIC = 0, 816 KSP_CG_HERMITIAN = 1 817 } KSPCGType; 818 PETSC_EXTERN const char *const KSPCGTypes[]; 819 820 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP, KSPCGType); 821 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP, PetscBool); 822 823 PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP, PetscReal); 824 PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP, PetscReal *); 825 PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP, PetscReal *); 826 827 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP, PetscReal *); 828 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP, PetscReal *); 829 PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp, PetscReal *x) 830 { 831 return KSPGLTRGetMinEig(ksp, x); 832 } 833 PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp, PetscReal *x) 834 { 835 return KSPGLTRGetLambda(ksp, x); 836 } 837 838 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP, const char[]); 839 PETSC_EXTERN PetscErrorCode KSPPythonGetType(KSP, const char *[]); 840 841 PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC, PetscErrorCode (*)(PC, KSP)); 842 PETSC_EXTERN PetscErrorCode PCPreSolve(PC, KSP); 843 PETSC_EXTERN PetscErrorCode PCPostSolve(PC, KSP); 844 845 #include <petscdrawtypes.h> 846 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP, PetscInt, PetscReal, void *); 847 848 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec)); 849 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC, PetscErrorCode (*)(PC, KSP, Vec, Vec)); 850 851 /*S 852 KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods. 853 854 Level: beginner 855 856 .seealso: `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType` 857 S*/ 858 typedef struct _p_KSPGuess *KSPGuess; 859 /*J 860 KSPGuessType - String with the name of a PETSc initial guess approach for Krylov methods. 861 862 Level: beginner 863 864 .seealso: `KSPGuess` 865 J*/ 866 typedef const char *KSPGuessType; 867 #define KSPGUESSFISCHER "fischer" 868 #define KSPGUESSPOD "pod" 869 PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[], PetscErrorCode (*)(KSPGuess)); 870 PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP, KSPGuess); 871 PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP, KSPGuess *); 872 PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess, PetscViewer); 873 PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess *); 874 PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm, KSPGuess *); 875 PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess, KSPGuessType); 876 PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess, KSPGuessType *); 877 PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess, PetscReal); 878 PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 879 PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess, Vec, Vec); 880 PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess, Vec, Vec); 881 PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 882 PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess, PetscInt, PetscInt); 883 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP, PetscInt, PetscInt); 884 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP, PetscBool); 885 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP, PetscBool *); 886 887 /*E 888 MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 889 890 Level: intermediate 891 892 .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`, `MatCreateSchurComplementPmat()` 893 E*/ 894 typedef enum { 895 MAT_SCHUR_COMPLEMENT_AINV_DIAG, 896 MAT_SCHUR_COMPLEMENT_AINV_LUMP, 897 MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG, 898 MAT_SCHUR_COMPLEMENT_AINV_FULL 899 } MatSchurComplementAinvType; 900 PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 901 902 typedef enum { 903 MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0, 904 MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1, 905 MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2, 906 MAT_LMVM_SYMBROYDEN_SCALE_USER = 3 907 } MatLMVMSymBroydenScaleType; 908 PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[]; 909 910 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat, Mat, Mat, Mat, Mat, Mat *); 911 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat, KSP *); 912 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat, KSP); 913 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 914 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat, Mat, Mat, Mat, Mat, Mat); 915 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat, Mat *, Mat *, Mat *, Mat *, Mat *); 916 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat, MatSchurComplementAinvType); 917 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat, MatSchurComplementAinvType *); 918 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat, MatReuse, Mat *); 919 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat, Mat *); 920 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat, IS, IS, IS, IS, MatReuse, Mat *, MatSchurComplementAinvType, MatReuse, Mat *); 921 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat, Mat, Mat, Mat, MatSchurComplementAinvType, MatReuse, Mat *); 922 923 PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm, PetscInt, PetscInt, Mat *); 924 PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm, PetscInt, PetscInt, Mat *); 925 PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm, PetscInt, PetscInt, Mat *); 926 PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 927 PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 928 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 929 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 930 PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm, PetscInt, PetscInt, Mat *); 931 932 PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 933 PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool *); 934 PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 935 PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 936 PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat); 937 PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat); 938 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 939 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 940 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 941 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 942 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 943 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec); 944 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 945 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat *); 946 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC *); 947 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP *); 948 PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt); 949 PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt *); 950 PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt *); 951 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar); 952 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType); 953 954 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP, DM); 955 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP, PetscBool); 956 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP, DM *); 957 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP, void *); 958 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP, void *); 959 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP, PetscErrorCode (*func)(KSP, Vec, void *), void *); 960 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *); 961 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP, PetscErrorCode (*)(KSP, Vec, void *), void *); 962 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM, PetscErrorCode (*)(KSP, Mat, Mat, void *), void *); 963 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM, PetscErrorCode (**)(KSP, Mat, Mat, void *), void *); 964 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM, PetscErrorCode (*)(KSP, Vec, void *), void *); 965 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM, PetscErrorCode (**)(KSP, Vec, void *), void *); 966 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM, PetscErrorCode (*)(KSP, Vec, void *), void *); 967 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM, PetscErrorCode (**)(KSP, Vec, void *), void *); 968 969 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM, Vec, Vec); 970 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); 971 972 PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *); 973 PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal); 974 #endif 975