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