126bd1501SBarry Smith #if !defined(PETSCPCTYPES_H) 226bd1501SBarry Smith #define PETSCPCTYPES_H 3b0753f9dSMatthew G. Knepley 4ac09b921SBarry Smith /* SUBMANSEC = PC */ 5ac09b921SBarry Smith 6b0753f9dSMatthew G. Knepley /*S 7b0753f9dSMatthew G. Knepley 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 2095452b02SPatrick Sanan Notes: 2195452b02SPatrick Sanan Click on the links above to see details on a particular solver 22b0753f9dSMatthew G. Knepley 23b0753f9dSMatthew G. Knepley PCRegister() is used to register preconditioners that are then accessible via PCSetType() 24b0753f9dSMatthew G. Knepley 25db781477SPatrick Sanan .seealso: `PCSetType()`, `PC`, `PCCreate()`, `PCRegister()`, `PCSetFromOptions()` 26b0753f9dSMatthew G. Knepley J*/ 27b0753f9dSMatthew G. Knepley typedef const char* PCType; 28b0753f9dSMatthew G. Knepley #define PCNONE "none" 29b0753f9dSMatthew G. Knepley #define PCJACOBI "jacobi" 30b0753f9dSMatthew G. Knepley #define PCSOR "sor" 31b0753f9dSMatthew G. Knepley #define PCLU "lu" 32a2fc1e05SToby Isaac #define PCQR "qr" 33b0753f9dSMatthew G. Knepley #define PCSHELL "shell" 34b0753f9dSMatthew G. Knepley #define PCBJACOBI "bjacobi" 35b0753f9dSMatthew G. Knepley #define PCMG "mg" 36b0753f9dSMatthew G. Knepley #define PCEISENSTAT "eisenstat" 37b0753f9dSMatthew G. Knepley #define PCILU "ilu" 38b0753f9dSMatthew G. Knepley #define PCICC "icc" 39b0753f9dSMatthew G. Knepley #define PCASM "asm" 40b0753f9dSMatthew G. Knepley #define PCGASM "gasm" 41b0753f9dSMatthew G. Knepley #define PCKSP "ksp" 42e607c864SMark Adams #define PCBJKOKKOS "bjkokkos" 43b0753f9dSMatthew G. Knepley #define PCCOMPOSITE "composite" 44b0753f9dSMatthew G. Knepley #define PCREDUNDANT "redundant" 45b0753f9dSMatthew G. Knepley #define PCSPAI "spai" 46b0753f9dSMatthew G. Knepley #define PCNN "nn" 47b0753f9dSMatthew G. Knepley #define PCCHOLESKY "cholesky" 48b0753f9dSMatthew G. Knepley #define PCPBJACOBI "pbjacobi" 490da83c2eSBarry Smith #define PCVPBJACOBI "vpbjacobi" 50b0753f9dSMatthew G. Knepley #define PCMAT "mat" 51b0753f9dSMatthew G. Knepley #define PCHYPRE "hypre" 52b0753f9dSMatthew G. Knepley #define PCPARMS "parms" 53b0753f9dSMatthew G. Knepley #define PCFIELDSPLIT "fieldsplit" 54b0753f9dSMatthew G. Knepley #define PCTFS "tfs" 55b0753f9dSMatthew G. Knepley #define PCML "ml" 56b0753f9dSMatthew G. Knepley #define PCGALERKIN "galerkin" 57b0753f9dSMatthew G. Knepley #define PCEXOTIC "exotic" 58b0753f9dSMatthew G. Knepley #define PCCP "cp" 59b0753f9dSMatthew G. Knepley #define PCBFBT "bfbt" 60b0753f9dSMatthew G. Knepley #define PCLSC "lsc" 61b0753f9dSMatthew G. Knepley #define PCPYTHON "python" 62b0753f9dSMatthew G. Knepley #define PCPFMG "pfmg" 63*1c188c59Sftrigaux #define PCSMG "smg" 64b0753f9dSMatthew G. Knepley #define PCSYSPFMG "syspfmg" 65b0753f9dSMatthew G. Knepley #define PCREDISTRIBUTE "redistribute" 66b0753f9dSMatthew G. Knepley #define PCSVD "svd" 67b0753f9dSMatthew G. Knepley #define PCGAMG "gamg" 684b3f184cSKarl Rupp #define PCCHOWILUVIENNACL "chowiluviennacl" 6970baa948SKarl Rupp #define PCROWSCALINGVIENNACL "rowscalingviennacl" 7007598726SKarl Rupp #define PCSAVIENNACL "saviennacl" 71b0753f9dSMatthew G. Knepley #define PCBDDC "bddc" 72b0753f9dSMatthew G. Knepley #define PCKACZMARZ "kaczmarz" 7368ddcbeaSBarry Smith #define PCTELESCOPE "telescope" 744bbf5ea8SMatthew G. Knepley #define PCPATCH "patch" 75b9ac7092SAlp Dener #define PCLMVM "lmvm" 76360ee056SFande Kong #define PCHMG "hmg" 7737eeb815SJakub Kruzik #define PCDEFLATION "deflation" 7838cfc46eSPierre Jolivet #define PCHPDDM "hpddm" 7953022affSStefano Zampini #define PCH2OPUS "h2opus" 80f1f2ae84SBarry Smith #define PCMPI "mpi" 81b0753f9dSMatthew G. Knepley 82b0753f9dSMatthew G. Knepley /*E 83b0753f9dSMatthew G. Knepley PCSide - If the preconditioner is to be applied to the left, right 84b0753f9dSMatthew G. Knepley or symmetrically around the operator. 85b0753f9dSMatthew G. Knepley 86b0753f9dSMatthew G. Knepley Level: beginner 87b0753f9dSMatthew G. Knepley 88b0753f9dSMatthew G. Knepley .seealso: 89b0753f9dSMatthew G. Knepley E*/ 90b0753f9dSMatthew G. Knepley typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide; 91b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 92b0753f9dSMatthew G. Knepley 93b0753f9dSMatthew G. Knepley /*E 94b0753f9dSMatthew G. Knepley PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates 95b0753f9dSMatthew G. Knepley 96b0753f9dSMatthew G. Knepley Level: advanced 97b0753f9dSMatthew G. Knepley 9895452b02SPatrick Sanan Notes: 9995452b02SPatrick Sanan this must match petsc/finclude/petscpc.h and the KSPConvergedReason values in petscksp.h 100b0753f9dSMatthew G. Knepley 101db781477SPatrick Sanan .seealso: `PCApplyRichardson()` 102b0753f9dSMatthew G. Knepley E*/ 103b0753f9dSMatthew G. Knepley typedef enum { 104b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_RTOL = 2, 105b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ATOL = 3, 106b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ITS = 4, 107b0753f9dSMatthew G. Knepley PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason; 108b0753f9dSMatthew G. Knepley 109b0753f9dSMatthew G. Knepley /*E 110b0753f9dSMatthew G. Knepley PCJacobiType - What elements are used to form the Jacobi preconditioner 111b0753f9dSMatthew G. Knepley 112b0753f9dSMatthew G. Knepley Level: intermediate 113b0753f9dSMatthew G. Knepley 114b0753f9dSMatthew G. Knepley .seealso: 115b0753f9dSMatthew G. Knepley E*/ 116b0753f9dSMatthew G. Knepley typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType; 117b0753f9dSMatthew G. Knepley 118b0753f9dSMatthew G. Knepley /*E 119b0753f9dSMatthew G. Knepley PCASMType - Type of additive Schwarz method to use 120b0753f9dSMatthew G. Knepley 121b0753f9dSMatthew G. Knepley $ PC_ASM_BASIC - Symmetric version where residuals from the ghost points are used 122b0753f9dSMatthew G. Knepley $ and computed values in ghost regions are added together. 123b0753f9dSMatthew G. Knepley $ Classical standard additive Schwarz. 124b0753f9dSMatthew G. Knepley $ PC_ASM_RESTRICT - Residuals from ghost points are used but computed values in ghost 125b0753f9dSMatthew G. Knepley $ region are discarded. 126b0753f9dSMatthew G. Knepley $ Default. 127b0753f9dSMatthew G. Knepley $ PC_ASM_INTERPOLATE - Residuals from ghost points are not used, computed values in ghost 128b0753f9dSMatthew G. Knepley $ region are added back in. 129b0753f9dSMatthew G. Knepley $ PC_ASM_NONE - Residuals from ghost points are not used, computed ghost values are 130b0753f9dSMatthew G. Knepley $ discarded. 131b0753f9dSMatthew G. Knepley $ Not very good. 132b0753f9dSMatthew G. Knepley 133b0753f9dSMatthew G. Knepley Level: beginner 134b0753f9dSMatthew G. Knepley 135db781477SPatrick Sanan .seealso: `PCASMSetType()` 136b0753f9dSMatthew G. Knepley E*/ 137b0753f9dSMatthew G. Knepley typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType; 138b0753f9dSMatthew G. Knepley 139b0753f9dSMatthew G. Knepley /*E 140b0753f9dSMatthew G. Knepley PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain). 141b0753f9dSMatthew G. Knepley 142b0753f9dSMatthew G. Knepley Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 143b0753f9dSMatthew G. Knepley domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 144b0753f9dSMatthew G. Knepley a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 145b0753f9dSMatthew G. Knepley (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 146b0753f9dSMatthew G. Knepley 147b0753f9dSMatthew G. Knepley $ PC_GASM_BASIC - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 148b0753f9dSMatthew G. Knepley $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 149b0753f9dSMatthew G. Knepley $ from neighboring subdomains. 150b0753f9dSMatthew G. Knepley $ Classical standard additive Schwarz. 151b0753f9dSMatthew G. Knepley $ PC_GASM_RESTRICT - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 152b0753f9dSMatthew G. Knepley $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 153b0753f9dSMatthew G. Knepley $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 154b0753f9dSMatthew G. Knepley $ assumption). 155b0753f9dSMatthew G. Knepley $ Default. 156b0753f9dSMatthew G. Knepley $ PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 157b0753f9dSMatthew G. Knepley $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 158b0753f9dSMatthew G. Knepley $ from neighboring subdomains. 159b0753f9dSMatthew G. Knepley $ 160b0753f9dSMatthew G. Knepley $ PC_GASM_NONE - Residuals and corrections are zeroed out outside the local subdomains. 161b0753f9dSMatthew G. Knepley $ Not very good. 162b0753f9dSMatthew G. Knepley 163b0753f9dSMatthew G. Knepley Level: beginner 164b0753f9dSMatthew G. Knepley 165db781477SPatrick Sanan .seealso: `PCGASMSetType()` 166b0753f9dSMatthew G. Knepley E*/ 167b0753f9dSMatthew G. Knepley typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType; 168b0753f9dSMatthew G. Knepley 169b0753f9dSMatthew G. Knepley /*E 170b0753f9dSMatthew G. Knepley PCCompositeType - Determines how two or more preconditioner are composed 171b0753f9dSMatthew G. Knepley 172b0753f9dSMatthew G. Knepley $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together 173b0753f9dSMatthew G. Knepley $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 174b0753f9dSMatthew G. Knepley $ computed after the previous preconditioner application 175b0753f9dSMatthew G. Knepley $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 176d0ecd4dfSBarry Smith $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners) 177b0753f9dSMatthew G. Knepley $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S 178b0753f9dSMatthew G. Knepley $ where first preconditioner is built from alpha I + S and second from 179b0753f9dSMatthew G. Knepley $ alpha I + R 180b0753f9dSMatthew G. Knepley 181b0753f9dSMatthew G. Knepley Level: beginner 182b0753f9dSMatthew G. Knepley 183db781477SPatrick Sanan .seealso: `PCCompositeSetType()` 184b0753f9dSMatthew G. Knepley E*/ 185a51937d4SCarola Kruse typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR,PC_COMPOSITE_GKB} PCCompositeType; 186b0753f9dSMatthew G. Knepley 187b0753f9dSMatthew G. Knepley /*E 188b0753f9dSMatthew G. Knepley PCFieldSplitSchurPreType - Determines how to precondition Schur complement 189b0753f9dSMatthew G. Knepley 190b0753f9dSMatthew G. Knepley Level: intermediate 191b0753f9dSMatthew G. Knepley 192db781477SPatrick Sanan .seealso: `PCFieldSplitSetSchurPre()` 193b0753f9dSMatthew G. Knepley E*/ 194b0753f9dSMatthew 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; 195b0753f9dSMatthew G. Knepley 196b0753f9dSMatthew G. Knepley /*E 197b0753f9dSMatthew G. Knepley PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 198b0753f9dSMatthew G. Knepley 199b0753f9dSMatthew G. Knepley Level: intermediate 200b0753f9dSMatthew G. Knepley 201db781477SPatrick Sanan .seealso: `PCFieldSplitSetSchurFactType()` 202b0753f9dSMatthew G. Knepley E*/ 203b0753f9dSMatthew G. Knepley typedef enum { 204b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_DIAG, 205b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_LOWER, 206b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_UPPER, 207b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_FULL 208b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType; 209b0753f9dSMatthew G. Knepley 210b0753f9dSMatthew G. Knepley /*E 211b0753f9dSMatthew G. Knepley PCPARMSGlobalType - Determines the global preconditioner method in PARMS 212b0753f9dSMatthew G. Knepley 213b0753f9dSMatthew G. Knepley Level: intermediate 214b0753f9dSMatthew G. Knepley 215db781477SPatrick Sanan .seealso: `PCPARMSSetGlobal()` 216b0753f9dSMatthew G. Knepley E*/ 217b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType; 2189d8081ecSMatthew G. Knepley 219b0753f9dSMatthew G. Knepley /*E 220b0753f9dSMatthew G. Knepley PCPARMSLocalType - Determines the local preconditioner method in PARMS 221b0753f9dSMatthew G. Knepley 222b0753f9dSMatthew G. Knepley Level: intermediate 223b0753f9dSMatthew G. Knepley 224db781477SPatrick Sanan .seealso: `PCPARMSSetLocal()` 225b0753f9dSMatthew G. Knepley E*/ 226b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType; 227b0753f9dSMatthew G. Knepley 228f0fc11ceSJed Brown /*J 229b0753f9dSMatthew G. Knepley PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method 230b0753f9dSMatthew G. Knepley 231b0753f9dSMatthew G. Knepley Level: intermediate 232b0753f9dSMatthew G. Knepley 233f6130c64SBarry Smith $ PCGAMGAGG - (the default) smoothed aggregation algorithm, robust, very well tested 234f6130c64SBarry Smith $ PCGAMGGEO - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not well tested 235f6130c64SBarry Smith $ PCGAMGCLASSICAL - classical algebraic multigrid preconditioner, incomplete, poorly tested 236f6130c64SBarry Smith 237db781477SPatrick Sanan .seealso: `PCMG`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()` 238f0fc11ceSJed Brown J*/ 239b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType; 240b0753f9dSMatthew G. Knepley #define PCGAMGAGG "agg" 241b0753f9dSMatthew G. Knepley #define PCGAMGGEO "geo" 242b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL "classical" 243b0753f9dSMatthew G. Knepley 244b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType; 245b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT "direct" 246b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard" 247b0753f9dSMatthew G. Knepley 248b0753f9dSMatthew G. Knepley /*E 249b0753f9dSMatthew G. Knepley PCMGType - Determines the type of multigrid method that is run. 250b0753f9dSMatthew G. Knepley 251b0753f9dSMatthew G. Knepley Level: beginner 252b0753f9dSMatthew G. Knepley 253b0753f9dSMatthew G. Knepley Values: 254c1cbb1deSBarry Smith + PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycleType() 255b0753f9dSMatthew G. Knepley . PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are 256b0753f9dSMatthew G. Knepley smoothed before updating the residual. This only uses the 257b0753f9dSMatthew G. Knepley down smoother, in the preconditioner the upper smoother is ignored 258b0753f9dSMatthew G. Knepley . PC_MG_FULL - same as multiplicative except one also performs grid sequencing, 259b0753f9dSMatthew G. Knepley that is starts on the coarsest grid, performs a cycle, interpolates 260b0753f9dSMatthew 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 261b0753f9dSMatthew G. Knepley algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 262b0753f9dSMatthew G. Knepley calls the V-cycle only on the coarser level and has a post-smoother instead. 263b0753f9dSMatthew G. Knepley - PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level 264b0753f9dSMatthew G. Knepley from a finer 265b0753f9dSMatthew G. Knepley 266db781477SPatrick Sanan .seealso: `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()` 267b0753f9dSMatthew G. Knepley 268b0753f9dSMatthew G. Knepley E*/ 269b0753f9dSMatthew G. Knepley typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType; 270b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE; 271b0753f9dSMatthew G. Knepley 272b0753f9dSMatthew G. Knepley /*E 273b0753f9dSMatthew G. Knepley PCMGCycleType - Use V-cycle or W-cycle 274b0753f9dSMatthew G. Knepley 275b0753f9dSMatthew G. Knepley Level: beginner 276b0753f9dSMatthew G. Knepley 277b0753f9dSMatthew G. Knepley Values: 278147403d9SBarry Smith + PC_MG_V_CYCLE - use the v cycle 279147403d9SBarry Smith - PC_MG_W_CYCLE - use the w cycle 280b0753f9dSMatthew G. Knepley 281db781477SPatrick Sanan .seealso: `PCMGSetCycleType()` 282b0753f9dSMatthew G. Knepley 283b0753f9dSMatthew G. Knepley E*/ 284b0753f9dSMatthew G. Knepley typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType; 285b0753f9dSMatthew G. Knepley 286b0753f9dSMatthew G. Knepley /*E 2872134b1e4SBarry Smith PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process 2882134b1e4SBarry Smith 2892134b1e4SBarry Smith Level: beginner 2902134b1e4SBarry Smith 2912134b1e4SBarry Smith Values: 2922134b1e4SBarry Smith + PC_MG_GALERKIN_PMAT - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid 2932134b1e4SBarry Smith . PC_MG_GALERKIN_MAT - computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid 2942134b1e4SBarry Smith . PC_MG_GALERKIN_BOTH - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once 2952134b1e4SBarry Smith - PC_MG_GALERKIN_NONE - neither operator is computed via the Galerkin process, the user must provide the operator 2962134b1e4SBarry Smith 2972134b1e4SBarry Smith Users should never set PC_MG_GALERKIN_EXTERNAL, it is used by GAMG and ML 2982134b1e4SBarry Smith 299db781477SPatrick Sanan .seealso: `PCMGSetCycleType()` 3002134b1e4SBarry Smith 3012134b1e4SBarry Smith E*/ 3022134b1e4SBarry Smith typedef enum { PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, PC_MG_GALERKIN_EXTERNAL} PCMGGalerkinType; 3032134b1e4SBarry Smith 3042134b1e4SBarry Smith /*E 305b0753f9dSMatthew G. Knepley PCExoticType - Face based or wirebasket based coarse grid space 306b0753f9dSMatthew G. Knepley 307b0753f9dSMatthew G. Knepley Level: beginner 308b0753f9dSMatthew G. Knepley 309db781477SPatrick Sanan .seealso: `PCExoticSetType()`, `PCEXOTIC` 310b0753f9dSMatthew G. Knepley E*/ 311b0753f9dSMatthew G. Knepley typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType; 312b0753f9dSMatthew G. Knepley 3138c1cd74cSHong Zhang /*E 314bc960bbfSJed Brown PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains. 315bc960bbfSJed Brown 316bc960bbfSJed Brown Level: intermediate 317bc960bbfSJed Brown 318bc960bbfSJed Brown Values: 319bc960bbfSJed Brown + PC_BDDC_INTERFACE_EXT_DIRICHLET - solves Dirichlet interior problem; this is the standard BDDC algorithm 320bc960bbfSJed Brown - PC_BDDC_INTERFACE_EXT_LUMP - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP" 321bc960bbfSJed Brown 322bc960bbfSJed Brown E*/ 323bc960bbfSJed Brown typedef enum { 324bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_DIRICHLET, 325bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_LUMP 326bc960bbfSJed Brown } PCBDDCInterfaceExtType; 327bc960bbfSJed Brown 328bc960bbfSJed Brown /*E 329f3b08a26SMatthew G. Knepley PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation 330f3b08a26SMatthew G. Knepley 331f3b08a26SMatthew G. Knepley Level: beginner 332f3b08a26SMatthew G. Knepley 333db781477SPatrick Sanan .seealso: `PCMGSetAdaptCoarseSpaceType()`, `PCMG` 334f3b08a26SMatthew G. Knepley E*/ 3352b3cbbdaSStefano Zampini typedef enum { PCMG_ADAPT_NONE, PCMG_ADAPT_POLYNOMIAL, PCMG_ADAPT_HARMONIC, PCMG_ADAPT_EIGENVECTOR, PCMG_ADAPT_GENERALIZED_EIGENVECTOR, PCMG_ADAPT_GDSW } PCMGCoarseSpaceType; 336f3b08a26SMatthew G. Knepley 337f3b08a26SMatthew G. Knepley /*E 338ee6acaa3SMatthew G. Knepley PCPatchConstructType - The algorithm used to construct patches for the preconditioner 3394bbf5ea8SMatthew G. Knepley 3404bbf5ea8SMatthew G. Knepley Level: beginner 3414bbf5ea8SMatthew G. Knepley 342db781477SPatrick Sanan .seealso: `PCPatchSetConstructType()`, `PCEXOTIC` 3434bbf5ea8SMatthew G. Knepley E*/ 344e5b9877fSPatrick Farrell typedef enum {PC_PATCH_STAR, PC_PATCH_VANKA, PC_PATCH_PARDECOMP, PC_PATCH_USER, PC_PATCH_PYTHON} PCPatchConstructType; 3454bbf5ea8SMatthew G. Knepley 3464bbf5ea8SMatthew G. Knepley /*E 347e26ad46dSJakub Kruzik PCDeflationSpaceType - Type of deflation 348e662bc50SJakub Kruzik 349e662bc50SJakub Kruzik Values: 3507e9012c5SJakub Kruzik + PC_DEFLATION_SPACE_HAAR - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off 351e26ad46dSJakub Kruzik . PC_DEFLATION_SPACE_DB2 - MATCOMPOSITE of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet) 352e26ad46dSJakub Kruzik . PC_DEFLATION_SPACE_DB4 - same as above, but with db4 (4 coefficient Daubechies) 353e26ad46dSJakub Kruzik . PC_DEFLATION_SPACE_DB8 - same as above, but with db8 (8 coefficient Daubechies) 354e26ad46dSJakub Kruzik . PC_DEFLATION_SPACE_DB16 - same as above, but with db16 (16 coefficient Daubechies) 355e26ad46dSJakub Kruzik . PC_DEFLATION_SPACE_BIORTH22 - same as above, but with biorthogonal 2.2 (6 coefficients) 356e26ad46dSJakub Kruzik . PC_DEFLATION_SPACE_MEYER - same as above, but with Meyer/FIR (62 coefficients) 357e26ad46dSJakub Kruzik . PC_DEFLATION_SPACE_AGGREGATION - aggregates local indices (given by operator matix distribution) into a subdomain 358e662bc50SJakub Kruzik - PC_DEFLATION_SPACE_USER - indicates space set by user 359e662bc50SJakub Kruzik 360e662bc50SJakub Kruzik Notes: 3617e9012c5SJakub Kruzik Wavelet-based space (except Haar) can be used in multilevel deflation. 362e662bc50SJakub Kruzik 3631fdb00f9SJakub Kruzik Level: intermediate 3641fdb00f9SJakub Kruzik 365db781477SPatrick Sanan .seealso: `PCDeflationSetSpaceToCompute()`, `PCDEFLATION` 366e662bc50SJakub Kruzik E*/ 367e662bc50SJakub Kruzik typedef enum { 368e662bc50SJakub Kruzik PC_DEFLATION_SPACE_HAAR, 369e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB2, 370e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB4, 371e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB8, 372e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB16, 373e662bc50SJakub Kruzik PC_DEFLATION_SPACE_BIORTH22, 374e662bc50SJakub Kruzik PC_DEFLATION_SPACE_MEYER, 375e662bc50SJakub Kruzik PC_DEFLATION_SPACE_AGGREGATION, 376e662bc50SJakub Kruzik PC_DEFLATION_SPACE_USER 377e662bc50SJakub Kruzik } PCDeflationSpaceType; 378e662bc50SJakub Kruzik 379e662bc50SJakub Kruzik /*E 38038cfc46eSPierre Jolivet PCHPDDMCoarseCorrectionType - Type of coarse correction used by PCHPDDM 38138cfc46eSPierre Jolivet 38238cfc46eSPierre Jolivet Level: intermediate 38338cfc46eSPierre Jolivet 38438cfc46eSPierre Jolivet Values: 385e5cd072eSprj- + PC_HPDDM_COARSE_CORRECTION_DEFLATED (default) - eq. (1) in PCHPDDMShellApply() 386e5cd072eSprj- . PC_HPDDM_COARSE_CORRECTION_ADDITIVE - eq. (2) 387e5cd072eSprj- - PC_HPDDM_COARSE_CORRECTION_BALANCED - eq. (3) 38838cfc46eSPierre Jolivet 389db781477SPatrick Sanan .seealso: `PCHPDDM`, `PCSetType()`, `PCHPDDMShellApply()` 39038cfc46eSPierre Jolivet E*/ 39138cfc46eSPierre Jolivet typedef enum { PC_HPDDM_COARSE_CORRECTION_DEFLATED, PC_HPDDM_COARSE_CORRECTION_ADDITIVE, PC_HPDDM_COARSE_CORRECTION_BALANCED } PCHPDDMCoarseCorrectionType; 39238cfc46eSPierre Jolivet 39338cfc46eSPierre Jolivet /*E 3948c1cd74cSHong Zhang PCFailedReason - indicates type of PC failure 3958c1cd74cSHong Zhang 3968c1cd74cSHong Zhang Level: beginner 3978c1cd74cSHong Zhang 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