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 Concepts: Krylov methods 17 18 Notes: When a direct solver is used but no Krylov solver is used the KSP object is still used by with a 19 KSPType of KSPPREONLY (meaning application of the preconditioner is only used as the linear solver). 20 21 .seealso: KSPCreate(), KSPSetType(), KSPType, SNES, TS, PC, KSP, KSPDestroy() 22 S*/ 23 typedef struct _p_KSP* KSP; 24 25 /*J 26 KSPType - String with the name of a PETSc Krylov method. 27 28 Level: beginner 29 30 .seealso: KSPSetType(), KSP, KSPRegister(), KSPCreate(), KSPSetFromOptions() 31 J*/ 32 typedef const char* KSPType; 33 #define KSPRICHARDSON "richardson" 34 #define KSPCHEBYSHEV "chebyshev" 35 #define KSPCG "cg" 36 #define KSPGROPPCG "groppcg" 37 #define KSPPIPECG "pipecg" 38 #define KSPCGNE "cgne" 39 #define KSPNASH "nash" 40 #define KSPSTCG "stcg" 41 #define KSPGLTR "gltr" 42 #define KSPFCG "fcg" 43 #define KSPGMRES "gmres" 44 #define KSPFGMRES "fgmres" 45 #define KSPLGMRES "lgmres" 46 #define KSPDGMRES "dgmres" 47 #define KSPPGMRES "pgmres" 48 #define KSPTCQMR "tcqmr" 49 #define KSPBCGS "bcgs" 50 #define KSPIBCGS "ibcgs" 51 #define KSPFBCGS "fbcgs" 52 #define KSPFBCGSR "fbcgsr" 53 #define KSPBCGSL "bcgsl" 54 #define KSPCGS "cgs" 55 #define KSPTFQMR "tfqmr" 56 #define KSPCR "cr" 57 #define KSPPIPECR "pipecr" 58 #define KSPLSQR "lsqr" 59 #define KSPPREONLY "preonly" 60 #define KSPQCG "qcg" 61 #define KSPBICG "bicg" 62 #define KSPMINRES "minres" 63 #define KSPSYMMLQ "symmlq" 64 #define KSPLCD "lcd" 65 #define KSPPYTHON "python" 66 #define KSPGCR "gcr" 67 68 /* Logging support */ 69 PETSC_EXTERN PetscClassId KSP_CLASSID; 70 PETSC_EXTERN PetscClassId DMKSP_CLASSID; 71 72 PETSC_EXTERN PetscErrorCode KSPCreate(MPI_Comm,KSP *); 73 PETSC_EXTERN PetscErrorCode KSPSetType(KSP,KSPType); 74 PETSC_EXTERN PetscErrorCode KSPGetType(KSP,KSPType *); 75 PETSC_EXTERN PetscErrorCode KSPSetUp(KSP); 76 PETSC_EXTERN PetscErrorCode KSPSetUpOnBlocks(KSP); 77 PETSC_EXTERN PetscErrorCode KSPSolve(KSP,Vec,Vec); 78 PETSC_EXTERN PetscErrorCode KSPSolveTranspose(KSP,Vec,Vec); 79 PETSC_EXTERN PetscErrorCode KSPReset(KSP); 80 PETSC_EXTERN PetscErrorCode KSPDestroy(KSP*); 81 PETSC_EXTERN PetscErrorCode KSPSetReusePreconditioner(KSP,PetscBool); 82 PETSC_EXTERN PetscErrorCode KSPSetSkipPCSetFromOptions(KSP,PetscBool); 83 84 PETSC_EXTERN PetscFunctionList KSPList; 85 PETSC_EXTERN PetscErrorCode KSPRegister(const char[],PetscErrorCode (*)(KSP)); 86 87 PETSC_EXTERN PetscErrorCode KSPSetPCSide(KSP,PCSide); 88 PETSC_EXTERN PetscErrorCode KSPGetPCSide(KSP,PCSide*); 89 PETSC_EXTERN PetscErrorCode KSPSetTolerances(KSP,PetscReal,PetscReal,PetscReal,PetscInt); 90 PETSC_EXTERN PetscErrorCode KSPGetTolerances(KSP,PetscReal*,PetscReal*,PetscReal*,PetscInt*); 91 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessNonzero(KSP,PetscBool ); 92 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessNonzero(KSP,PetscBool *); 93 PETSC_EXTERN PetscErrorCode KSPSetInitialGuessKnoll(KSP,PetscBool ); 94 PETSC_EXTERN PetscErrorCode KSPGetInitialGuessKnoll(KSP,PetscBool *); 95 PETSC_EXTERN PetscErrorCode KSPSetErrorIfNotConverged(KSP,PetscBool ); 96 PETSC_EXTERN PetscErrorCode KSPGetErrorIfNotConverged(KSP,PetscBool *); 97 PETSC_EXTERN PetscErrorCode KSPSetComputeEigenvalues(KSP,PetscBool ); 98 PETSC_EXTERN PetscErrorCode KSPGetComputeEigenvalues(KSP,PetscBool *); 99 PETSC_EXTERN PetscErrorCode KSPSetComputeSingularValues(KSP,PetscBool ); 100 PETSC_EXTERN PetscErrorCode KSPGetComputeSingularValues(KSP,PetscBool *); 101 PETSC_EXTERN PetscErrorCode KSPGetRhs(KSP,Vec *); 102 PETSC_EXTERN PetscErrorCode KSPGetSolution(KSP,Vec *); 103 PETSC_EXTERN PetscErrorCode KSPGetResidualNorm(KSP,PetscReal*); 104 PETSC_EXTERN PetscErrorCode KSPGetIterationNumber(KSP,PetscInt*); 105 PETSC_EXTERN PetscErrorCode KSPGetTotalIterations(KSP,PetscInt*); 106 PETSC_EXTERN PetscErrorCode KSPCreateVecs(KSP,PetscInt,Vec**,PetscInt,Vec**); 107 PETSC_DEPRECATED("Use KSPCreateVecs()") PETSC_STATIC_INLINE PetscErrorCode KSPGetVecs(KSP ksp,PetscInt n,Vec **x,PetscInt m,Vec **y) {return KSPCreateVecs(ksp,n,x,m,y);} 108 109 PETSC_EXTERN PetscErrorCode KSPSetPreSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 110 PETSC_EXTERN PetscErrorCode KSPSetPostSolve(KSP,PetscErrorCode (*)(KSP,Vec,Vec,void*),void*); 111 112 PETSC_EXTERN PetscErrorCode KSPSetPC(KSP,PC); 113 PETSC_EXTERN PetscErrorCode KSPGetPC(KSP,PC*); 114 115 PETSC_EXTERN PetscErrorCode KSPMonitor(KSP,PetscInt,PetscReal); 116 PETSC_EXTERN PetscErrorCode KSPMonitorSet(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *,PetscErrorCode (*)(void**)); 117 PETSC_EXTERN PetscErrorCode KSPMonitorCancel(KSP); 118 PETSC_EXTERN PetscErrorCode KSPGetMonitorContext(KSP,void **); 119 PETSC_EXTERN PetscErrorCode KSPGetResidualHistory(KSP,PetscReal*[],PetscInt *); 120 PETSC_EXTERN PetscErrorCode KSPSetResidualHistory(KSP,PetscReal[],PetscInt,PetscBool ); 121 122 PETSC_EXTERN PetscErrorCode KSPBuildSolutionDefault(KSP,Vec,Vec*); 123 PETSC_EXTERN PetscErrorCode KSPBuildResidualDefault(KSP,Vec,Vec,Vec *); 124 PETSC_EXTERN PetscErrorCode KSPDestroyDefault(KSP); 125 PETSC_EXTERN PetscErrorCode KSPSetWorkVecs(KSP,PetscInt); 126 127 PETSC_EXTERN PetscErrorCode PCKSPGetKSP(PC,KSP*); 128 PETSC_EXTERN PetscErrorCode PCBJacobiGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 129 PETSC_EXTERN PetscErrorCode PCASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 130 PETSC_EXTERN PetscErrorCode PCGASMGetSubKSP(PC,PetscInt*,PetscInt*,KSP*[]); 131 PETSC_EXTERN PetscErrorCode PCFieldSplitGetSubKSP(PC,PetscInt*,KSP*[]); 132 PETSC_EXTERN PetscErrorCode PCMGGetSmoother(PC,PetscInt,KSP*); 133 PETSC_EXTERN PetscErrorCode PCMGGetSmootherDown(PC,PetscInt,KSP*); 134 PETSC_EXTERN PetscErrorCode PCMGGetSmootherUp(PC,PetscInt,KSP*); 135 PETSC_EXTERN PetscErrorCode PCMGGetCoarseSolve(PC,KSP*); 136 PETSC_EXTERN PetscErrorCode PCGalerkinGetKSP(PC,KSP *); 137 138 PETSC_EXTERN PetscErrorCode KSPBuildSolution(KSP,Vec,Vec *); 139 PETSC_EXTERN PetscErrorCode KSPBuildResidual(KSP,Vec,Vec,Vec *); 140 141 PETSC_EXTERN PetscErrorCode KSPRichardsonSetScale(KSP,PetscReal); 142 PETSC_EXTERN PetscErrorCode KSPRichardsonSetSelfScale(KSP,PetscBool ); 143 PETSC_EXTERN PetscErrorCode KSPChebyshevSetEigenvalues(KSP,PetscReal,PetscReal); 144 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSet(KSP,PetscReal,PetscReal,PetscReal,PetscReal); 145 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetUseRandom(KSP,PetscBool); 146 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigSetRandom(KSP,PetscRandom); 147 PETSC_EXTERN PetscErrorCode KSPChebyshevEstEigGetKSP(KSP,KSP*); 148 PETSC_EXTERN PetscErrorCode KSPComputeExtremeSingularValues(KSP,PetscReal*,PetscReal*); 149 PETSC_EXTERN PetscErrorCode KSPComputeEigenvalues(KSP,PetscInt,PetscReal[],PetscReal[],PetscInt *); 150 PETSC_EXTERN PetscErrorCode KSPComputeEigenvaluesExplicitly(KSP,PetscInt,PetscReal[],PetscReal[]); 151 152 /*E 153 154 KSPFCGTruncationType - Define how stored directions are used to orthogonalize in FCG 155 156 KSP_FCG_TRUNC_TYPE_STANDARD uses all (up to mmax) stored directions 157 KSP_FCG_TRUNC_TYPE_NOTAY uses the last max(1,mod(i,mmax)) stored directions at iteration i=0,1.. 158 159 Level: intermediate 160 161 .seealso : KSPFCG,KSPFCGSetTruncationType(),KSPFCGGetTruncationType() 162 163 E*/ 164 typedef enum {KSP_FCG_TRUNC_TYPE_STANDARD,KSP_FCG_TRUNC_TYPE_NOTAY} KSPFCGTruncationType; 165 PETSC_EXTERN const char *const KSPFCGTruncationTypes[]; 166 167 PETSC_EXTERN PetscErrorCode KSPFCGSetMmax(KSP,PetscInt); 168 PETSC_EXTERN PetscErrorCode KSPFCGGetMmax(KSP,PetscInt*); 169 PETSC_EXTERN PetscErrorCode KSPFCGSetNprealloc(KSP,PetscInt); 170 PETSC_EXTERN PetscErrorCode KSPFCGGetNprealloc(KSP,PetscInt*); 171 PETSC_EXTERN PetscErrorCode KSPFCGSetTruncationType(KSP,KSPFCGTruncationType); 172 PETSC_EXTERN PetscErrorCode KSPFCGGetTruncationType(KSP,KSPFCGTruncationType*); 173 174 PETSC_EXTERN PetscErrorCode KSPGMRESSetRestart(KSP, PetscInt); 175 PETSC_EXTERN PetscErrorCode KSPGMRESGetRestart(KSP, PetscInt*); 176 PETSC_EXTERN PetscErrorCode KSPGMRESSetHapTol(KSP,PetscReal); 177 178 PETSC_EXTERN PetscErrorCode KSPGMRESSetPreAllocateVectors(KSP); 179 PETSC_EXTERN PetscErrorCode KSPGMRESSetOrthogonalization(KSP,PetscErrorCode (*)(KSP,PetscInt)); 180 PETSC_EXTERN PetscErrorCode KSPGMRESGetOrthogonalization(KSP,PetscErrorCode (**)(KSP,PetscInt)); 181 PETSC_EXTERN PetscErrorCode KSPGMRESModifiedGramSchmidtOrthogonalization(KSP,PetscInt); 182 PETSC_EXTERN PetscErrorCode KSPGMRESClassicalGramSchmidtOrthogonalization(KSP,PetscInt); 183 184 PETSC_EXTERN PetscErrorCode KSPLGMRESSetAugDim(KSP,PetscInt); 185 PETSC_EXTERN PetscErrorCode KSPLGMRESSetConstant(KSP); 186 187 PETSC_EXTERN PetscErrorCode KSPGCRSetRestart(KSP,PetscInt); 188 PETSC_EXTERN PetscErrorCode KSPGCRGetRestart(KSP,PetscInt*); 189 PETSC_EXTERN PetscErrorCode KSPGCRSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 190 191 /*E 192 KSPGMRESCGSRefinementType - How the classical (unmodified) Gram-Schmidt is performed. 193 194 Level: advanced 195 196 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 197 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSPGMRESModifiedGramSchmidtOrthogonalization() 198 199 E*/ 200 typedef enum {KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS} KSPGMRESCGSRefinementType; 201 PETSC_EXTERN const char *const KSPGMRESCGSRefinementTypes[]; 202 /*MC 203 KSP_GMRES_CGS_REFINE_NEVER - Do the classical (unmodified) Gram-Schmidt process 204 205 Level: advanced 206 207 Note: Possible unstable, but the fastest to compute 208 209 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 210 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 211 KSPGMRESModifiedGramSchmidtOrthogonalization() 212 M*/ 213 214 /*MC 215 KSP_GMRES_CGS_REFINE_IFNEEDED - Do the classical (unmodified) Gram-Schmidt process and one step of 216 iterative refinement if an estimate of the orthogonality of the resulting vectors indicates 217 poor orthogonality. 218 219 Level: advanced 220 221 Note: This is slower than KSP_GMRES_CGS_REFINE_NEVER because it requires an extra norm computation to 222 estimate the orthogonality but is more stable. 223 224 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 225 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_NEVER, KSP_GMRES_CGS_REFINE_ALWAYS, 226 KSPGMRESModifiedGramSchmidtOrthogonalization() 227 M*/ 228 229 /*MC 230 KSP_GMRES_CGS_REFINE_NEVER - Do two steps of the classical (unmodified) Gram-Schmidt process. 231 232 Level: advanced 233 234 Note: This is roughly twice the cost of KSP_GMRES_CGS_REFINE_NEVER because it performs the process twice 235 but it saves the extra norm calculation needed by KSP_GMRES_CGS_REFINE_IFNEEDED. 236 237 You should only use this if you absolutely know that the iterative refinement is needed. 238 239 .seealso: KSPGMRESClassicalGramSchmidtOrthogonalization(), KSPGMRESSetOrthogonalization(), KSPGMRESGetOrthogonalization(), 240 KSPGMRESSetCGSRefinementType(), KSPGMRESGetCGSRefinementType(), KSP_GMRES_CGS_REFINE_IFNEEDED, KSP_GMRES_CGS_REFINE_ALWAYS, 241 KSPGMRESModifiedGramSchmidtOrthogonalization() 242 M*/ 243 244 PETSC_EXTERN PetscErrorCode KSPGMRESSetCGSRefinementType(KSP,KSPGMRESCGSRefinementType); 245 PETSC_EXTERN PetscErrorCode KSPGMRESGetCGSRefinementType(KSP,KSPGMRESCGSRefinementType*); 246 247 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCNoChange(KSP,PetscInt,PetscInt,PetscReal,void*); 248 PETSC_EXTERN PetscErrorCode KSPFGMRESModifyPCKSP(KSP,PetscInt,PetscInt,PetscReal,void*); 249 PETSC_EXTERN PetscErrorCode KSPFGMRESSetModifyPC(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscInt,PetscReal,void*),void*,PetscErrorCode(*)(void*)); 250 251 PETSC_EXTERN PetscErrorCode KSPQCGSetTrustRegionRadius(KSP,PetscReal); 252 PETSC_EXTERN PetscErrorCode KSPQCGGetQuadratic(KSP,PetscReal*); 253 PETSC_EXTERN PetscErrorCode KSPQCGGetTrialStepNorm(KSP,PetscReal*); 254 255 PETSC_EXTERN PetscErrorCode KSPBCGSLSetXRes(KSP,PetscReal); 256 PETSC_EXTERN PetscErrorCode KSPBCGSLSetPol(KSP,PetscBool ); 257 PETSC_EXTERN PetscErrorCode KSPBCGSLSetEll(KSP,PetscInt); 258 PETSC_EXTERN PetscErrorCode KSPBCGSLSetUsePseudoinverse(KSP,PetscBool); 259 260 PETSC_EXTERN PetscErrorCode KSPSetFromOptions(KSP); 261 PETSC_EXTERN PetscErrorCode KSPAddOptionsChecker(PetscErrorCode (*)(KSP)); 262 263 PETSC_EXTERN PetscErrorCode KSPMonitorSingularValue(KSP,PetscInt,PetscReal,void *); 264 PETSC_EXTERN PetscErrorCode KSPMonitorDefault(KSP,PetscInt,PetscReal,void *); 265 PETSC_EXTERN PetscErrorCode KSPLSQRMonitorDefault(KSP,PetscInt,PetscReal,void *); 266 PETSC_EXTERN PetscErrorCode KSPMonitorRange(KSP,PetscInt,PetscReal,void *); 267 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicTolerance(KSP ksp,PetscInt its,PetscReal fnorm,void *dummy); 268 PETSC_EXTERN PetscErrorCode KSPMonitorDynamicToleranceDestroy(void **dummy); 269 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualNorm(KSP,PetscInt,PetscReal,void *); 270 PETSC_EXTERN PetscErrorCode KSPMonitorTrueResidualMaxNorm(KSP,PetscInt,PetscReal,void *); 271 PETSC_EXTERN PetscErrorCode KSPMonitorDefaultShort(KSP,PetscInt,PetscReal,void *); 272 PETSC_EXTERN PetscErrorCode KSPMonitorSolution(KSP,PetscInt,PetscReal,void *); 273 PETSC_EXTERN PetscErrorCode KSPMonitorSAWs(KSP,PetscInt,PetscReal,void*); 274 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsCreate(KSP,void**); 275 PETSC_EXTERN PetscErrorCode KSPMonitorSAWsDestroy(void**); 276 PETSC_EXTERN PetscErrorCode KSPGMRESMonitorKrylov(KSP,PetscInt,PetscReal,void *); 277 278 PETSC_EXTERN PetscErrorCode KSPUnwindPreconditioner(KSP,Vec,Vec); 279 PETSC_EXTERN PetscErrorCode KSPInitialResidual(KSP,Vec,Vec,Vec,Vec,Vec); 280 281 PETSC_EXTERN PetscErrorCode KSPSetOperators(KSP,Mat,Mat); 282 PETSC_EXTERN PetscErrorCode KSPGetOperators(KSP,Mat*,Mat*); 283 PETSC_EXTERN PetscErrorCode KSPGetOperatorsSet(KSP,PetscBool *,PetscBool *); 284 PETSC_EXTERN PetscErrorCode KSPSetOptionsPrefix(KSP,const char[]); 285 PETSC_EXTERN PetscErrorCode KSPAppendOptionsPrefix(KSP,const char[]); 286 PETSC_EXTERN PetscErrorCode KSPGetOptionsPrefix(KSP,const char*[]); 287 PETSC_EXTERN PetscErrorCode KSPSetTabLevel(KSP,PetscInt); 288 PETSC_EXTERN PetscErrorCode KSPGetTabLevel(KSP,PetscInt*); 289 290 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScale(KSP,PetscBool ); 291 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScale(KSP,PetscBool *); 292 PETSC_EXTERN PetscErrorCode KSPSetDiagonalScaleFix(KSP,PetscBool ); 293 PETSC_EXTERN PetscErrorCode KSPGetDiagonalScaleFix(KSP,PetscBool *); 294 295 PETSC_EXTERN PetscErrorCode KSPView(KSP,PetscViewer); 296 PETSC_EXTERN PetscErrorCode KSPLoad(KSP,PetscViewer); 297 PETSC_STATIC_INLINE PetscErrorCode KSPViewFromOptions(KSP A,PetscObject obj,const char name[]) {return PetscObjectViewFromOptions((PetscObject)A,obj,name);} 298 PETSC_EXTERN PetscErrorCode KSPReasonView(KSP,PetscViewer); 299 PETSC_EXTERN PetscErrorCode KSPReasonViewFromOptions(KSP); 300 301 #define KSP_FILE_CLASSID 1211223 302 303 PETSC_EXTERN PetscErrorCode KSPLSQRSetStandardErrorVec(KSP,Vec); 304 PETSC_EXTERN PetscErrorCode KSPLSQRGetStandardErrorVec(KSP,Vec*); 305 306 PETSC_EXTERN PetscErrorCode PCRedundantGetKSP(PC,KSP*); 307 PETSC_EXTERN PetscErrorCode PCRedistributeGetKSP(PC,KSP*); 308 PETSC_EXTERN PetscErrorCode PCTelescopeGetKSP(PC,KSP*); 309 310 /*E 311 KSPNormType - Norm that is passed in the Krylov convergence 312 test routines. 313 314 Level: advanced 315 316 Each solver only supports a subset of these and some may support different ones 317 depending on left or right preconditioning, see KSPSetPCSide() 318 319 Notes: this must match petsc/finclude/petscksp.h 320 321 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetNormType(), 322 KSPSetConvergenceTest(), KSPSetPCSide() 323 E*/ 324 typedef enum {KSP_NORM_DEFAULT = -1,KSP_NORM_NONE = 0,KSP_NORM_PRECONDITIONED = 1,KSP_NORM_UNPRECONDITIONED = 2,KSP_NORM_NATURAL = 3} KSPNormType; 325 #define KSP_NORM_MAX (KSP_NORM_NATURAL + 1) 326 PETSC_EXTERN const char *const*const KSPNormTypes; 327 328 /*MC 329 KSP_NORM_NONE - Do not compute a norm during the Krylov process. This will 330 possibly save some computation but means the convergence test cannot 331 be based on a norm of a residual etc. 332 333 Level: advanced 334 335 Note: Some Krylov methods need to compute a residual norm (such as GMRES) and then this option is ignored 336 337 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL 338 M*/ 339 340 /*MC 341 KSP_NORM_PRECONDITIONED - Compute the norm of the preconditioned residual B*(b - A*x), if left preconditioning, and pass that to the 342 convergence test routine. 343 344 Level: advanced 345 346 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_UNPRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 347 M*/ 348 349 /*MC 350 KSP_NORM_UNPRECONDITIONED - Compute the norm of the true residual (b - A*x) and pass that to the 351 convergence test routine. 352 353 Level: advanced 354 355 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_NATURAL, KSPSetConvergenceTest() 356 M*/ 357 358 /*MC 359 KSP_NORM_NATURAL - Compute the 'natural norm' of residual sqrt((b - A*x)*B*(b - A*x)) and pass that to the 360 convergence test routine. This is only supported by KSPCG, KSPCR, KSPCGNE, KSPCGS, KSPFCG 361 362 Level: advanced 363 364 .seealso: KSPNormType, KSPSetNormType(), KSP_NORM_NONE, KSP_NORM_PRECONDITIONED, KSP_NORM_UNPRECONDITIONED, KSPSetConvergenceTest() 365 M*/ 366 367 PETSC_EXTERN PetscErrorCode KSPSetNormType(KSP,KSPNormType); 368 PETSC_EXTERN PetscErrorCode KSPGetNormType(KSP,KSPNormType*); 369 PETSC_EXTERN PetscErrorCode KSPSetSupportedNorm(KSP ksp,KSPNormType,PCSide,PetscInt); 370 PETSC_EXTERN PetscErrorCode KSPSetCheckNormIteration(KSP,PetscInt); 371 PETSC_EXTERN PetscErrorCode KSPSetLagNorm(KSP,PetscBool); 372 373 /*E 374 KSPConvergedReason - reason a Krylov method was said to have converged or diverged 375 376 Level: beginner 377 378 Notes: See KSPGetConvergedReason() for explanation of each value 379 380 Developer notes: this must match petsc/finclude/petscksp.h 381 382 The string versions of these are KSPConvergedReasons; if you change 383 any of the values here also change them that array of names. 384 385 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPSetTolerances() 386 E*/ 387 typedef enum {/* converged */ 388 KSP_CONVERGED_RTOL_NORMAL = 1, 389 KSP_CONVERGED_ATOL_NORMAL = 9, 390 KSP_CONVERGED_RTOL = 2, 391 KSP_CONVERGED_ATOL = 3, 392 KSP_CONVERGED_ITS = 4, 393 KSP_CONVERGED_CG_NEG_CURVE = 5, 394 KSP_CONVERGED_CG_CONSTRAINED = 6, 395 KSP_CONVERGED_STEP_LENGTH = 7, 396 KSP_CONVERGED_HAPPY_BREAKDOWN = 8, 397 /* diverged */ 398 KSP_DIVERGED_NULL = -2, 399 KSP_DIVERGED_ITS = -3, 400 KSP_DIVERGED_DTOL = -4, 401 KSP_DIVERGED_BREAKDOWN = -5, 402 KSP_DIVERGED_BREAKDOWN_BICG = -6, 403 KSP_DIVERGED_NONSYMMETRIC = -7, 404 KSP_DIVERGED_INDEFINITE_PC = -8, 405 KSP_DIVERGED_NANORINF = -9, 406 KSP_DIVERGED_INDEFINITE_MAT = -10, 407 KSP_DIVERGED_PCSETUP_FAILED = -11, 408 409 KSP_CONVERGED_ITERATING = 0} KSPConvergedReason; 410 PETSC_EXTERN const char *const*KSPConvergedReasons; 411 412 /*MC 413 KSP_CONVERGED_RTOL - norm(r) <= rtol*norm(b) 414 415 Level: beginner 416 417 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 418 for left preconditioning it is the 2-norm of the preconditioned residual, and the 419 2-norm of the residual for right preconditioning 420 421 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 422 423 M*/ 424 425 /*MC 426 KSP_CONVERGED_ATOL - norm(r) <= atol 427 428 Level: beginner 429 430 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 431 for left preconditioning it is the 2-norm of the preconditioned residual, and the 432 2-norm of the residual for right preconditioning 433 434 Level: beginner 435 436 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 437 438 M*/ 439 440 /*MC 441 KSP_DIVERGED_DTOL - norm(r) >= dtol*norm(b) 442 443 Level: beginner 444 445 See KSPNormType and KSPSetNormType() for possible norms that may be used. By default 446 for left preconditioning it is the 2-norm of the preconditioned residual, and the 447 2-norm of the residual for right preconditioning 448 449 Level: beginner 450 451 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 452 453 M*/ 454 455 /*MC 456 KSP_DIVERGED_ITS - Ran out of iterations before any convergence criteria was 457 reached 458 459 Level: beginner 460 461 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 462 463 M*/ 464 465 /*MC 466 KSP_CONVERGED_ITS - Used by the KSPPREONLY solver after the single iteration of 467 the preconditioner is applied. Also used when the KSPConvergedSkip() convergence 468 test routine is set in KSP. 469 470 Level: beginner 471 472 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 473 474 M*/ 475 476 /*MC 477 KSP_DIVERGED_BREAKDOWN - A breakdown in the Krylov method was detected so the 478 method could not continue to enlarge the Krylov space. Could be due to a singlular matrix or 479 preconditioner. 480 481 Level: beginner 482 483 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 484 485 M*/ 486 487 /*MC 488 KSP_DIVERGED_BREAKDOWN_BICG - A breakdown in the KSPBICG method was detected so the 489 method could not continue to enlarge the Krylov space. 490 491 Level: beginner 492 493 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 494 495 M*/ 496 497 /*MC 498 KSP_DIVERGED_NONSYMMETRIC - It appears the operator or preconditioner is not 499 symmetric and this Krylov method (KSPCG, KSPMINRES, KSPCR) requires symmetry 500 501 Level: beginner 502 503 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 504 505 M*/ 506 507 /*MC 508 KSP_DIVERGED_INDEFINITE_PC - It appears the preconditioner is indefinite (has both 509 positive and negative eigenvalues) and this Krylov method (KSPCG) requires it to 510 be positive definite 511 512 Level: beginner 513 514 Notes: This can happen with the PCICC preconditioner, use -pc_factor_shift_positive_definite to force 515 the PCICC preconditioner to generate a positive definite preconditioner 516 517 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 518 519 M*/ 520 521 /*MC 522 KSP_CONVERGED_ITERATING - This flag is returned if you call KSPGetConvergedReason() 523 while the KSPSolve() is still running. 524 525 Level: beginner 526 527 .seealso: KSPSolve(), KSPGetConvergedReason(), KSPConvergedReason, KSPSetTolerances() 528 529 M*/ 530 531 PETSC_EXTERN PetscErrorCode KSPSetConvergenceTest(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,KSPConvergedReason*,void*),void *,PetscErrorCode (*)(void*)); 532 PETSC_EXTERN PetscErrorCode KSPGetConvergenceContext(KSP,void **); 533 PETSC_EXTERN PetscErrorCode KSPConvergedDefault(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 534 PETSC_EXTERN PetscErrorCode KSPConvergedLSQR(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 535 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultDestroy(void *); 536 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultCreate(void **); 537 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUIRNorm(KSP); 538 PETSC_EXTERN PetscErrorCode KSPConvergedDefaultSetUMIRNorm(KSP); 539 PETSC_EXTERN PetscErrorCode KSPConvergedSkip(KSP,PetscInt,PetscReal,KSPConvergedReason*,void *); 540 PETSC_EXTERN PetscErrorCode KSPGetConvergedReason(KSP,KSPConvergedReason *); 541 542 PETSC_DEPRECATED("Use KSPConvergedDefault()") PETSC_STATIC_INLINE void KSPDefaultConverged(void) { /* never called */ } 543 #define KSPDefaultConverged (KSPDefaultConverged, KSPConvergedDefault) 544 PETSC_DEPRECATED("Use KSPConvergedDefaultDestroy()") PETSC_STATIC_INLINE void KSPDefaultConvergedDestroy(void) { /* never called */ } 545 #define KSPDefaultConvergedDestroy (KSPDefaultConvergedDestroy, KSPConvergedDefaultDestroy) 546 PETSC_DEPRECATED("Use KSPConvergedDefaultCreate()") PETSC_STATIC_INLINE void KSPDefaultConvergedCreate(void) { /* never called */ } 547 #define KSPDefaultConvergedCreate (KSPDefaultConvergedCreate, KSPConvergedDefaultCreate) 548 PETSC_DEPRECATED("Use KSPConvergedDefaultSetUIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUIRNorm(void) { /* never called */ } 549 #define KSPDefaultConvergedSetUIRNorm (KSPDefaultConvergedSetUIRNorm, KSPConvergedDefaultSetUIRNorm) 550 PETSC_DEPRECATED("Use KSPConvergedDefaultSetUMIRNorm()") PETSC_STATIC_INLINE void KSPDefaultConvergedSetUMIRNorm(void) { /* never called */ } 551 #define KSPDefaultConvergedSetUMIRNorm (KSPDefaultConvergedSetUMIRNorm, KSPConvergedDefaultSetUMIRNorm) 552 PETSC_DEPRECATED("Use KSPConvergedSkip()") PETSC_STATIC_INLINE void KSPSkipConverged(void) { /* never called */ } 553 #define KSPSkipConverged (KSPSkipConverged, KSPConvergedSkip) 554 555 PETSC_EXTERN PetscErrorCode KSPComputeExplicitOperator(KSP,Mat *); 556 557 /*E 558 KSPCGType - Determines what type of CG to use 559 560 Level: beginner 561 562 .seealso: KSPCGSetType() 563 E*/ 564 typedef enum {KSP_CG_SYMMETRIC=0,KSP_CG_HERMITIAN=1} KSPCGType; 565 PETSC_EXTERN const char *const KSPCGTypes[]; 566 567 PETSC_EXTERN PetscErrorCode KSPCGSetType(KSP,KSPCGType); 568 PETSC_EXTERN PetscErrorCode KSPCGUseSingleReduction(KSP,PetscBool ); 569 570 PETSC_EXTERN PetscErrorCode KSPNASHSetRadius(KSP,PetscReal); 571 PETSC_EXTERN PetscErrorCode KSPNASHGetNormD(KSP,PetscReal *); 572 PETSC_EXTERN PetscErrorCode KSPNASHGetObjFcn(KSP,PetscReal *); 573 574 PETSC_EXTERN PetscErrorCode KSPSTCGSetRadius(KSP,PetscReal); 575 PETSC_EXTERN PetscErrorCode KSPSTCGGetNormD(KSP,PetscReal *); 576 PETSC_EXTERN PetscErrorCode KSPSTCGGetObjFcn(KSP,PetscReal *); 577 578 PETSC_EXTERN PetscErrorCode KSPGLTRSetRadius(KSP,PetscReal); 579 PETSC_EXTERN PetscErrorCode KSPGLTRGetNormD(KSP,PetscReal *); 580 PETSC_EXTERN PetscErrorCode KSPGLTRGetObjFcn(KSP,PetscReal *); 581 PETSC_EXTERN PetscErrorCode KSPGLTRGetMinEig(KSP,PetscReal *); 582 PETSC_EXTERN PetscErrorCode KSPGLTRGetLambda(KSP,PetscReal *); 583 584 PETSC_EXTERN PetscErrorCode KSPPythonSetType(KSP,const char[]); 585 586 PETSC_EXTERN PetscErrorCode PCPreSolve(PC,KSP); 587 PETSC_EXTERN PetscErrorCode PCPostSolve(PC,KSP); 588 589 #include <petscdrawtypes.h> 590 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 591 PETSC_EXTERN PetscErrorCode KSPMonitorLGResidualNorm(KSP,PetscInt,PetscReal,void*); 592 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNormCreate(MPI_Comm,const char[],const char[],int,int,int,int,PetscDrawLG*); 593 PETSC_EXTERN PetscErrorCode KSPMonitorLGTrueResidualNorm(KSP,PetscInt,PetscReal,void*); 594 PETSC_EXTERN PetscErrorCode KSPMonitorLGRange(KSP,PetscInt,PetscReal,void*); 595 596 PETSC_EXTERN PetscErrorCode PCShellSetPreSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 597 PETSC_EXTERN PetscErrorCode PCShellSetPostSolve(PC,PetscErrorCode (*)(PC,KSP,Vec,Vec)); 598 599 /* see src/ksp/ksp/interface/iguess.c */ 600 typedef struct _p_KSPFischerGuess {PetscInt method,curl,maxl,refcnt;PetscBool monitor;Mat mat; KSP ksp;}* KSPFischerGuess; 601 602 PETSC_EXTERN PetscErrorCode KSPFischerGuessCreate(KSP,PetscInt,PetscInt,KSPFischerGuess*); 603 PETSC_EXTERN PetscErrorCode KSPFischerGuessDestroy(KSPFischerGuess*); 604 PETSC_EXTERN PetscErrorCode KSPFischerGuessReset(KSPFischerGuess); 605 PETSC_EXTERN PetscErrorCode KSPFischerGuessUpdate(KSPFischerGuess,Vec); 606 PETSC_EXTERN PetscErrorCode KSPFischerGuessFormGuess(KSPFischerGuess,Vec,Vec); 607 PETSC_EXTERN PetscErrorCode KSPFischerGuessSetFromOptions(KSPFischerGuess); 608 609 PETSC_EXTERN PetscErrorCode KSPSetUseFischerGuess(KSP,PetscInt,PetscInt); 610 PETSC_EXTERN PetscErrorCode KSPSetFischerGuess(KSP,KSPFischerGuess); 611 PETSC_EXTERN PetscErrorCode KSPGetFischerGuess(KSP,KSPFischerGuess*); 612 613 /*E 614 MatSchurComplementAinvType - Determines how to approximate the inverse of the (0,0) block in Schur complement preconditioning matrix assembly routines 615 616 Level: intermediate 617 618 .seealso: MatSchurComplementGetAinvType(), MatSchurComplementSetAinvType(), MatSchurComplementGetPmat(), MatGetSchurComplement(), MatCreateSchurComplementPmat() 619 E*/ 620 typedef enum {MAT_SCHUR_COMPLEMENT_AINV_DIAG, MAT_SCHUR_COMPLEMENT_AINV_LUMP} MatSchurComplementAinvType; 621 PETSC_EXTERN const char *const MatSchurComplementAinvTypes[]; 622 623 PETSC_EXTERN PetscErrorCode MatCreateSchurComplement(Mat,Mat,Mat,Mat,Mat,Mat*); 624 PETSC_EXTERN PetscErrorCode MatSchurComplementGetKSP(Mat,KSP*); 625 PETSC_EXTERN PetscErrorCode MatSchurComplementSetKSP(Mat,KSP); 626 PETSC_EXTERN PetscErrorCode MatSchurComplementSetSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 627 PETSC_EXTERN PetscErrorCode MatSchurComplementUpdateSubMatrices(Mat,Mat,Mat,Mat,Mat,Mat); 628 PETSC_EXTERN PetscErrorCode MatSchurComplementGetSubMatrices(Mat,Mat*,Mat*,Mat*,Mat*,Mat*); 629 PETSC_EXTERN PetscErrorCode MatSchurComplementSetAinvType(Mat,MatSchurComplementAinvType); 630 PETSC_EXTERN PetscErrorCode MatSchurComplementGetAinvType(Mat,MatSchurComplementAinvType*); 631 PETSC_EXTERN PetscErrorCode MatSchurComplementGetPmat(Mat,MatReuse,Mat*); 632 PETSC_EXTERN PetscErrorCode MatSchurComplementComputeExplicitOperator(Mat,Mat*); 633 PETSC_EXTERN PetscErrorCode MatGetSchurComplement(Mat,IS,IS,IS,IS,MatReuse,Mat *,MatSchurComplementAinvType,MatReuse,Mat *); 634 PETSC_EXTERN PetscErrorCode MatCreateSchurComplementPmat(Mat,Mat,Mat,Mat,MatSchurComplementAinvType,MatReuse,Mat*); 635 636 PETSC_EXTERN PetscErrorCode KSPSetDM(KSP,DM); 637 PETSC_EXTERN PetscErrorCode KSPSetDMActive(KSP,PetscBool ); 638 PETSC_EXTERN PetscErrorCode KSPGetDM(KSP,DM*); 639 PETSC_EXTERN PetscErrorCode KSPSetApplicationContext(KSP,void*); 640 PETSC_EXTERN PetscErrorCode KSPGetApplicationContext(KSP,void*); 641 PETSC_EXTERN PetscErrorCode KSPSetComputeRHS(KSP,PetscErrorCode (*func)(KSP,Vec,void*),void *); 642 PETSC_EXTERN PetscErrorCode KSPSetComputeOperators(KSP,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 643 PETSC_EXTERN PetscErrorCode KSPSetComputeInitialGuess(KSP,PetscErrorCode(*)(KSP,Vec,void*),void*); 644 PETSC_EXTERN PetscErrorCode DMKSPSetComputeOperators(DM,PetscErrorCode(*)(KSP,Mat,Mat,void*),void*); 645 PETSC_EXTERN PetscErrorCode DMKSPGetComputeOperators(DM,PetscErrorCode(**)(KSP,Mat,Mat,void*),void*); 646 PETSC_EXTERN PetscErrorCode DMKSPSetComputeRHS(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 647 PETSC_EXTERN PetscErrorCode DMKSPGetComputeRHS(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 648 PETSC_EXTERN PetscErrorCode DMKSPSetComputeInitialGuess(DM,PetscErrorCode(*)(KSP,Vec,void*),void*); 649 PETSC_EXTERN PetscErrorCode DMKSPGetComputeInitialGuess(DM,PetscErrorCode(**)(KSP,Vec,void*),void*); 650 651 PETSC_EXTERN PetscErrorCode DMGlobalToLocalSolve(DM,Vec,Vec); 652 PETSC_EXTERN PetscErrorCode DMPlexProjectField(DM, Vec, void (**)(const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscScalar[], const PetscReal [], PetscScalar []), InsertMode, Vec); 653 #endif 654