16524c165SJacob Faibussowitsch #ifndef 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 11*09e53591SBarry Smith .seealso: [](sec_pc), `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 23*09e53591SBarry Smith .seealso: [](sec_pc), `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" 32e6f8f311SMark 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 87*09e53591SBarry Smith .seealso: [](sec_pc), `PC` 88b0753f9dSMatthew G. Knepley E*/ 899371c9d4SSatish Balay typedef enum { 909371c9d4SSatish Balay PC_SIDE_DEFAULT = -1, 919371c9d4SSatish Balay PC_LEFT, 929371c9d4SSatish Balay PC_RIGHT, 939371c9d4SSatish Balay PC_SYMMETRIC 949371c9d4SSatish Balay } PCSide; 95b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 96b0753f9dSMatthew G. Knepley 97b0753f9dSMatthew G. Knepley /*E 98*09e53591SBarry Smith PCRichardsonConvergedReason - reason a `PCRICHARDSON` `PCApplyRichardson()` method terminated 99b0753f9dSMatthew G. Knepley 100b0753f9dSMatthew G. Knepley Level: advanced 101b0753f9dSMatthew G. Knepley 10287497f52SBarry Smith Developer Note: 10387497f52SBarry Smith this must match petsc/finclude/petscpc.h and the `KSPConvergedReason` values in petscksp.h 104b0753f9dSMatthew G. Knepley 105*09e53591SBarry Smith .seealso: [](sec_pc), `PCRICHARDSON`, `PC`, `PCApplyRichardson()` 106b0753f9dSMatthew G. Knepley E*/ 107b0753f9dSMatthew G. Knepley typedef enum { 108b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_RTOL = 2, 109b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ATOL = 3, 110b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ITS = 4, 1119371c9d4SSatish Balay PCRICHARDSON_DIVERGED_DTOL = -4 1129371c9d4SSatish Balay } PCRichardsonConvergedReason; 113b0753f9dSMatthew G. Knepley 114b0753f9dSMatthew G. Knepley /*E 115b0753f9dSMatthew G. Knepley PCJacobiType - What elements are used to form the Jacobi preconditioner 116b0753f9dSMatthew G. Knepley 117*09e53591SBarry Smith Values: 118*09e53591SBarry Smith + `PC_JACOBI_DIAGONAL` - use the diagonal entry, if it is zero use one 119*09e53591SBarry Smith . `PC_JACOBI_ROWMAX` - use the maximum absolute value in the row 120*09e53591SBarry Smith - `PC_JACOBI_ROWSUM` - use the sum of the values in the row (not the absolute values) 121*09e53591SBarry Smith 122b0753f9dSMatthew G. Knepley Level: intermediate 123b0753f9dSMatthew G. Knepley 124*09e53591SBarry Smith .seealso: [](sec_pc), `PCJACOBI`, `PC` 125b0753f9dSMatthew G. Knepley E*/ 1269371c9d4SSatish Balay typedef enum { 1279371c9d4SSatish Balay PC_JACOBI_DIAGONAL, 1289371c9d4SSatish Balay PC_JACOBI_ROWMAX, 1299371c9d4SSatish Balay PC_JACOBI_ROWSUM 1309371c9d4SSatish Balay } PCJacobiType; 131b0753f9dSMatthew G. Knepley 132b0753f9dSMatthew G. Knepley /*E 133b0753f9dSMatthew G. Knepley PCASMType - Type of additive Schwarz method to use 134b0753f9dSMatthew G. Knepley 135*09e53591SBarry Smith Values: 136*09e53591SBarry Smith + `PC_ASM_BASIC` - Symmetric version where residuals from the ghost points are used 137*09e53591SBarry Smith and computed values in ghost regions are added together. 138*09e53591SBarry Smith Classical standard additive Schwarz. 139*09e53591SBarry Smith . `PC_ASM_RESTRICT` - Residuals from ghost points are used but computed values in ghost 140*09e53591SBarry Smith region are discarded. 141*09e53591SBarry Smith Default. 142*09e53591SBarry Smith . `PC_ASM_INTERPOLATE` - Residuals from ghost points are not used, computed values in ghost 143*09e53591SBarry Smith region are added back in. 144*09e53591SBarry Smith - `PC_ASM_NONE` - Residuals from ghost points are not used, computed ghost values are 145*09e53591SBarry Smith discarded. 146*09e53591SBarry Smith Not very good. 147b0753f9dSMatthew G. Knepley 148b0753f9dSMatthew G. Knepley Level: beginner 149b0753f9dSMatthew G. Knepley 150*09e53591SBarry Smith .seealso: [](sec_pc), `PC`, `PCASM`, `PCASMSetType()` 151b0753f9dSMatthew G. Knepley E*/ 1529371c9d4SSatish Balay typedef enum { 1539371c9d4SSatish Balay PC_ASM_BASIC = 3, 1549371c9d4SSatish Balay PC_ASM_RESTRICT = 1, 1559371c9d4SSatish Balay PC_ASM_INTERPOLATE = 2, 1569371c9d4SSatish Balay PC_ASM_NONE = 0 1579371c9d4SSatish Balay } PCASMType; 158b0753f9dSMatthew G. Knepley 159b0753f9dSMatthew G. Knepley /*E 16087497f52SBarry Smith PCGASMType - Type of generalized additive Schwarz method to use (differs from `PCASM` in allowing multiple processors per subdomain). 161b0753f9dSMatthew G. Knepley 162b0753f9dSMatthew G. Knepley Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 163b0753f9dSMatthew G. Knepley domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 164b0753f9dSMatthew G. Knepley a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 165b0753f9dSMatthew G. Knepley (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 166b0753f9dSMatthew G. Knepley 167*09e53591SBarry Smith Values: 168*09e53591SBarry Smith + `PC_GASM_BASIC` - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 169*09e53591SBarry Smith over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 170*09e53591SBarry Smith from neighboring subdomains. 171*09e53591SBarry Smith Classical standard additive Schwarz. 172*09e53591SBarry Smith . `PC_GASM_RESTRICT` - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 173*09e53591SBarry Smith (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 174*09e53591SBarry Smith each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 175*09e53591SBarry Smith assumption). 176*09e53591SBarry Smith Default. 177*09e53591SBarry Smith . `PC_GASM_INTERPOLATE` - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 178*09e53591SBarry Smith applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 179*09e53591SBarry Smith from neighboring subdomains. 180*09e53591SBarry Smith 181*09e53591SBarry Smith - `PC_GASM_NONE` - Residuals and corrections are zeroed out outside the local subdomains. 182*09e53591SBarry Smith Not very good. 183b0753f9dSMatthew G. Knepley 184b0753f9dSMatthew G. Knepley Level: beginner 185b0753f9dSMatthew G. Knepley 186*09e53591SBarry Smith .seealso: [](sec_pc), `PCGASM`, `PCASM`, `PC`, `PCGASMSetType()` 187b0753f9dSMatthew G. Knepley E*/ 1889371c9d4SSatish Balay typedef enum { 1899371c9d4SSatish Balay PC_GASM_BASIC = 3, 1909371c9d4SSatish Balay PC_GASM_RESTRICT = 1, 1919371c9d4SSatish Balay PC_GASM_INTERPOLATE = 2, 1929371c9d4SSatish Balay PC_GASM_NONE = 0 1939371c9d4SSatish Balay } PCGASMType; 194b0753f9dSMatthew G. Knepley 195b0753f9dSMatthew G. Knepley /*E 196b0753f9dSMatthew G. Knepley PCCompositeType - Determines how two or more preconditioner are composed 197b0753f9dSMatthew G. Knepley 198*09e53591SBarry Smith Values: 199*09e53591SBarry Smith + `PC_COMPOSITE_ADDITIVE` - results from application of all preconditioners are added together 200*09e53591SBarry Smith . `PC_COMPOSITE_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly 201*09e53591SBarry Smith computed after the previous preconditioner application 202*09e53591SBarry Smith . `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly 203*09e53591SBarry Smith computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners) 204*09e53591SBarry Smith . `PC_COMPOSITE_SPECIAL` - This is very special for a matrix of the form alpha I + R + S 205*09e53591SBarry Smith where first preconditioner is built from alpha I + S and second from 206*09e53591SBarry Smith alpha I + R 207*09e53591SBarry Smith . `PC_COMPOSITE_SCHUR` - composes the Schur complement of the matrix from two blocks, see `PCFIELDSPLIT` 208*09e53591SBarry Smith - `PC_COMPOSITE_GKB` - the generalized Golub-Kahan bidiagonalization preconditioner, see `PCFIELDSPLIT` 209b0753f9dSMatthew G. Knepley 210b0753f9dSMatthew G. Knepley Level: beginner 211b0753f9dSMatthew G. Knepley 212*09e53591SBarry Smith .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()` 213b0753f9dSMatthew G. Knepley E*/ 2149371c9d4SSatish Balay typedef enum { 2159371c9d4SSatish Balay PC_COMPOSITE_ADDITIVE, 2169371c9d4SSatish Balay PC_COMPOSITE_MULTIPLICATIVE, 2179371c9d4SSatish Balay PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, 2189371c9d4SSatish Balay PC_COMPOSITE_SPECIAL, 2199371c9d4SSatish Balay PC_COMPOSITE_SCHUR, 2209371c9d4SSatish Balay PC_COMPOSITE_GKB 2219371c9d4SSatish Balay } PCCompositeType; 222b0753f9dSMatthew G. Knepley 223b0753f9dSMatthew G. Knepley /*E 22487497f52SBarry Smith PCFieldSplitSchurPreType - Determines how to precondition a Schur complement 225b0753f9dSMatthew G. Knepley 226b0753f9dSMatthew G. Knepley Level: intermediate 227b0753f9dSMatthew G. Knepley 228*09e53591SBarry Smith .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurPre()`, `PC` 229b0753f9dSMatthew G. Knepley E*/ 2309371c9d4SSatish Balay typedef enum { 2319371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_SELF, 2329371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_SELFP, 2339371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_A11, 2349371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_USER, 2359371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_FULL 2369371c9d4SSatish Balay } PCFieldSplitSchurPreType; 237b0753f9dSMatthew G. Knepley 238b0753f9dSMatthew G. Knepley /*E 239b0753f9dSMatthew G. Knepley PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 240b0753f9dSMatthew G. Knepley 241b0753f9dSMatthew G. Knepley Level: intermediate 242b0753f9dSMatthew G. Knepley 243*09e53591SBarry Smith .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurFactType()`, `PC` 244b0753f9dSMatthew G. Knepley E*/ 245b0753f9dSMatthew G. Knepley typedef enum { 246b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_DIAG, 247b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_LOWER, 248b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_UPPER, 249b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_FULL 250b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType; 251b0753f9dSMatthew G. Knepley 252b0753f9dSMatthew G. Knepley /*E 25387497f52SBarry Smith PCPARMSGlobalType - Determines the global preconditioner method in `PCPARMS` 254b0753f9dSMatthew G. Knepley 255b0753f9dSMatthew G. Knepley Level: intermediate 256b0753f9dSMatthew G. Knepley 257*09e53591SBarry Smith .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetGlobal()`, `PC` 258b0753f9dSMatthew G. Knepley E*/ 2599371c9d4SSatish Balay typedef enum { 2609371c9d4SSatish Balay PC_PARMS_GLOBAL_RAS, 2619371c9d4SSatish Balay PC_PARMS_GLOBAL_SCHUR, 2629371c9d4SSatish Balay PC_PARMS_GLOBAL_BJ 2639371c9d4SSatish Balay } PCPARMSGlobalType; 2649d8081ecSMatthew G. Knepley 265b0753f9dSMatthew G. Knepley /*E 26687497f52SBarry Smith PCPARMSLocalType - Determines the local preconditioner method in `PCPARMS` 267b0753f9dSMatthew G. Knepley 268b0753f9dSMatthew G. Knepley Level: intermediate 269b0753f9dSMatthew G. Knepley 270*09e53591SBarry Smith .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetLocal()`, `PC` 271b0753f9dSMatthew G. Knepley E*/ 2729371c9d4SSatish Balay typedef enum { 2739371c9d4SSatish Balay PC_PARMS_LOCAL_ILU0, 2749371c9d4SSatish Balay PC_PARMS_LOCAL_ILUK, 2759371c9d4SSatish Balay PC_PARMS_LOCAL_ILUT, 2769371c9d4SSatish Balay PC_PARMS_LOCAL_ARMS 2779371c9d4SSatish Balay } PCPARMSLocalType; 278b0753f9dSMatthew G. Knepley 279f0fc11ceSJed Brown /*J 28087497f52SBarry Smith PCGAMGType - type of generalized algebraic multigrid `PCGAMG` method 281b0753f9dSMatthew G. Knepley 282*09e53591SBarry Smith Values: 283*09e53591SBarry Smith + `PCGAMGAGG` - (the default) smoothed aggregation algorithm, robust, very well tested 284*09e53591SBarry Smith . `PCGAMGGEO` - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not well tested 285*09e53591SBarry Smith - `PCGAMGCLASSICAL` - classical algebraic multigrid preconditioner, incomplete, poorly tested 286*09e53591SBarry Smith 287b0753f9dSMatthew G. Knepley Level: intermediate 288b0753f9dSMatthew G. Knepley 289*09e53591SBarry Smith .seealso: [](sec_pc), `PCGAMG`, `PCMG`, `PC`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()` 290f0fc11ceSJed Brown J*/ 291b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType; 292b0753f9dSMatthew G. Knepley #define PCGAMGAGG "agg" 293b0753f9dSMatthew G. Knepley #define PCGAMGGEO "geo" 294b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL "classical" 295b0753f9dSMatthew G. Knepley 296b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType; 297b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT "direct" 298b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard" 299b0753f9dSMatthew G. Knepley 300b0753f9dSMatthew G. Knepley /*E 301b0753f9dSMatthew G. Knepley PCMGType - Determines the type of multigrid method that is run. 302b0753f9dSMatthew G. Knepley 303b0753f9dSMatthew G. Knepley Level: beginner 304b0753f9dSMatthew G. Knepley 305b0753f9dSMatthew G. Knepley Values: 30687497f52SBarry Smith + `PC_MG_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `PCMGSetCycleType()` 30787497f52SBarry Smith . `PC_MG_ADDITIVE` - the additive multigrid preconditioner where all levels are 308b0753f9dSMatthew G. Knepley smoothed before updating the residual. This only uses the 309b0753f9dSMatthew G. Knepley down smoother, in the preconditioner the upper smoother is ignored 31087497f52SBarry Smith . `PC_MG_FULL` - same as multiplicative except one also performs grid sequencing, 311b0753f9dSMatthew G. Knepley that is starts on the coarsest grid, performs a cycle, interpolates 312b0753f9dSMatthew 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 313b0753f9dSMatthew G. Knepley algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 314b0753f9dSMatthew G. Knepley calls the V-cycle only on the coarser level and has a post-smoother instead. 31587497f52SBarry Smith - `PC_MG_KASKADE` - like full multigrid except one never goes back to a coarser level 316b0753f9dSMatthew G. Knepley from a finer 317b0753f9dSMatthew G. Knepley 318*09e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()` 319b0753f9dSMatthew G. Knepley E*/ 3209371c9d4SSatish Balay typedef enum { 3219371c9d4SSatish Balay PC_MG_MULTIPLICATIVE, 3229371c9d4SSatish Balay PC_MG_ADDITIVE, 3239371c9d4SSatish Balay PC_MG_FULL, 3249371c9d4SSatish Balay PC_MG_KASKADE 3259371c9d4SSatish Balay } PCMGType; 326b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE; 327b0753f9dSMatthew G. Knepley 328b0753f9dSMatthew G. Knepley /*E 329b0753f9dSMatthew G. Knepley PCMGCycleType - Use V-cycle or W-cycle 330b0753f9dSMatthew G. Knepley 331b0753f9dSMatthew G. Knepley Level: beginner 332b0753f9dSMatthew G. Knepley 333b0753f9dSMatthew G. Knepley Values: 33487497f52SBarry Smith + `PC_MG_V_CYCLE` - use the v cycle 33587497f52SBarry Smith - `PC_MG_W_CYCLE` - use the w cycle 336b0753f9dSMatthew G. Knepley 337*09e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()` 338b0753f9dSMatthew G. Knepley E*/ 3399371c9d4SSatish Balay typedef enum { 3409371c9d4SSatish Balay PC_MG_CYCLE_V = 1, 3419371c9d4SSatish Balay PC_MG_CYCLE_W = 2 3429371c9d4SSatish Balay } PCMGCycleType; 343b0753f9dSMatthew G. Knepley 344b0753f9dSMatthew G. Knepley /*E 3452134b1e4SBarry Smith PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process 3462134b1e4SBarry Smith 3472134b1e4SBarry Smith Values: 34887497f52SBarry Smith + `PC_MG_GALERKIN_PMAT` - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid 34987497f52SBarry Smith . `PC_MG_GALERKIN_MAT` - computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid 35087497f52SBarry Smith . `PC_MG_GALERKIN_BOTH` - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once 35187497f52SBarry Smith - `PC_MG_GALERKIN_NONE` - neither operator is computed via the Galerkin process, the user must provide the operator 3522134b1e4SBarry Smith 353*09e53591SBarry Smith Level: beginner 354*09e53591SBarry Smith 355*09e53591SBarry Smith Note: 35687497f52SBarry Smith Users should never set `PC_MG_GALERKIN_EXTERNAL`, it is used by `PCHYPRE` and `PCML` 3572134b1e4SBarry Smith 358*09e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()` 3592134b1e4SBarry Smith E*/ 3609371c9d4SSatish Balay typedef enum { 3619371c9d4SSatish Balay PC_MG_GALERKIN_BOTH, 3629371c9d4SSatish Balay PC_MG_GALERKIN_PMAT, 3639371c9d4SSatish Balay PC_MG_GALERKIN_MAT, 3649371c9d4SSatish Balay PC_MG_GALERKIN_NONE, 3659371c9d4SSatish Balay PC_MG_GALERKIN_EXTERNAL 3669371c9d4SSatish Balay } PCMGGalerkinType; 3672134b1e4SBarry Smith 3682134b1e4SBarry Smith /*E 369b0753f9dSMatthew G. Knepley PCExoticType - Face based or wirebasket based coarse grid space 370b0753f9dSMatthew G. Knepley 371b0753f9dSMatthew G. Knepley Level: beginner 372b0753f9dSMatthew G. Knepley 373*09e53591SBarry Smith .seealso: [](sec_pc), `PCExoticSetType()`, `PCEXOTIC` 374b0753f9dSMatthew G. Knepley E*/ 3759371c9d4SSatish Balay typedef enum { 3769371c9d4SSatish Balay PC_EXOTIC_FACE, 3779371c9d4SSatish Balay PC_EXOTIC_WIREBASKET 3789371c9d4SSatish Balay } PCExoticType; 379b0753f9dSMatthew G. Knepley 3808c1cd74cSHong Zhang /*E 381bc960bbfSJed Brown PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains. 382bc960bbfSJed Brown 383bc960bbfSJed Brown Level: intermediate 384bc960bbfSJed Brown 385bc960bbfSJed Brown Values: 38687497f52SBarry Smith + `PC_BDDC_INTERFACE_EXT_DIRICHLET` - solves Dirichlet interior problem; this is the standard BDDC algorithm 38787497f52SBarry Smith - `PC_BDDC_INTERFACE_EXT_LUMP` - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP" 388bc960bbfSJed Brown 389*09e53591SBarry Smith .seealso: [](sec_pc), `PCBDDC`, `PC` 390bc960bbfSJed Brown E*/ 391bc960bbfSJed Brown typedef enum { 392bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_DIRICHLET, 393bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_LUMP 394bc960bbfSJed Brown } PCBDDCInterfaceExtType; 395bc960bbfSJed Brown 396bc960bbfSJed Brown /*E 397f3b08a26SMatthew G. Knepley PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation 398f3b08a26SMatthew G. Knepley 399f3b08a26SMatthew G. Knepley Level: beginner 400f3b08a26SMatthew G. Knepley 401*09e53591SBarry Smith .seealso: [](sec_pc), `PCMGSetAdaptCoarseSpaceType()`, `PCMG`, `PC` 402f3b08a26SMatthew G. Knepley E*/ 4039371c9d4SSatish Balay typedef enum { 4049371c9d4SSatish Balay PCMG_ADAPT_NONE, 4059371c9d4SSatish Balay PCMG_ADAPT_POLYNOMIAL, 4069371c9d4SSatish Balay PCMG_ADAPT_HARMONIC, 4079371c9d4SSatish Balay PCMG_ADAPT_EIGENVECTOR, 4089371c9d4SSatish Balay PCMG_ADAPT_GENERALIZED_EIGENVECTOR, 4099371c9d4SSatish Balay PCMG_ADAPT_GDSW 4109371c9d4SSatish Balay } PCMGCoarseSpaceType; 411f3b08a26SMatthew G. Knepley 412f3b08a26SMatthew G. Knepley /*E 413ee6acaa3SMatthew G. Knepley PCPatchConstructType - The algorithm used to construct patches for the preconditioner 4144bbf5ea8SMatthew G. Knepley 4154bbf5ea8SMatthew G. Knepley Level: beginner 4164bbf5ea8SMatthew G. Knepley 417*09e53591SBarry Smith .seealso: [](sec_pc), `PCPatchSetConstructType()`, `PCPATCH`, `PC` 4184bbf5ea8SMatthew G. Knepley E*/ 4199371c9d4SSatish Balay typedef enum { 4209371c9d4SSatish Balay PC_PATCH_STAR, 4219371c9d4SSatish Balay PC_PATCH_VANKA, 4229371c9d4SSatish Balay PC_PATCH_PARDECOMP, 4239371c9d4SSatish Balay PC_PATCH_USER, 4249371c9d4SSatish Balay PC_PATCH_PYTHON 4259371c9d4SSatish Balay } PCPatchConstructType; 4264bbf5ea8SMatthew G. Knepley 4274bbf5ea8SMatthew G. Knepley /*E 428e26ad46dSJakub Kruzik PCDeflationSpaceType - Type of deflation 429e662bc50SJakub Kruzik 430e662bc50SJakub Kruzik Values: 43187497f52SBarry Smith + `PC_DEFLATION_SPACE_HAAR` - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off 432*09e53591SBarry Smith . `PC_DEFLATION_SPACE_DB2` - `MATCOMPOSITE` of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet) 43387497f52SBarry Smith . `PC_DEFLATION_SPACE_DB4` - same as above, but with db4 (4 coefficient Daubechies) 43487497f52SBarry Smith . `PC_DEFLATION_SPACE_DB8` - same as above, but with db8 (8 coefficient Daubechies) 43587497f52SBarry Smith . `PC_DEFLATION_SPACE_DB16` - same as above, but with db16 (16 coefficient Daubechies) 43687497f52SBarry Smith . `PC_DEFLATION_SPACE_BIORTH22` - same as above, but with biorthogonal 2.2 (6 coefficients) 43787497f52SBarry Smith . `PC_DEFLATION_SPACE_MEYER` - same as above, but with Meyer/FIR (62 coefficients) 438d5b43468SJose E. Roman . `PC_DEFLATION_SPACE_AGGREGATION` - aggregates local indices (given by operator matrix distribution) into a subdomain 43987497f52SBarry Smith - `PC_DEFLATION_SPACE_USER` - indicates space set by user 440e662bc50SJakub Kruzik 4411fdb00f9SJakub Kruzik Level: intermediate 4421fdb00f9SJakub Kruzik 443*09e53591SBarry Smith Note: 444*09e53591SBarry Smith Wavelet-based space (except Haar) can be used in multilevel deflation. 445*09e53591SBarry Smith 446*09e53591SBarry Smith .seealso: [](sec_pc), `PCDeflationSetSpaceToCompute()`, `PCDEFLATION`, `PC` 447e662bc50SJakub Kruzik E*/ 448e662bc50SJakub Kruzik typedef enum { 449e662bc50SJakub Kruzik PC_DEFLATION_SPACE_HAAR, 450e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB2, 451e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB4, 452e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB8, 453e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB16, 454e662bc50SJakub Kruzik PC_DEFLATION_SPACE_BIORTH22, 455e662bc50SJakub Kruzik PC_DEFLATION_SPACE_MEYER, 456e662bc50SJakub Kruzik PC_DEFLATION_SPACE_AGGREGATION, 457e662bc50SJakub Kruzik PC_DEFLATION_SPACE_USER 458e662bc50SJakub Kruzik } PCDeflationSpaceType; 459e662bc50SJakub Kruzik 460e662bc50SJakub Kruzik /*E 46187497f52SBarry Smith PCHPDDMCoarseCorrectionType - Type of coarse correction used by `PCHPDDM` 46238cfc46eSPierre Jolivet 46338cfc46eSPierre Jolivet Level: intermediate 46438cfc46eSPierre Jolivet 46538cfc46eSPierre Jolivet Values: 466*09e53591SBarry Smith + `PC_HPDDM_COARSE_CORRECTION_DEFLATED` (default) - eq. (1) in `PCHPDDMShellApply()` 46787497f52SBarry Smith . `PC_HPDDM_COARSE_CORRECTION_ADDITIVE` - eq. (2) 46887497f52SBarry Smith - `PC_HPDDM_COARSE_CORRECTION_BALANCED` - eq. (3) 46938cfc46eSPierre Jolivet 470*09e53591SBarry Smith .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCSetType()`, `PCHPDDMShellApply()` 47138cfc46eSPierre Jolivet E*/ 4729371c9d4SSatish Balay typedef enum { 4739371c9d4SSatish Balay PC_HPDDM_COARSE_CORRECTION_DEFLATED, 4749371c9d4SSatish Balay PC_HPDDM_COARSE_CORRECTION_ADDITIVE, 4759371c9d4SSatish Balay PC_HPDDM_COARSE_CORRECTION_BALANCED 4769371c9d4SSatish Balay } PCHPDDMCoarseCorrectionType; 47738cfc46eSPierre Jolivet 47838cfc46eSPierre Jolivet /*E 47987497f52SBarry Smith PCFailedReason - indicates type of `PC` failure 4808c1cd74cSHong Zhang 4818c1cd74cSHong Zhang Level: beginner 4828c1cd74cSHong Zhang 48387497f52SBarry Smith Developer Note: 4848c1cd74cSHong Zhang Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h 485*09e53591SBarry Smith 486*09e53591SBarry Smith .seealso: [](sec_pc), `PC` 4878c1cd74cSHong Zhang E*/ 4889371c9d4SSatish Balay typedef enum { 4899371c9d4SSatish Balay PC_SETUP_ERROR = -1, 4909371c9d4SSatish Balay PC_NOERROR, 4919371c9d4SSatish Balay PC_FACTOR_STRUCT_ZEROPIVOT, 4929371c9d4SSatish Balay PC_FACTOR_NUMERIC_ZEROPIVOT, 4939371c9d4SSatish Balay PC_FACTOR_OUTMEMORY, 4949371c9d4SSatish Balay PC_FACTOR_OTHER, 49587b47708SBarry Smith PC_INCONSISTENT_RHS, 4969371c9d4SSatish Balay PC_SUBPC_ERROR 4979371c9d4SSatish Balay } PCFailedReason; 498ce7c7f2fSMark Adams 499ce7c7f2fSMark Adams /*E 500ce7c7f2fSMark Adams PCGAMGLayoutType - Layout for reduced grids 501ce7c7f2fSMark Adams 502ce7c7f2fSMark Adams Level: intermediate 503ce7c7f2fSMark Adams 504*09e53591SBarry Smith Developer Note: 505ce7c7f2fSMark Adams Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h 506*09e53591SBarry Smith 507*09e53591SBarry Smith .seealso: [](sec_pc), `PCGAMG`, `PC`, `PCGAMGSetCoarseGridLayoutType()` 508ce7c7f2fSMark Adams E*/ 5099371c9d4SSatish Balay typedef enum { 5109371c9d4SSatish Balay PCGAMG_LAYOUT_COMPACT, 5119371c9d4SSatish Balay PCGAMG_LAYOUT_SPREAD 5129371c9d4SSatish Balay } PCGAMGLayoutType; 513ce7c7f2fSMark Adams 514b0753f9dSMatthew G. Knepley #endif 515