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