1a4963045SJacob Faibussowitsch #pragma once 2b0753f9dSMatthew G. Knepley 3*ce78bad3SBarry Smith /* MANSEC = KSP */ 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 1116a05f60SBarry Smith .seealso: [](doc_linsolve), [](sec_pc), `PCCreate()`, `PCSetType()`, `PCType` 12b0753f9dSMatthew G. Knepley S*/ 13b0753f9dSMatthew G. Knepley typedef struct _p_PC *PC; 14b0753f9dSMatthew G. Knepley 15b0753f9dSMatthew G. Knepley /*J 160b4b7b1cSBarry Smith PCType - String with the name of a PETSc preconditioner. These are all the preconditioners and direct solvers that PETSc provides. 17b0753f9dSMatthew G. Knepley 18b0753f9dSMatthew G. Knepley Level: beginner 19b0753f9dSMatthew G. Knepley 200b4b7b1cSBarry Smith Notes: 210b4b7b1cSBarry Smith Use `PCSetType()` or the options database key `-pc_type` to set the preconditioner to use with a given `PC` object 220b4b7b1cSBarry Smith 2387497f52SBarry Smith `PCRegister()` is used to register preconditioners that are then accessible via `PCSetType()` 24b0753f9dSMatthew G. Knepley 2516a05f60SBarry Smith .seealso: [](doc_linsolve), [](sec_pc), `PCSetType()`, `PC`, `PCCreate()`, `PCRegister()`, `PCSetFromOptions()`, `PCLU`, `PCJACOBI`, `PCBJACOBI` 26b0753f9dSMatthew G. Knepley J*/ 27b0753f9dSMatthew G. Knepley typedef const char *PCType; 28b0753f9dSMatthew G. Knepley #define PCNONE "none" 29b0753f9dSMatthew G. Knepley #define PCJACOBI "jacobi" 30b0753f9dSMatthew G. Knepley #define PCSOR "sor" 31b0753f9dSMatthew G. Knepley #define PCLU "lu" 32a2fc1e05SToby Isaac #define PCQR "qr" 33b0753f9dSMatthew G. Knepley #define PCSHELL "shell" 34e6f8f311SMark Adams #define PCAMGX "amgx" 35b0753f9dSMatthew G. Knepley #define PCBJACOBI "bjacobi" 36b0753f9dSMatthew G. Knepley #define PCMG "mg" 37b0753f9dSMatthew G. Knepley #define PCEISENSTAT "eisenstat" 38b0753f9dSMatthew G. Knepley #define PCILU "ilu" 39b0753f9dSMatthew G. Knepley #define PCICC "icc" 40b0753f9dSMatthew G. Knepley #define PCASM "asm" 41b0753f9dSMatthew G. Knepley #define PCGASM "gasm" 42b0753f9dSMatthew G. Knepley #define PCKSP "ksp" 43e607c864SMark Adams #define PCBJKOKKOS "bjkokkos" 44b0753f9dSMatthew G. Knepley #define PCCOMPOSITE "composite" 45b0753f9dSMatthew G. Knepley #define PCREDUNDANT "redundant" 46b0753f9dSMatthew G. Knepley #define PCSPAI "spai" 47b0753f9dSMatthew G. Knepley #define PCNN "nn" 48b0753f9dSMatthew G. Knepley #define PCCHOLESKY "cholesky" 49b0753f9dSMatthew G. Knepley #define PCPBJACOBI "pbjacobi" 500da83c2eSBarry Smith #define PCVPBJACOBI "vpbjacobi" 51b0753f9dSMatthew G. Knepley #define PCMAT "mat" 52b0753f9dSMatthew G. Knepley #define PCHYPRE "hypre" 53b0753f9dSMatthew G. Knepley #define PCPARMS "parms" 54b0753f9dSMatthew G. Knepley #define PCFIELDSPLIT "fieldsplit" 55b0753f9dSMatthew G. Knepley #define PCTFS "tfs" 56b0753f9dSMatthew G. Knepley #define PCML "ml" 57b0753f9dSMatthew G. Knepley #define PCGALERKIN "galerkin" 58b0753f9dSMatthew G. Knepley #define PCEXOTIC "exotic" 59b0753f9dSMatthew G. Knepley #define PCCP "cp" 60b0753f9dSMatthew G. Knepley #define PCLSC "lsc" 61b0753f9dSMatthew G. Knepley #define PCPYTHON "python" 62b0753f9dSMatthew G. Knepley #define PCPFMG "pfmg" 631c188c59Sftrigaux #define PCSMG "smg" 64b0753f9dSMatthew G. Knepley #define PCSYSPFMG "syspfmg" 65b0753f9dSMatthew G. Knepley #define PCREDISTRIBUTE "redistribute" 66b0753f9dSMatthew G. Knepley #define PCSVD "svd" 67b0753f9dSMatthew G. Knepley #define PCGAMG "gamg" 684b3f184cSKarl Rupp #define PCCHOWILUVIENNACL "chowiluviennacl" 6970baa948SKarl Rupp #define PCROWSCALINGVIENNACL "rowscalingviennacl" 7007598726SKarl Rupp #define PCSAVIENNACL "saviennacl" 71b0753f9dSMatthew G. Knepley #define PCBDDC "bddc" 72b0753f9dSMatthew G. Knepley #define PCKACZMARZ "kaczmarz" 7368ddcbeaSBarry Smith #define PCTELESCOPE "telescope" 744bbf5ea8SMatthew G. Knepley #define PCPATCH "patch" 75b9ac7092SAlp Dener #define PCLMVM "lmvm" 76360ee056SFande Kong #define PCHMG "hmg" 7737eeb815SJakub Kruzik #define PCDEFLATION "deflation" 7838cfc46eSPierre Jolivet #define PCHPDDM "hpddm" 7953022affSStefano Zampini #define PCH2OPUS "h2opus" 80f1f2ae84SBarry Smith #define PCMPI "mpi" 81b0753f9dSMatthew G. Knepley 82b0753f9dSMatthew G. Knepley /*E 830b4b7b1cSBarry Smith PCSide - Determines if the preconditioner is to be applied to the left, right 840b4b7b1cSBarry Smith or symmetrically around the operator in `KSPSolve()`. 85b0753f9dSMatthew G. Knepley 8616a05f60SBarry Smith Values: 8716a05f60SBarry Smith + `PC_LEFT` - applied after the operator is applied 8816a05f60SBarry Smith . `PC_RIGHT` - applied before the operator is applied 8916a05f60SBarry 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. 9016a05f60SBarry Smith 91b0753f9dSMatthew G. Knepley Level: beginner 92b0753f9dSMatthew G. Knepley 9316a05f60SBarry Smith Note: 9416a05f60SBarry Smith Certain `KSPType` support only a subset of `PCSide` values 9516a05f60SBarry Smith 960b4b7b1cSBarry Smith .seealso: [](sec_pc), `PC`, `KSPSetPCSide()`, `KSP`, `KSPType`, `KSPGetPCSide()`, `KSPSolve()` 97b0753f9dSMatthew G. Knepley E*/ 989371c9d4SSatish Balay typedef enum { 999371c9d4SSatish Balay PC_SIDE_DEFAULT = -1, 100*ce78bad3SBarry Smith PC_LEFT = 0, 101*ce78bad3SBarry Smith PC_RIGHT = 1, 102*ce78bad3SBarry Smith PC_SYMMETRIC = 2 1039371c9d4SSatish Balay } PCSide; 104b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 105b0753f9dSMatthew G. Knepley 106b0753f9dSMatthew G. Knepley /*E 10739b1ba1bSPierre Jolivet PCRichardsonConvergedReason - reason a `PCApplyRichardson()` method terminated 108b0753f9dSMatthew G. Knepley 109b0753f9dSMatthew G. Knepley Level: advanced 110b0753f9dSMatthew G. Knepley 11139b1ba1bSPierre Jolivet .seealso: [](sec_pc), `KSPRICHARDSON`, `PC`, `PCApplyRichardson()` 112b0753f9dSMatthew G. Knepley E*/ 113b0753f9dSMatthew G. Knepley typedef enum { 114dd460d27SBarry Smith PCRICHARDSON_NOT_SET = 0, 115b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_RTOL = 2, 116b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ATOL = 3, 117b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ITS = 4, 1189371c9d4SSatish Balay PCRICHARDSON_DIVERGED_DTOL = -4 1199371c9d4SSatish Balay } PCRichardsonConvergedReason; 120b0753f9dSMatthew G. Knepley 121b0753f9dSMatthew G. Knepley /*E 1220b4b7b1cSBarry Smith PCJacobiType - Determines what elements of the matrix are used to form the Jacobi preconditioner, that is with the `PCType` of `PCJACOBI` 123b0753f9dSMatthew G. Knepley 12409e53591SBarry Smith Values: 12509e53591SBarry Smith + `PC_JACOBI_DIAGONAL` - use the diagonal entry, if it is zero use one 126eede4a3fSMark Adams . `PC_JACOBI_ROWL1` - add sum of absolute values in row i, j != i, to diag_ii 12709e53591SBarry Smith . `PC_JACOBI_ROWMAX` - use the maximum absolute value in the row 12809e53591SBarry Smith - `PC_JACOBI_ROWSUM` - use the sum of the values in the row (not the absolute values) 12909e53591SBarry Smith 130b0753f9dSMatthew G. Knepley Level: intermediate 131b0753f9dSMatthew G. Knepley 13209e53591SBarry Smith .seealso: [](sec_pc), `PCJACOBI`, `PC` 133b0753f9dSMatthew G. Knepley E*/ 1349371c9d4SSatish Balay typedef enum { 1359371c9d4SSatish Balay PC_JACOBI_DIAGONAL, 136eede4a3fSMark Adams PC_JACOBI_ROWL1, 1379371c9d4SSatish Balay PC_JACOBI_ROWMAX, 1389371c9d4SSatish Balay PC_JACOBI_ROWSUM 1399371c9d4SSatish Balay } PCJacobiType; 140b0753f9dSMatthew G. Knepley 141b0753f9dSMatthew G. Knepley /*E 1420b4b7b1cSBarry Smith PCASMType - Determines the type of additive Schwarz method, `PCASM`, to use 143b0753f9dSMatthew G. Knepley 14409e53591SBarry Smith Values: 14509e53591SBarry Smith + `PC_ASM_BASIC` - Symmetric version where residuals from the ghost points are used 14609e53591SBarry Smith and computed values in ghost regions are added together. 147af27ebaaSBarry Smith Classical standard additive Schwarz as introduced in {cite}`dryja1987additive`. 14809e53591SBarry Smith . `PC_ASM_RESTRICT` - Residuals from ghost points are used but computed values in ghost 149af27ebaaSBarry Smith region are discarded {cite}`cs99`. Default. 15009e53591SBarry Smith . `PC_ASM_INTERPOLATE` - Residuals from ghost points are not used, computed values in ghost 15109e53591SBarry Smith region are added back in. 15209e53591SBarry Smith - `PC_ASM_NONE` - Residuals from ghost points are not used, computed ghost values are 153af27ebaaSBarry Smith discarded. Not very good. 154b0753f9dSMatthew G. Knepley 155b0753f9dSMatthew G. Knepley Level: beginner 156b0753f9dSMatthew G. Knepley 15716a05f60SBarry Smith .seealso: [](sec_pc), `PC`, `PCASM`, `PCASMSetType()`, `PCGASMType` 158b0753f9dSMatthew G. Knepley E*/ 1599371c9d4SSatish Balay typedef enum { 1609371c9d4SSatish Balay PC_ASM_BASIC = 3, 1619371c9d4SSatish Balay PC_ASM_RESTRICT = 1, 1629371c9d4SSatish Balay PC_ASM_INTERPOLATE = 2, 1639371c9d4SSatish Balay PC_ASM_NONE = 0 1649371c9d4SSatish Balay } PCASMType; 165b0753f9dSMatthew G. Knepley 166b0753f9dSMatthew G. Knepley /*E 1670b4b7b1cSBarry Smith PCGASMType - Determines the type of generalized additive Schwarz method to use (differs from `PCASM` in allowing multiple processors per subdomain) with the `PCType` of `PCGASM` 168b0753f9dSMatthew G. Knepley 16909e53591SBarry Smith Values: 17009e53591SBarry Smith + `PC_GASM_BASIC` - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 17109e53591SBarry Smith over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 172af27ebaaSBarry Smith from neighboring subdomains. Classical standard additive Schwarz {cite}`dryja1987additive`. 17309e53591SBarry Smith . `PC_GASM_RESTRICT` - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 17409e53591SBarry Smith (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 17509e53591SBarry Smith each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 176af27ebaaSBarry Smith assumption) {cite}`cs99`. Default. 17709e53591SBarry Smith . `PC_GASM_INTERPOLATE` - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 17809e53591SBarry Smith applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 17909e53591SBarry Smith from neighboring subdomains. 180af27ebaaSBarry Smith - `PC_GASM_NONE` - Residuals and corrections are zeroed out outside the local subdomains. Not very good. 181b0753f9dSMatthew G. Knepley 182b0753f9dSMatthew G. Knepley Level: beginner 183b0753f9dSMatthew G. Knepley 18416a05f60SBarry Smith Note: 18516a05f60SBarry Smith Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 1860b4b7b1cSBarry Smith domain, while the outer subdomains contain the inner subdomains and overlap with each other. The `PCGASM` preconditioner will compute 18716a05f60SBarry Smith a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 18816a05f60SBarry Smith (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 18916a05f60SBarry Smith 190af27ebaaSBarry Smith Developer Note: 191af27ebaaSBarry Smith Perhaps better to remove this since it matches `PCASMType` 192af27ebaaSBarry Smith 19316a05f60SBarry Smith .seealso: [](sec_pc), `PCGASM`, `PCASM`, `PC`, `PCGASMSetType()`, `PCASMType` 194b0753f9dSMatthew G. Knepley E*/ 1959371c9d4SSatish Balay typedef enum { 1969371c9d4SSatish Balay PC_GASM_BASIC = 3, 1979371c9d4SSatish Balay PC_GASM_RESTRICT = 1, 1989371c9d4SSatish Balay PC_GASM_INTERPOLATE = 2, 1999371c9d4SSatish Balay PC_GASM_NONE = 0 2009371c9d4SSatish Balay } PCGASMType; 201b0753f9dSMatthew G. Knepley 202b0753f9dSMatthew G. Knepley /*E 20316a05f60SBarry Smith PCCompositeType - Determines how two or more preconditioner are composed with the `PCType` of `PCCOMPOSITE` 204b0753f9dSMatthew G. Knepley 20509e53591SBarry Smith Values: 20609e53591SBarry Smith + `PC_COMPOSITE_ADDITIVE` - results from application of all preconditioners are added together 20709e53591SBarry Smith . `PC_COMPOSITE_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly 20809e53591SBarry Smith computed after the previous preconditioner application 20909e53591SBarry Smith . `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly 21009e53591SBarry Smith computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners) 211af27ebaaSBarry Smith . `PC_COMPOSITE_SPECIAL` - This is very special for a matrix of the form $ \alpha I + R + S$ 212af27ebaaSBarry Smith where the first preconditioner is built from $\alpha I + S$ and second from $\alpha I + R$ 21309e53591SBarry Smith . `PC_COMPOSITE_SCHUR` - composes the Schur complement of the matrix from two blocks, see `PCFIELDSPLIT` 21409e53591SBarry Smith - `PC_COMPOSITE_GKB` - the generalized Golub-Kahan bidiagonalization preconditioner, see `PCFIELDSPLIT` 215b0753f9dSMatthew G. Knepley 216b0753f9dSMatthew G. Knepley Level: beginner 217b0753f9dSMatthew G. Knepley 218af27ebaaSBarry Smith .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()`, `SNESCompositeType`, `PCCompositeSpecialSetAlpha()` 219b0753f9dSMatthew G. Knepley E*/ 2209371c9d4SSatish Balay typedef enum { 2219371c9d4SSatish Balay PC_COMPOSITE_ADDITIVE, 2229371c9d4SSatish Balay PC_COMPOSITE_MULTIPLICATIVE, 2239371c9d4SSatish Balay PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, 2249371c9d4SSatish Balay PC_COMPOSITE_SPECIAL, 2259371c9d4SSatish Balay PC_COMPOSITE_SCHUR, 2269371c9d4SSatish Balay PC_COMPOSITE_GKB 2279371c9d4SSatish Balay } PCCompositeType; 228b0753f9dSMatthew G. Knepley 229b0753f9dSMatthew G. Knepley /*E 2300b4b7b1cSBarry Smith PCFieldSplitSchurPreType - Determines how to precondition a Schur complement arising with the `PCType` of `PCFIELDSPLIT` 231b0753f9dSMatthew G. Knepley 23216a05f60SBarry Smith Values: 23313044ca3SPierre Jolivet + `PC_FIELDSPLIT_SCHUR_PRE_SELF` - the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix. 23413044ca3SPierre Jolivet The only preconditioners that currently work with this symbolic representation matrix object are `PCLSC` and `PCHPDDM` 235a077d33dSBarry 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$. 236a077d33dSBarry Smith This is only a good preconditioner when $diag(A00)$ is a good preconditioner for $A00$. Optionally, $A00$ can be 23716a05f60SBarry Smith lumped before extracting the diagonal using the additional option `-fieldsplit_1_mat_schur_complement_ainv_type lump` 238a077d33dSBarry Smith . `PC_FIELDSPLIT_SCHUR_PRE_A11` - the preconditioner for the Schur complement is generated from $A11$, not the Schur complement matrix 23916a05f60SBarry Smith . `PC_FIELDSPLIT_SCHUR_PRE_USER` - the preconditioner for the Schur complement is generated from the user provided matrix (pre argument 24016a05f60SBarry Smith to this function). 24116a05f60SBarry Smith - `PC_FIELDSPLIT_SCHUR_PRE_FULL` - the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation 24216a05f60SBarry Smith computed internally by `PCFIELDSPLIT` (this is expensive) useful mostly as a test that the Schur complement approach can work for your problem 24316a05f60SBarry Smith 244b0753f9dSMatthew G. Knepley Level: intermediate 245b0753f9dSMatthew G. Knepley 24609e53591SBarry Smith .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurPre()`, `PC` 247b0753f9dSMatthew G. Knepley E*/ 2489371c9d4SSatish Balay typedef enum { 2499371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_SELF, 2509371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_SELFP, 2519371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_A11, 2529371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_USER, 2539371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_FULL 2549371c9d4SSatish Balay } PCFieldSplitSchurPreType; 255b0753f9dSMatthew G. Knepley 256b0753f9dSMatthew G. Knepley /*E 2570b4b7b1cSBarry Smith PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use with the `PCType` of `PCFIELDSPLIT` 258b0753f9dSMatthew G. Knepley 25916a05f60SBarry Smith Values: 26016a05f60SBarry Smith + `PC_FIELDSPLIT_SCHUR_FACT_DIAG` - the preconditioner is solving `D` 26116a05f60SBarry Smith . `PC_FIELDSPLIT_SCHUR_FACT_LOWER` - the preconditioner is solving `L D` 26216a05f60SBarry Smith . `PC_FIELDSPLIT_SCHUR_FACT_UPPER` - the preconditioner is solving `D U` 26316a05f60SBarry Smith - `PC_FIELDSPLIT_SCHUR_FACT_FULL` - the preconditioner is solving `L(D U)` 26416a05f60SBarry Smith 26516a05f60SBarry Smith where the matrix is factorized as 26616a05f60SBarry Smith .vb 26716a05f60SBarry Smith (A B) = (1 0) (A 0) (1 Ainv*B) = L D U 26816a05f60SBarry Smith (C E) (C*Ainv 1) (0 S) (0 1) 26916a05f60SBarry Smith .ve 27016a05f60SBarry Smith 271b0753f9dSMatthew G. Knepley Level: intermediate 272b0753f9dSMatthew G. Knepley 27309e53591SBarry Smith .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurFactType()`, `PC` 274b0753f9dSMatthew G. Knepley E*/ 275b0753f9dSMatthew G. Knepley typedef enum { 276b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_DIAG, 277b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_LOWER, 278b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_UPPER, 279b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_FULL 280b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType; 281b0753f9dSMatthew G. Knepley 282b0753f9dSMatthew G. Knepley /*E 28387497f52SBarry Smith PCPARMSGlobalType - Determines the global preconditioner method in `PCPARMS` 284b0753f9dSMatthew G. Knepley 285b0753f9dSMatthew G. Knepley Level: intermediate 286b0753f9dSMatthew G. Knepley 28709e53591SBarry Smith .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetGlobal()`, `PC` 288b0753f9dSMatthew G. Knepley E*/ 2899371c9d4SSatish Balay typedef enum { 2909371c9d4SSatish Balay PC_PARMS_GLOBAL_RAS, 2919371c9d4SSatish Balay PC_PARMS_GLOBAL_SCHUR, 2929371c9d4SSatish Balay PC_PARMS_GLOBAL_BJ 2939371c9d4SSatish Balay } PCPARMSGlobalType; 2949d8081ecSMatthew G. Knepley 295b0753f9dSMatthew G. Knepley /*E 29687497f52SBarry Smith PCPARMSLocalType - Determines the local preconditioner method in `PCPARMS` 297b0753f9dSMatthew G. Knepley 298b0753f9dSMatthew G. Knepley Level: intermediate 299b0753f9dSMatthew G. Knepley 30009e53591SBarry Smith .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetLocal()`, `PC` 301b0753f9dSMatthew G. Knepley E*/ 3029371c9d4SSatish Balay typedef enum { 3039371c9d4SSatish Balay PC_PARMS_LOCAL_ILU0, 3049371c9d4SSatish Balay PC_PARMS_LOCAL_ILUK, 3059371c9d4SSatish Balay PC_PARMS_LOCAL_ILUT, 3069371c9d4SSatish Balay PC_PARMS_LOCAL_ARMS 3079371c9d4SSatish Balay } PCPARMSLocalType; 308b0753f9dSMatthew G. Knepley 309f0fc11ceSJed Brown /*J 31087497f52SBarry Smith PCGAMGType - type of generalized algebraic multigrid `PCGAMG` method 311b0753f9dSMatthew G. Knepley 31209e53591SBarry Smith Values: 31309e53591SBarry Smith + `PCGAMGAGG` - (the default) smoothed aggregation algorithm, robust, very well tested 314baca6076SPierre Jolivet . `PCGAMGGEO` - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not supported, reference implementation (2D) 315baca6076SPierre Jolivet - `PCGAMGCLASSICAL` - classical algebraic multigrid preconditioner, incomplete, not supported, reference implementation 31609e53591SBarry Smith 317b0753f9dSMatthew G. Knepley Level: intermediate 318b0753f9dSMatthew G. Knepley 31909e53591SBarry Smith .seealso: [](sec_pc), `PCGAMG`, `PCMG`, `PC`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()` 320f0fc11ceSJed Brown J*/ 321b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType; 322b0753f9dSMatthew G. Knepley #define PCGAMGAGG "agg" 323b0753f9dSMatthew G. Knepley #define PCGAMGGEO "geo" 324b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL "classical" 325b0753f9dSMatthew G. Knepley 326b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType; 327b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT "direct" 328b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard" 329b0753f9dSMatthew G. Knepley 330b0753f9dSMatthew G. Knepley /*E 3310b4b7b1cSBarry Smith PCMGType - Determines the type of multigrid method that is run with the `PCType` of `PCMG` 332b0753f9dSMatthew G. Knepley 333b0753f9dSMatthew G. Knepley Values: 33487497f52SBarry Smith + `PC_MG_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `PCMGSetCycleType()` 33587497f52SBarry Smith . `PC_MG_ADDITIVE` - the additive multigrid preconditioner where all levels are 336b0753f9dSMatthew G. Knepley smoothed before updating the residual. This only uses the 337b0753f9dSMatthew G. Knepley down smoother, in the preconditioner the upper smoother is ignored 33887497f52SBarry Smith . `PC_MG_FULL` - same as multiplicative except one also performs grid sequencing, 339b0753f9dSMatthew G. Knepley that is starts on the coarsest grid, performs a cycle, interpolates 340b0753f9dSMatthew 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 341b0753f9dSMatthew G. Knepley algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 342b0753f9dSMatthew G. Knepley calls the V-cycle only on the coarser level and has a post-smoother instead. 3430b4b7b1cSBarry Smith - `PC_MG_KASKADE` - Cascadic or Kaskadic multigrid, like full multigrid except one never goes back to a coarser level from a finer 34416a05f60SBarry Smith 34516a05f60SBarry Smith Level: beginner 346b0753f9dSMatthew G. Knepley 34709e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()` 348b0753f9dSMatthew G. Knepley E*/ 3499371c9d4SSatish Balay typedef enum { 3509371c9d4SSatish Balay PC_MG_MULTIPLICATIVE, 3519371c9d4SSatish Balay PC_MG_ADDITIVE, 3529371c9d4SSatish Balay PC_MG_FULL, 3539371c9d4SSatish Balay PC_MG_KASKADE 3549371c9d4SSatish Balay } PCMGType; 355b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE; 356b0753f9dSMatthew G. Knepley 357b0753f9dSMatthew G. Knepley /*E 3580b4b7b1cSBarry Smith PCMGCycleType - Determines which of V-cycle or W-cycle to use with the `PCType` of `PCMG` or `PCGAMG` 359b0753f9dSMatthew G. Knepley 360b0753f9dSMatthew G. Knepley Values: 361af27ebaaSBarry Smith + `PC_MG_V_CYCLE` - use the V cycle 362af27ebaaSBarry Smith - `PC_MG_W_CYCLE` - use the W cycle 363b0753f9dSMatthew G. Knepley 36416a05f60SBarry Smith Level: beginner 36516a05f60SBarry Smith 36609e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()` 367b0753f9dSMatthew G. Knepley E*/ 3689371c9d4SSatish Balay typedef enum { 3699371c9d4SSatish Balay PC_MG_CYCLE_V = 1, 3709371c9d4SSatish Balay PC_MG_CYCLE_W = 2 3719371c9d4SSatish Balay } PCMGCycleType; 372b0753f9dSMatthew G. Knepley 373b0753f9dSMatthew G. Knepley /*E 3740b4b7b1cSBarry Smith PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process with the `PCType` of `PCMG` 3752134b1e4SBarry Smith 3762134b1e4SBarry Smith Values: 3770b4b7b1cSBarry Smith + `PC_MG_GALERKIN_PMAT` - computes the `pmat` (matrix from which the preconditioner is built) via the Galerkin process from the finest grid 3780b4b7b1cSBarry Smith . `PC_MG_GALERKIN_MAT` - computes the `mat` (matrix used to apply the operator) via the Galerkin process from the finest grid 3790b4b7b1cSBarry Smith . `PC_MG_GALERKIN_BOTH` - computes both the `mat` and `pmat` via the Galerkin process (if pmat == mat the construction is only done once 38087497f52SBarry Smith - `PC_MG_GALERKIN_NONE` - neither operator is computed via the Galerkin process, the user must provide the operator 3812134b1e4SBarry Smith 38209e53591SBarry Smith Level: beginner 38309e53591SBarry Smith 38409e53591SBarry Smith Note: 38587497f52SBarry Smith Users should never set `PC_MG_GALERKIN_EXTERNAL`, it is used by `PCHYPRE` and `PCML` 3862134b1e4SBarry Smith 38709e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()` 3882134b1e4SBarry Smith E*/ 3899371c9d4SSatish Balay typedef enum { 3909371c9d4SSatish Balay PC_MG_GALERKIN_BOTH, 3919371c9d4SSatish Balay PC_MG_GALERKIN_PMAT, 3929371c9d4SSatish Balay PC_MG_GALERKIN_MAT, 3939371c9d4SSatish Balay PC_MG_GALERKIN_NONE, 3949371c9d4SSatish Balay PC_MG_GALERKIN_EXTERNAL 3959371c9d4SSatish Balay } PCMGGalerkinType; 3962134b1e4SBarry Smith 3972134b1e4SBarry Smith /*E 3980b4b7b1cSBarry Smith PCExoticType - Determines which of the face-based or wirebasket-based coarse grid space to use with the `PCType` of `PCEXOTIC` 399b0753f9dSMatthew G. Knepley 400b0753f9dSMatthew G. Knepley Level: beginner 401b0753f9dSMatthew G. Knepley 40209e53591SBarry Smith .seealso: [](sec_pc), `PCExoticSetType()`, `PCEXOTIC` 403b0753f9dSMatthew G. Knepley E*/ 4049371c9d4SSatish Balay typedef enum { 4059371c9d4SSatish Balay PC_EXOTIC_FACE, 4069371c9d4SSatish Balay PC_EXOTIC_WIREBASKET 4079371c9d4SSatish Balay } PCExoticType; 408b0753f9dSMatthew G. Knepley 4098c1cd74cSHong Zhang /*E 4100b4b7b1cSBarry Smith PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains with the `PCType` of `PCBDDC` 411bc960bbfSJed Brown 412bc960bbfSJed Brown Values: 41387497f52SBarry Smith + `PC_BDDC_INTERFACE_EXT_DIRICHLET` - solves Dirichlet interior problem; this is the standard BDDC algorithm 4140b4b7b1cSBarry Smith - `PC_BDDC_INTERFACE_EXT_LUMP` - skips interior solve; sometimes called $M_1$ and associated with "lumped FETI-DP" 415bc960bbfSJed Brown 41616a05f60SBarry Smith Level: intermediate 41716a05f60SBarry Smith 41809e53591SBarry Smith .seealso: [](sec_pc), `PCBDDC`, `PC` 419bc960bbfSJed Brown E*/ 420bc960bbfSJed Brown typedef enum { 421bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_DIRICHLET, 422bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_LUMP 423bc960bbfSJed Brown } PCBDDCInterfaceExtType; 424bc960bbfSJed Brown 425bc960bbfSJed Brown /*E 426f3b08a26SMatthew G. Knepley PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation 427f3b08a26SMatthew G. Knepley 428f3b08a26SMatthew G. Knepley Level: beginner 429f3b08a26SMatthew G. Knepley 43009e53591SBarry Smith .seealso: [](sec_pc), `PCMGSetAdaptCoarseSpaceType()`, `PCMG`, `PC` 431f3b08a26SMatthew G. Knepley E*/ 4329371c9d4SSatish Balay typedef enum { 4339371c9d4SSatish Balay PCMG_ADAPT_NONE, 4349371c9d4SSatish Balay PCMG_ADAPT_POLYNOMIAL, 4359371c9d4SSatish Balay PCMG_ADAPT_HARMONIC, 4369371c9d4SSatish Balay PCMG_ADAPT_EIGENVECTOR, 4379371c9d4SSatish Balay PCMG_ADAPT_GENERALIZED_EIGENVECTOR, 4389371c9d4SSatish Balay PCMG_ADAPT_GDSW 4399371c9d4SSatish Balay } PCMGCoarseSpaceType; 440f3b08a26SMatthew G. Knepley 441f3b08a26SMatthew G. Knepley /*E 4420b4b7b1cSBarry Smith PCPatchConstructType - Determines the algorithm used to construct patches for the `PCPATCH` preconditioner 4434bbf5ea8SMatthew G. Knepley 4444bbf5ea8SMatthew G. Knepley Level: beginner 4454bbf5ea8SMatthew G. Knepley 44609e53591SBarry Smith .seealso: [](sec_pc), `PCPatchSetConstructType()`, `PCPATCH`, `PC` 4474bbf5ea8SMatthew G. Knepley E*/ 4489371c9d4SSatish Balay typedef enum { 4499371c9d4SSatish Balay PC_PATCH_STAR, 4509371c9d4SSatish Balay PC_PATCH_VANKA, 4519371c9d4SSatish Balay PC_PATCH_PARDECOMP, 4529371c9d4SSatish Balay PC_PATCH_USER, 4539371c9d4SSatish Balay PC_PATCH_PYTHON 4549371c9d4SSatish Balay } PCPatchConstructType; 4554bbf5ea8SMatthew G. Knepley 4564bbf5ea8SMatthew G. Knepley /*E 4570b4b7b1cSBarry Smith PCDeflationSpaceType - Type of deflation used by `PCType` `PCDEFLATION` 458e662bc50SJakub Kruzik 459e662bc50SJakub Kruzik Values: 46087497f52SBarry Smith + `PC_DEFLATION_SPACE_HAAR` - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off 46109e53591SBarry Smith . `PC_DEFLATION_SPACE_DB2` - `MATCOMPOSITE` of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet) 46287497f52SBarry Smith . `PC_DEFLATION_SPACE_DB4` - same as above, but with db4 (4 coefficient Daubechies) 46387497f52SBarry Smith . `PC_DEFLATION_SPACE_DB8` - same as above, but with db8 (8 coefficient Daubechies) 46487497f52SBarry Smith . `PC_DEFLATION_SPACE_DB16` - same as above, but with db16 (16 coefficient Daubechies) 46587497f52SBarry Smith . `PC_DEFLATION_SPACE_BIORTH22` - same as above, but with biorthogonal 2.2 (6 coefficients) 46687497f52SBarry Smith . `PC_DEFLATION_SPACE_MEYER` - same as above, but with Meyer/FIR (62 coefficients) 467d5b43468SJose E. Roman . `PC_DEFLATION_SPACE_AGGREGATION` - aggregates local indices (given by operator matrix distribution) into a subdomain 46887497f52SBarry Smith - `PC_DEFLATION_SPACE_USER` - indicates space set by user 469e662bc50SJakub Kruzik 4701fdb00f9SJakub Kruzik Level: intermediate 4711fdb00f9SJakub Kruzik 47209e53591SBarry Smith Note: 47309e53591SBarry Smith Wavelet-based space (except Haar) can be used in multilevel deflation. 47409e53591SBarry Smith 47509e53591SBarry Smith .seealso: [](sec_pc), `PCDeflationSetSpaceToCompute()`, `PCDEFLATION`, `PC` 476e662bc50SJakub Kruzik E*/ 477e662bc50SJakub Kruzik typedef enum { 478e662bc50SJakub Kruzik PC_DEFLATION_SPACE_HAAR, 479e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB2, 480e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB4, 481e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB8, 482e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB16, 483e662bc50SJakub Kruzik PC_DEFLATION_SPACE_BIORTH22, 484e662bc50SJakub Kruzik PC_DEFLATION_SPACE_MEYER, 485e662bc50SJakub Kruzik PC_DEFLATION_SPACE_AGGREGATION, 486e662bc50SJakub Kruzik PC_DEFLATION_SPACE_USER 487e662bc50SJakub Kruzik } PCDeflationSpaceType; 488e662bc50SJakub Kruzik 489e662bc50SJakub Kruzik /*E 4900b4b7b1cSBarry Smith PCHPDDMCoarseCorrectionType - Type of coarse correction used by `PCType` `PCHPDDM` 49138cfc46eSPierre Jolivet 49238cfc46eSPierre Jolivet Values: 49309e53591SBarry Smith + `PC_HPDDM_COARSE_CORRECTION_DEFLATED` (default) - eq. (1) in `PCHPDDMShellApply()` 49487497f52SBarry Smith . `PC_HPDDM_COARSE_CORRECTION_ADDITIVE` - eq. (2) 495aa1539e9SPierre Jolivet . `PC_HPDDM_COARSE_CORRECTION_BALANCED` - eq. (3) 496aa1539e9SPierre Jolivet - `PC_HPDDM_COARSE_CORRECTION_NONE` - no coarse correction (mostly useful for debugging) 49738cfc46eSPierre Jolivet 49816a05f60SBarry Smith Level: intermediate 49916a05f60SBarry Smith 50009e53591SBarry Smith .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCSetType()`, `PCHPDDMShellApply()` 50138cfc46eSPierre Jolivet E*/ 5029371c9d4SSatish Balay typedef enum { 5039371c9d4SSatish Balay PC_HPDDM_COARSE_CORRECTION_DEFLATED, 5049371c9d4SSatish Balay PC_HPDDM_COARSE_CORRECTION_ADDITIVE, 505aa1539e9SPierre Jolivet PC_HPDDM_COARSE_CORRECTION_BALANCED, 506aa1539e9SPierre Jolivet PC_HPDDM_COARSE_CORRECTION_NONE 5079371c9d4SSatish Balay } PCHPDDMCoarseCorrectionType; 50838cfc46eSPierre Jolivet 50938cfc46eSPierre Jolivet /*E 51013044ca3SPierre Jolivet PCHPDDMSchurPreType - Type of `PCHPDDM` preconditioner for a `MATSCHURCOMPLEMENT` generated by `PCFIELDSPLIT` with `PCFieldSplitSchurPreType` set to `PC_FIELDSPLIT_SCHUR_PRE_SELF` 51113044ca3SPierre Jolivet 51213044ca3SPierre Jolivet Values: 5130b4b7b1cSBarry Smith + `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 5140b4b7b1cSBarry Smith is built by approximating the Schur complement with (inv(sqrt(diag(A00))) A01)^T (inv(sqrt(diag(A00))) A01) 5150b4b7b1cSBarry Smith and by considering the associated linear least squares problem 5160b4b7b1cSBarry Smith - `PC_HPDDM_SCHUR_PRE_GENEO` - only with A10 = A01^T, `PCHPDDMSetAuxiliaryMat()` called on the `PC` of the A00 block, and if A11 is nonzero, 5170b4b7b1cSBarry Smith then `PCHPDDMSetAuxiliaryMat()` must be called on the associated `PC` as well (it is built automatically for the 5180b4b7b1cSBarry Smith user otherwise); the Schur complement `PC` is set internally to `PCKSP`, with the prefix `-fieldsplit_1_pc_hpddm_`; 5190b4b7b1cSBarry Smith the operator associated to the `PC` is spectrally equivalent to the original Schur complement 52013044ca3SPierre Jolivet 52113044ca3SPierre Jolivet Level: advanced 52213044ca3SPierre Jolivet 52313044ca3SPierre Jolivet .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCFIELDSPLIT`, `PC_FIELDSPLIT_SCHUR_PRE_SELF`, `PCFieldSplitSetSchurPre()`, `PCHPDDMSetAuxiliaryMat()` 52413044ca3SPierre Jolivet E*/ 52513044ca3SPierre Jolivet typedef enum { 52613044ca3SPierre Jolivet PC_HPDDM_SCHUR_PRE_LEAST_SQUARES, 527*ce78bad3SBarry Smith PC_HPDDM_SCHUR_PRE_GENEO 52813044ca3SPierre Jolivet } PCHPDDMSchurPreType; 52913044ca3SPierre Jolivet 53013044ca3SPierre Jolivet /*E 5310b4b7b1cSBarry Smith PCFailedReason - indicates the type of `PC` failure. That is why the construction of the preconditioner, `PCSetUp()`, or its use, `PCApply()`, failed 5328c1cd74cSHong Zhang 5338c1cd74cSHong Zhang Level: beginner 5348c1cd74cSHong Zhang 5350b4b7b1cSBarry Smith .seealso: [](sec_pc), `PC`, `PCGetFailedReason()`, `PCSetUp()` 5368c1cd74cSHong Zhang E*/ 5379371c9d4SSatish Balay typedef enum { 5389371c9d4SSatish Balay PC_SETUP_ERROR = -1, 539*ce78bad3SBarry Smith PC_NOERROR = 0, 540*ce78bad3SBarry Smith PC_FACTOR_STRUCT_ZEROPIVOT = 1, 541*ce78bad3SBarry Smith PC_FACTOR_NUMERIC_ZEROPIVOT = 2, 542*ce78bad3SBarry Smith PC_FACTOR_OUTMEMORY = 3, 543*ce78bad3SBarry Smith PC_FACTOR_OTHER = 4, 544*ce78bad3SBarry Smith PC_INCONSISTENT_RHS = 5, 545*ce78bad3SBarry Smith PC_SUBPC_ERROR = 6 5469371c9d4SSatish Balay } PCFailedReason; 547ce7c7f2fSMark Adams 548ce7c7f2fSMark Adams /*E 5490b4b7b1cSBarry Smith PCGAMGLayoutType - Layout for reduced grids for `PCType` `PCGAMG` 550ce7c7f2fSMark Adams 551ce7c7f2fSMark Adams Level: intermediate 552ce7c7f2fSMark Adams 55309e53591SBarry Smith .seealso: [](sec_pc), `PCGAMG`, `PC`, `PCGAMGSetCoarseGridLayoutType()` 554ce7c7f2fSMark Adams E*/ 5559371c9d4SSatish Balay typedef enum { 5569371c9d4SSatish Balay PCGAMG_LAYOUT_COMPACT, 5579371c9d4SSatish Balay PCGAMG_LAYOUT_SPREAD 5589371c9d4SSatish Balay } PCGAMGLayoutType; 559