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" 30b0753f9dSMatthew G. Knepley #define PCSHELL "shell" 31b0753f9dSMatthew G. Knepley #define PCBJACOBI "bjacobi" 32b0753f9dSMatthew G. Knepley #define PCMG "mg" 33b0753f9dSMatthew G. Knepley #define PCEISENSTAT "eisenstat" 34b0753f9dSMatthew G. Knepley #define PCILU "ilu" 35b0753f9dSMatthew G. Knepley #define PCICC "icc" 36b0753f9dSMatthew G. Knepley #define PCASM "asm" 37b0753f9dSMatthew G. Knepley #define PCGASM "gasm" 38b0753f9dSMatthew G. Knepley #define PCKSP "ksp" 39b0753f9dSMatthew G. Knepley #define PCCOMPOSITE "composite" 40b0753f9dSMatthew G. Knepley #define PCREDUNDANT "redundant" 41b0753f9dSMatthew G. Knepley #define PCSPAI "spai" 42b0753f9dSMatthew G. Knepley #define PCNN "nn" 43b0753f9dSMatthew G. Knepley #define PCCHOLESKY "cholesky" 44b0753f9dSMatthew G. Knepley #define PCPBJACOBI "pbjacobi" 450da83c2eSBarry Smith #define PCVPBJACOBI "vpbjacobi" 46b0753f9dSMatthew G. Knepley #define PCMAT "mat" 47b0753f9dSMatthew G. Knepley #define PCHYPRE "hypre" 48b0753f9dSMatthew G. Knepley #define PCPARMS "parms" 49b0753f9dSMatthew G. Knepley #define PCFIELDSPLIT "fieldsplit" 50b0753f9dSMatthew G. Knepley #define PCTFS "tfs" 51b0753f9dSMatthew G. Knepley #define PCML "ml" 52b0753f9dSMatthew G. Knepley #define PCGALERKIN "galerkin" 53b0753f9dSMatthew G. Knepley #define PCEXOTIC "exotic" 54b0753f9dSMatthew G. Knepley #define PCCP "cp" 55b0753f9dSMatthew G. Knepley #define PCBFBT "bfbt" 56b0753f9dSMatthew G. Knepley #define PCLSC "lsc" 57b0753f9dSMatthew G. Knepley #define PCPYTHON "python" 58b0753f9dSMatthew G. Knepley #define PCPFMG "pfmg" 59b0753f9dSMatthew G. Knepley #define PCSYSPFMG "syspfmg" 60b0753f9dSMatthew G. Knepley #define PCREDISTRIBUTE "redistribute" 61b0753f9dSMatthew G. Knepley #define PCSVD "svd" 62b0753f9dSMatthew G. Knepley #define PCGAMG "gamg" 634b3f184cSKarl Rupp #define PCCHOWILUVIENNACL "chowiluviennacl" 6470baa948SKarl Rupp #define PCROWSCALINGVIENNACL "rowscalingviennacl" 6507598726SKarl Rupp #define PCSAVIENNACL "saviennacl" 66b0753f9dSMatthew G. Knepley #define PCBDDC "bddc" 67b0753f9dSMatthew G. Knepley #define PCKACZMARZ "kaczmarz" 6868ddcbeaSBarry Smith #define PCTELESCOPE "telescope" 694bbf5ea8SMatthew G. Knepley #define PCPATCH "patch" 70b9ac7092SAlp Dener #define PCLMVM "lmvm" 71360ee056SFande Kong #define PCHMG "hmg" 72*37eeb815SJakub Kruzik #define PCDEFLATION "deflation" 73b0753f9dSMatthew G. Knepley 74b0753f9dSMatthew G. Knepley /*E 75b0753f9dSMatthew G. Knepley PCSide - If the preconditioner is to be applied to the left, right 76b0753f9dSMatthew G. Knepley or symmetrically around the operator. 77b0753f9dSMatthew G. Knepley 78b0753f9dSMatthew G. Knepley Level: beginner 79b0753f9dSMatthew G. Knepley 80b0753f9dSMatthew G. Knepley .seealso: 81b0753f9dSMatthew G. Knepley E*/ 82b0753f9dSMatthew G. Knepley typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide; 83b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 84b0753f9dSMatthew G. Knepley 85b0753f9dSMatthew G. Knepley /*E 86b0753f9dSMatthew G. Knepley PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates 87b0753f9dSMatthew G. Knepley 88b0753f9dSMatthew G. Knepley Level: advanced 89b0753f9dSMatthew G. Knepley 9095452b02SPatrick Sanan Notes: 9195452b02SPatrick Sanan this must match petsc/finclude/petscpc.h and the KSPConvergedReason values in petscksp.h 92b0753f9dSMatthew G. Knepley 93b0753f9dSMatthew G. Knepley .seealso: PCApplyRichardson() 94b0753f9dSMatthew G. Knepley E*/ 95b0753f9dSMatthew G. Knepley typedef enum { 96b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_RTOL = 2, 97b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ATOL = 3, 98b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ITS = 4, 99b0753f9dSMatthew G. Knepley PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason; 100b0753f9dSMatthew G. Knepley 101b0753f9dSMatthew G. Knepley /*E 102b0753f9dSMatthew G. Knepley PCJacobiType - What elements are used to form the Jacobi preconditioner 103b0753f9dSMatthew G. Knepley 104b0753f9dSMatthew G. Knepley Level: intermediate 105b0753f9dSMatthew G. Knepley 106b0753f9dSMatthew G. Knepley .seealso: 107b0753f9dSMatthew G. Knepley E*/ 108b0753f9dSMatthew G. Knepley typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType; 109b0753f9dSMatthew G. Knepley 110b0753f9dSMatthew G. Knepley /*E 111b0753f9dSMatthew G. Knepley PCASMType - Type of additive Schwarz method to use 112b0753f9dSMatthew G. Knepley 113b0753f9dSMatthew G. Knepley $ PC_ASM_BASIC - Symmetric version where residuals from the ghost points are used 114b0753f9dSMatthew G. Knepley $ and computed values in ghost regions are added together. 115b0753f9dSMatthew G. Knepley $ Classical standard additive Schwarz. 116b0753f9dSMatthew G. Knepley $ PC_ASM_RESTRICT - Residuals from ghost points are used but computed values in ghost 117b0753f9dSMatthew G. Knepley $ region are discarded. 118b0753f9dSMatthew G. Knepley $ Default. 119b0753f9dSMatthew G. Knepley $ PC_ASM_INTERPOLATE - Residuals from ghost points are not used, computed values in ghost 120b0753f9dSMatthew G. Knepley $ region are added back in. 121b0753f9dSMatthew G. Knepley $ PC_ASM_NONE - Residuals from ghost points are not used, computed ghost values are 122b0753f9dSMatthew G. Knepley $ discarded. 123b0753f9dSMatthew G. Knepley $ Not very good. 124b0753f9dSMatthew G. Knepley 125b0753f9dSMatthew G. Knepley Level: beginner 126b0753f9dSMatthew G. Knepley 127b0753f9dSMatthew G. Knepley .seealso: PCASMSetType() 128b0753f9dSMatthew G. Knepley E*/ 129b0753f9dSMatthew G. Knepley typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType; 130b0753f9dSMatthew G. Knepley 131b0753f9dSMatthew G. Knepley /*E 132b0753f9dSMatthew G. Knepley PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain). 133b0753f9dSMatthew G. Knepley 134b0753f9dSMatthew G. Knepley Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 135b0753f9dSMatthew G. Knepley domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 136b0753f9dSMatthew G. Knepley a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 137b0753f9dSMatthew G. Knepley (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 138b0753f9dSMatthew G. Knepley 139b0753f9dSMatthew G. Knepley $ PC_GASM_BASIC - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 140b0753f9dSMatthew G. Knepley $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 141b0753f9dSMatthew G. Knepley $ from neighboring subdomains. 142b0753f9dSMatthew G. Knepley $ Classical standard additive Schwarz. 143b0753f9dSMatthew G. Knepley $ PC_GASM_RESTRICT - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 144b0753f9dSMatthew G. Knepley $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 145b0753f9dSMatthew G. Knepley $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 146b0753f9dSMatthew G. Knepley $ assumption). 147b0753f9dSMatthew G. Knepley $ Default. 148b0753f9dSMatthew G. Knepley $ PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 149b0753f9dSMatthew G. Knepley $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 150b0753f9dSMatthew G. Knepley $ from neighboring subdomains. 151b0753f9dSMatthew G. Knepley $ 152b0753f9dSMatthew G. Knepley $ PC_GASM_NONE - Residuals and corrections are zeroed out outside the local subdomains. 153b0753f9dSMatthew G. Knepley $ Not very good. 154b0753f9dSMatthew G. Knepley 155b0753f9dSMatthew G. Knepley Level: beginner 156b0753f9dSMatthew G. Knepley 157b0753f9dSMatthew G. Knepley .seealso: PCGASMSetType() 158b0753f9dSMatthew G. Knepley E*/ 159b0753f9dSMatthew G. Knepley typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType; 160b0753f9dSMatthew G. Knepley 161b0753f9dSMatthew G. Knepley /*E 162b0753f9dSMatthew G. Knepley PCCompositeType - Determines how two or more preconditioner are composed 163b0753f9dSMatthew G. Knepley 164b0753f9dSMatthew G. Knepley $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together 165b0753f9dSMatthew G. Knepley $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 166b0753f9dSMatthew G. Knepley $ computed after the previous preconditioner application 167b0753f9dSMatthew G. Knepley $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly 168d0ecd4dfSBarry Smith $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners) 169b0753f9dSMatthew G. Knepley $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S 170b0753f9dSMatthew G. Knepley $ where first preconditioner is built from alpha I + S and second from 171b0753f9dSMatthew G. Knepley $ alpha I + R 172b0753f9dSMatthew G. Knepley 173b0753f9dSMatthew G. Knepley Level: beginner 174b0753f9dSMatthew G. Knepley 175b0753f9dSMatthew G. Knepley .seealso: PCCompositeSetType() 176b0753f9dSMatthew G. Knepley E*/ 177a51937d4SCarola Kruse typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR,PC_COMPOSITE_GKB} PCCompositeType; 178b0753f9dSMatthew G. Knepley 179b0753f9dSMatthew G. Knepley /*E 180b0753f9dSMatthew G. Knepley PCFieldSplitSchurPreType - Determines how to precondition Schur complement 181b0753f9dSMatthew G. Knepley 182b0753f9dSMatthew G. Knepley Level: intermediate 183b0753f9dSMatthew G. Knepley 184b0753f9dSMatthew G. Knepley .seealso: PCFieldSplitSetSchurPre() 185b0753f9dSMatthew G. Knepley E*/ 186b0753f9dSMatthew 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; 187b0753f9dSMatthew G. Knepley 188b0753f9dSMatthew G. Knepley /*E 189b0753f9dSMatthew G. Knepley PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 190b0753f9dSMatthew G. Knepley 191b0753f9dSMatthew G. Knepley Level: intermediate 192b0753f9dSMatthew G. Knepley 193b0753f9dSMatthew G. Knepley .seealso: PCFieldSplitSetSchurFactType() 194b0753f9dSMatthew G. Knepley E*/ 195b0753f9dSMatthew G. Knepley typedef enum { 196b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_DIAG, 197b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_LOWER, 198b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_UPPER, 199b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_FULL 200b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType; 201b0753f9dSMatthew G. Knepley 202b0753f9dSMatthew G. Knepley /*E 203b0753f9dSMatthew G. Knepley PCPARMSGlobalType - Determines the global preconditioner method in PARMS 204b0753f9dSMatthew G. Knepley 205b0753f9dSMatthew G. Knepley Level: intermediate 206b0753f9dSMatthew G. Knepley 207b0753f9dSMatthew G. Knepley .seealso: PCPARMSSetGlobal() 208b0753f9dSMatthew G. Knepley E*/ 209b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType; 2109d8081ecSMatthew G. Knepley 211b0753f9dSMatthew G. Knepley /*E 212b0753f9dSMatthew G. Knepley PCPARMSLocalType - Determines the local preconditioner method in PARMS 213b0753f9dSMatthew G. Knepley 214b0753f9dSMatthew G. Knepley Level: intermediate 215b0753f9dSMatthew G. Knepley 216b0753f9dSMatthew G. Knepley .seealso: PCPARMSSetLocal() 217b0753f9dSMatthew G. Knepley E*/ 218b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType; 219b0753f9dSMatthew G. Knepley 220b0753f9dSMatthew G. Knepley /*E 221b0753f9dSMatthew G. Knepley PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method 222b0753f9dSMatthew G. Knepley 223b0753f9dSMatthew G. Knepley Level: intermediate 224b0753f9dSMatthew G. Knepley 225b0753f9dSMatthew G. Knepley .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseInterpolation() 226b0753f9dSMatthew G. Knepley E*/ 227b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType; 228b0753f9dSMatthew G. Knepley #define PCGAMGAGG "agg" 229b0753f9dSMatthew G. Knepley #define PCGAMGGEO "geo" 230b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL "classical" 231b0753f9dSMatthew G. Knepley 232b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType; 233b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT "direct" 234b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard" 235b0753f9dSMatthew G. Knepley 236b0753f9dSMatthew G. Knepley /*E 237b0753f9dSMatthew G. Knepley PCMGType - Determines the type of multigrid method that is run. 238b0753f9dSMatthew G. Knepley 239b0753f9dSMatthew G. Knepley Level: beginner 240b0753f9dSMatthew G. Knepley 241b0753f9dSMatthew G. Knepley Values: 242c1cbb1deSBarry Smith + PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycleType() 243b0753f9dSMatthew G. Knepley . PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are 244b0753f9dSMatthew G. Knepley smoothed before updating the residual. This only uses the 245b0753f9dSMatthew G. Knepley down smoother, in the preconditioner the upper smoother is ignored 246b0753f9dSMatthew G. Knepley . PC_MG_FULL - same as multiplicative except one also performs grid sequencing, 247b0753f9dSMatthew G. Knepley that is starts on the coarsest grid, performs a cycle, interpolates 248b0753f9dSMatthew 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 249b0753f9dSMatthew G. Knepley algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 250b0753f9dSMatthew G. Knepley calls the V-cycle only on the coarser level and has a post-smoother instead. 251b0753f9dSMatthew G. Knepley - PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level 252b0753f9dSMatthew G. Knepley from a finer 253b0753f9dSMatthew G. Knepley 254c1cbb1deSBarry Smith .seealso: PCMGSetType(), PCMGSetCycleType(), PCMGSetCycleTypeOnLevel() 255b0753f9dSMatthew G. Knepley 256b0753f9dSMatthew G. Knepley E*/ 257b0753f9dSMatthew G. Knepley typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType; 258b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE; 259b0753f9dSMatthew G. Knepley 260b0753f9dSMatthew G. Knepley /*E 261b0753f9dSMatthew G. Knepley PCMGCycleType - Use V-cycle or W-cycle 262b0753f9dSMatthew G. Knepley 263b0753f9dSMatthew G. Knepley Level: beginner 264b0753f9dSMatthew G. Knepley 265b0753f9dSMatthew G. Knepley Values: 266b0753f9dSMatthew G. Knepley + PC_MG_V_CYCLE 267b0753f9dSMatthew G. Knepley - PC_MG_W_CYCLE 268b0753f9dSMatthew G. Knepley 269b0753f9dSMatthew G. Knepley .seealso: PCMGSetCycleType() 270b0753f9dSMatthew G. Knepley 271b0753f9dSMatthew G. Knepley E*/ 272b0753f9dSMatthew G. Knepley typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType; 273b0753f9dSMatthew G. Knepley 274b0753f9dSMatthew G. Knepley /*E 2752134b1e4SBarry Smith PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process 2762134b1e4SBarry Smith 2772134b1e4SBarry Smith Level: beginner 2782134b1e4SBarry Smith 2792134b1e4SBarry Smith Values: 2802134b1e4SBarry Smith + PC_MG_GALERKIN_PMAT - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid 2812134b1e4SBarry Smith . PC_MG_GALERKIN_MAT - computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid 2822134b1e4SBarry Smith . PC_MG_GALERKIN_BOTH - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once 2832134b1e4SBarry Smith - PC_MG_GALERKIN_NONE - neither operator is computed via the Galerkin process, the user must provide the operator 2842134b1e4SBarry Smith 2852134b1e4SBarry Smith Users should never set PC_MG_GALERKIN_EXTERNAL, it is used by GAMG and ML 2862134b1e4SBarry Smith 2872134b1e4SBarry Smith .seealso: PCMGSetCycleType() 2882134b1e4SBarry Smith 2892134b1e4SBarry Smith E*/ 2902134b1e4SBarry Smith typedef enum { PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, PC_MG_GALERKIN_EXTERNAL} PCMGGalerkinType; 2912134b1e4SBarry Smith 2922134b1e4SBarry Smith /*E 293b0753f9dSMatthew G. Knepley PCExoticType - Face based or wirebasket based coarse grid space 294b0753f9dSMatthew G. Knepley 295b0753f9dSMatthew G. Knepley Level: beginner 296b0753f9dSMatthew G. Knepley 297b0753f9dSMatthew G. Knepley .seealso: PCExoticSetType(), PCEXOTIC 298b0753f9dSMatthew G. Knepley E*/ 299b0753f9dSMatthew G. Knepley typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType; 300b0753f9dSMatthew G. Knepley 3018c1cd74cSHong Zhang /*E 302ee6acaa3SMatthew G. Knepley PCPatchConstructType - The algorithm used to construct patches for the preconditioner 3034bbf5ea8SMatthew G. Knepley 3044bbf5ea8SMatthew G. Knepley Level: beginner 3054bbf5ea8SMatthew G. Knepley 3064bbf5ea8SMatthew G. Knepley .seealso: PCPatchSetConstructType(), PCEXOTIC 3074bbf5ea8SMatthew G. Knepley E*/ 308e5b9877fSPatrick Farrell typedef enum {PC_PATCH_STAR, PC_PATCH_VANKA, PC_PATCH_PARDECOMP, PC_PATCH_USER, PC_PATCH_PYTHON} PCPatchConstructType; 3094bbf5ea8SMatthew G. Knepley 3104bbf5ea8SMatthew G. Knepley /*E 311*37eeb815SJakub Kruzik PCDeflationType - Type of deflation to use 312*37eeb815SJakub Kruzik Level: beginner 313*37eeb815SJakub Kruzik 314*37eeb815SJakub Kruzik Values: 315*37eeb815SJakub Kruzik + PC_DEFLATION_PRE - use special initial guess and keep the iterations in the correct space (default) 316*37eeb815SJakub Kruzik . PC_DEFLATION_INIT - use special initial guess (good for deflation space consiting of eigenvectors approximated to high precision) 317*37eeb815SJakub Kruzik - PC_DEFLATION_POST - keep the iterations in the correct space and correct the solution after it was found 318*37eeb815SJakub Kruzik 319*37eeb815SJakub Kruzik .seealso: PCDeflationSetType() 320*37eeb815SJakub Kruzik 321*37eeb815SJakub Kruzik E*/ 322*37eeb815SJakub Kruzik typedef enum {PC_DEFLATION_INIT,PC_DEFLATION_PRE,PC_DEFLATION_POST} PCDeflationType; 323*37eeb815SJakub Kruzik 324*37eeb815SJakub Kruzik /*E 3258c1cd74cSHong Zhang PCFailedReason - indicates type of PC failure 3268c1cd74cSHong Zhang 3278c1cd74cSHong Zhang Level: beginner 3288c1cd74cSHong Zhang 3298c1cd74cSHong Zhang Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h 3308c1cd74cSHong Zhang E*/ 331284f0211SHong Zhang typedef enum {PC_NOERROR,PC_FACTOR_STRUCT_ZEROPIVOT,PC_FACTOR_NUMERIC_ZEROPIVOT,PC_FACTOR_OUTMEMORY,PC_FACTOR_OTHER,PC_SUBPC_ERROR} PCFailedReason; 332b0753f9dSMatthew G. Knepley #endif 333