1 #if !defined(PETSCPCTYPES_H) 2 #define PETSCPCTYPES_H 3 4 /*S 5 PC - Abstract PETSc object that manages all preconditioners including direct solvers such as PCLU 6 7 Level: beginner 8 9 .seealso: PCCreate(), PCSetType(), PCType (for list of available types) 10 S*/ 11 typedef struct _p_PC* PC; 12 13 /*J 14 PCType - String with the name of a PETSc preconditioner method. 15 16 Level: beginner 17 18 Notes: 19 Click on the links above to see details on a particular solver 20 21 PCRegister() is used to register preconditioners that are then accessible via PCSetType() 22 23 .seealso: PCSetType(), PC, PCCreate(), PCRegister(), PCSetFromOptions() 24 J*/ 25 typedef const char* PCType; 26 #define PCNONE "none" 27 #define PCJACOBI "jacobi" 28 #define PCSOR "sor" 29 #define PCLU "lu" 30 #define PCSHELL "shell" 31 #define PCBJACOBI "bjacobi" 32 #define PCMG "mg" 33 #define PCEISENSTAT "eisenstat" 34 #define PCILU "ilu" 35 #define PCICC "icc" 36 #define PCASM "asm" 37 #define PCGASM "gasm" 38 #define PCKSP "ksp" 39 #define PCCOMPOSITE "composite" 40 #define PCREDUNDANT "redundant" 41 #define PCSPAI "spai" 42 #define PCNN "nn" 43 #define PCCHOLESKY "cholesky" 44 #define PCPBJACOBI "pbjacobi" 45 #define PCVPBJACOBI "vpbjacobi" 46 #define PCMAT "mat" 47 #define PCHYPRE "hypre" 48 #define PCPARMS "parms" 49 #define PCFIELDSPLIT "fieldsplit" 50 #define PCTFS "tfs" 51 #define PCML "ml" 52 #define PCGALERKIN "galerkin" 53 #define PCEXOTIC "exotic" 54 #define PCCP "cp" 55 #define PCBFBT "bfbt" 56 #define PCLSC "lsc" 57 #define PCPYTHON "python" 58 #define PCPFMG "pfmg" 59 #define PCSYSPFMG "syspfmg" 60 #define PCREDISTRIBUTE "redistribute" 61 #define PCSVD "svd" 62 #define PCGAMG "gamg" 63 #define PCCHOWILUVIENNACL "chowiluviennacl" 64 #define PCROWSCALINGVIENNACL "rowscalingviennacl" 65 #define PCSAVIENNACL "saviennacl" 66 #define PCBDDC "bddc" 67 #define PCKACZMARZ "kaczmarz" 68 #define PCTELESCOPE "telescope" 69 #define PCPATCH "patch" 70 #define PCLMVM "lmvm" 71 #define PCHMG "hmg" 72 73 /*E 74 PCSide - If the preconditioner is to be applied to the left, right 75 or symmetrically around the operator. 76 77 Level: beginner 78 79 .seealso: 80 E*/ 81 typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide; 82 #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 83 84 /*E 85 PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates 86 87 Level: advanced 88 89 Notes: 90 this must match petsc/finclude/petscpc.h and the KSPConvergedReason values in petscksp.h 91 92 .seealso: PCApplyRichardson() 93 E*/ 94 typedef enum { 95 PCRICHARDSON_CONVERGED_RTOL = 2, 96 PCRICHARDSON_CONVERGED_ATOL = 3, 97 PCRICHARDSON_CONVERGED_ITS = 4, 98 PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason; 99 100 /*E 101 PCJacobiType - What elements are used to form the Jacobi preconditioner 102 103 Level: intermediate 104 105 .seealso: 106 E*/ 107 typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType; 108 109 /*E 110 PCASMType - Type of additive Schwarz method to use 111 112 $ PC_ASM_BASIC - Symmetric version where residuals from the ghost points are used 113 $ and computed values in ghost regions are added together. 114 $ Classical standard additive Schwarz. 115 $ PC_ASM_RESTRICT - Residuals from ghost points are used but computed values in ghost 116 $ region are discarded. 117 $ Default. 118 $ PC_ASM_INTERPOLATE - Residuals from ghost points are not used, computed values in ghost 119 $ region are added back in. 120 $ PC_ASM_NONE - Residuals from ghost points are not used, computed ghost values are 121 $ discarded. 122 $ Not very good. 123 124 Level: beginner 125 126 .seealso: PCASMSetType() 127 E*/ 128 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType; 129 130 /*E 131 PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain). 132 133 Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 134 domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 135 a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 136 (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 137 138 $ PC_GASM_BASIC - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 139 $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 140 $ from neighboring subdomains. 141 $ Classical standard additive Schwarz. 142 $ PC_GASM_RESTRICT - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 143 $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 144 $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 145 $ assumption). 146 $ Default. 147 $ PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 148 $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 149 $ from neighboring subdomains. 150 $ 151 $ PC_GASM_NONE - Residuals and corrections are zeroed out outside the local subdomains. 152 $ Not very good. 153 154 Level: beginner 155 156 .seealso: PCGASMSetType() 157 E*/ 158 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType; 159 160 /*E 161 PCCompositeType - Determines how two or more preconditioner are composed 162 163 $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together 164 $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 165 $ computed after the previous preconditioner application 166 $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 167 $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners) 168 $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S 169 $ where first preconditioner is built from alpha I + S and second from 170 $ alpha I + R 171 172 Level: beginner 173 174 .seealso: PCCompositeSetType() 175 E*/ 176 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR,PC_COMPOSITE_GKB} PCCompositeType; 177 178 /*E 179 PCFieldSplitSchurPreType - Determines how to precondition Schur complement 180 181 Level: intermediate 182 183 .seealso: PCFieldSplitSetSchurPre() 184 E*/ 185 typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_SELFP,PC_FIELDSPLIT_SCHUR_PRE_A11,PC_FIELDSPLIT_SCHUR_PRE_USER,PC_FIELDSPLIT_SCHUR_PRE_FULL} PCFieldSplitSchurPreType; 186 187 /*E 188 PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 189 190 Level: intermediate 191 192 .seealso: PCFieldSplitSetSchurFactType() 193 E*/ 194 typedef enum { 195 PC_FIELDSPLIT_SCHUR_FACT_DIAG, 196 PC_FIELDSPLIT_SCHUR_FACT_LOWER, 197 PC_FIELDSPLIT_SCHUR_FACT_UPPER, 198 PC_FIELDSPLIT_SCHUR_FACT_FULL 199 } PCFieldSplitSchurFactType; 200 201 /*E 202 PCPARMSGlobalType - Determines the global preconditioner method in PARMS 203 204 Level: intermediate 205 206 .seealso: PCPARMSSetGlobal() 207 E*/ 208 typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType; 209 210 /*E 211 PCPARMSLocalType - Determines the local preconditioner method in PARMS 212 213 Level: intermediate 214 215 .seealso: PCPARMSSetLocal() 216 E*/ 217 typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType; 218 219 /*E 220 PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method 221 222 Level: intermediate 223 224 .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseInterpolation() 225 E*/ 226 typedef const char *PCGAMGType; 227 #define PCGAMGAGG "agg" 228 #define PCGAMGGEO "geo" 229 #define PCGAMGCLASSICAL "classical" 230 231 typedef const char *PCGAMGClassicalType; 232 #define PCGAMGCLASSICALDIRECT "direct" 233 #define PCGAMGCLASSICALSTANDARD "standard" 234 235 /*E 236 PCMGType - Determines the type of multigrid method that is run. 237 238 Level: beginner 239 240 Values: 241 + PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycleType() 242 . PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are 243 smoothed before updating the residual. This only uses the 244 down smoother, in the preconditioner the upper smoother is ignored 245 . PC_MG_FULL - same as multiplicative except one also performs grid sequencing, 246 that is starts on the coarsest grid, performs a cycle, interpolates 247 to the next, performs a cycle etc. This is much like the F-cycle presented in "Multigrid" by Trottenberg, Oosterlee, Schuller page 49, but that 248 algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 249 calls the V-cycle only on the coarser level and has a post-smoother instead. 250 - PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level 251 from a finer 252 253 .seealso: PCMGSetType(), PCMGSetCycleType(), PCMGSetCycleTypeOnLevel() 254 255 E*/ 256 typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType; 257 #define PC_MG_CASCADE PC_MG_KASKADE; 258 259 /*E 260 PCMGCycleType - Use V-cycle or W-cycle 261 262 Level: beginner 263 264 Values: 265 + PC_MG_V_CYCLE 266 - PC_MG_W_CYCLE 267 268 .seealso: PCMGSetCycleType() 269 270 E*/ 271 typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType; 272 273 /*E 274 PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process 275 276 Level: beginner 277 278 Values: 279 + PC_MG_GALERKIN_PMAT - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid 280 . PC_MG_GALERKIN_MAT - computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid 281 . PC_MG_GALERKIN_BOTH - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once 282 - PC_MG_GALERKIN_NONE - neither operator is computed via the Galerkin process, the user must provide the operator 283 284 Users should never set PC_MG_GALERKIN_EXTERNAL, it is used by GAMG and ML 285 286 .seealso: PCMGSetCycleType() 287 288 E*/ 289 typedef enum { PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, PC_MG_GALERKIN_EXTERNAL} PCMGGalerkinType; 290 291 /*E 292 PCExoticType - Face based or wirebasket based coarse grid space 293 294 Level: beginner 295 296 .seealso: PCExoticSetType(), PCEXOTIC 297 E*/ 298 typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType; 299 300 /*E 301 PCPatchConstructType - The algorithm used to construct patches for the preconditioner 302 303 Level: beginner 304 305 .seealso: PCPatchSetConstructType(), PCEXOTIC 306 E*/ 307 typedef enum {PC_PATCH_STAR, PC_PATCH_VANKA, PC_PATCH_PARDECOMP, PC_PATCH_USER, PC_PATCH_PYTHON} PCPatchConstructType; 308 309 /*E 310 PCFailedReason - indicates type of PC failure 311 312 Level: beginner 313 314 Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h 315 E*/ 316 typedef enum {PC_NOERROR,PC_FACTOR_STRUCT_ZEROPIVOT,PC_FACTOR_NUMERIC_ZEROPIVOT,PC_FACTOR_OUTMEMORY,PC_FACTOR_OTHER,PC_SUBPC_ERROR} PCFailedReason; 317 #endif 318