126bd1501SBarry Smith #if !defined(PETSCPCTYPES_H) 226bd1501SBarry Smith #define PETSCPCTYPES_H 3b0753f9dSMatthew G. Knepley 4b0753f9dSMatthew G. Knepley /*S 5b0753f9dSMatthew G. Knepley PC - Abstract PETSc object that manages all preconditioners including direct solvers such as PCLU 6b0753f9dSMatthew G. Knepley 7b0753f9dSMatthew G. Knepley Level: beginner 8b0753f9dSMatthew G. Knepley 9b0753f9dSMatthew G. Knepley .seealso: PCCreate(), PCSetType(), PCType (for list of available types) 10b0753f9dSMatthew G. Knepley S*/ 11b0753f9dSMatthew G. Knepley typedef struct _p_PC* PC; 12b0753f9dSMatthew G. Knepley 13b0753f9dSMatthew G. Knepley /*J 14b0753f9dSMatthew G. Knepley PCType - String with the name of a PETSc preconditioner method. 15b0753f9dSMatthew G. Knepley 16b0753f9dSMatthew G. Knepley Level: beginner 17b0753f9dSMatthew G. Knepley 1895452b02SPatrick Sanan Notes: 1995452b02SPatrick Sanan Click on the links above to see details on a particular solver 20b0753f9dSMatthew G. Knepley 21b0753f9dSMatthew G. Knepley PCRegister() is used to register preconditioners that are then accessible via PCSetType() 22b0753f9dSMatthew G. Knepley 23b0753f9dSMatthew G. Knepley .seealso: PCSetType(), PC, PCCreate(), PCRegister(), PCSetFromOptions() 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" 32b0753f9dSMatthew G. Knepley #define PCBJACOBI "bjacobi" 33b0753f9dSMatthew G. Knepley #define PCMG "mg" 34b0753f9dSMatthew G. Knepley #define PCEISENSTAT "eisenstat" 35b0753f9dSMatthew G. Knepley #define PCILU "ilu" 36b0753f9dSMatthew G. Knepley #define PCICC "icc" 37b0753f9dSMatthew G. Knepley #define PCASM "asm" 38b0753f9dSMatthew G. Knepley #define PCGASM "gasm" 39b0753f9dSMatthew G. Knepley #define PCKSP "ksp" 40b0753f9dSMatthew G. Knepley #define PCCOMPOSITE "composite" 41b0753f9dSMatthew G. Knepley #define PCREDUNDANT "redundant" 42b0753f9dSMatthew G. Knepley #define PCSPAI "spai" 43b0753f9dSMatthew G. Knepley #define PCNN "nn" 44b0753f9dSMatthew G. Knepley #define PCCHOLESKY "cholesky" 45b0753f9dSMatthew G. Knepley #define PCPBJACOBI "pbjacobi" 460da83c2eSBarry Smith #define PCVPBJACOBI "vpbjacobi" 47b0753f9dSMatthew G. Knepley #define PCMAT "mat" 48b0753f9dSMatthew G. Knepley #define PCHYPRE "hypre" 49b0753f9dSMatthew G. Knepley #define PCPARMS "parms" 50b0753f9dSMatthew G. Knepley #define PCFIELDSPLIT "fieldsplit" 51b0753f9dSMatthew G. Knepley #define PCTFS "tfs" 52b0753f9dSMatthew G. Knepley #define PCML "ml" 53b0753f9dSMatthew G. Knepley #define PCGALERKIN "galerkin" 54b0753f9dSMatthew G. Knepley #define PCEXOTIC "exotic" 55b0753f9dSMatthew G. Knepley #define PCCP "cp" 56b0753f9dSMatthew G. Knepley #define PCBFBT "bfbt" 57b0753f9dSMatthew G. Knepley #define PCLSC "lsc" 58b0753f9dSMatthew G. Knepley #define PCPYTHON "python" 59b0753f9dSMatthew G. Knepley #define PCPFMG "pfmg" 60b0753f9dSMatthew G. Knepley #define PCSYSPFMG "syspfmg" 61b0753f9dSMatthew G. Knepley #define PCREDISTRIBUTE "redistribute" 62b0753f9dSMatthew G. Knepley #define PCSVD "svd" 63b0753f9dSMatthew G. Knepley #define PCGAMG "gamg" 644b3f184cSKarl Rupp #define PCCHOWILUVIENNACL "chowiluviennacl" 6570baa948SKarl Rupp #define PCROWSCALINGVIENNACL "rowscalingviennacl" 6607598726SKarl Rupp #define PCSAVIENNACL "saviennacl" 67b0753f9dSMatthew G. Knepley #define PCBDDC "bddc" 68b0753f9dSMatthew G. Knepley #define PCKACZMARZ "kaczmarz" 6968ddcbeaSBarry Smith #define PCTELESCOPE "telescope" 704bbf5ea8SMatthew G. Knepley #define PCPATCH "patch" 71b9ac7092SAlp Dener #define PCLMVM "lmvm" 72360ee056SFande Kong #define PCHMG "hmg" 7337eeb815SJakub Kruzik #define PCDEFLATION "deflation" 7438cfc46eSPierre Jolivet #define PCHPDDM "hpddm" 75*53022affSStefano Zampini #define PCH2OPUS "h2opus" 76b0753f9dSMatthew G. Knepley 77b0753f9dSMatthew G. Knepley /*E 78b0753f9dSMatthew G. Knepley PCSide - If the preconditioner is to be applied to the left, right 79b0753f9dSMatthew G. Knepley or symmetrically around the operator. 80b0753f9dSMatthew G. Knepley 81b0753f9dSMatthew G. Knepley Level: beginner 82b0753f9dSMatthew G. Knepley 83b0753f9dSMatthew G. Knepley .seealso: 84b0753f9dSMatthew G. Knepley E*/ 85b0753f9dSMatthew G. Knepley typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide; 86b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 87b0753f9dSMatthew G. Knepley 88b0753f9dSMatthew G. Knepley /*E 89b0753f9dSMatthew G. Knepley PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates 90b0753f9dSMatthew G. Knepley 91b0753f9dSMatthew G. Knepley Level: advanced 92b0753f9dSMatthew G. Knepley 9395452b02SPatrick Sanan Notes: 9495452b02SPatrick Sanan this must match petsc/finclude/petscpc.h and the KSPConvergedReason values in petscksp.h 95b0753f9dSMatthew G. Knepley 96b0753f9dSMatthew G. Knepley .seealso: PCApplyRichardson() 97b0753f9dSMatthew G. Knepley E*/ 98b0753f9dSMatthew G. Knepley typedef enum { 99b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_RTOL = 2, 100b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ATOL = 3, 101b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ITS = 4, 102b0753f9dSMatthew G. Knepley PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason; 103b0753f9dSMatthew G. Knepley 104b0753f9dSMatthew G. Knepley /*E 105b0753f9dSMatthew G. Knepley PCJacobiType - What elements are used to form the Jacobi preconditioner 106b0753f9dSMatthew G. Knepley 107b0753f9dSMatthew G. Knepley Level: intermediate 108b0753f9dSMatthew G. Knepley 109b0753f9dSMatthew G. Knepley .seealso: 110b0753f9dSMatthew G. Knepley E*/ 111b0753f9dSMatthew G. Knepley typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType; 112b0753f9dSMatthew G. Knepley 113b0753f9dSMatthew G. Knepley /*E 114b0753f9dSMatthew G. Knepley PCASMType - Type of additive Schwarz method to use 115b0753f9dSMatthew G. Knepley 116b0753f9dSMatthew G. Knepley $ PC_ASM_BASIC - Symmetric version where residuals from the ghost points are used 117b0753f9dSMatthew G. Knepley $ and computed values in ghost regions are added together. 118b0753f9dSMatthew G. Knepley $ Classical standard additive Schwarz. 119b0753f9dSMatthew G. Knepley $ PC_ASM_RESTRICT - Residuals from ghost points are used but computed values in ghost 120b0753f9dSMatthew G. Knepley $ region are discarded. 121b0753f9dSMatthew G. Knepley $ Default. 122b0753f9dSMatthew G. Knepley $ PC_ASM_INTERPOLATE - Residuals from ghost points are not used, computed values in ghost 123b0753f9dSMatthew G. Knepley $ region are added back in. 124b0753f9dSMatthew G. Knepley $ PC_ASM_NONE - Residuals from ghost points are not used, computed ghost values are 125b0753f9dSMatthew G. Knepley $ discarded. 126b0753f9dSMatthew G. Knepley $ Not very good. 127b0753f9dSMatthew G. Knepley 128b0753f9dSMatthew G. Knepley Level: beginner 129b0753f9dSMatthew G. Knepley 130b0753f9dSMatthew G. Knepley .seealso: PCASMSetType() 131b0753f9dSMatthew G. Knepley E*/ 132b0753f9dSMatthew G. Knepley typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType; 133b0753f9dSMatthew G. Knepley 134b0753f9dSMatthew G. Knepley /*E 135b0753f9dSMatthew G. Knepley PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain). 136b0753f9dSMatthew G. Knepley 137b0753f9dSMatthew G. Knepley Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 138b0753f9dSMatthew G. Knepley domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 139b0753f9dSMatthew G. Knepley a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 140b0753f9dSMatthew G. Knepley (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 141b0753f9dSMatthew G. Knepley 142b0753f9dSMatthew G. Knepley $ PC_GASM_BASIC - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 143b0753f9dSMatthew G. Knepley $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 144b0753f9dSMatthew G. Knepley $ from neighboring subdomains. 145b0753f9dSMatthew G. Knepley $ Classical standard additive Schwarz. 146b0753f9dSMatthew G. Knepley $ PC_GASM_RESTRICT - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 147b0753f9dSMatthew G. Knepley $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 148b0753f9dSMatthew G. Knepley $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 149b0753f9dSMatthew G. Knepley $ assumption). 150b0753f9dSMatthew G. Knepley $ Default. 151b0753f9dSMatthew G. Knepley $ PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 152b0753f9dSMatthew G. Knepley $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 153b0753f9dSMatthew G. Knepley $ from neighboring subdomains. 154b0753f9dSMatthew G. Knepley $ 155b0753f9dSMatthew G. Knepley $ PC_GASM_NONE - Residuals and corrections are zeroed out outside the local subdomains. 156b0753f9dSMatthew G. Knepley $ Not very good. 157b0753f9dSMatthew G. Knepley 158b0753f9dSMatthew G. Knepley Level: beginner 159b0753f9dSMatthew G. Knepley 160b0753f9dSMatthew G. Knepley .seealso: PCGASMSetType() 161b0753f9dSMatthew G. Knepley E*/ 162b0753f9dSMatthew G. Knepley typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType; 163b0753f9dSMatthew G. Knepley 164b0753f9dSMatthew G. Knepley /*E 165b0753f9dSMatthew G. Knepley PCCompositeType - Determines how two or more preconditioner are composed 166b0753f9dSMatthew G. Knepley 167b0753f9dSMatthew G. Knepley $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together 168b0753f9dSMatthew G. Knepley $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 169b0753f9dSMatthew G. Knepley $ computed after the previous preconditioner application 170b0753f9dSMatthew G. Knepley $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 171d0ecd4dfSBarry Smith $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners) 172b0753f9dSMatthew G. Knepley $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S 173b0753f9dSMatthew G. Knepley $ where first preconditioner is built from alpha I + S and second from 174b0753f9dSMatthew G. Knepley $ alpha I + R 175b0753f9dSMatthew G. Knepley 176b0753f9dSMatthew G. Knepley Level: beginner 177b0753f9dSMatthew G. Knepley 178b0753f9dSMatthew G. Knepley .seealso: PCCompositeSetType() 179b0753f9dSMatthew G. Knepley E*/ 180a51937d4SCarola Kruse typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR,PC_COMPOSITE_GKB} PCCompositeType; 181b0753f9dSMatthew G. Knepley 182b0753f9dSMatthew G. Knepley /*E 183b0753f9dSMatthew G. Knepley PCFieldSplitSchurPreType - Determines how to precondition Schur complement 184b0753f9dSMatthew G. Knepley 185b0753f9dSMatthew G. Knepley Level: intermediate 186b0753f9dSMatthew G. Knepley 187b0753f9dSMatthew G. Knepley .seealso: PCFieldSplitSetSchurPre() 188b0753f9dSMatthew G. Knepley E*/ 189b0753f9dSMatthew 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; 190b0753f9dSMatthew G. Knepley 191b0753f9dSMatthew G. Knepley /*E 192b0753f9dSMatthew G. Knepley PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 193b0753f9dSMatthew G. Knepley 194b0753f9dSMatthew G. Knepley Level: intermediate 195b0753f9dSMatthew G. Knepley 196b0753f9dSMatthew G. Knepley .seealso: PCFieldSplitSetSchurFactType() 197b0753f9dSMatthew G. Knepley E*/ 198b0753f9dSMatthew G. Knepley typedef enum { 199b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_DIAG, 200b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_LOWER, 201b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_UPPER, 202b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_FULL 203b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType; 204b0753f9dSMatthew G. Knepley 205b0753f9dSMatthew G. Knepley /*E 206b0753f9dSMatthew G. Knepley PCPARMSGlobalType - Determines the global preconditioner method in PARMS 207b0753f9dSMatthew G. Knepley 208b0753f9dSMatthew G. Knepley Level: intermediate 209b0753f9dSMatthew G. Knepley 210b0753f9dSMatthew G. Knepley .seealso: PCPARMSSetGlobal() 211b0753f9dSMatthew G. Knepley E*/ 212b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType; 2139d8081ecSMatthew G. Knepley 214b0753f9dSMatthew G. Knepley /*E 215b0753f9dSMatthew G. Knepley PCPARMSLocalType - Determines the local preconditioner method in PARMS 216b0753f9dSMatthew G. Knepley 217b0753f9dSMatthew G. Knepley Level: intermediate 218b0753f9dSMatthew G. Knepley 219b0753f9dSMatthew G. Knepley .seealso: PCPARMSSetLocal() 220b0753f9dSMatthew G. Knepley E*/ 221b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType; 222b0753f9dSMatthew G. Knepley 223f0fc11ceSJed Brown /*J 224b0753f9dSMatthew G. Knepley PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method 225b0753f9dSMatthew G. Knepley 226b0753f9dSMatthew G. Knepley Level: intermediate 227b0753f9dSMatthew G. Knepley 228f6130c64SBarry Smith $ PCGAMGAGG - (the default) smoothed aggregation algorithm, robust, very well tested 229f6130c64SBarry Smith $ PCGAMGGEO - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not well tested 230f6130c64SBarry Smith $ PCGAMGCLASSICAL - classical algebraic multigrid preconditioner, incomplete, poorly tested 231f6130c64SBarry Smith 232b0753f9dSMatthew G. Knepley .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseInterpolation() 233f0fc11ceSJed Brown J*/ 234b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType; 235b0753f9dSMatthew G. Knepley #define PCGAMGAGG "agg" 236b0753f9dSMatthew G. Knepley #define PCGAMGGEO "geo" 237b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL "classical" 238b0753f9dSMatthew G. Knepley 239b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType; 240b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT "direct" 241b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard" 242b0753f9dSMatthew G. Knepley 243b0753f9dSMatthew G. Knepley /*E 244b0753f9dSMatthew G. Knepley PCMGType - Determines the type of multigrid method that is run. 245b0753f9dSMatthew G. Knepley 246b0753f9dSMatthew G. Knepley Level: beginner 247b0753f9dSMatthew G. Knepley 248b0753f9dSMatthew G. Knepley Values: 249c1cbb1deSBarry Smith + PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycleType() 250b0753f9dSMatthew G. Knepley . PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are 251b0753f9dSMatthew G. Knepley smoothed before updating the residual. This only uses the 252b0753f9dSMatthew G. Knepley down smoother, in the preconditioner the upper smoother is ignored 253b0753f9dSMatthew G. Knepley . PC_MG_FULL - same as multiplicative except one also performs grid sequencing, 254b0753f9dSMatthew G. Knepley that is starts on the coarsest grid, performs a cycle, interpolates 255b0753f9dSMatthew 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 256b0753f9dSMatthew G. Knepley algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 257b0753f9dSMatthew G. Knepley calls the V-cycle only on the coarser level and has a post-smoother instead. 258b0753f9dSMatthew G. Knepley - PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level 259b0753f9dSMatthew G. Knepley from a finer 260b0753f9dSMatthew G. Knepley 261c1cbb1deSBarry Smith .seealso: PCMGSetType(), PCMGSetCycleType(), PCMGSetCycleTypeOnLevel() 262b0753f9dSMatthew G. Knepley 263b0753f9dSMatthew G. Knepley E*/ 264b0753f9dSMatthew G. Knepley typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType; 265b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE; 266b0753f9dSMatthew G. Knepley 267b0753f9dSMatthew G. Knepley /*E 268b0753f9dSMatthew G. Knepley PCMGCycleType - Use V-cycle or W-cycle 269b0753f9dSMatthew G. Knepley 270b0753f9dSMatthew G. Knepley Level: beginner 271b0753f9dSMatthew G. Knepley 272b0753f9dSMatthew G. Knepley Values: 273b0753f9dSMatthew G. Knepley + PC_MG_V_CYCLE 274b0753f9dSMatthew G. Knepley - PC_MG_W_CYCLE 275b0753f9dSMatthew G. Knepley 276b0753f9dSMatthew G. Knepley .seealso: PCMGSetCycleType() 277b0753f9dSMatthew G. Knepley 278b0753f9dSMatthew G. Knepley E*/ 279b0753f9dSMatthew G. Knepley typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType; 280b0753f9dSMatthew G. Knepley 281b0753f9dSMatthew G. Knepley /*E 2822134b1e4SBarry Smith PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process 2832134b1e4SBarry Smith 2842134b1e4SBarry Smith Level: beginner 2852134b1e4SBarry Smith 2862134b1e4SBarry Smith Values: 2872134b1e4SBarry Smith + PC_MG_GALERKIN_PMAT - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid 2882134b1e4SBarry Smith . PC_MG_GALERKIN_MAT - computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid 2892134b1e4SBarry Smith . PC_MG_GALERKIN_BOTH - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once 2902134b1e4SBarry Smith - PC_MG_GALERKIN_NONE - neither operator is computed via the Galerkin process, the user must provide the operator 2912134b1e4SBarry Smith 2922134b1e4SBarry Smith Users should never set PC_MG_GALERKIN_EXTERNAL, it is used by GAMG and ML 2932134b1e4SBarry Smith 2942134b1e4SBarry Smith .seealso: PCMGSetCycleType() 2952134b1e4SBarry Smith 2962134b1e4SBarry Smith E*/ 2972134b1e4SBarry Smith typedef enum { PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, PC_MG_GALERKIN_EXTERNAL} PCMGGalerkinType; 2982134b1e4SBarry Smith 2992134b1e4SBarry Smith /*E 300b0753f9dSMatthew G. Knepley PCExoticType - Face based or wirebasket based coarse grid space 301b0753f9dSMatthew G. Knepley 302b0753f9dSMatthew G. Knepley Level: beginner 303b0753f9dSMatthew G. Knepley 304b0753f9dSMatthew G. Knepley .seealso: PCExoticSetType(), PCEXOTIC 305b0753f9dSMatthew G. Knepley E*/ 306b0753f9dSMatthew G. Knepley typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType; 307b0753f9dSMatthew G. Knepley 3088c1cd74cSHong Zhang /*E 309bc960bbfSJed Brown PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains. 310bc960bbfSJed Brown 311bc960bbfSJed Brown Level: intermediate 312bc960bbfSJed Brown 313bc960bbfSJed Brown Values: 314bc960bbfSJed Brown + PC_BDDC_INTERFACE_EXT_DIRICHLET - solves Dirichlet interior problem; this is the standard BDDC algorithm 315bc960bbfSJed Brown - PC_BDDC_INTERFACE_EXT_LUMP - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP" 316bc960bbfSJed Brown 317bc960bbfSJed Brown E*/ 318bc960bbfSJed Brown typedef enum { 319bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_DIRICHLET, 320bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_LUMP 321bc960bbfSJed Brown } PCBDDCInterfaceExtType; 322bc960bbfSJed Brown 323bc960bbfSJed Brown /*E 324f3b08a26SMatthew G. Knepley PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation 325f3b08a26SMatthew G. Knepley 326f3b08a26SMatthew G. Knepley Level: beginner 327f3b08a26SMatthew G. Knepley 328f3b08a26SMatthew G. Knepley .seealso: PCMGSetAdaptCoarseSpaceType(), PCMG 329f3b08a26SMatthew G. Knepley E*/ 330f3b08a26SMatthew G. Knepley typedef enum { PCMG_POLYNOMIAL, PCMG_HARMONIC, PCMG_EIGENVECTOR, PCMG_GENERALIZED_EIGENVECTOR } PCMGCoarseSpaceType; 331f3b08a26SMatthew G. Knepley 332f3b08a26SMatthew G. Knepley /*E 333ee6acaa3SMatthew G. Knepley PCPatchConstructType - The algorithm used to construct patches for the preconditioner 3344bbf5ea8SMatthew G. Knepley 3354bbf5ea8SMatthew G. Knepley Level: beginner 3364bbf5ea8SMatthew G. Knepley 3374bbf5ea8SMatthew G. Knepley .seealso: PCPatchSetConstructType(), PCEXOTIC 3384bbf5ea8SMatthew G. Knepley E*/ 339e5b9877fSPatrick Farrell typedef enum {PC_PATCH_STAR, PC_PATCH_VANKA, PC_PATCH_PARDECOMP, PC_PATCH_USER, PC_PATCH_PYTHON} PCPatchConstructType; 3404bbf5ea8SMatthew G. Knepley 3414bbf5ea8SMatthew G. Knepley /*E 342e26ad46dSJakub Kruzik PCDeflationSpaceType - Type of deflation 343e662bc50SJakub Kruzik 344e662bc50SJakub Kruzik Values: 3457e9012c5SJakub Kruzik + PC_DEFLATION_SPACE_HAAR - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off 346e26ad46dSJakub Kruzik . PC_DEFLATION_SPACE_DB2 - MATCOMPOSITE of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet) 347e26ad46dSJakub Kruzik . PC_DEFLATION_SPACE_DB4 - same as above, but with db4 (4 coefficient Daubechies) 348e26ad46dSJakub Kruzik . PC_DEFLATION_SPACE_DB8 - same as above, but with db8 (8 coefficient Daubechies) 349e26ad46dSJakub Kruzik . PC_DEFLATION_SPACE_DB16 - same as above, but with db16 (16 coefficient Daubechies) 350e26ad46dSJakub Kruzik . PC_DEFLATION_SPACE_BIORTH22 - same as above, but with biorthogonal 2.2 (6 coefficients) 351e26ad46dSJakub Kruzik . PC_DEFLATION_SPACE_MEYER - same as above, but with Meyer/FIR (62 coefficients) 352e26ad46dSJakub Kruzik . PC_DEFLATION_SPACE_AGGREGATION - aggregates local indices (given by operator matix distribution) into a subdomain 353e662bc50SJakub Kruzik - PC_DEFLATION_SPACE_USER - indicates space set by user 354e662bc50SJakub Kruzik 355e662bc50SJakub Kruzik Notes: 3567e9012c5SJakub Kruzik Wavelet-based space (except Haar) can be used in multilevel deflation. 357e662bc50SJakub Kruzik 3581fdb00f9SJakub Kruzik Level: intermediate 3591fdb00f9SJakub Kruzik 3601fdb00f9SJakub Kruzik .seealso: PCDeflationSetSpaceToCompute(), PCDEFLATION 361e662bc50SJakub Kruzik E*/ 362e662bc50SJakub Kruzik typedef enum { 363e662bc50SJakub Kruzik PC_DEFLATION_SPACE_HAAR, 364e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB2, 365e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB4, 366e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB8, 367e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB16, 368e662bc50SJakub Kruzik PC_DEFLATION_SPACE_BIORTH22, 369e662bc50SJakub Kruzik PC_DEFLATION_SPACE_MEYER, 370e662bc50SJakub Kruzik PC_DEFLATION_SPACE_AGGREGATION, 371e662bc50SJakub Kruzik PC_DEFLATION_SPACE_USER 372e662bc50SJakub Kruzik } PCDeflationSpaceType; 373e662bc50SJakub Kruzik 374e662bc50SJakub Kruzik /*E 37538cfc46eSPierre Jolivet PCHPDDMCoarseCorrectionType - Type of coarse correction used by PCHPDDM 37638cfc46eSPierre Jolivet 37738cfc46eSPierre Jolivet Level: intermediate 37838cfc46eSPierre Jolivet 37938cfc46eSPierre Jolivet Values: 380e5cd072eSprj- + PC_HPDDM_COARSE_CORRECTION_DEFLATED (default) - eq. (1) in PCHPDDMShellApply() 381e5cd072eSprj- . PC_HPDDM_COARSE_CORRECTION_ADDITIVE - eq. (2) 382e5cd072eSprj- - PC_HPDDM_COARSE_CORRECTION_BALANCED - eq. (3) 38338cfc46eSPierre Jolivet 38438cfc46eSPierre Jolivet .seealso: PCHPDDM, PCSetType(), PCHPDDMShellApply() 38538cfc46eSPierre Jolivet E*/ 38638cfc46eSPierre Jolivet typedef enum { PC_HPDDM_COARSE_CORRECTION_DEFLATED, PC_HPDDM_COARSE_CORRECTION_ADDITIVE, PC_HPDDM_COARSE_CORRECTION_BALANCED } PCHPDDMCoarseCorrectionType; 38738cfc46eSPierre Jolivet 38838cfc46eSPierre Jolivet /*E 3898c1cd74cSHong Zhang PCFailedReason - indicates type of PC failure 3908c1cd74cSHong Zhang 3918c1cd74cSHong Zhang Level: beginner 3928c1cd74cSHong Zhang 3938c1cd74cSHong Zhang Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h 3948c1cd74cSHong Zhang E*/ 3951b2b9847SBarry 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; 396ce7c7f2fSMark Adams 397ce7c7f2fSMark Adams /*E 398ce7c7f2fSMark Adams PCGAMGLayoutType - Layout for reduced grids 399ce7c7f2fSMark Adams 400ce7c7f2fSMark Adams Level: intermediate 401ce7c7f2fSMark Adams 402f73ac439SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType() 403ce7c7f2fSMark Adams Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h 404ce7c7f2fSMark Adams E*/ 405a0095786SMark typedef enum {PCGAMG_LAYOUT_COMPACT,PCGAMG_LAYOUT_SPREAD} PCGAMGLayoutType; 406ce7c7f2fSMark Adams 407b0753f9dSMatthew G. Knepley #endif 408