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