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