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 Notes: 20 When a direct solver is used, but no Krylov solver is used, the KSP object is still used but with a 21 KSPType of KSPPREONLY, 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) {return KSPSetMatSolveBatchSize(ksp,n);} 103 PETSC_EXTERN PetscErrorCode KSPGetMatSolveBatchSize(KSP,PetscInt*); 104 PETSC_DEPRECATED_FUNCTION("Use KSPGetMatSolveBatchSize() (since version 3.15)") static inline PetscErrorCode KSPGetMatSolveBlockSize(KSP ksp,PetscInt *n) {return KSPGetMatSolveBatchSize(ksp,n);} 105 PETSC_EXTERN PetscErrorCode KSPReset(KSP); 106 PETSC_EXTERN PetscErrorCode KSPResetViewers(KSP); 107 PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 108 PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 109 PETSC_EXTERN PetscErrorCode KSPGetReusePreconditioner(KSP,PetscBool*); 110 PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool); 111 PETSC_EXTERN PetscErrorCode KSPCheckSolve(KSP,PC,Vec); 112 113 PETSC_EXTERN PetscFunctionList KSPList; 114 PETSC_EXTERN PetscFunctionList KSPGuessList; 115 PETSC_EXTERN PetscFunctionList KSPMonitorList; 116 PETSC_EXTERN PetscFunctionList KSPMonitorCreateList; 117 PETSC_EXTERN PetscFunctionList KSPMonitorDestroyList; 118 PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 119 PETSC_EXTERN PetscErrorCode KSPMonitorRegister(const char[], PetscViewerType, PetscViewerFormat, PetscErrorCode (*)(KSP, PetscInt, PetscReal, PetscViewerAndFormat *), PetscErrorCode (*)(PetscViewer, PetscViewerFormat, void *, PetscViewerAndFormat **), PetscErrorCode (*)(PetscViewerAndFormat **)); 120 121 PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 122 PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 123 PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 124 PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 125 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool); 126 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool*); 127 PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool); 128 PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool*); 129 PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool); 130 PETSC_EXTERN PetscErrorCode KSPSetComputeRitz(KSP,PetscBool); 131 PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool*); 132 PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool); 133 PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool*); 134 PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec*); 135 PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec*); 136 PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 137 PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 138 PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*); 139 PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 140 PETSC_DEPRECATED_FUNCTION("Use KSPCreateVecs() (since version 3.6)") static inline PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);} 141 142 PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 143 PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 144 145 PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 146 PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 147 148 PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 149 PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode (*)(void**)); 150 PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 151 PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void*); 152 PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,const PetscReal*[],PetscInt*); 153 PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool); 154 PETSC_EXTERN PetscErrorCode KSPGetErrorHistory(KSP,const PetscReal*[],PetscInt*); 155 PETSC_EXTERN PetscErrorCode KSPSetErrorHistory(KSP,PetscReal[],PetscInt,PetscBool); 156 157 PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 158 PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec*); 159 PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 160 PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 161 162 PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 163 PETSC_EXTERN PetscErrorCode PCKSPSetKSP(PC,KSP); 164 PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 165 PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 166 PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 167 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 168 PETSC_EXTERN PetscErrorCode PCFieldSplitSchurGetSubKSP(PC,PetscInt*,KSP*[]); 169 PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 170 PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 171 PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 172 PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 173 PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP*); 174 PETSC_EXTERN PetscErrorCode PCDeflationGetCoarseKSP(PC,KSP*); 175 /* 176 PCMGCoarseList contains the list of coarse space constructor currently registered 177 These are added with PCMGRegisterCoarseSpaceConstructor() 178 */ 179 PETSC_EXTERN PetscFunctionList PCMGCoarseList; 180 PETSC_EXTERN PetscErrorCode PCMGRegisterCoarseSpaceConstructor(const char[], PetscErrorCode (*)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat*)); 181 PETSC_EXTERN PetscErrorCode PCMGGetCoarseSpaceConstructor(const char[], PetscErrorCode (**)(PC, PetscInt, DM, KSP, PetscInt, Mat, Mat*)); 182 183 PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec*); 184 PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec*); 185 186 PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 187 PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool); 188 PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 189 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 190 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseNoisy(KSP,PetscBool); 191 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*); 192 PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 193 PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt*); 194 PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]); 195 PETSC_EXTERN PetscErrorCode KSPComputeRitz(KSP,PetscBool,PetscBool,PetscInt*,Vec[],PetscReal[],PetscReal[]); 196 197 /*E 198 199 KSPFCDTruncationType - Define how stored directions are used to orthogonalize in flexible conjugate directions (FCD) methods 200 201 Values: 202 $ KSP_FCD_TRUNC_TYPE_STANDARD - uses all (up to mmax) stored directions 203 $ KSP_FCD_TRUNC_TYPE_NOTAY - uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 204 205 Level: intermediate 206 .seealso `:` `KSPFCG`, `KSPPIPEFCG`, `KSPPIPEGCR`, `KSPFCGSetTruncationType()`, `KSPFCGGetTruncationType()` 207 208 E*/ 209 typedef enum {KSP_FCD_TRUNC_TYPE_STANDARD,KSP_FCD_TRUNC_TYPE_NOTAY} KSPFCDTruncationType; 210 PETSC_EXTERN const char *const KSPFCDTruncationTypes[]; 211 212 PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 213 PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 214 PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 215 PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 216 PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCDTruncationType); 217 PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCDTruncationType*); 218 219 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetMmax(KSP,PetscInt); 220 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetMmax(KSP,PetscInt*); 221 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetNprealloc(KSP,PetscInt); 222 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetNprealloc(KSP,PetscInt*); 223 PETSC_EXTERN PetscErrorCode KSPPIPEFCGSetTruncationType(KSP,KSPFCDTruncationType); 224 PETSC_EXTERN PetscErrorCode KSPPIPEFCGGetTruncationType(KSP,KSPFCDTruncationType*); 225 226 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetMmax(KSP,PetscInt); 227 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetMmax(KSP,PetscInt*); 228 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetNprealloc(KSP,PetscInt); 229 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetNprealloc(KSP,PetscInt*); 230 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetTruncationType(KSP,KSPFCDTruncationType); 231 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetTruncationType(KSP,KSPFCDTruncationType*); 232 PETSC_EXTERN PetscErrorCode KSPPIPEGCRSetUnrollW(KSP,PetscBool); 233 PETSC_EXTERN PetscErrorCode KSPPIPEGCRGetUnrollW(KSP,PetscBool*); 234 235 PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 236 PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 237 PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 238 PETSC_EXTERN PetscErrorCode KSPGMRESSetBreakdownTolerance(KSP,PetscReal); 239 240 PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 241 PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 242 PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 243 PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 244 PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 245 246 PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 247 PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 248 249 PETSC_EXTERN PetscErrorCode KSPPIPEFGMRESSetShift(KSP,PetscScalar); 250 251 PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 252 PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 253 PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 254 255 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerBDDC(KSP,PC*); 256 PETSC_EXTERN PetscErrorCode KSPFETIDPSetInnerBDDC(KSP,PC); 257 PETSC_EXTERN PetscErrorCode KSPFETIDPGetInnerKSP(KSP,KSP*); 258 PETSC_EXTERN PetscErrorCode KSPFETIDPSetPressureOperator(KSP,Mat); 259 260 PETSC_EXTERN PetscErrorCode KSPHPDDMSetDeflationMat(KSP,Mat); 261 PETSC_EXTERN PetscErrorCode KSPHPDDMGetDeflationMat(KSP,Mat*); 262 #if PetscDefined(HAVE_HPDDM) 263 PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMSetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMSetDeflationSpace(KSP ksp, Mat U) { return KSPHPDDMSetDeflationMat(ksp, U); } 264 PETSC_DEPRECATED_FUNCTION("Use KSPHPDDMGetDeflationMat() (since version 3.18)") static inline PetscErrorCode KSPHPDDMGetDeflationSpace(KSP ksp, Mat *U) { return KSPHPDDMGetDeflationMat(ksp, U); } 265 #endif 266 PETSC_DEPRECATED_FUNCTION("Use KSPMatSolve() (since version 3.14)") static inline PetscErrorCode KSPHPDDMMatSolve(KSP ksp, Mat B, Mat X) { return KSPMatSolve(ksp, B, X); } 267 /*E 268 KSPHPDDMType - Type of Krylov method used by KSPHPDDM 269 270 Level: intermediate 271 272 Values: 273 $ KSP_HPDDM_TYPE_GMRES (default) 274 $ KSP_HPDDM_TYPE_BGMRES 275 $ KSP_HPDDM_TYPE_CG 276 $ KSP_HPDDM_TYPE_BCG 277 $ KSP_HPDDM_TYPE_GCRODR 278 $ KSP_HPDDM_TYPE_BGCRODR 279 $ KSP_HPDDM_TYPE_BFBCG 280 $ KSP_HPDDM_TYPE_PREONLY 281 282 .seealso: `KSPHPDDM`, `KSPHPDDMSetType()` 283 E*/ 284 typedef enum { 285 KSP_HPDDM_TYPE_GMRES = 0, 286 KSP_HPDDM_TYPE_BGMRES = 1, 287 KSP_HPDDM_TYPE_CG = 2, 288 KSP_HPDDM_TYPE_BCG = 3, 289 KSP_HPDDM_TYPE_GCRODR = 4, 290 KSP_HPDDM_TYPE_BGCRODR = 5, 291 KSP_HPDDM_TYPE_BFBCG = 6, 292 KSP_HPDDM_TYPE_PREONLY = 7 293 } KSPHPDDMType; 294 PETSC_EXTERN const char *const KSPHPDDMTypes[]; 295 /*E 296 KSPHPDDMPrecision - Precision of Krylov bases used by KSPHPDDM 297 298 Level: intermediate 299 300 Values: 301 $ KSP_HPDDM_PRECISION_HALF (currently unsupported) 302 $ KSP_HPDDM_PRECISION_SINGLE (default when PETSc is configured --with-precision=single) 303 $ KSP_HPDDM_PRECISION_DOUBLE (default when PETSc is configured --with-precision=double) 304 $ KSP_HPDDM_PRECISION_QUADRUPLE (currently unsupported) 305 306 .seealso: `KSPHPDDM` 307 E*/ 308 typedef enum { 309 KSP_HPDDM_PRECISION_HALF = 0, 310 KSP_HPDDM_PRECISION_SINGLE = 1, 311 KSP_HPDDM_PRECISION_DOUBLE = 2, 312 KSP_HPDDM_PRECISION_QUADRUPLE = 3 313 } KSPHPDDMPrecision; 314 PETSC_EXTERN PetscErrorCode KSPHPDDMSetType(KSP,KSPHPDDMType); 315 PETSC_EXTERN PetscErrorCode KSPHPDDMGetType(KSP,KSPHPDDMType*); 316 317 /*E 318 KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 319 320 Level: advanced 321 322 .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 323 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSPGMRESModifiedGramSchmidtOrthogonalization()` 324 325 E*/ 326 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 327 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 328 /*MC 329 KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 330 331 Level: advanced 332 333 Note: Possible unstable, but the fastest to compute 334 335 .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 336 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 337 `KSPGMRESModifiedGramSchmidtOrthogonalization()` 338 M*/ 339 340 /*MC 341 KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 342 iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 343 poor orthogonality. 344 345 Level: advanced 346 347 Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 348 estimate the orthogonality but is more stable. 349 350 .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 351 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_NEVER`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 352 `KSPGMRESModifiedGramSchmidtOrthogonalization()` 353 M*/ 354 355 /*MC 356 KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 357 358 Level: advanced 359 360 Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 361 but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 362 363 You should only use this if you absolutely know that the iterative refinement is needed. 364 365 .seealso: `KSPGMRESClassicalGramSchmidtOrthogonalization()`, `KSPGMRESSetOrthogonalization()`, `KSPGMRESGetOrthogonalization()`, 366 `KSPGMRESSetCGSRefinementType()`, `KSPGMRESGetCGSRefinementType()`, `KSP_GMRES_CGS_REFINE_IFNEEDED`, `KSP_GMRES_CGS_REFINE_ALWAYS`, 367 `KSPGMRESModifiedGramSchmidtOrthogonalization()` 368 M*/ 369 370 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 371 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 372 373 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 374 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 375 PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 376 377 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 378 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 379 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 380 381 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 382 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool); 383 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 384 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 385 386 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 387 PETSC_EXTERN PetscErrorCode KSPResetFromOptions(KSP); 388 PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 389 390 PETSC_EXTERN PetscErrorCode KSPMonitorSetFromOptions(KSP,const char[],const char[],void*); 391 PETSC_EXTERN PetscErrorCode KSPMonitorLGCreate(MPI_Comm,const char[],const char[],const char[],PetscInt,const char *[],int,int,int,int,PetscDrawLG *); 392 PETSC_EXTERN PetscErrorCode KSPMonitorResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 393 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 394 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 395 PETSC_EXTERN PetscErrorCode KSPMonitorResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 396 PETSC_EXTERN PetscErrorCode KSPMonitorResidualShort(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 397 PETSC_EXTERN PetscErrorCode KSPMonitorResidualRange(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 398 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 399 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 400 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 401 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 402 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMax(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 403 PETSC_EXTERN PetscErrorCode KSPMonitorError(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 404 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 405 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 406 PETSC_EXTERN PetscErrorCode KSPMonitorErrorDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 407 PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 408 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDraw(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 409 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 410 PETSC_EXTERN PetscErrorCode KSPMonitorSolutionDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 411 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 412 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValueCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 413 PETSC_DEPRECATED_FUNCTION("Use KSPMonitorResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorDefault(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorResidual(ksp,n,rnorm,vf);} 414 PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidual() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorTrueResidual(ksp,n,rnorm,vf);} 415 PETSC_DEPRECATED_FUNCTION("Use KSPMonitorTrueResidualMax() (since version 3.15)") static inline PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP ksp, PetscInt n, PetscReal rnorm, PetscViewerAndFormat *vf) {return KSPMonitorTrueResidualMax(ksp,n,rnorm,vf);} 416 417 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void*); 418 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 419 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 420 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 421 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 422 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 423 424 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 425 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 426 427 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 428 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 429 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool*,PetscBool*); 430 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 431 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 432 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 433 434 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool); 435 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool*); 436 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool); 437 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool*); 438 439 PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 440 PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 441 PETSC_EXTERN PetscErrorCode KSPViewFromOptions(KSP,PetscObject,const char[]); 442 PETSC_EXTERN PetscErrorCode KSPConvergedReasonView(KSP,PetscViewer); 443 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewSet(KSP,PetscErrorCode (*)(KSP,void*),void *vctx,PetscErrorCode (*)(void**)); 444 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewFromOptions(KSP); 445 PETSC_EXTERN PetscErrorCode KSPConvergedReasonViewCancel(KSP); 446 PETSC_EXTERN PetscErrorCode KSPConvergedRateView(KSP,PetscViewer); 447 448 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonView() (since version 3.14)") static inline PetscErrorCode KSPReasonView(KSP ksp,PetscViewer v) {return KSPConvergedReasonView(ksp,v);} 449 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedReasonViewFromOptions() (since version 3.14)") static inline PetscErrorCode KSPReasonViewFromOptions(KSP ksp) {return KSPConvergedReasonViewFromOptions(ksp);} 450 451 #define KSP_FILE_CLASSID 1211223 452 453 PETSC_EXTERN PetscErrorCode KSPLSQRSetExactMatNorm(KSP,PetscBool); 454 PETSC_EXTERN PetscErrorCode KSPLSQRSetComputeStandardErrorVec(KSP,PetscBool); 455 PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 456 PETSC_EXTERN PetscErrorCode KSPLSQRGetNorms(KSP,PetscReal*,PetscReal*); 457 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidual(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 458 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLG(KSP,PetscInt,PetscReal,PetscViewerAndFormat*); 459 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorResidualDrawLGCreate(PetscViewer,PetscViewerFormat,void *,PetscViewerAndFormat **); 460 461 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 462 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 463 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*); 464 465 /*E 466 KSPNormType - Norm that is passed in the Krylov convergence 467 test routines. 468 469 Level: advanced 470 471 Each solver only supports a subset of these and some may support different ones 472 depending on left or right preconditioning, see KSPSetPCSide() 473 474 Notes: 475 this must match petsc/finclude/petscksp.h 476 477 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetNormType()`, 478 `KSPSetConvergenceTest()`, `KSPSetPCSide()` 479 E*/ 480 typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 481 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 482 PETSC_EXTERN const char *const*const KSPNormTypes; 483 484 /*MC 485 KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 486 possibly save some computation but means the convergence test cannot 487 be based on a norm of a residual etc. 488 489 Level: advanced 490 491 Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 492 493 .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL` 494 M*/ 495 496 /*MC 497 KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 498 convergence test routine. 499 500 Level: advanced 501 502 .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_UNPRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 503 M*/ 504 505 /*MC 506 KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 507 convergence test routine. 508 509 Level: advanced 510 511 .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_NATURAL`, `KSPSetConvergenceTest()` 512 M*/ 513 514 /*MC 515 KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 516 convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG, KSPPIPEFCG, KSPPIPEGCR 517 518 Level: advanced 519 520 .seealso: `KSPNormType`, `KSPSetNormType()`, `KSP_NORM_NONE`, `KSP_NORM_PRECONDITIONED`, `KSP_NORM_UNPRECONDITIONED`, `KSPSetConvergenceTest()` 521 M*/ 522 523 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 524 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 525 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 526 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 527 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 528 529 #define KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED KSP_DIVERGED_PCSETUP_FAILED PETSC_DEPRECATED_ENUM("Use KSP_DIVERGED_PC_FAILED (since version 3.11)") 530 /*E 531 KSPConvergedReason - reason a Krylov method was said to have converged or diverged 532 533 Level: beginner 534 535 Notes: 536 See KSPGetConvergedReason() for explanation of each value 537 538 Developer Notes: 539 this must match petsc/finclude/petscksp.h 540 541 The string versions of these are KSPConvergedReasons; if you change 542 any of the values here also change them that array of names. 543 544 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPSetTolerances()`, `KSPConvergedReasonView()` 545 E*/ 546 typedef enum {/* converged */ 547 KSP_CONVERGED_RTOL_NORMAL = 1, 548 KSP_CONVERGED_ATOL_NORMAL = 9, 549 KSP_CONVERGED_RTOL = 2, 550 KSP_CONVERGED_ATOL = 3, 551 KSP_CONVERGED_ITS = 4, 552 KSP_CONVERGED_CG_NEG_CURVE = 5, 553 KSP_CONVERGED_CG_CONSTRAINED = 6, 554 KSP_CONVERGED_STEP_LENGTH = 7, 555 KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 556 /* diverged */ 557 KSP_DIVERGED_NULL = -2, 558 KSP_DIVERGED_ITS = -3, 559 KSP_DIVERGED_DTOL = -4, 560 KSP_DIVERGED_BREAKDOWN = -5, 561 KSP_DIVERGED_BREAKDOWN_BICG = -6, 562 KSP_DIVERGED_NONSYMMETRIC = -7, 563 KSP_DIVERGED_INDEFINITE_PC = -8, 564 KSP_DIVERGED_NANORINF = -9, 565 KSP_DIVERGED_INDEFINITE_MAT = -10, 566 KSP_DIVERGED_PC_FAILED = -11, 567 KSP_DIVERGED_PCSETUP_FAILED_DEPRECATED = -11, 568 569 KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 570 PETSC_EXTERN const char *const*KSPConvergedReasons; 571 572 /*MC 573 KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) or rtol*norm(b - A*x_0) if KSPConvergedDefaultSetUIRNorm() was called 574 575 Level: beginner 576 577 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 578 for left preconditioning it is the 2-norm of the preconditioned residual, and the 579 2-norm of the residual for right preconditioning 580 581 See also KSP_CONVERGED_ATOL which may apply before this tolerance. 582 583 .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 584 585 M*/ 586 587 /*MC 588 KSP_CONVERGED_ATOL - norm(r) <= atol 589 590 Level: beginner 591 592 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 593 for left preconditioning it is the 2-norm of the preconditioned residual, and the 594 2-norm of the residual for right preconditioning 595 596 See also KSP_CONVERGED_RTOL which may apply before this tolerance. 597 598 .seealso: `KSP_CONVERGED_RTOL`, `KSP_DIVERGED_DTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 599 600 M*/ 601 602 /*MC 603 KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 604 605 Level: beginner 606 607 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 608 for left preconditioning it is the 2-norm of the preconditioned residual, and the 609 2-norm of the residual for right preconditioning 610 611 Level: beginner 612 613 .seealso: `KSP_CONVERGED_ATOL`, `KSP_DIVERGED_RTOL`, `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 614 615 M*/ 616 617 /*MC 618 KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 619 reached 620 621 Level: beginner 622 623 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 624 625 M*/ 626 627 /*MC 628 KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 629 the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 630 test routine is set in KSP. 631 632 Level: beginner 633 634 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 635 636 M*/ 637 638 /*MC 639 KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 640 method could not continue to enlarge the Krylov space. Could be due to a singular matrix or 641 preconditioner. In KSPHPDDM, this is also returned when some search directions within a block 642 are colinear. 643 644 Level: beginner 645 646 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 647 648 M*/ 649 650 /*MC 651 KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 652 method could not continue to enlarge the Krylov space. 653 654 Level: beginner 655 656 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 657 658 M*/ 659 660 /*MC 661 KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 662 symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 663 664 Level: beginner 665 666 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 667 668 M*/ 669 670 /*MC 671 KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 672 positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 673 be positive definite 674 675 Level: beginner 676 677 Notes: 678 This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 679 the PCICC preconditioner to generate a positive definite preconditioner 680 681 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 682 683 M*/ 684 685 /*MC 686 KSP_DIVERGED_PC_FAILED - It was not possible to build or use the requested preconditioner. This is usually due to a 687 zero pivot in a factorization. It can also result from a failure in a subpreconditioner inside a nested preconditioner 688 such as PCFIELDSPLIT. 689 690 Level: beginner 691 692 Notes: 693 Run with -ksp_error_if_not_converged to stop the program when the error is detected and print an error message with details. 694 695 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 696 697 M*/ 698 699 /*MC 700 KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 701 while the KSPSolve() is still running. 702 703 Level: beginner 704 705 .seealso: `KSPSolve()`, `KSPGetConvergedReason()`, `KSPConvergedReason`, `KSPSetTolerances()` 706 707 M*/ 708 709 PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void*,PetscErrorCode (*)(void*)); 710 PETSC_EXTERN PetscErrorCode KSPGetConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*)); 711 PETSC_EXTERN PetscErrorCode KSPGetAndClearConvergenceTest(KSP,PetscErrorCode (**)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void**,PetscErrorCode (**)(void*)); 712 PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void*); 713 PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 714 PETSC_EXTERN PetscErrorCode KSPLSQRConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 715 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void*); 716 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void**); 717 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 718 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 719 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetConvergedMaxits(KSP,PetscBool); 720 PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*); 721 PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason*); 722 PETSC_EXTERN PetscErrorCode KSPGetConvergedReasonString(KSP,const char**); 723 PETSC_EXTERN PetscErrorCode KSPComputeConvergenceRate(KSP,PetscReal*,PetscReal*,PetscReal*,PetscReal*); 724 725 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefault() (since version 3.5)") static inline void KSPDefaultConverged(void) { /* never called */ } 726 #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 727 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultDestroy() (since version 3.5)") static inline void KSPDefaultConvergedDestroy(void) { /* never called */ } 728 #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 729 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultCreate() (since version 3.5)") static inline void KSPDefaultConvergedCreate(void) { /* never called */ } 730 #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 731 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 732 #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 733 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedDefaultSetUMIRNorm() (since version 3.5)") static inline void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 734 #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 735 PETSC_DEPRECATED_FUNCTION("Use KSPConvergedSkip() (since version 3.5)") static inline void KSPSkipConverged(void) { /* never called */ } 736 #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 737 738 PETSC_EXTERN PetscErrorCode KSPComputeOperator(KSP,MatType,Mat*); 739 PETSC_DEPRECATED_FUNCTION("Use KSPComputeOperator() (since version 3.12)") static inline PetscErrorCode KSPComputeExplicitOperator(KSP A,Mat* B) { return KSPComputeOperator(A,NULL,B); } 740 741 /*E 742 KSPCGType - Determines what type of CG to use 743 744 Level: beginner 745 746 .seealso: `KSPCGSetType()` 747 E*/ 748 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 749 PETSC_EXTERN const char *const KSPCGTypes[]; 750 751 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 752 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool); 753 754 PETSC_EXTERN PetscErrorCode KSPCGSetRadius(KSP,PetscReal); 755 PETSC_EXTERN PetscErrorCode KSPCGGetNormD(KSP,PetscReal*); 756 PETSC_EXTERN PetscErrorCode KSPCGGetObjFcn(KSP,PetscReal*); 757 758 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal*); 759 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal*); 760 PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetMinEig (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetMinEig(KSP ksp,PetscReal *x) {return KSPGLTRGetMinEig(ksp,x);} 761 PETSC_DEPRECATED_FUNCTION("Use KSPGLTRGetLambda (since v3.12)") static inline PetscErrorCode KSPCGGLTRGetLambda(KSP ksp,PetscReal *x) {return KSPGLTRGetLambda(ksp,x);} 762 763 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 764 765 PETSC_EXTERN PetscErrorCode PCSetPreSolve(PC,PetscErrorCode (*)(PC,KSP)); 766 PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 767 PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 768 769 #include <petscdrawtypes.h> 770 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 771 772 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 773 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 774 775 /*S 776 KSPGuess - Abstract PETSc object that manages all initial guess methods in Krylov methods. 777 778 Level: beginner 779 780 .seealso: `KSPCreate()`, `KSPSetGuessType()`, `KSPGuessType` 781 S*/ 782 typedef struct _p_KSPGuess* KSPGuess; 783 /*J 784 KSPGuessType - String with the name of a PETSc initial guess for Krylov methods. 785 786 Level: beginner 787 788 .seealso: `KSPGuess` 789 J*/ 790 typedef const char* KSPGuessType; 791 #define KSPGUESSFISCHER "fischer" 792 #define KSPGUESSPOD "pod" 793 PETSC_EXTERN PetscErrorCode KSPGuessRegister(const char[],PetscErrorCode (*)(KSPGuess)); 794 PETSC_EXTERN PetscErrorCode KSPSetGuess(KSP,KSPGuess); 795 PETSC_EXTERN PetscErrorCode KSPGetGuess(KSP,KSPGuess*); 796 PETSC_EXTERN PetscErrorCode KSPGuessView(KSPGuess,PetscViewer); 797 PETSC_EXTERN PetscErrorCode KSPGuessDestroy(KSPGuess*); 798 PETSC_EXTERN PetscErrorCode KSPGuessCreate(MPI_Comm,KSPGuess*); 799 PETSC_EXTERN PetscErrorCode KSPGuessSetType(KSPGuess,KSPGuessType); 800 PETSC_EXTERN PetscErrorCode KSPGuessGetType(KSPGuess,KSPGuessType*); 801 PETSC_EXTERN PetscErrorCode KSPGuessSetTolerance(KSPGuess,PetscReal); 802 PETSC_EXTERN PetscErrorCode KSPGuessSetUp(KSPGuess); 803 PETSC_EXTERN PetscErrorCode KSPGuessUpdate(KSPGuess,Vec,Vec); 804 PETSC_EXTERN PetscErrorCode KSPGuessFormGuess(KSPGuess,Vec,Vec); 805 PETSC_EXTERN PetscErrorCode KSPGuessSetFromOptions(KSPGuess); 806 PETSC_EXTERN PetscErrorCode KSPGuessFischerSetModel(KSPGuess,PetscInt,PetscInt); 807 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 808 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool); 809 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool*); 810 811 /*E 812 MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 813 814 Level: intermediate 815 816 .seealso: `MatSchurComplementGetAinvType()`, `MatSchurComplementSetAinvType()`, `MatSchurComplementGetPmat()`, `MatGetSchurComplement()`, `MatCreateSchurComplementPmat()` 817 E*/ 818 typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP, MAT_SCHUR_COMPLEMENT_AINV_BLOCK_DIAG} MatSchurComplementAinvType; 819 PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 820 821 typedef enum {MAT_LMVM_SYMBROYDEN_SCALE_NONE = 0, 822 MAT_LMVM_SYMBROYDEN_SCALE_SCALAR = 1, 823 MAT_LMVM_SYMBROYDEN_SCALE_DIAGONAL = 2, 824 MAT_LMVM_SYMBROYDEN_SCALE_USER = 3} MatLMVMSymBroydenScaleType; 825 PETSC_EXTERN const char *const MatLMVMSymBroydenScaleTypes[]; 826 827 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 828 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 829 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 830 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 831 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 832 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 833 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 834 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 835 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 836 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 837 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat*,MatSchurComplementAinvType,MatReuse,Mat*); 838 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 839 840 PETSC_EXTERN PetscErrorCode MatCreateLMVMDFP(MPI_Comm,PetscInt,PetscInt,Mat*); 841 PETSC_EXTERN PetscErrorCode MatCreateLMVMBFGS(MPI_Comm,PetscInt,PetscInt,Mat*); 842 PETSC_EXTERN PetscErrorCode MatCreateLMVMSR1(MPI_Comm,PetscInt,PetscInt,Mat*); 843 PETSC_EXTERN PetscErrorCode MatCreateLMVMBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 844 PETSC_EXTERN PetscErrorCode MatCreateLMVMBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 845 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 846 PETSC_EXTERN PetscErrorCode MatCreateLMVMSymBadBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 847 PETSC_EXTERN PetscErrorCode MatCreateLMVMDiagBroyden(MPI_Comm,PetscInt,PetscInt,Mat*); 848 849 PETSC_EXTERN PetscErrorCode MatLMVMUpdate(Mat, Vec, Vec); 850 PETSC_EXTERN PetscErrorCode MatLMVMIsAllocated(Mat, PetscBool*); 851 PETSC_EXTERN PetscErrorCode MatLMVMAllocate(Mat, Vec, Vec); 852 PETSC_EXTERN PetscErrorCode MatLMVMReset(Mat, PetscBool); 853 PETSC_EXTERN PetscErrorCode MatLMVMResetShift(Mat); 854 PETSC_EXTERN PetscErrorCode MatLMVMClearJ0(Mat); 855 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0(Mat, Mat); 856 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Scale(Mat, PetscReal); 857 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0Diag(Mat, Vec); 858 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0PC(Mat, PC); 859 PETSC_EXTERN PetscErrorCode MatLMVMSetJ0KSP(Mat, KSP); 860 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Fwd(Mat, Vec, Vec); 861 PETSC_EXTERN PetscErrorCode MatLMVMApplyJ0Inv(Mat, Vec, Vec); 862 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0(Mat, Mat*); 863 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0PC(Mat, PC*); 864 PETSC_EXTERN PetscErrorCode MatLMVMGetJ0KSP(Mat, KSP*); 865 PETSC_EXTERN PetscErrorCode MatLMVMSetHistorySize(Mat, PetscInt); 866 PETSC_EXTERN PetscErrorCode MatLMVMGetUpdateCount(Mat, PetscInt*); 867 PETSC_EXTERN PetscErrorCode MatLMVMGetRejectCount(Mat, PetscInt*); 868 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetDelta(Mat, PetscScalar); 869 PETSC_EXTERN PetscErrorCode MatLMVMSymBroydenSetScaleType(Mat, MatLMVMSymBroydenScaleType); 870 871 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 872 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool); 873 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 874 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 875 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 876 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void*); 877 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 878 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 879 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 880 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 881 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 882 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 883 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 884 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 885 886 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 887 PETSC_EXTERN PetscErrorCode DMProjectField(DM, PetscReal, Vec, 888 void (**)(PetscInt, PetscInt, PetscInt, 889 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 890 const PetscInt[], const PetscInt[], const PetscScalar[], const PetscScalar[], const PetscScalar[], 891 PetscReal, const PetscReal [], PetscInt, const PetscScalar[], PetscScalar []), InsertMode, Vec); 892 893 PETSC_EXTERN PetscErrorCode DMAdaptInterpolator(DM, DM, Mat, KSP, Mat, Mat, Mat *, void *); 894 PETSC_EXTERN PetscErrorCode DMCheckInterpolator(DM, Mat, Mat, Mat, PetscReal); 895 #endif 896