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 PCBFBT "bfbt" 61b0753f9dSMatthew G. Knepley #define PCLSC "lsc" 62b0753f9dSMatthew G. Knepley #define PCPYTHON "python" 63b0753f9dSMatthew G. Knepley #define PCPFMG "pfmg" 641c188c59Sftrigaux #define PCSMG "smg" 65b0753f9dSMatthew G. Knepley #define PCSYSPFMG "syspfmg" 66b0753f9dSMatthew G. Knepley #define PCREDISTRIBUTE "redistribute" 67b0753f9dSMatthew G. Knepley #define PCSVD "svd" 68b0753f9dSMatthew G. Knepley #define PCGAMG "gamg" 694b3f184cSKarl Rupp #define PCCHOWILUVIENNACL "chowiluviennacl" 7070baa948SKarl Rupp #define PCROWSCALINGVIENNACL "rowscalingviennacl" 7107598726SKarl Rupp #define PCSAVIENNACL "saviennacl" 72b0753f9dSMatthew G. Knepley #define PCBDDC "bddc" 73b0753f9dSMatthew G. Knepley #define PCKACZMARZ "kaczmarz" 7468ddcbeaSBarry Smith #define PCTELESCOPE "telescope" 754bbf5ea8SMatthew G. Knepley #define PCPATCH "patch" 76b9ac7092SAlp Dener #define PCLMVM "lmvm" 77360ee056SFande Kong #define PCHMG "hmg" 7837eeb815SJakub Kruzik #define PCDEFLATION "deflation" 7938cfc46eSPierre Jolivet #define PCHPDDM "hpddm" 8053022affSStefano Zampini #define PCH2OPUS "h2opus" 81f1f2ae84SBarry Smith #define PCMPI "mpi" 82b0753f9dSMatthew G. Knepley 83b0753f9dSMatthew G. Knepley /*E 840b4b7b1cSBarry Smith PCSide - Determines if the preconditioner is to be applied to the left, right 850b4b7b1cSBarry Smith or symmetrically around the operator in `KSPSolve()`. 86b0753f9dSMatthew G. Knepley 8716a05f60SBarry Smith Values: 8816a05f60SBarry Smith + `PC_LEFT` - applied after the operator is applied 8916a05f60SBarry Smith . `PC_RIGHT` - applied before the operator is applied 9016a05f60SBarry 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. 9116a05f60SBarry Smith 92b0753f9dSMatthew G. Knepley Level: beginner 93b0753f9dSMatthew G. Knepley 9416a05f60SBarry Smith Note: 9516a05f60SBarry Smith Certain `KSPType` support only a subset of `PCSide` values 9616a05f60SBarry Smith 970b4b7b1cSBarry Smith .seealso: [](sec_pc), `PC`, `KSPSetPCSide()`, `KSP`, `KSPType`, `KSPGetPCSide()`, `KSPSolve()` 98b0753f9dSMatthew G. Knepley E*/ 999371c9d4SSatish Balay typedef enum { 1009371c9d4SSatish Balay PC_SIDE_DEFAULT = -1, 101*ce78bad3SBarry Smith PC_LEFT = 0, 102*ce78bad3SBarry Smith PC_RIGHT = 1, 103*ce78bad3SBarry Smith PC_SYMMETRIC = 2 1049371c9d4SSatish Balay } PCSide; 105b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 106b0753f9dSMatthew G. Knepley 107b0753f9dSMatthew G. Knepley /*E 10839b1ba1bSPierre Jolivet PCRichardsonConvergedReason - reason a `PCApplyRichardson()` method terminated 109b0753f9dSMatthew G. Knepley 110b0753f9dSMatthew G. Knepley Level: advanced 111b0753f9dSMatthew G. Knepley 11239b1ba1bSPierre Jolivet .seealso: [](sec_pc), `KSPRICHARDSON`, `PC`, `PCApplyRichardson()` 113b0753f9dSMatthew G. Knepley E*/ 114b0753f9dSMatthew G. Knepley typedef enum { 115dd460d27SBarry Smith PCRICHARDSON_NOT_SET = 0, 116b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_RTOL = 2, 117b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ATOL = 3, 118b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ITS = 4, 1199371c9d4SSatish Balay PCRICHARDSON_DIVERGED_DTOL = -4 1209371c9d4SSatish Balay } PCRichardsonConvergedReason; 121b0753f9dSMatthew G. Knepley 122b0753f9dSMatthew G. Knepley /*E 1230b4b7b1cSBarry Smith PCJacobiType - Determines what elements of the matrix are used to form the Jacobi preconditioner, that is with the `PCType` of `PCJACOBI` 124b0753f9dSMatthew G. Knepley 12509e53591SBarry Smith Values: 12609e53591SBarry Smith + `PC_JACOBI_DIAGONAL` - use the diagonal entry, if it is zero use one 127eede4a3fSMark Adams . `PC_JACOBI_ROWL1` - add sum of absolute values in row i, j != i, to diag_ii 12809e53591SBarry Smith . `PC_JACOBI_ROWMAX` - use the maximum absolute value in the row 12909e53591SBarry Smith - `PC_JACOBI_ROWSUM` - use the sum of the values in the row (not the absolute values) 13009e53591SBarry Smith 131b0753f9dSMatthew G. Knepley Level: intermediate 132b0753f9dSMatthew G. Knepley 13309e53591SBarry Smith .seealso: [](sec_pc), `PCJACOBI`, `PC` 134b0753f9dSMatthew G. Knepley E*/ 1359371c9d4SSatish Balay typedef enum { 1369371c9d4SSatish Balay PC_JACOBI_DIAGONAL, 137eede4a3fSMark Adams PC_JACOBI_ROWL1, 1389371c9d4SSatish Balay PC_JACOBI_ROWMAX, 1399371c9d4SSatish Balay PC_JACOBI_ROWSUM 1409371c9d4SSatish Balay } PCJacobiType; 141b0753f9dSMatthew G. Knepley 142b0753f9dSMatthew G. Knepley /*E 1430b4b7b1cSBarry Smith PCASMType - Determines the type of additive Schwarz method, `PCASM`, to use 144b0753f9dSMatthew G. Knepley 14509e53591SBarry Smith Values: 14609e53591SBarry Smith + `PC_ASM_BASIC` - Symmetric version where residuals from the ghost points are used 14709e53591SBarry Smith and computed values in ghost regions are added together. 148af27ebaaSBarry Smith Classical standard additive Schwarz as introduced in {cite}`dryja1987additive`. 14909e53591SBarry Smith . `PC_ASM_RESTRICT` - Residuals from ghost points are used but computed values in ghost 150af27ebaaSBarry Smith region are discarded {cite}`cs99`. Default. 15109e53591SBarry Smith . `PC_ASM_INTERPOLATE` - Residuals from ghost points are not used, computed values in ghost 15209e53591SBarry Smith region are added back in. 15309e53591SBarry Smith - `PC_ASM_NONE` - Residuals from ghost points are not used, computed ghost values are 154af27ebaaSBarry Smith discarded. Not very good. 155b0753f9dSMatthew G. Knepley 156b0753f9dSMatthew G. Knepley Level: beginner 157b0753f9dSMatthew G. Knepley 15816a05f60SBarry Smith .seealso: [](sec_pc), `PC`, `PCASM`, `PCASMSetType()`, `PCGASMType` 159b0753f9dSMatthew G. Knepley E*/ 1609371c9d4SSatish Balay typedef enum { 1619371c9d4SSatish Balay PC_ASM_BASIC = 3, 1629371c9d4SSatish Balay PC_ASM_RESTRICT = 1, 1639371c9d4SSatish Balay PC_ASM_INTERPOLATE = 2, 1649371c9d4SSatish Balay PC_ASM_NONE = 0 1659371c9d4SSatish Balay } PCASMType; 166b0753f9dSMatthew G. Knepley 167b0753f9dSMatthew G. Knepley /*E 1680b4b7b1cSBarry 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` 169b0753f9dSMatthew G. Knepley 17009e53591SBarry Smith Values: 17109e53591SBarry Smith + `PC_GASM_BASIC` - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 17209e53591SBarry Smith over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 173af27ebaaSBarry Smith from neighboring subdomains. Classical standard additive Schwarz {cite}`dryja1987additive`. 17409e53591SBarry Smith . `PC_GASM_RESTRICT` - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 17509e53591SBarry Smith (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 17609e53591SBarry Smith each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 177af27ebaaSBarry Smith assumption) {cite}`cs99`. Default. 17809e53591SBarry Smith . `PC_GASM_INTERPOLATE` - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 17909e53591SBarry Smith applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 18009e53591SBarry Smith from neighboring subdomains. 181af27ebaaSBarry Smith - `PC_GASM_NONE` - Residuals and corrections are zeroed out outside the local subdomains. Not very good. 182b0753f9dSMatthew G. Knepley 183b0753f9dSMatthew G. Knepley Level: beginner 184b0753f9dSMatthew G. Knepley 18516a05f60SBarry Smith Note: 18616a05f60SBarry Smith Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 1870b4b7b1cSBarry Smith domain, while the outer subdomains contain the inner subdomains and overlap with each other. The `PCGASM` preconditioner will compute 18816a05f60SBarry Smith a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 18916a05f60SBarry Smith (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 19016a05f60SBarry Smith 191af27ebaaSBarry Smith Developer Note: 192af27ebaaSBarry Smith Perhaps better to remove this since it matches `PCASMType` 193af27ebaaSBarry Smith 19416a05f60SBarry Smith .seealso: [](sec_pc), `PCGASM`, `PCASM`, `PC`, `PCGASMSetType()`, `PCASMType` 195b0753f9dSMatthew G. Knepley E*/ 1969371c9d4SSatish Balay typedef enum { 1979371c9d4SSatish Balay PC_GASM_BASIC = 3, 1989371c9d4SSatish Balay PC_GASM_RESTRICT = 1, 1999371c9d4SSatish Balay PC_GASM_INTERPOLATE = 2, 2009371c9d4SSatish Balay PC_GASM_NONE = 0 2019371c9d4SSatish Balay } PCGASMType; 202b0753f9dSMatthew G. Knepley 203b0753f9dSMatthew G. Knepley /*E 20416a05f60SBarry Smith PCCompositeType - Determines how two or more preconditioner are composed with the `PCType` of `PCCOMPOSITE` 205b0753f9dSMatthew G. Knepley 20609e53591SBarry Smith Values: 20709e53591SBarry Smith + `PC_COMPOSITE_ADDITIVE` - results from application of all preconditioners are added together 20809e53591SBarry Smith . `PC_COMPOSITE_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly 20909e53591SBarry Smith computed after the previous preconditioner application 21009e53591SBarry Smith . `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly 21109e53591SBarry Smith computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners) 212af27ebaaSBarry Smith . `PC_COMPOSITE_SPECIAL` - This is very special for a matrix of the form $ \alpha I + R + S$ 213af27ebaaSBarry Smith where the first preconditioner is built from $\alpha I + S$ and second from $\alpha I + R$ 21409e53591SBarry Smith . `PC_COMPOSITE_SCHUR` - composes the Schur complement of the matrix from two blocks, see `PCFIELDSPLIT` 21509e53591SBarry Smith - `PC_COMPOSITE_GKB` - the generalized Golub-Kahan bidiagonalization preconditioner, see `PCFIELDSPLIT` 216b0753f9dSMatthew G. Knepley 217b0753f9dSMatthew G. Knepley Level: beginner 218b0753f9dSMatthew G. Knepley 219af27ebaaSBarry Smith .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()`, `SNESCompositeType`, `PCCompositeSpecialSetAlpha()` 220b0753f9dSMatthew G. Knepley E*/ 2219371c9d4SSatish Balay typedef enum { 2229371c9d4SSatish Balay PC_COMPOSITE_ADDITIVE, 2239371c9d4SSatish Balay PC_COMPOSITE_MULTIPLICATIVE, 2249371c9d4SSatish Balay PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, 2259371c9d4SSatish Balay PC_COMPOSITE_SPECIAL, 2269371c9d4SSatish Balay PC_COMPOSITE_SCHUR, 2279371c9d4SSatish Balay PC_COMPOSITE_GKB 2289371c9d4SSatish Balay } PCCompositeType; 229b0753f9dSMatthew G. Knepley 230b0753f9dSMatthew G. Knepley /*E 2310b4b7b1cSBarry Smith PCFieldSplitSchurPreType - Determines how to precondition a Schur complement arising with the `PCType` of `PCFIELDSPLIT` 232b0753f9dSMatthew G. Knepley 23316a05f60SBarry Smith Values: 23413044ca3SPierre Jolivet + `PC_FIELDSPLIT_SCHUR_PRE_SELF` - the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix. 23513044ca3SPierre Jolivet The only preconditioners that currently work with this symbolic representation matrix object are `PCLSC` and `PCHPDDM` 236a077d33dSBarry 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$. 237a077d33dSBarry Smith This is only a good preconditioner when $diag(A00)$ is a good preconditioner for $A00$. Optionally, $A00$ can be 23816a05f60SBarry Smith lumped before extracting the diagonal using the additional option `-fieldsplit_1_mat_schur_complement_ainv_type lump` 239a077d33dSBarry Smith . `PC_FIELDSPLIT_SCHUR_PRE_A11` - the preconditioner for the Schur complement is generated from $A11$, not the Schur complement matrix 24016a05f60SBarry Smith . `PC_FIELDSPLIT_SCHUR_PRE_USER` - the preconditioner for the Schur complement is generated from the user provided matrix (pre argument 24116a05f60SBarry Smith to this function). 24216a05f60SBarry Smith - `PC_FIELDSPLIT_SCHUR_PRE_FULL` - the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation 24316a05f60SBarry Smith computed internally by `PCFIELDSPLIT` (this is expensive) useful mostly as a test that the Schur complement approach can work for your problem 24416a05f60SBarry Smith 245b0753f9dSMatthew G. Knepley Level: intermediate 246b0753f9dSMatthew G. Knepley 24709e53591SBarry Smith .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurPre()`, `PC` 248b0753f9dSMatthew G. Knepley E*/ 2499371c9d4SSatish Balay typedef enum { 2509371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_SELF, 2519371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_SELFP, 2529371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_A11, 2539371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_USER, 2549371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_FULL 2559371c9d4SSatish Balay } PCFieldSplitSchurPreType; 256b0753f9dSMatthew G. Knepley 257b0753f9dSMatthew G. Knepley /*E 2580b4b7b1cSBarry Smith PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use with the `PCType` of `PCFIELDSPLIT` 259b0753f9dSMatthew G. Knepley 26016a05f60SBarry Smith Values: 26116a05f60SBarry Smith + `PC_FIELDSPLIT_SCHUR_FACT_DIAG` - the preconditioner is solving `D` 26216a05f60SBarry Smith . `PC_FIELDSPLIT_SCHUR_FACT_LOWER` - the preconditioner is solving `L D` 26316a05f60SBarry Smith . `PC_FIELDSPLIT_SCHUR_FACT_UPPER` - the preconditioner is solving `D U` 26416a05f60SBarry Smith - `PC_FIELDSPLIT_SCHUR_FACT_FULL` - the preconditioner is solving `L(D U)` 26516a05f60SBarry Smith 26616a05f60SBarry Smith where the matrix is factorized as 26716a05f60SBarry Smith .vb 26816a05f60SBarry Smith (A B) = (1 0) (A 0) (1 Ainv*B) = L D U 26916a05f60SBarry Smith (C E) (C*Ainv 1) (0 S) (0 1) 27016a05f60SBarry Smith .ve 27116a05f60SBarry Smith 272b0753f9dSMatthew G. Knepley Level: intermediate 273b0753f9dSMatthew G. Knepley 27409e53591SBarry Smith .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurFactType()`, `PC` 275b0753f9dSMatthew G. Knepley E*/ 276b0753f9dSMatthew G. Knepley typedef enum { 277b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_DIAG, 278b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_LOWER, 279b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_UPPER, 280b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_FULL 281b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType; 282b0753f9dSMatthew G. Knepley 283b0753f9dSMatthew G. Knepley /*E 28487497f52SBarry Smith PCPARMSGlobalType - Determines the global preconditioner method in `PCPARMS` 285b0753f9dSMatthew G. Knepley 286b0753f9dSMatthew G. Knepley Level: intermediate 287b0753f9dSMatthew G. Knepley 28809e53591SBarry Smith .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetGlobal()`, `PC` 289b0753f9dSMatthew G. Knepley E*/ 2909371c9d4SSatish Balay typedef enum { 2919371c9d4SSatish Balay PC_PARMS_GLOBAL_RAS, 2929371c9d4SSatish Balay PC_PARMS_GLOBAL_SCHUR, 2939371c9d4SSatish Balay PC_PARMS_GLOBAL_BJ 2949371c9d4SSatish Balay } PCPARMSGlobalType; 2959d8081ecSMatthew G. Knepley 296b0753f9dSMatthew G. Knepley /*E 29787497f52SBarry Smith PCPARMSLocalType - Determines the local preconditioner method in `PCPARMS` 298b0753f9dSMatthew G. Knepley 299b0753f9dSMatthew G. Knepley Level: intermediate 300b0753f9dSMatthew G. Knepley 30109e53591SBarry Smith .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetLocal()`, `PC` 302b0753f9dSMatthew G. Knepley E*/ 3039371c9d4SSatish Balay typedef enum { 3049371c9d4SSatish Balay PC_PARMS_LOCAL_ILU0, 3059371c9d4SSatish Balay PC_PARMS_LOCAL_ILUK, 3069371c9d4SSatish Balay PC_PARMS_LOCAL_ILUT, 3079371c9d4SSatish Balay PC_PARMS_LOCAL_ARMS 3089371c9d4SSatish Balay } PCPARMSLocalType; 309b0753f9dSMatthew G. Knepley 310f0fc11ceSJed Brown /*J 31187497f52SBarry Smith PCGAMGType - type of generalized algebraic multigrid `PCGAMG` method 312b0753f9dSMatthew G. Knepley 31309e53591SBarry Smith Values: 31409e53591SBarry Smith + `PCGAMGAGG` - (the default) smoothed aggregation algorithm, robust, very well tested 315baca6076SPierre Jolivet . `PCGAMGGEO` - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not supported, reference implementation (2D) 316baca6076SPierre Jolivet - `PCGAMGCLASSICAL` - classical algebraic multigrid preconditioner, incomplete, not supported, reference implementation 31709e53591SBarry Smith 318b0753f9dSMatthew G. Knepley Level: intermediate 319b0753f9dSMatthew G. Knepley 32009e53591SBarry Smith .seealso: [](sec_pc), `PCGAMG`, `PCMG`, `PC`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()` 321f0fc11ceSJed Brown J*/ 322b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType; 323b0753f9dSMatthew G. Knepley #define PCGAMGAGG "agg" 324b0753f9dSMatthew G. Knepley #define PCGAMGGEO "geo" 325b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL "classical" 326b0753f9dSMatthew G. Knepley 327b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType; 328b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT "direct" 329b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard" 330b0753f9dSMatthew G. Knepley 331b0753f9dSMatthew G. Knepley /*E 3320b4b7b1cSBarry Smith PCMGType - Determines the type of multigrid method that is run with the `PCType` of `PCMG` 333b0753f9dSMatthew G. Knepley 334b0753f9dSMatthew G. Knepley Values: 33587497f52SBarry Smith + `PC_MG_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `PCMGSetCycleType()` 33687497f52SBarry Smith . `PC_MG_ADDITIVE` - the additive multigrid preconditioner where all levels are 337b0753f9dSMatthew G. Knepley smoothed before updating the residual. This only uses the 338b0753f9dSMatthew G. Knepley down smoother, in the preconditioner the upper smoother is ignored 33987497f52SBarry Smith . `PC_MG_FULL` - same as multiplicative except one also performs grid sequencing, 340b0753f9dSMatthew G. Knepley that is starts on the coarsest grid, performs a cycle, interpolates 341b0753f9dSMatthew 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 342b0753f9dSMatthew G. Knepley algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 343b0753f9dSMatthew G. Knepley calls the V-cycle only on the coarser level and has a post-smoother instead. 3440b4b7b1cSBarry Smith - `PC_MG_KASKADE` - Cascadic or Kaskadic multigrid, like full multigrid except one never goes back to a coarser level from a finer 34516a05f60SBarry Smith 34616a05f60SBarry Smith Level: beginner 347b0753f9dSMatthew G. Knepley 34809e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()` 349b0753f9dSMatthew G. Knepley E*/ 3509371c9d4SSatish Balay typedef enum { 3519371c9d4SSatish Balay PC_MG_MULTIPLICATIVE, 3529371c9d4SSatish Balay PC_MG_ADDITIVE, 3539371c9d4SSatish Balay PC_MG_FULL, 3549371c9d4SSatish Balay PC_MG_KASKADE 3559371c9d4SSatish Balay } PCMGType; 356b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE; 357b0753f9dSMatthew G. Knepley 358b0753f9dSMatthew G. Knepley /*E 3590b4b7b1cSBarry Smith PCMGCycleType - Determines which of V-cycle or W-cycle to use with the `PCType` of `PCMG` or `PCGAMG` 360b0753f9dSMatthew G. Knepley 361b0753f9dSMatthew G. Knepley Values: 362af27ebaaSBarry Smith + `PC_MG_V_CYCLE` - use the V cycle 363af27ebaaSBarry Smith - `PC_MG_W_CYCLE` - use the W cycle 364b0753f9dSMatthew G. Knepley 36516a05f60SBarry Smith Level: beginner 36616a05f60SBarry Smith 36709e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()` 368b0753f9dSMatthew G. Knepley E*/ 3699371c9d4SSatish Balay typedef enum { 3709371c9d4SSatish Balay PC_MG_CYCLE_V = 1, 3719371c9d4SSatish Balay PC_MG_CYCLE_W = 2 3729371c9d4SSatish Balay } PCMGCycleType; 373b0753f9dSMatthew G. Knepley 374b0753f9dSMatthew G. Knepley /*E 3750b4b7b1cSBarry Smith PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process with the `PCType` of `PCMG` 3762134b1e4SBarry Smith 3772134b1e4SBarry Smith Values: 3780b4b7b1cSBarry Smith + `PC_MG_GALERKIN_PMAT` - computes the `pmat` (matrix from which the preconditioner is built) via the Galerkin process from the finest grid 3790b4b7b1cSBarry Smith . `PC_MG_GALERKIN_MAT` - computes the `mat` (matrix used to apply the operator) via the Galerkin process from the finest grid 3800b4b7b1cSBarry Smith . `PC_MG_GALERKIN_BOTH` - computes both the `mat` and `pmat` via the Galerkin process (if pmat == mat the construction is only done once 38187497f52SBarry Smith - `PC_MG_GALERKIN_NONE` - neither operator is computed via the Galerkin process, the user must provide the operator 3822134b1e4SBarry Smith 38309e53591SBarry Smith Level: beginner 38409e53591SBarry Smith 38509e53591SBarry Smith Note: 38687497f52SBarry Smith Users should never set `PC_MG_GALERKIN_EXTERNAL`, it is used by `PCHYPRE` and `PCML` 3872134b1e4SBarry Smith 38809e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()` 3892134b1e4SBarry Smith E*/ 3909371c9d4SSatish Balay typedef enum { 3919371c9d4SSatish Balay PC_MG_GALERKIN_BOTH, 3929371c9d4SSatish Balay PC_MG_GALERKIN_PMAT, 3939371c9d4SSatish Balay PC_MG_GALERKIN_MAT, 3949371c9d4SSatish Balay PC_MG_GALERKIN_NONE, 3959371c9d4SSatish Balay PC_MG_GALERKIN_EXTERNAL 3969371c9d4SSatish Balay } PCMGGalerkinType; 3972134b1e4SBarry Smith 3982134b1e4SBarry Smith /*E 3990b4b7b1cSBarry Smith PCExoticType - Determines which of the face-based or wirebasket-based coarse grid space to use with the `PCType` of `PCEXOTIC` 400b0753f9dSMatthew G. Knepley 401b0753f9dSMatthew G. Knepley Level: beginner 402b0753f9dSMatthew G. Knepley 40309e53591SBarry Smith .seealso: [](sec_pc), `PCExoticSetType()`, `PCEXOTIC` 404b0753f9dSMatthew G. Knepley E*/ 4059371c9d4SSatish Balay typedef enum { 4069371c9d4SSatish Balay PC_EXOTIC_FACE, 4079371c9d4SSatish Balay PC_EXOTIC_WIREBASKET 4089371c9d4SSatish Balay } PCExoticType; 409b0753f9dSMatthew G. Knepley 4108c1cd74cSHong Zhang /*E 4110b4b7b1cSBarry Smith PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains with the `PCType` of `PCBDDC` 412bc960bbfSJed Brown 413bc960bbfSJed Brown Values: 41487497f52SBarry Smith + `PC_BDDC_INTERFACE_EXT_DIRICHLET` - solves Dirichlet interior problem; this is the standard BDDC algorithm 4150b4b7b1cSBarry Smith - `PC_BDDC_INTERFACE_EXT_LUMP` - skips interior solve; sometimes called $M_1$ and associated with "lumped FETI-DP" 416bc960bbfSJed Brown 41716a05f60SBarry Smith Level: intermediate 41816a05f60SBarry Smith 41909e53591SBarry Smith .seealso: [](sec_pc), `PCBDDC`, `PC` 420bc960bbfSJed Brown E*/ 421bc960bbfSJed Brown typedef enum { 422bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_DIRICHLET, 423bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_LUMP 424bc960bbfSJed Brown } PCBDDCInterfaceExtType; 425bc960bbfSJed Brown 426bc960bbfSJed Brown /*E 427f3b08a26SMatthew G. Knepley PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation 428f3b08a26SMatthew G. Knepley 429f3b08a26SMatthew G. Knepley Level: beginner 430f3b08a26SMatthew G. Knepley 43109e53591SBarry Smith .seealso: [](sec_pc), `PCMGSetAdaptCoarseSpaceType()`, `PCMG`, `PC` 432f3b08a26SMatthew G. Knepley E*/ 4339371c9d4SSatish Balay typedef enum { 4349371c9d4SSatish Balay PCMG_ADAPT_NONE, 4359371c9d4SSatish Balay PCMG_ADAPT_POLYNOMIAL, 4369371c9d4SSatish Balay PCMG_ADAPT_HARMONIC, 4379371c9d4SSatish Balay PCMG_ADAPT_EIGENVECTOR, 4389371c9d4SSatish Balay PCMG_ADAPT_GENERALIZED_EIGENVECTOR, 4399371c9d4SSatish Balay PCMG_ADAPT_GDSW 4409371c9d4SSatish Balay } PCMGCoarseSpaceType; 441f3b08a26SMatthew G. Knepley 442f3b08a26SMatthew G. Knepley /*E 4430b4b7b1cSBarry Smith PCPatchConstructType - Determines the algorithm used to construct patches for the `PCPATCH` preconditioner 4444bbf5ea8SMatthew G. Knepley 4454bbf5ea8SMatthew G. Knepley Level: beginner 4464bbf5ea8SMatthew G. Knepley 44709e53591SBarry Smith .seealso: [](sec_pc), `PCPatchSetConstructType()`, `PCPATCH`, `PC` 4484bbf5ea8SMatthew G. Knepley E*/ 4499371c9d4SSatish Balay typedef enum { 4509371c9d4SSatish Balay PC_PATCH_STAR, 4519371c9d4SSatish Balay PC_PATCH_VANKA, 4529371c9d4SSatish Balay PC_PATCH_PARDECOMP, 4539371c9d4SSatish Balay PC_PATCH_USER, 4549371c9d4SSatish Balay PC_PATCH_PYTHON 4559371c9d4SSatish Balay } PCPatchConstructType; 4564bbf5ea8SMatthew G. Knepley 4574bbf5ea8SMatthew G. Knepley /*E 4580b4b7b1cSBarry Smith PCDeflationSpaceType - Type of deflation used by `PCType` `PCDEFLATION` 459e662bc50SJakub Kruzik 460e662bc50SJakub Kruzik Values: 46187497f52SBarry Smith + `PC_DEFLATION_SPACE_HAAR` - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off 46209e53591SBarry Smith . `PC_DEFLATION_SPACE_DB2` - `MATCOMPOSITE` of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet) 46387497f52SBarry Smith . `PC_DEFLATION_SPACE_DB4` - same as above, but with db4 (4 coefficient Daubechies) 46487497f52SBarry Smith . `PC_DEFLATION_SPACE_DB8` - same as above, but with db8 (8 coefficient Daubechies) 46587497f52SBarry Smith . `PC_DEFLATION_SPACE_DB16` - same as above, but with db16 (16 coefficient Daubechies) 46687497f52SBarry Smith . `PC_DEFLATION_SPACE_BIORTH22` - same as above, but with biorthogonal 2.2 (6 coefficients) 46787497f52SBarry Smith . `PC_DEFLATION_SPACE_MEYER` - same as above, but with Meyer/FIR (62 coefficients) 468d5b43468SJose E. Roman . `PC_DEFLATION_SPACE_AGGREGATION` - aggregates local indices (given by operator matrix distribution) into a subdomain 46987497f52SBarry Smith - `PC_DEFLATION_SPACE_USER` - indicates space set by user 470e662bc50SJakub Kruzik 4711fdb00f9SJakub Kruzik Level: intermediate 4721fdb00f9SJakub Kruzik 47309e53591SBarry Smith Note: 47409e53591SBarry Smith Wavelet-based space (except Haar) can be used in multilevel deflation. 47509e53591SBarry Smith 47609e53591SBarry Smith .seealso: [](sec_pc), `PCDeflationSetSpaceToCompute()`, `PCDEFLATION`, `PC` 477e662bc50SJakub Kruzik E*/ 478e662bc50SJakub Kruzik typedef enum { 479e662bc50SJakub Kruzik PC_DEFLATION_SPACE_HAAR, 480e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB2, 481e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB4, 482e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB8, 483e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB16, 484e662bc50SJakub Kruzik PC_DEFLATION_SPACE_BIORTH22, 485e662bc50SJakub Kruzik PC_DEFLATION_SPACE_MEYER, 486e662bc50SJakub Kruzik PC_DEFLATION_SPACE_AGGREGATION, 487e662bc50SJakub Kruzik PC_DEFLATION_SPACE_USER 488e662bc50SJakub Kruzik } PCDeflationSpaceType; 489e662bc50SJakub Kruzik 490e662bc50SJakub Kruzik /*E 4910b4b7b1cSBarry Smith PCHPDDMCoarseCorrectionType - Type of coarse correction used by `PCType` `PCHPDDM` 49238cfc46eSPierre Jolivet 49338cfc46eSPierre Jolivet Values: 49409e53591SBarry Smith + `PC_HPDDM_COARSE_CORRECTION_DEFLATED` (default) - eq. (1) in `PCHPDDMShellApply()` 49587497f52SBarry Smith . `PC_HPDDM_COARSE_CORRECTION_ADDITIVE` - eq. (2) 496aa1539e9SPierre Jolivet . `PC_HPDDM_COARSE_CORRECTION_BALANCED` - eq. (3) 497aa1539e9SPierre Jolivet - `PC_HPDDM_COARSE_CORRECTION_NONE` - no coarse correction (mostly useful for debugging) 49838cfc46eSPierre Jolivet 49916a05f60SBarry Smith Level: intermediate 50016a05f60SBarry Smith 50109e53591SBarry Smith .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCSetType()`, `PCHPDDMShellApply()` 50238cfc46eSPierre Jolivet E*/ 5039371c9d4SSatish Balay typedef enum { 5049371c9d4SSatish Balay PC_HPDDM_COARSE_CORRECTION_DEFLATED, 5059371c9d4SSatish Balay PC_HPDDM_COARSE_CORRECTION_ADDITIVE, 506aa1539e9SPierre Jolivet PC_HPDDM_COARSE_CORRECTION_BALANCED, 507aa1539e9SPierre Jolivet PC_HPDDM_COARSE_CORRECTION_NONE 5089371c9d4SSatish Balay } PCHPDDMCoarseCorrectionType; 50938cfc46eSPierre Jolivet 51038cfc46eSPierre Jolivet /*E 51113044ca3SPierre Jolivet PCHPDDMSchurPreType - Type of `PCHPDDM` preconditioner for a `MATSCHURCOMPLEMENT` generated by `PCFIELDSPLIT` with `PCFieldSplitSchurPreType` set to `PC_FIELDSPLIT_SCHUR_PRE_SELF` 51213044ca3SPierre Jolivet 51313044ca3SPierre Jolivet Values: 5140b4b7b1cSBarry 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 5150b4b7b1cSBarry Smith is built by approximating the Schur complement with (inv(sqrt(diag(A00))) A01)^T (inv(sqrt(diag(A00))) A01) 5160b4b7b1cSBarry Smith and by considering the associated linear least squares problem 5170b4b7b1cSBarry 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, 5180b4b7b1cSBarry Smith then `PCHPDDMSetAuxiliaryMat()` must be called on the associated `PC` as well (it is built automatically for the 5190b4b7b1cSBarry Smith user otherwise); the Schur complement `PC` is set internally to `PCKSP`, with the prefix `-fieldsplit_1_pc_hpddm_`; 5200b4b7b1cSBarry Smith the operator associated to the `PC` is spectrally equivalent to the original Schur complement 52113044ca3SPierre Jolivet 52213044ca3SPierre Jolivet Level: advanced 52313044ca3SPierre Jolivet 52413044ca3SPierre Jolivet .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCFIELDSPLIT`, `PC_FIELDSPLIT_SCHUR_PRE_SELF`, `PCFieldSplitSetSchurPre()`, `PCHPDDMSetAuxiliaryMat()` 52513044ca3SPierre Jolivet E*/ 52613044ca3SPierre Jolivet typedef enum { 52713044ca3SPierre Jolivet PC_HPDDM_SCHUR_PRE_LEAST_SQUARES, 528*ce78bad3SBarry Smith PC_HPDDM_SCHUR_PRE_GENEO 52913044ca3SPierre Jolivet } PCHPDDMSchurPreType; 53013044ca3SPierre Jolivet 53113044ca3SPierre Jolivet /*E 5320b4b7b1cSBarry Smith PCFailedReason - indicates the type of `PC` failure. That is why the construction of the preconditioner, `PCSetUp()`, or its use, `PCApply()`, failed 5338c1cd74cSHong Zhang 5348c1cd74cSHong Zhang Level: beginner 5358c1cd74cSHong Zhang 5360b4b7b1cSBarry Smith .seealso: [](sec_pc), `PC`, `PCGetFailedReason()`, `PCSetUp()` 5378c1cd74cSHong Zhang E*/ 5389371c9d4SSatish Balay typedef enum { 5399371c9d4SSatish Balay PC_SETUP_ERROR = -1, 540*ce78bad3SBarry Smith PC_NOERROR = 0, 541*ce78bad3SBarry Smith PC_FACTOR_STRUCT_ZEROPIVOT = 1, 542*ce78bad3SBarry Smith PC_FACTOR_NUMERIC_ZEROPIVOT = 2, 543*ce78bad3SBarry Smith PC_FACTOR_OUTMEMORY = 3, 544*ce78bad3SBarry Smith PC_FACTOR_OTHER = 4, 545*ce78bad3SBarry Smith PC_INCONSISTENT_RHS = 5, 546*ce78bad3SBarry Smith PC_SUBPC_ERROR = 6 5479371c9d4SSatish Balay } PCFailedReason; 548ce7c7f2fSMark Adams 549ce7c7f2fSMark Adams /*E 5500b4b7b1cSBarry Smith PCGAMGLayoutType - Layout for reduced grids for `PCType` `PCGAMG` 551ce7c7f2fSMark Adams 552ce7c7f2fSMark Adams Level: intermediate 553ce7c7f2fSMark Adams 55409e53591SBarry Smith .seealso: [](sec_pc), `PCGAMG`, `PC`, `PCGAMGSetCoarseGridLayoutType()` 555ce7c7f2fSMark Adams E*/ 5569371c9d4SSatish Balay typedef enum { 5579371c9d4SSatish Balay PCGAMG_LAYOUT_COMPACT, 5589371c9d4SSatish Balay PCGAMG_LAYOUT_SPREAD 5599371c9d4SSatish Balay } PCGAMGLayoutType; 560