126bd1501SBarry Smith #if !defined(PETSCPCTYPES_H) 226bd1501SBarry Smith #define PETSCPCTYPES_H 3b0753f9dSMatthew G. Knepley 4ac09b921SBarry Smith /* SUBMANSEC = PC */ 5ac09b921SBarry Smith 6b0753f9dSMatthew G. Knepley /*S 787497f52SBarry Smith PC - Abstract PETSc object that manages all preconditioners including direct solvers such as `PCLU` 8b0753f9dSMatthew G. Knepley 9b0753f9dSMatthew G. Knepley Level: beginner 10b0753f9dSMatthew G. Knepley 11db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType` 12b0753f9dSMatthew G. Knepley S*/ 13b0753f9dSMatthew G. Knepley typedef struct _p_PC* PC; 14b0753f9dSMatthew G. Knepley 15b0753f9dSMatthew G. Knepley /*J 16b0753f9dSMatthew G. Knepley PCType - String with the name of a PETSc preconditioner method. 17b0753f9dSMatthew G. Knepley 18b0753f9dSMatthew G. Knepley Level: beginner 19b0753f9dSMatthew G. Knepley 2087497f52SBarry Smith Note: 2187497f52SBarry Smith `PCRegister()` is used to register preconditioners that are then accessible via `PCSetType()` 22b0753f9dSMatthew G. Knepley 2387497f52SBarry Smith .seealso: `PCSetType()`, `PC`, `PCCreate()`, `PCRegister()`, `PCSetFromOptions()`, `PCLU`, `PCJACOBI`, `PCBJACOBI` 24b0753f9dSMatthew G. Knepley J*/ 25b0753f9dSMatthew G. Knepley typedef const char* PCType; 26b0753f9dSMatthew G. Knepley #define PCNONE "none" 27b0753f9dSMatthew G. Knepley #define PCJACOBI "jacobi" 28b0753f9dSMatthew G. Knepley #define PCSOR "sor" 29b0753f9dSMatthew G. Knepley #define PCLU "lu" 30a2fc1e05SToby Isaac #define PCQR "qr" 31b0753f9dSMatthew G. Knepley #define PCSHELL "shell" 32*e6f8f311SMark Adams #define PCAMGX "amgx" 33b0753f9dSMatthew G. Knepley #define PCBJACOBI "bjacobi" 34b0753f9dSMatthew G. Knepley #define PCMG "mg" 35b0753f9dSMatthew G. Knepley #define PCEISENSTAT "eisenstat" 36b0753f9dSMatthew G. Knepley #define PCILU "ilu" 37b0753f9dSMatthew G. Knepley #define PCICC "icc" 38b0753f9dSMatthew G. Knepley #define PCASM "asm" 39b0753f9dSMatthew G. Knepley #define PCGASM "gasm" 40b0753f9dSMatthew G. Knepley #define PCKSP "ksp" 41e607c864SMark Adams #define PCBJKOKKOS "bjkokkos" 42b0753f9dSMatthew G. Knepley #define PCCOMPOSITE "composite" 43b0753f9dSMatthew G. Knepley #define PCREDUNDANT "redundant" 44b0753f9dSMatthew G. Knepley #define PCSPAI "spai" 45b0753f9dSMatthew G. Knepley #define PCNN "nn" 46b0753f9dSMatthew G. Knepley #define PCCHOLESKY "cholesky" 47b0753f9dSMatthew G. Knepley #define PCPBJACOBI "pbjacobi" 480da83c2eSBarry Smith #define PCVPBJACOBI "vpbjacobi" 49b0753f9dSMatthew G. Knepley #define PCMAT "mat" 50b0753f9dSMatthew G. Knepley #define PCHYPRE "hypre" 51b0753f9dSMatthew G. Knepley #define PCPARMS "parms" 52b0753f9dSMatthew G. Knepley #define PCFIELDSPLIT "fieldsplit" 53b0753f9dSMatthew G. Knepley #define PCTFS "tfs" 54b0753f9dSMatthew G. Knepley #define PCML "ml" 55b0753f9dSMatthew G. Knepley #define PCGALERKIN "galerkin" 56b0753f9dSMatthew G. Knepley #define PCEXOTIC "exotic" 57b0753f9dSMatthew G. Knepley #define PCCP "cp" 58b0753f9dSMatthew G. Knepley #define PCBFBT "bfbt" 59b0753f9dSMatthew G. Knepley #define PCLSC "lsc" 60b0753f9dSMatthew G. Knepley #define PCPYTHON "python" 61b0753f9dSMatthew G. Knepley #define PCPFMG "pfmg" 621c188c59Sftrigaux #define PCSMG "smg" 63b0753f9dSMatthew G. Knepley #define PCSYSPFMG "syspfmg" 64b0753f9dSMatthew G. Knepley #define PCREDISTRIBUTE "redistribute" 65b0753f9dSMatthew G. Knepley #define PCSVD "svd" 66b0753f9dSMatthew G. Knepley #define PCGAMG "gamg" 674b3f184cSKarl Rupp #define PCCHOWILUVIENNACL "chowiluviennacl" 6870baa948SKarl Rupp #define PCROWSCALINGVIENNACL "rowscalingviennacl" 6907598726SKarl Rupp #define PCSAVIENNACL "saviennacl" 70b0753f9dSMatthew G. Knepley #define PCBDDC "bddc" 71b0753f9dSMatthew G. Knepley #define PCKACZMARZ "kaczmarz" 7268ddcbeaSBarry Smith #define PCTELESCOPE "telescope" 734bbf5ea8SMatthew G. Knepley #define PCPATCH "patch" 74b9ac7092SAlp Dener #define PCLMVM "lmvm" 75360ee056SFande Kong #define PCHMG "hmg" 7637eeb815SJakub Kruzik #define PCDEFLATION "deflation" 7738cfc46eSPierre Jolivet #define PCHPDDM "hpddm" 7853022affSStefano Zampini #define PCH2OPUS "h2opus" 79f1f2ae84SBarry Smith #define PCMPI "mpi" 80b0753f9dSMatthew G. Knepley 81b0753f9dSMatthew G. Knepley /*E 82b0753f9dSMatthew G. Knepley PCSide - If the preconditioner is to be applied to the left, right 83b0753f9dSMatthew G. Knepley or symmetrically around the operator. 84b0753f9dSMatthew G. Knepley 85b0753f9dSMatthew G. Knepley Level: beginner 86b0753f9dSMatthew G. Knepley 87b0753f9dSMatthew G. Knepley .seealso: 88b0753f9dSMatthew G. Knepley E*/ 89b0753f9dSMatthew G. Knepley typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide; 90b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 91b0753f9dSMatthew G. Knepley 92b0753f9dSMatthew G. Knepley /*E 9387497f52SBarry Smith PCRichardsonConvergedReason - reason a `PCApplyRichardson() method terminated 94b0753f9dSMatthew G. Knepley 95b0753f9dSMatthew G. Knepley Level: advanced 96b0753f9dSMatthew G. Knepley 9787497f52SBarry Smith Developer Note: 9887497f52SBarry Smith this must match petsc/finclude/petscpc.h and the `KSPConvergedReason` values in petscksp.h 99b0753f9dSMatthew G. Knepley 100db781477SPatrick Sanan .seealso: `PCApplyRichardson()` 101b0753f9dSMatthew G. Knepley E*/ 102b0753f9dSMatthew G. Knepley typedef enum { 103b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_RTOL = 2, 104b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ATOL = 3, 105b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ITS = 4, 106b0753f9dSMatthew G. Knepley PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason; 107b0753f9dSMatthew G. Knepley 108b0753f9dSMatthew G. Knepley /*E 109b0753f9dSMatthew G. Knepley PCJacobiType - What elements are used to form the Jacobi preconditioner 110b0753f9dSMatthew G. Knepley 111b0753f9dSMatthew G. Knepley Level: intermediate 112b0753f9dSMatthew G. Knepley 113b0753f9dSMatthew G. Knepley .seealso: 114b0753f9dSMatthew G. Knepley E*/ 115b0753f9dSMatthew G. Knepley typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType; 116b0753f9dSMatthew G. Knepley 117b0753f9dSMatthew G. Knepley /*E 118b0753f9dSMatthew G. Knepley PCASMType - Type of additive Schwarz method to use 119b0753f9dSMatthew G. Knepley 12087497f52SBarry Smith $ `PC_ASM_BASIC` - Symmetric version where residuals from the ghost points are used 121b0753f9dSMatthew G. Knepley $ and computed values in ghost regions are added together. 122b0753f9dSMatthew G. Knepley $ Classical standard additive Schwarz. 12387497f52SBarry Smith $ `PC_ASM_RESTRICT` - Residuals from ghost points are used but computed values in ghost 124b0753f9dSMatthew G. Knepley $ region are discarded. 125b0753f9dSMatthew G. Knepley $ Default. 12687497f52SBarry Smith $ `PC_ASM_INTERPOLATE` - Residuals from ghost points are not used, computed values in ghost 127b0753f9dSMatthew G. Knepley $ region are added back in. 12887497f52SBarry Smith $ `PC_ASM_NONE` - Residuals from ghost points are not used, computed ghost values are 129b0753f9dSMatthew G. Knepley $ discarded. 130b0753f9dSMatthew G. Knepley $ Not very good. 131b0753f9dSMatthew G. Knepley 132b0753f9dSMatthew G. Knepley Level: beginner 133b0753f9dSMatthew G. Knepley 134db781477SPatrick Sanan .seealso: `PCASMSetType()` 135b0753f9dSMatthew G. Knepley E*/ 136b0753f9dSMatthew G. Knepley typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType; 137b0753f9dSMatthew G. Knepley 138b0753f9dSMatthew G. Knepley /*E 13987497f52SBarry Smith PCGASMType - Type of generalized additive Schwarz method to use (differs from `PCASM` in allowing multiple processors per subdomain). 140b0753f9dSMatthew G. Knepley 141b0753f9dSMatthew G. Knepley Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 142b0753f9dSMatthew G. Knepley domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 143b0753f9dSMatthew G. Knepley a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 144b0753f9dSMatthew G. Knepley (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 145b0753f9dSMatthew G. Knepley 14687497f52SBarry Smith $ `PC_GASM_BASIC` - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 147b0753f9dSMatthew G. Knepley $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 148b0753f9dSMatthew G. Knepley $ from neighboring subdomains. 149b0753f9dSMatthew G. Knepley $ Classical standard additive Schwarz. 15087497f52SBarry Smith $ `PC_GASM_RESTRICT` - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 151b0753f9dSMatthew G. Knepley $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 152b0753f9dSMatthew G. Knepley $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 153b0753f9dSMatthew G. Knepley $ assumption). 154b0753f9dSMatthew G. Knepley $ Default. 15587497f52SBarry Smith $ `PC_GASM_INTERPOLATE` - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 156b0753f9dSMatthew G. Knepley $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 157b0753f9dSMatthew G. Knepley $ from neighboring subdomains. 158b0753f9dSMatthew G. Knepley $ 15987497f52SBarry Smith $ `PC_GASM_NONE` - Residuals and corrections are zeroed out outside the local subdomains. 160b0753f9dSMatthew G. Knepley $ Not very good. 161b0753f9dSMatthew G. Knepley 162b0753f9dSMatthew G. Knepley Level: beginner 163b0753f9dSMatthew G. Knepley 164db781477SPatrick Sanan .seealso: `PCGASMSetType()` 165b0753f9dSMatthew G. Knepley E*/ 166b0753f9dSMatthew G. Knepley typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType; 167b0753f9dSMatthew G. Knepley 168b0753f9dSMatthew G. Knepley /*E 169b0753f9dSMatthew G. Knepley PCCompositeType - Determines how two or more preconditioner are composed 170b0753f9dSMatthew G. Knepley 17187497f52SBarry Smith $ `PC_COMPOSITE_ADDITIVE` - results from application of all preconditioners are added together 17287497f52SBarry Smith $ `PC_COMPOSITE_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly 173b0753f9dSMatthew G. Knepley $ computed after the previous preconditioner application 17487497f52SBarry Smith $ `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly 175d0ecd4dfSBarry Smith $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners) 17687497f52SBarry Smith $ `PC_COMPOSITE_SPECIAL` - This is very special for a matrix of the form alpha I + R + S 177b0753f9dSMatthew G. Knepley $ where first preconditioner is built from alpha I + S and second from 178b0753f9dSMatthew G. Knepley $ alpha I + R 179b0753f9dSMatthew G. Knepley 180b0753f9dSMatthew G. Knepley Level: beginner 181b0753f9dSMatthew G. Knepley 182db781477SPatrick Sanan .seealso: `PCCompositeSetType()` 183b0753f9dSMatthew G. Knepley E*/ 184a51937d4SCarola Kruse typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR,PC_COMPOSITE_GKB} PCCompositeType; 185b0753f9dSMatthew G. Knepley 186b0753f9dSMatthew G. Knepley /*E 18787497f52SBarry Smith PCFieldSplitSchurPreType - Determines how to precondition a Schur complement 188b0753f9dSMatthew G. Knepley 189b0753f9dSMatthew G. Knepley Level: intermediate 190b0753f9dSMatthew G. Knepley 191db781477SPatrick Sanan .seealso: `PCFieldSplitSetSchurPre()` 192b0753f9dSMatthew G. Knepley E*/ 193b0753f9dSMatthew G. Knepley 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; 194b0753f9dSMatthew G. Knepley 195b0753f9dSMatthew G. Knepley /*E 196b0753f9dSMatthew G. Knepley PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 197b0753f9dSMatthew G. Knepley 198b0753f9dSMatthew G. Knepley Level: intermediate 199b0753f9dSMatthew G. Knepley 200db781477SPatrick Sanan .seealso: `PCFieldSplitSetSchurFactType()` 201b0753f9dSMatthew G. Knepley E*/ 202b0753f9dSMatthew G. Knepley typedef enum { 203b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_DIAG, 204b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_LOWER, 205b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_UPPER, 206b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_FULL 207b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType; 208b0753f9dSMatthew G. Knepley 209b0753f9dSMatthew G. Knepley /*E 21087497f52SBarry Smith PCPARMSGlobalType - Determines the global preconditioner method in `PCPARMS` 211b0753f9dSMatthew G. Knepley 212b0753f9dSMatthew G. Knepley Level: intermediate 213b0753f9dSMatthew G. Knepley 214db781477SPatrick Sanan .seealso: `PCPARMSSetGlobal()` 215b0753f9dSMatthew G. Knepley E*/ 216b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType; 2179d8081ecSMatthew G. Knepley 218b0753f9dSMatthew G. Knepley /*E 21987497f52SBarry Smith PCPARMSLocalType - Determines the local preconditioner method in `PCPARMS` 220b0753f9dSMatthew G. Knepley 221b0753f9dSMatthew G. Knepley Level: intermediate 222b0753f9dSMatthew G. Knepley 223db781477SPatrick Sanan .seealso: `PCPARMSSetLocal()` 224b0753f9dSMatthew G. Knepley E*/ 225b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType; 226b0753f9dSMatthew G. Knepley 227f0fc11ceSJed Brown /*J 22887497f52SBarry Smith PCGAMGType - type of generalized algebraic multigrid `PCGAMG` method 229b0753f9dSMatthew G. Knepley 230b0753f9dSMatthew G. Knepley Level: intermediate 231b0753f9dSMatthew G. Knepley 23287497f52SBarry Smith $ `PCGAMGAGG` - (the default) smoothed aggregation algorithm, robust, very well tested 23387497f52SBarry Smith $ `PCGAMGGEO` - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not well tested 23487497f52SBarry Smith $ `PCGAMGCLASSICAL` - classical algebraic multigrid preconditioner, incomplete, poorly tested 235f6130c64SBarry Smith 236db781477SPatrick Sanan .seealso: `PCMG`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()` 237f0fc11ceSJed Brown J*/ 238b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType; 239b0753f9dSMatthew G. Knepley #define PCGAMGAGG "agg" 240b0753f9dSMatthew G. Knepley #define PCGAMGGEO "geo" 241b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL "classical" 242b0753f9dSMatthew G. Knepley 243b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType; 244b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT "direct" 245b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard" 246b0753f9dSMatthew G. Knepley 247b0753f9dSMatthew G. Knepley /*E 248b0753f9dSMatthew G. Knepley PCMGType - Determines the type of multigrid method that is run. 249b0753f9dSMatthew G. Knepley 250b0753f9dSMatthew G. Knepley Level: beginner 251b0753f9dSMatthew G. Knepley 252b0753f9dSMatthew G. Knepley Values: 25387497f52SBarry Smith + `PC_MG_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `PCMGSetCycleType()` 25487497f52SBarry Smith . `PC_MG_ADDITIVE` - the additive multigrid preconditioner where all levels are 255b0753f9dSMatthew G. Knepley smoothed before updating the residual. This only uses the 256b0753f9dSMatthew G. Knepley down smoother, in the preconditioner the upper smoother is ignored 25787497f52SBarry Smith . `PC_MG_FULL` - same as multiplicative except one also performs grid sequencing, 258b0753f9dSMatthew G. Knepley that is starts on the coarsest grid, performs a cycle, interpolates 259b0753f9dSMatthew G. Knepley 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 260b0753f9dSMatthew G. Knepley algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 261b0753f9dSMatthew G. Knepley calls the V-cycle only on the coarser level and has a post-smoother instead. 26287497f52SBarry Smith - `PC_MG_KASKADE` - like full multigrid except one never goes back to a coarser level 263b0753f9dSMatthew G. Knepley from a finer 264b0753f9dSMatthew G. Knepley 265db781477SPatrick Sanan .seealso: `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()` 266b0753f9dSMatthew G. Knepley 267b0753f9dSMatthew G. Knepley E*/ 268b0753f9dSMatthew G. Knepley typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType; 269b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE; 270b0753f9dSMatthew G. Knepley 271b0753f9dSMatthew G. Knepley /*E 272b0753f9dSMatthew G. Knepley PCMGCycleType - Use V-cycle or W-cycle 273b0753f9dSMatthew G. Knepley 274b0753f9dSMatthew G. Knepley Level: beginner 275b0753f9dSMatthew G. Knepley 276b0753f9dSMatthew G. Knepley Values: 27787497f52SBarry Smith + `PC_MG_V_CYCLE` - use the v cycle 27887497f52SBarry Smith - `PC_MG_W_CYCLE` - use the w cycle 279b0753f9dSMatthew G. Knepley 280db781477SPatrick Sanan .seealso: `PCMGSetCycleType()` 281b0753f9dSMatthew G. Knepley 282b0753f9dSMatthew G. Knepley E*/ 283b0753f9dSMatthew G. Knepley typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType; 284b0753f9dSMatthew G. Knepley 285b0753f9dSMatthew G. Knepley /*E 2862134b1e4SBarry Smith PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process 2872134b1e4SBarry Smith 2882134b1e4SBarry Smith Level: beginner 2892134b1e4SBarry Smith 2902134b1e4SBarry Smith Values: 29187497f52SBarry Smith + `PC_MG_GALERKIN_PMAT` - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid 29287497f52SBarry Smith . `PC_MG_GALERKIN_MAT` - computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid 29387497f52SBarry Smith . `PC_MG_GALERKIN_BOTH` - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once 29487497f52SBarry Smith - `PC_MG_GALERKIN_NONE` - neither operator is computed via the Galerkin process, the user must provide the operator 2952134b1e4SBarry Smith 29687497f52SBarry Smith Users should never set `PC_MG_GALERKIN_EXTERNAL`, it is used by `PCHYPRE` and `PCML` 2972134b1e4SBarry Smith 298db781477SPatrick Sanan .seealso: `PCMGSetCycleType()` 2992134b1e4SBarry Smith 3002134b1e4SBarry Smith E*/ 3012134b1e4SBarry Smith typedef enum { PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, PC_MG_GALERKIN_EXTERNAL} PCMGGalerkinType; 3022134b1e4SBarry Smith 3032134b1e4SBarry Smith /*E 304b0753f9dSMatthew G. Knepley PCExoticType - Face based or wirebasket based coarse grid space 305b0753f9dSMatthew G. Knepley 306b0753f9dSMatthew G. Knepley Level: beginner 307b0753f9dSMatthew G. Knepley 308db781477SPatrick Sanan .seealso: `PCExoticSetType()`, `PCEXOTIC` 309b0753f9dSMatthew G. Knepley E*/ 310b0753f9dSMatthew G. Knepley typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType; 311b0753f9dSMatthew G. Knepley 3128c1cd74cSHong Zhang /*E 313bc960bbfSJed Brown PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains. 314bc960bbfSJed Brown 315bc960bbfSJed Brown Level: intermediate 316bc960bbfSJed Brown 317bc960bbfSJed Brown Values: 31887497f52SBarry Smith + `PC_BDDC_INTERFACE_EXT_DIRICHLET` - solves Dirichlet interior problem; this is the standard BDDC algorithm 31987497f52SBarry Smith - `PC_BDDC_INTERFACE_EXT_LUMP` - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP" 320bc960bbfSJed Brown 321bc960bbfSJed Brown E*/ 322bc960bbfSJed Brown typedef enum { 323bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_DIRICHLET, 324bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_LUMP 325bc960bbfSJed Brown } PCBDDCInterfaceExtType; 326bc960bbfSJed Brown 327bc960bbfSJed Brown /*E 328f3b08a26SMatthew G. Knepley PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation 329f3b08a26SMatthew G. Knepley 330f3b08a26SMatthew G. Knepley Level: beginner 331f3b08a26SMatthew G. Knepley 332db781477SPatrick Sanan .seealso: `PCMGSetAdaptCoarseSpaceType()`, `PCMG` 333f3b08a26SMatthew G. Knepley E*/ 3342b3cbbdaSStefano Zampini typedef enum { PCMG_ADAPT_NONE, PCMG_ADAPT_POLYNOMIAL, PCMG_ADAPT_HARMONIC, PCMG_ADAPT_EIGENVECTOR, PCMG_ADAPT_GENERALIZED_EIGENVECTOR, PCMG_ADAPT_GDSW } PCMGCoarseSpaceType; 335f3b08a26SMatthew G. Knepley 336f3b08a26SMatthew G. Knepley /*E 337ee6acaa3SMatthew G. Knepley PCPatchConstructType - The algorithm used to construct patches for the preconditioner 3384bbf5ea8SMatthew G. Knepley 3394bbf5ea8SMatthew G. Knepley Level: beginner 3404bbf5ea8SMatthew G. Knepley 341db781477SPatrick Sanan .seealso: `PCPatchSetConstructType()`, `PCEXOTIC` 3424bbf5ea8SMatthew G. Knepley E*/ 343e5b9877fSPatrick Farrell typedef enum {PC_PATCH_STAR, PC_PATCH_VANKA, PC_PATCH_PARDECOMP, PC_PATCH_USER, PC_PATCH_PYTHON} PCPatchConstructType; 3444bbf5ea8SMatthew G. Knepley 3454bbf5ea8SMatthew G. Knepley /*E 346e26ad46dSJakub Kruzik PCDeflationSpaceType - Type of deflation 347e662bc50SJakub Kruzik 348e662bc50SJakub Kruzik Values: 34987497f52SBarry Smith + `PC_DEFLATION_SPACE_HAAR` - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off 35087497f52SBarry Smith . `PC_DEFLATION_SPACE_DB2` - MATCOMPOSITE of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet) 35187497f52SBarry Smith . `PC_DEFLATION_SPACE_DB4` - same as above, but with db4 (4 coefficient Daubechies) 35287497f52SBarry Smith . `PC_DEFLATION_SPACE_DB8` - same as above, but with db8 (8 coefficient Daubechies) 35387497f52SBarry Smith . `PC_DEFLATION_SPACE_DB16` - same as above, but with db16 (16 coefficient Daubechies) 35487497f52SBarry Smith . `PC_DEFLATION_SPACE_BIORTH22` - same as above, but with biorthogonal 2.2 (6 coefficients) 35587497f52SBarry Smith . `PC_DEFLATION_SPACE_MEYER` - same as above, but with Meyer/FIR (62 coefficients) 35687497f52SBarry Smith . `PC_DEFLATION_SPACE_AGGREGATION` - aggregates local indices (given by operator matix distribution) into a subdomain 35787497f52SBarry Smith - `PC_DEFLATION_SPACE_USER` - indicates space set by user 358e662bc50SJakub Kruzik 359e662bc50SJakub Kruzik Notes: 3607e9012c5SJakub Kruzik Wavelet-based space (except Haar) can be used in multilevel deflation. 361e662bc50SJakub Kruzik 3621fdb00f9SJakub Kruzik Level: intermediate 3631fdb00f9SJakub Kruzik 364db781477SPatrick Sanan .seealso: `PCDeflationSetSpaceToCompute()`, `PCDEFLATION` 365e662bc50SJakub Kruzik E*/ 366e662bc50SJakub Kruzik typedef enum { 367e662bc50SJakub Kruzik PC_DEFLATION_SPACE_HAAR, 368e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB2, 369e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB4, 370e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB8, 371e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB16, 372e662bc50SJakub Kruzik PC_DEFLATION_SPACE_BIORTH22, 373e662bc50SJakub Kruzik PC_DEFLATION_SPACE_MEYER, 374e662bc50SJakub Kruzik PC_DEFLATION_SPACE_AGGREGATION, 375e662bc50SJakub Kruzik PC_DEFLATION_SPACE_USER 376e662bc50SJakub Kruzik } PCDeflationSpaceType; 377e662bc50SJakub Kruzik 378e662bc50SJakub Kruzik /*E 37987497f52SBarry Smith PCHPDDMCoarseCorrectionType - Type of coarse correction used by `PCHPDDM` 38038cfc46eSPierre Jolivet 38138cfc46eSPierre Jolivet Level: intermediate 38238cfc46eSPierre Jolivet 38338cfc46eSPierre Jolivet Values: 38487497f52SBarry Smith + `PC_HPDDM_COARSE_CORRECTION_DEFLATED` (default) - eq. (1) in PCHPDDMShellApply() 38587497f52SBarry Smith . `PC_HPDDM_COARSE_CORRECTION_ADDITIVE` - eq. (2) 38687497f52SBarry Smith - `PC_HPDDM_COARSE_CORRECTION_BALANCED` - eq. (3) 38738cfc46eSPierre Jolivet 388db781477SPatrick Sanan .seealso: `PCHPDDM`, `PCSetType()`, `PCHPDDMShellApply()` 38938cfc46eSPierre Jolivet E*/ 39038cfc46eSPierre Jolivet typedef enum { PC_HPDDM_COARSE_CORRECTION_DEFLATED, PC_HPDDM_COARSE_CORRECTION_ADDITIVE, PC_HPDDM_COARSE_CORRECTION_BALANCED } PCHPDDMCoarseCorrectionType; 39138cfc46eSPierre Jolivet 39238cfc46eSPierre Jolivet /*E 39387497f52SBarry Smith PCFailedReason - indicates type of `PC` failure 3948c1cd74cSHong Zhang 3958c1cd74cSHong Zhang Level: beginner 3968c1cd74cSHong Zhang 39787497f52SBarry Smith Developer Note: 3988c1cd74cSHong Zhang Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h 3998c1cd74cSHong Zhang E*/ 4001b2b9847SBarry Smith typedef enum {PC_SETUP_ERROR = -1,PC_NOERROR,PC_FACTOR_STRUCT_ZEROPIVOT,PC_FACTOR_NUMERIC_ZEROPIVOT,PC_FACTOR_OUTMEMORY,PC_FACTOR_OTHER,PC_SUBPC_ERROR} PCFailedReason; 401ce7c7f2fSMark Adams 402ce7c7f2fSMark Adams /*E 403ce7c7f2fSMark Adams PCGAMGLayoutType - Layout for reduced grids 404ce7c7f2fSMark Adams 405ce7c7f2fSMark Adams Level: intermediate 406ce7c7f2fSMark Adams 407db781477SPatrick Sanan .seealso: `PCGAMGSetCoarseGridLayoutType()` 4089969d477SPatrick Sanan 409ce7c7f2fSMark Adams Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h 410ce7c7f2fSMark Adams E*/ 411a0095786SMark typedef enum {PCGAMG_LAYOUT_COMPACT,PCGAMG_LAYOUT_SPREAD} PCGAMGLayoutType; 412ce7c7f2fSMark Adams 413b0753f9dSMatthew G. Knepley #endif 414