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