1a4963045SJacob Faibussowitsch #pragma once 2b0753f9dSMatthew G. Knepley 3ac09b921SBarry Smith /* SUBMANSEC = PC */ 4ac09b921SBarry Smith 5b0753f9dSMatthew G. Knepley /*S 687497f52SBarry Smith PC - Abstract PETSc object that manages all preconditioners including direct solvers such as `PCLU` 7b0753f9dSMatthew G. Knepley 8b0753f9dSMatthew G. Knepley Level: beginner 9b0753f9dSMatthew G. Knepley 1016a05f60SBarry Smith .seealso: [](doc_linsolve), [](sec_pc), `PCCreate()`, `PCSetType()`, `PCType` 11b0753f9dSMatthew G. Knepley S*/ 12b0753f9dSMatthew G. Knepley typedef struct _p_PC *PC; 13b0753f9dSMatthew G. Knepley 14b0753f9dSMatthew G. Knepley /*J 1516a05f60SBarry Smith PCType - String with the name of a PETSc preconditioner 16b0753f9dSMatthew G. Knepley 17b0753f9dSMatthew G. Knepley Level: beginner 18b0753f9dSMatthew G. Knepley 1987497f52SBarry Smith Note: 2087497f52SBarry Smith `PCRegister()` is used to register preconditioners that are then accessible via `PCSetType()` 21b0753f9dSMatthew G. Knepley 2216a05f60SBarry Smith .seealso: [](doc_linsolve), [](sec_pc), `PCSetType()`, `PC`, `PCCreate()`, `PCRegister()`, `PCSetFromOptions()`, `PCLU`, `PCJACOBI`, `PCBJACOBI` 23b0753f9dSMatthew G. Knepley J*/ 24b0753f9dSMatthew G. Knepley typedef const char *PCType; 25b0753f9dSMatthew G. Knepley #define PCNONE "none" 26b0753f9dSMatthew G. Knepley #define PCJACOBI "jacobi" 27b0753f9dSMatthew G. Knepley #define PCSOR "sor" 28b0753f9dSMatthew G. Knepley #define PCLU "lu" 29a2fc1e05SToby Isaac #define PCQR "qr" 30b0753f9dSMatthew G. Knepley #define PCSHELL "shell" 31e6f8f311SMark Adams #define PCAMGX "amgx" 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" 40e607c864SMark Adams #define PCBJKOKKOS "bjkokkos" 41b0753f9dSMatthew G. Knepley #define PCCOMPOSITE "composite" 42b0753f9dSMatthew G. Knepley #define PCREDUNDANT "redundant" 43b0753f9dSMatthew G. Knepley #define PCSPAI "spai" 44b0753f9dSMatthew G. Knepley #define PCNN "nn" 45b0753f9dSMatthew G. Knepley #define PCCHOLESKY "cholesky" 46b0753f9dSMatthew G. Knepley #define PCPBJACOBI "pbjacobi" 470da83c2eSBarry Smith #define PCVPBJACOBI "vpbjacobi" 48b0753f9dSMatthew G. Knepley #define PCMAT "mat" 49b0753f9dSMatthew G. Knepley #define PCHYPRE "hypre" 50b0753f9dSMatthew G. Knepley #define PCPARMS "parms" 51b0753f9dSMatthew G. Knepley #define PCFIELDSPLIT "fieldsplit" 52b0753f9dSMatthew G. Knepley #define PCTFS "tfs" 53b0753f9dSMatthew G. Knepley #define PCML "ml" 54b0753f9dSMatthew G. Knepley #define PCGALERKIN "galerkin" 55b0753f9dSMatthew G. Knepley #define PCEXOTIC "exotic" 56b0753f9dSMatthew G. Knepley #define PCCP "cp" 57b0753f9dSMatthew G. Knepley #define PCBFBT "bfbt" 58b0753f9dSMatthew G. Knepley #define PCLSC "lsc" 59b0753f9dSMatthew G. Knepley #define PCPYTHON "python" 60b0753f9dSMatthew G. Knepley #define PCPFMG "pfmg" 611c188c59Sftrigaux #define PCSMG "smg" 62b0753f9dSMatthew G. Knepley #define PCSYSPFMG "syspfmg" 63b0753f9dSMatthew G. Knepley #define PCREDISTRIBUTE "redistribute" 64b0753f9dSMatthew G. Knepley #define PCSVD "svd" 65b0753f9dSMatthew G. Knepley #define PCGAMG "gamg" 664b3f184cSKarl Rupp #define PCCHOWILUVIENNACL "chowiluviennacl" 6770baa948SKarl Rupp #define PCROWSCALINGVIENNACL "rowscalingviennacl" 6807598726SKarl Rupp #define PCSAVIENNACL "saviennacl" 69b0753f9dSMatthew G. Knepley #define PCBDDC "bddc" 70b0753f9dSMatthew G. Knepley #define PCKACZMARZ "kaczmarz" 7168ddcbeaSBarry Smith #define PCTELESCOPE "telescope" 724bbf5ea8SMatthew G. Knepley #define PCPATCH "patch" 73b9ac7092SAlp Dener #define PCLMVM "lmvm" 74360ee056SFande Kong #define PCHMG "hmg" 7537eeb815SJakub Kruzik #define PCDEFLATION "deflation" 7638cfc46eSPierre Jolivet #define PCHPDDM "hpddm" 7753022affSStefano Zampini #define PCH2OPUS "h2opus" 78f1f2ae84SBarry Smith #define PCMPI "mpi" 79b0753f9dSMatthew G. Knepley 80b0753f9dSMatthew G. Knepley /*E 81b0753f9dSMatthew G. Knepley PCSide - If the preconditioner is to be applied to the left, right 82b0753f9dSMatthew G. Knepley or symmetrically around the operator. 83b0753f9dSMatthew G. Knepley 8416a05f60SBarry Smith Values: 8516a05f60SBarry Smith + `PC_LEFT` - applied after the operator is applied 8616a05f60SBarry Smith . `PC_RIGHT` - applied before the operator is applied 8716a05f60SBarry Smith - `PC_SYMMETRIC` - a portion of the preconditioner is applied before the operator and the transpose of this portion is applied after the operator is applied. 8816a05f60SBarry Smith 89b0753f9dSMatthew G. Knepley Level: beginner 90b0753f9dSMatthew G. Knepley 9116a05f60SBarry Smith Note: 9216a05f60SBarry Smith Certain `KSPType` support only a subset of `PCSide` values 9316a05f60SBarry Smith 94af27ebaaSBarry Smith .seealso: [](sec_pc), `PC`, `KSPSetPCSide()`, `KSP`, `KSPType` 95b0753f9dSMatthew G. Knepley E*/ 969371c9d4SSatish Balay typedef enum { 979371c9d4SSatish Balay PC_SIDE_DEFAULT = -1, 989371c9d4SSatish Balay PC_LEFT, 999371c9d4SSatish Balay PC_RIGHT, 1009371c9d4SSatish Balay PC_SYMMETRIC 1019371c9d4SSatish Balay } PCSide; 102b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 103b0753f9dSMatthew G. Knepley 104b0753f9dSMatthew G. Knepley /*E 10539b1ba1bSPierre Jolivet PCRichardsonConvergedReason - reason a `PCApplyRichardson()` method terminated 106b0753f9dSMatthew G. Knepley 107b0753f9dSMatthew G. Knepley Level: advanced 108b0753f9dSMatthew G. Knepley 10939b1ba1bSPierre Jolivet .seealso: [](sec_pc), `KSPRICHARDSON`, `PC`, `PCApplyRichardson()` 110b0753f9dSMatthew G. Knepley E*/ 111b0753f9dSMatthew G. Knepley typedef enum { 112b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_RTOL = 2, 113b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ATOL = 3, 114b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ITS = 4, 1159371c9d4SSatish Balay PCRICHARDSON_DIVERGED_DTOL = -4 1169371c9d4SSatish Balay } PCRichardsonConvergedReason; 117b0753f9dSMatthew G. Knepley 118b0753f9dSMatthew G. Knepley /*E 11916a05f60SBarry Smith PCJacobiType - What elements of the matrix are used to form the Jacobi preconditioner 120b0753f9dSMatthew G. Knepley 12109e53591SBarry Smith Values: 12209e53591SBarry Smith + `PC_JACOBI_DIAGONAL` - use the diagonal entry, if it is zero use one 123eede4a3fSMark Adams . `PC_JACOBI_ROWL1` - add sum of absolute values in row i, j != i, to diag_ii 12409e53591SBarry Smith . `PC_JACOBI_ROWMAX` - use the maximum absolute value in the row 12509e53591SBarry Smith - `PC_JACOBI_ROWSUM` - use the sum of the values in the row (not the absolute values) 12609e53591SBarry Smith 127b0753f9dSMatthew G. Knepley Level: intermediate 128b0753f9dSMatthew G. Knepley 12909e53591SBarry Smith .seealso: [](sec_pc), `PCJACOBI`, `PC` 130b0753f9dSMatthew G. Knepley E*/ 1319371c9d4SSatish Balay typedef enum { 1329371c9d4SSatish Balay PC_JACOBI_DIAGONAL, 133eede4a3fSMark Adams PC_JACOBI_ROWL1, 1349371c9d4SSatish Balay PC_JACOBI_ROWMAX, 1359371c9d4SSatish Balay PC_JACOBI_ROWSUM 1369371c9d4SSatish Balay } PCJacobiType; 137b0753f9dSMatthew G. Knepley 138b0753f9dSMatthew G. Knepley /*E 139b0753f9dSMatthew G. Knepley PCASMType - Type of additive Schwarz method to use 140b0753f9dSMatthew G. Knepley 14109e53591SBarry Smith Values: 14209e53591SBarry Smith + `PC_ASM_BASIC` - Symmetric version where residuals from the ghost points are used 14309e53591SBarry Smith and computed values in ghost regions are added together. 144af27ebaaSBarry Smith Classical standard additive Schwarz as introduced in {cite}`dryja1987additive`. 14509e53591SBarry Smith . `PC_ASM_RESTRICT` - Residuals from ghost points are used but computed values in ghost 146af27ebaaSBarry Smith region are discarded {cite}`cs99`. Default. 14709e53591SBarry Smith . `PC_ASM_INTERPOLATE` - Residuals from ghost points are not used, computed values in ghost 14809e53591SBarry Smith region are added back in. 14909e53591SBarry Smith - `PC_ASM_NONE` - Residuals from ghost points are not used, computed ghost values are 150af27ebaaSBarry Smith discarded. Not very good. 151b0753f9dSMatthew G. Knepley 152b0753f9dSMatthew G. Knepley Level: beginner 153b0753f9dSMatthew G. Knepley 15416a05f60SBarry Smith .seealso: [](sec_pc), `PC`, `PCASM`, `PCASMSetType()`, `PCGASMType` 155b0753f9dSMatthew G. Knepley E*/ 1569371c9d4SSatish Balay typedef enum { 1579371c9d4SSatish Balay PC_ASM_BASIC = 3, 1589371c9d4SSatish Balay PC_ASM_RESTRICT = 1, 1599371c9d4SSatish Balay PC_ASM_INTERPOLATE = 2, 1609371c9d4SSatish Balay PC_ASM_NONE = 0 1619371c9d4SSatish Balay } PCASMType; 162b0753f9dSMatthew G. Knepley 163b0753f9dSMatthew G. Knepley /*E 16487497f52SBarry Smith PCGASMType - Type of generalized additive Schwarz method to use (differs from `PCASM` in allowing multiple processors per subdomain). 165b0753f9dSMatthew G. Knepley 16609e53591SBarry Smith Values: 16709e53591SBarry Smith + `PC_GASM_BASIC` - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 16809e53591SBarry Smith over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 169af27ebaaSBarry Smith from neighboring subdomains. Classical standard additive Schwarz {cite}`dryja1987additive`. 17009e53591SBarry Smith . `PC_GASM_RESTRICT` - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 17109e53591SBarry Smith (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 17209e53591SBarry Smith each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 173af27ebaaSBarry Smith assumption) {cite}`cs99`. Default. 17409e53591SBarry Smith . `PC_GASM_INTERPOLATE` - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 17509e53591SBarry Smith applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 17609e53591SBarry Smith from neighboring subdomains. 177af27ebaaSBarry Smith - `PC_GASM_NONE` - Residuals and corrections are zeroed out outside the local subdomains. Not very good. 178b0753f9dSMatthew G. Knepley 179b0753f9dSMatthew G. Knepley Level: beginner 180b0753f9dSMatthew G. Knepley 18116a05f60SBarry Smith Note: 18216a05f60SBarry Smith Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 18316a05f60SBarry Smith domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 18416a05f60SBarry Smith a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 18516a05f60SBarry Smith (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 18616a05f60SBarry Smith 187af27ebaaSBarry Smith Developer Note: 188af27ebaaSBarry Smith Perhaps better to remove this since it matches `PCASMType` 189af27ebaaSBarry Smith 19016a05f60SBarry Smith .seealso: [](sec_pc), `PCGASM`, `PCASM`, `PC`, `PCGASMSetType()`, `PCASMType` 191b0753f9dSMatthew G. Knepley E*/ 1929371c9d4SSatish Balay typedef enum { 1939371c9d4SSatish Balay PC_GASM_BASIC = 3, 1949371c9d4SSatish Balay PC_GASM_RESTRICT = 1, 1959371c9d4SSatish Balay PC_GASM_INTERPOLATE = 2, 1969371c9d4SSatish Balay PC_GASM_NONE = 0 1979371c9d4SSatish Balay } PCGASMType; 198b0753f9dSMatthew G. Knepley 199b0753f9dSMatthew G. Knepley /*E 20016a05f60SBarry Smith PCCompositeType - Determines how two or more preconditioner are composed with the `PCType` of `PCCOMPOSITE` 201b0753f9dSMatthew G. Knepley 20209e53591SBarry Smith Values: 20309e53591SBarry Smith + `PC_COMPOSITE_ADDITIVE` - results from application of all preconditioners are added together 20409e53591SBarry Smith . `PC_COMPOSITE_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly 20509e53591SBarry Smith computed after the previous preconditioner application 20609e53591SBarry Smith . `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly 20709e53591SBarry Smith computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners) 208af27ebaaSBarry Smith . `PC_COMPOSITE_SPECIAL` - This is very special for a matrix of the form $ \alpha I + R + S$ 209af27ebaaSBarry Smith where the first preconditioner is built from $\alpha I + S$ and second from $\alpha I + R$ 21009e53591SBarry Smith . `PC_COMPOSITE_SCHUR` - composes the Schur complement of the matrix from two blocks, see `PCFIELDSPLIT` 21109e53591SBarry Smith - `PC_COMPOSITE_GKB` - the generalized Golub-Kahan bidiagonalization preconditioner, see `PCFIELDSPLIT` 212b0753f9dSMatthew G. Knepley 213b0753f9dSMatthew G. Knepley Level: beginner 214b0753f9dSMatthew G. Knepley 215af27ebaaSBarry Smith .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()`, `SNESCompositeType`, `PCCompositeSpecialSetAlpha()` 216b0753f9dSMatthew G. Knepley E*/ 2179371c9d4SSatish Balay typedef enum { 2189371c9d4SSatish Balay PC_COMPOSITE_ADDITIVE, 2199371c9d4SSatish Balay PC_COMPOSITE_MULTIPLICATIVE, 2209371c9d4SSatish Balay PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, 2219371c9d4SSatish Balay PC_COMPOSITE_SPECIAL, 2229371c9d4SSatish Balay PC_COMPOSITE_SCHUR, 2239371c9d4SSatish Balay PC_COMPOSITE_GKB 2249371c9d4SSatish Balay } PCCompositeType; 225b0753f9dSMatthew G. Knepley 226b0753f9dSMatthew G. Knepley /*E 22787497f52SBarry Smith PCFieldSplitSchurPreType - Determines how to precondition a Schur complement 228b0753f9dSMatthew G. Knepley 22916a05f60SBarry Smith Values: 23013044ca3SPierre Jolivet + `PC_FIELDSPLIT_SCHUR_PRE_SELF` - the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix. 23113044ca3SPierre Jolivet The only preconditioners that currently work with this symbolic representation matrix object are `PCLSC` and `PCHPDDM` 232*a077d33dSBarry Smith . `PC_FIELDSPLIT_SCHUR_PRE_SELFP` - the preconditioning for the Schur complement is generated from an explicitly-assembled approximation $Sp = A11 - A10 diag(A00)^{-1} A01$. 233*a077d33dSBarry Smith This is only a good preconditioner when $diag(A00)$ is a good preconditioner for $A00$. Optionally, $A00$ can be 23416a05f60SBarry Smith lumped before extracting the diagonal using the additional option `-fieldsplit_1_mat_schur_complement_ainv_type lump` 235*a077d33dSBarry Smith . `PC_FIELDSPLIT_SCHUR_PRE_A11` - the preconditioner for the Schur complement is generated from $A11$, not the Schur complement matrix 23616a05f60SBarry Smith . `PC_FIELDSPLIT_SCHUR_PRE_USER` - the preconditioner for the Schur complement is generated from the user provided matrix (pre argument 23716a05f60SBarry Smith to this function). 23816a05f60SBarry Smith - `PC_FIELDSPLIT_SCHUR_PRE_FULL` - the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation 23916a05f60SBarry Smith computed internally by `PCFIELDSPLIT` (this is expensive) useful mostly as a test that the Schur complement approach can work for your problem 24016a05f60SBarry Smith 241b0753f9dSMatthew G. Knepley Level: intermediate 242b0753f9dSMatthew G. Knepley 24309e53591SBarry Smith .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurPre()`, `PC` 244b0753f9dSMatthew G. Knepley E*/ 2459371c9d4SSatish Balay typedef enum { 2469371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_SELF, 2479371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_SELFP, 2489371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_A11, 2499371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_USER, 2509371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_FULL 2519371c9d4SSatish Balay } PCFieldSplitSchurPreType; 252b0753f9dSMatthew G. Knepley 253b0753f9dSMatthew G. Knepley /*E 254b0753f9dSMatthew G. Knepley PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 255b0753f9dSMatthew G. Knepley 25616a05f60SBarry Smith Values: 25716a05f60SBarry Smith + `PC_FIELDSPLIT_SCHUR_FACT_DIAG` - the preconditioner is solving `D` 25816a05f60SBarry Smith . `PC_FIELDSPLIT_SCHUR_FACT_LOWER` - the preconditioner is solving `L D` 25916a05f60SBarry Smith . `PC_FIELDSPLIT_SCHUR_FACT_UPPER` - the preconditioner is solving `D U` 26016a05f60SBarry Smith - `PC_FIELDSPLIT_SCHUR_FACT_FULL` - the preconditioner is solving `L(D U)` 26116a05f60SBarry Smith 26216a05f60SBarry Smith where the matrix is factorized as 26316a05f60SBarry Smith .vb 26416a05f60SBarry Smith (A B) = (1 0) (A 0) (1 Ainv*B) = L D U 26516a05f60SBarry Smith (C E) (C*Ainv 1) (0 S) (0 1) 26616a05f60SBarry Smith .ve 26716a05f60SBarry Smith 268b0753f9dSMatthew G. Knepley Level: intermediate 269b0753f9dSMatthew G. Knepley 27009e53591SBarry Smith .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurFactType()`, `PC` 271b0753f9dSMatthew G. Knepley E*/ 272b0753f9dSMatthew G. Knepley typedef enum { 273b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_DIAG, 274b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_LOWER, 275b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_UPPER, 276b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_FULL 277b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType; 278b0753f9dSMatthew G. Knepley 279b0753f9dSMatthew G. Knepley /*E 28087497f52SBarry Smith PCPARMSGlobalType - Determines the global preconditioner method in `PCPARMS` 281b0753f9dSMatthew G. Knepley 282b0753f9dSMatthew G. Knepley Level: intermediate 283b0753f9dSMatthew G. Knepley 28409e53591SBarry Smith .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetGlobal()`, `PC` 285b0753f9dSMatthew G. Knepley E*/ 2869371c9d4SSatish Balay typedef enum { 2879371c9d4SSatish Balay PC_PARMS_GLOBAL_RAS, 2889371c9d4SSatish Balay PC_PARMS_GLOBAL_SCHUR, 2899371c9d4SSatish Balay PC_PARMS_GLOBAL_BJ 2909371c9d4SSatish Balay } PCPARMSGlobalType; 2919d8081ecSMatthew G. Knepley 292b0753f9dSMatthew G. Knepley /*E 29387497f52SBarry Smith PCPARMSLocalType - Determines the local preconditioner method in `PCPARMS` 294b0753f9dSMatthew G. Knepley 295b0753f9dSMatthew G. Knepley Level: intermediate 296b0753f9dSMatthew G. Knepley 29709e53591SBarry Smith .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetLocal()`, `PC` 298b0753f9dSMatthew G. Knepley E*/ 2999371c9d4SSatish Balay typedef enum { 3009371c9d4SSatish Balay PC_PARMS_LOCAL_ILU0, 3019371c9d4SSatish Balay PC_PARMS_LOCAL_ILUK, 3029371c9d4SSatish Balay PC_PARMS_LOCAL_ILUT, 3039371c9d4SSatish Balay PC_PARMS_LOCAL_ARMS 3049371c9d4SSatish Balay } PCPARMSLocalType; 305b0753f9dSMatthew G. Knepley 306f0fc11ceSJed Brown /*J 30787497f52SBarry Smith PCGAMGType - type of generalized algebraic multigrid `PCGAMG` method 308b0753f9dSMatthew G. Knepley 30909e53591SBarry Smith Values: 31009e53591SBarry Smith + `PCGAMGAGG` - (the default) smoothed aggregation algorithm, robust, very well tested 311baca6076SPierre Jolivet . `PCGAMGGEO` - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not supported, reference implementation (2D) 312baca6076SPierre Jolivet - `PCGAMGCLASSICAL` - classical algebraic multigrid preconditioner, incomplete, not supported, reference implementation 31309e53591SBarry Smith 314b0753f9dSMatthew G. Knepley Level: intermediate 315b0753f9dSMatthew G. Knepley 31609e53591SBarry Smith .seealso: [](sec_pc), `PCGAMG`, `PCMG`, `PC`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()` 317f0fc11ceSJed Brown J*/ 318b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType; 319b0753f9dSMatthew G. Knepley #define PCGAMGAGG "agg" 320b0753f9dSMatthew G. Knepley #define PCGAMGGEO "geo" 321b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL "classical" 322b0753f9dSMatthew G. Knepley 323b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType; 324b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT "direct" 325b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard" 326b0753f9dSMatthew G. Knepley 327b0753f9dSMatthew G. Knepley /*E 328b0753f9dSMatthew G. Knepley PCMGType - Determines the type of multigrid method that is run. 329b0753f9dSMatthew G. Knepley 330b0753f9dSMatthew G. Knepley Values: 33187497f52SBarry Smith + `PC_MG_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `PCMGSetCycleType()` 33287497f52SBarry Smith . `PC_MG_ADDITIVE` - the additive multigrid preconditioner where all levels are 333b0753f9dSMatthew G. Knepley smoothed before updating the residual. This only uses the 334b0753f9dSMatthew G. Knepley down smoother, in the preconditioner the upper smoother is ignored 33587497f52SBarry Smith . `PC_MG_FULL` - same as multiplicative except one also performs grid sequencing, 336b0753f9dSMatthew G. Knepley that is starts on the coarsest grid, performs a cycle, interpolates 337b0753f9dSMatthew 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 338b0753f9dSMatthew G. Knepley algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 339b0753f9dSMatthew G. Knepley calls the V-cycle only on the coarser level and has a post-smoother instead. 34016a05f60SBarry Smith - `PC_MG_KASKADE` - like full multigrid except one never goes back to a coarser level from a finer 34116a05f60SBarry Smith 34216a05f60SBarry Smith Level: beginner 343b0753f9dSMatthew G. Knepley 34409e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()` 345b0753f9dSMatthew G. Knepley E*/ 3469371c9d4SSatish Balay typedef enum { 3479371c9d4SSatish Balay PC_MG_MULTIPLICATIVE, 3489371c9d4SSatish Balay PC_MG_ADDITIVE, 3499371c9d4SSatish Balay PC_MG_FULL, 3509371c9d4SSatish Balay PC_MG_KASKADE 3519371c9d4SSatish Balay } PCMGType; 352b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE; 353b0753f9dSMatthew G. Knepley 354b0753f9dSMatthew G. Knepley /*E 355b0753f9dSMatthew G. Knepley PCMGCycleType - Use V-cycle or W-cycle 356b0753f9dSMatthew G. Knepley 357b0753f9dSMatthew G. Knepley Values: 358af27ebaaSBarry Smith + `PC_MG_V_CYCLE` - use the V cycle 359af27ebaaSBarry Smith - `PC_MG_W_CYCLE` - use the W cycle 360b0753f9dSMatthew G. Knepley 36116a05f60SBarry Smith Level: beginner 36216a05f60SBarry Smith 36309e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()` 364b0753f9dSMatthew G. Knepley E*/ 3659371c9d4SSatish Balay typedef enum { 3669371c9d4SSatish Balay PC_MG_CYCLE_V = 1, 3679371c9d4SSatish Balay PC_MG_CYCLE_W = 2 3689371c9d4SSatish Balay } PCMGCycleType; 369b0753f9dSMatthew G. Knepley 370b0753f9dSMatthew G. Knepley /*E 3712134b1e4SBarry Smith PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process 3722134b1e4SBarry Smith 3732134b1e4SBarry Smith Values: 37487497f52SBarry Smith + `PC_MG_GALERKIN_PMAT` - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid 37587497f52SBarry Smith . `PC_MG_GALERKIN_MAT` - computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid 37687497f52SBarry Smith . `PC_MG_GALERKIN_BOTH` - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once 37787497f52SBarry Smith - `PC_MG_GALERKIN_NONE` - neither operator is computed via the Galerkin process, the user must provide the operator 3782134b1e4SBarry Smith 37909e53591SBarry Smith Level: beginner 38009e53591SBarry Smith 38109e53591SBarry Smith Note: 38287497f52SBarry Smith Users should never set `PC_MG_GALERKIN_EXTERNAL`, it is used by `PCHYPRE` and `PCML` 3832134b1e4SBarry Smith 38409e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()` 3852134b1e4SBarry Smith E*/ 3869371c9d4SSatish Balay typedef enum { 3879371c9d4SSatish Balay PC_MG_GALERKIN_BOTH, 3889371c9d4SSatish Balay PC_MG_GALERKIN_PMAT, 3899371c9d4SSatish Balay PC_MG_GALERKIN_MAT, 3909371c9d4SSatish Balay PC_MG_GALERKIN_NONE, 3919371c9d4SSatish Balay PC_MG_GALERKIN_EXTERNAL 3929371c9d4SSatish Balay } PCMGGalerkinType; 3932134b1e4SBarry Smith 3942134b1e4SBarry Smith /*E 395b0753f9dSMatthew G. Knepley PCExoticType - Face based or wirebasket based coarse grid space 396b0753f9dSMatthew G. Knepley 397b0753f9dSMatthew G. Knepley Level: beginner 398b0753f9dSMatthew G. Knepley 39909e53591SBarry Smith .seealso: [](sec_pc), `PCExoticSetType()`, `PCEXOTIC` 400b0753f9dSMatthew G. Knepley E*/ 4019371c9d4SSatish Balay typedef enum { 4029371c9d4SSatish Balay PC_EXOTIC_FACE, 4039371c9d4SSatish Balay PC_EXOTIC_WIREBASKET 4049371c9d4SSatish Balay } PCExoticType; 405b0753f9dSMatthew G. Knepley 4068c1cd74cSHong Zhang /*E 407bc960bbfSJed Brown PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains. 408bc960bbfSJed Brown 409bc960bbfSJed Brown Values: 41087497f52SBarry Smith + `PC_BDDC_INTERFACE_EXT_DIRICHLET` - solves Dirichlet interior problem; this is the standard BDDC algorithm 41187497f52SBarry Smith - `PC_BDDC_INTERFACE_EXT_LUMP` - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP" 412bc960bbfSJed Brown 41316a05f60SBarry Smith Level: intermediate 41416a05f60SBarry Smith 41509e53591SBarry Smith .seealso: [](sec_pc), `PCBDDC`, `PC` 416bc960bbfSJed Brown E*/ 417bc960bbfSJed Brown typedef enum { 418bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_DIRICHLET, 419bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_LUMP 420bc960bbfSJed Brown } PCBDDCInterfaceExtType; 421bc960bbfSJed Brown 422bc960bbfSJed Brown /*E 423f3b08a26SMatthew G. Knepley PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation 424f3b08a26SMatthew G. Knepley 425f3b08a26SMatthew G. Knepley Level: beginner 426f3b08a26SMatthew G. Knepley 42709e53591SBarry Smith .seealso: [](sec_pc), `PCMGSetAdaptCoarseSpaceType()`, `PCMG`, `PC` 428f3b08a26SMatthew G. Knepley E*/ 4299371c9d4SSatish Balay typedef enum { 4309371c9d4SSatish Balay PCMG_ADAPT_NONE, 4319371c9d4SSatish Balay PCMG_ADAPT_POLYNOMIAL, 4329371c9d4SSatish Balay PCMG_ADAPT_HARMONIC, 4339371c9d4SSatish Balay PCMG_ADAPT_EIGENVECTOR, 4349371c9d4SSatish Balay PCMG_ADAPT_GENERALIZED_EIGENVECTOR, 4359371c9d4SSatish Balay PCMG_ADAPT_GDSW 4369371c9d4SSatish Balay } PCMGCoarseSpaceType; 437f3b08a26SMatthew G. Knepley 438f3b08a26SMatthew G. Knepley /*E 43916a05f60SBarry Smith PCPatchConstructType - The algorithm used to construct patches for the `PCPATCH` preconditioner 4404bbf5ea8SMatthew G. Knepley 4414bbf5ea8SMatthew G. Knepley Level: beginner 4424bbf5ea8SMatthew G. Knepley 44309e53591SBarry Smith .seealso: [](sec_pc), `PCPatchSetConstructType()`, `PCPATCH`, `PC` 4444bbf5ea8SMatthew G. Knepley E*/ 4459371c9d4SSatish Balay typedef enum { 4469371c9d4SSatish Balay PC_PATCH_STAR, 4479371c9d4SSatish Balay PC_PATCH_VANKA, 4489371c9d4SSatish Balay PC_PATCH_PARDECOMP, 4499371c9d4SSatish Balay PC_PATCH_USER, 4509371c9d4SSatish Balay PC_PATCH_PYTHON 4519371c9d4SSatish Balay } PCPatchConstructType; 4524bbf5ea8SMatthew G. Knepley 4534bbf5ea8SMatthew G. Knepley /*E 454e26ad46dSJakub Kruzik PCDeflationSpaceType - Type of deflation 455e662bc50SJakub Kruzik 456e662bc50SJakub Kruzik Values: 45787497f52SBarry Smith + `PC_DEFLATION_SPACE_HAAR` - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off 45809e53591SBarry Smith . `PC_DEFLATION_SPACE_DB2` - `MATCOMPOSITE` of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet) 45987497f52SBarry Smith . `PC_DEFLATION_SPACE_DB4` - same as above, but with db4 (4 coefficient Daubechies) 46087497f52SBarry Smith . `PC_DEFLATION_SPACE_DB8` - same as above, but with db8 (8 coefficient Daubechies) 46187497f52SBarry Smith . `PC_DEFLATION_SPACE_DB16` - same as above, but with db16 (16 coefficient Daubechies) 46287497f52SBarry Smith . `PC_DEFLATION_SPACE_BIORTH22` - same as above, but with biorthogonal 2.2 (6 coefficients) 46387497f52SBarry Smith . `PC_DEFLATION_SPACE_MEYER` - same as above, but with Meyer/FIR (62 coefficients) 464d5b43468SJose E. Roman . `PC_DEFLATION_SPACE_AGGREGATION` - aggregates local indices (given by operator matrix distribution) into a subdomain 46587497f52SBarry Smith - `PC_DEFLATION_SPACE_USER` - indicates space set by user 466e662bc50SJakub Kruzik 4671fdb00f9SJakub Kruzik Level: intermediate 4681fdb00f9SJakub Kruzik 46909e53591SBarry Smith Note: 47009e53591SBarry Smith Wavelet-based space (except Haar) can be used in multilevel deflation. 47109e53591SBarry Smith 47209e53591SBarry Smith .seealso: [](sec_pc), `PCDeflationSetSpaceToCompute()`, `PCDEFLATION`, `PC` 473e662bc50SJakub Kruzik E*/ 474e662bc50SJakub Kruzik typedef enum { 475e662bc50SJakub Kruzik PC_DEFLATION_SPACE_HAAR, 476e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB2, 477e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB4, 478e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB8, 479e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB16, 480e662bc50SJakub Kruzik PC_DEFLATION_SPACE_BIORTH22, 481e662bc50SJakub Kruzik PC_DEFLATION_SPACE_MEYER, 482e662bc50SJakub Kruzik PC_DEFLATION_SPACE_AGGREGATION, 483e662bc50SJakub Kruzik PC_DEFLATION_SPACE_USER 484e662bc50SJakub Kruzik } PCDeflationSpaceType; 485e662bc50SJakub Kruzik 486e662bc50SJakub Kruzik /*E 48787497f52SBarry Smith PCHPDDMCoarseCorrectionType - Type of coarse correction used by `PCHPDDM` 48838cfc46eSPierre Jolivet 48938cfc46eSPierre Jolivet Values: 49009e53591SBarry Smith + `PC_HPDDM_COARSE_CORRECTION_DEFLATED` (default) - eq. (1) in `PCHPDDMShellApply()` 49187497f52SBarry Smith . `PC_HPDDM_COARSE_CORRECTION_ADDITIVE` - eq. (2) 492aa1539e9SPierre Jolivet . `PC_HPDDM_COARSE_CORRECTION_BALANCED` - eq. (3) 493aa1539e9SPierre Jolivet - `PC_HPDDM_COARSE_CORRECTION_NONE` - no coarse correction (mostly useful for debugging) 49438cfc46eSPierre Jolivet 49516a05f60SBarry Smith Level: intermediate 49616a05f60SBarry Smith 49709e53591SBarry Smith .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCSetType()`, `PCHPDDMShellApply()` 49838cfc46eSPierre Jolivet E*/ 4999371c9d4SSatish Balay typedef enum { 5009371c9d4SSatish Balay PC_HPDDM_COARSE_CORRECTION_DEFLATED, 5019371c9d4SSatish Balay PC_HPDDM_COARSE_CORRECTION_ADDITIVE, 502aa1539e9SPierre Jolivet PC_HPDDM_COARSE_CORRECTION_BALANCED, 503aa1539e9SPierre Jolivet PC_HPDDM_COARSE_CORRECTION_NONE 5049371c9d4SSatish Balay } PCHPDDMCoarseCorrectionType; 50538cfc46eSPierre Jolivet 50638cfc46eSPierre Jolivet /*E 50713044ca3SPierre Jolivet PCHPDDMSchurPreType - Type of `PCHPDDM` preconditioner for a `MATSCHURCOMPLEMENT` generated by `PCFIELDSPLIT` with `PCFieldSplitSchurPreType` set to `PC_FIELDSPLIT_SCHUR_PRE_SELF` 50813044ca3SPierre Jolivet 50913044ca3SPierre Jolivet Values: 51013044ca3SPierre Jolivet + `PC_HPDDM_SCHUR_PRE_LEAST_SQUARES` (default) - only with a near-zero A11 block and A10 = A01^T; a preconditioner for solving A01^T A00^-1 A01 x = b is built by approximating the Schur complement with (inv(sqrt(diag(A00))) A01)^T (inv(sqrt(diag(A00))) A01) and by considering the associated linear least squares problem 51113044ca3SPierre Jolivet - `PC_HPDDM_SCHUR_PRE_GENEO` - only with A10 = A01^T, `PCHPDDMSetAuxiliaryMat()` called on the `PC` of the A00 block, and if A11 is nonzero, then `PCHPDDMSetAuxiliaryMat()` must be called on the associated `PC` as well (it is built automatically for the user otherwise); the Schur complement `PC` is set internally to `PCKSP`, with the prefix `-fieldsplit_1_pc_hpddm_`; the operator associated to the `PC` is spectrally equivalent to the original Schur complement 51213044ca3SPierre Jolivet 51313044ca3SPierre Jolivet Level: advanced 51413044ca3SPierre Jolivet 51513044ca3SPierre Jolivet .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCFIELDSPLIT`, `PC_FIELDSPLIT_SCHUR_PRE_SELF`, `PCFieldSplitSetSchurPre()`, `PCHPDDMSetAuxiliaryMat()` 51613044ca3SPierre Jolivet E*/ 51713044ca3SPierre Jolivet typedef enum { 51813044ca3SPierre Jolivet PC_HPDDM_SCHUR_PRE_LEAST_SQUARES, 51913044ca3SPierre Jolivet PC_HPDDM_SCHUR_PRE_GENEO, 52013044ca3SPierre Jolivet } PCHPDDMSchurPreType; 52113044ca3SPierre Jolivet 52213044ca3SPierre Jolivet /*E 52387497f52SBarry Smith PCFailedReason - indicates type of `PC` failure 5248c1cd74cSHong Zhang 5258c1cd74cSHong Zhang Level: beginner 5268c1cd74cSHong Zhang 52709e53591SBarry Smith .seealso: [](sec_pc), `PC` 5288c1cd74cSHong Zhang E*/ 5299371c9d4SSatish Balay typedef enum { 5309371c9d4SSatish Balay PC_SETUP_ERROR = -1, 5319371c9d4SSatish Balay PC_NOERROR, 5329371c9d4SSatish Balay PC_FACTOR_STRUCT_ZEROPIVOT, 5339371c9d4SSatish Balay PC_FACTOR_NUMERIC_ZEROPIVOT, 5349371c9d4SSatish Balay PC_FACTOR_OUTMEMORY, 5359371c9d4SSatish Balay PC_FACTOR_OTHER, 53687b47708SBarry Smith PC_INCONSISTENT_RHS, 5379371c9d4SSatish Balay PC_SUBPC_ERROR 5389371c9d4SSatish Balay } PCFailedReason; 539ce7c7f2fSMark Adams 540ce7c7f2fSMark Adams /*E 541ce7c7f2fSMark Adams PCGAMGLayoutType - Layout for reduced grids 542ce7c7f2fSMark Adams 543ce7c7f2fSMark Adams Level: intermediate 544ce7c7f2fSMark Adams 54509e53591SBarry Smith .seealso: [](sec_pc), `PCGAMG`, `PC`, `PCGAMGSetCoarseGridLayoutType()` 546ce7c7f2fSMark Adams E*/ 5479371c9d4SSatish Balay typedef enum { 5489371c9d4SSatish Balay PCGAMG_LAYOUT_COMPACT, 5499371c9d4SSatish Balay PCGAMG_LAYOUT_SPREAD 5509371c9d4SSatish Balay } PCGAMGLayoutType; 551