xref: /petsc/include/petscpctypes.h (revision af27ebaa0199971c43fd2e2e162251afd1bcda49)
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 
94*af27ebaaSBarry 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
10509e53591SBarry Smith     PCRichardsonConvergedReason - reason a `PCRICHARDSON` `PCApplyRichardson()` method terminated
106b0753f9dSMatthew G. Knepley 
107b0753f9dSMatthew G. Knepley    Level: advanced
108b0753f9dSMatthew G. Knepley 
10909e53591SBarry Smith .seealso: [](sec_pc), `PCRICHARDSON`, `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
12309e53591SBarry Smith .  `PC_JACOBI_ROWMAX`   - use the maximum absolute value in the row
12409e53591SBarry Smith -  `PC_JACOBI_ROWSUM`   - use the sum of the values in the row (not the absolute values)
12509e53591SBarry Smith 
126b0753f9dSMatthew G. Knepley    Level: intermediate
127b0753f9dSMatthew G. Knepley 
12809e53591SBarry Smith .seealso: [](sec_pc), `PCJACOBI`, `PC`
129b0753f9dSMatthew G. Knepley E*/
1309371c9d4SSatish Balay typedef enum {
1319371c9d4SSatish Balay   PC_JACOBI_DIAGONAL,
1329371c9d4SSatish Balay   PC_JACOBI_ROWMAX,
1339371c9d4SSatish Balay   PC_JACOBI_ROWSUM
1349371c9d4SSatish Balay } PCJacobiType;
135b0753f9dSMatthew G. Knepley 
136b0753f9dSMatthew G. Knepley /*E
137b0753f9dSMatthew G. Knepley     PCASMType - Type of additive Schwarz method to use
138b0753f9dSMatthew G. Knepley 
13909e53591SBarry Smith    Values:
14009e53591SBarry Smith +  `PC_ASM_BASIC`        - Symmetric version where residuals from the ghost points are used
14109e53591SBarry Smith                            and computed values in ghost regions are added together.
142*af27ebaaSBarry Smith                            Classical standard additive Schwarz as introduced in {cite}`dryja1987additive`.
14309e53591SBarry Smith .  `PC_ASM_RESTRICT`     - Residuals from ghost points are used but computed values in ghost
144*af27ebaaSBarry Smith                            region are discarded {cite}`cs99`. Default.
14509e53591SBarry Smith .  `PC_ASM_INTERPOLATE`  - Residuals from ghost points are not used, computed values in ghost
14609e53591SBarry Smith                            region are added back in.
14709e53591SBarry Smith -  `PC_ASM_NONE`         - Residuals from ghost points are not used, computed ghost values are
148*af27ebaaSBarry Smith                            discarded. Not very good.
149b0753f9dSMatthew G. Knepley 
150b0753f9dSMatthew G. Knepley    Level: beginner
151b0753f9dSMatthew G. Knepley 
15216a05f60SBarry Smith .seealso: [](sec_pc), `PC`, `PCASM`, `PCASMSetType()`, `PCGASMType`
153b0753f9dSMatthew G. Knepley E*/
1549371c9d4SSatish Balay typedef enum {
1559371c9d4SSatish Balay   PC_ASM_BASIC       = 3,
1569371c9d4SSatish Balay   PC_ASM_RESTRICT    = 1,
1579371c9d4SSatish Balay   PC_ASM_INTERPOLATE = 2,
1589371c9d4SSatish Balay   PC_ASM_NONE        = 0
1599371c9d4SSatish Balay } PCASMType;
160b0753f9dSMatthew G. Knepley 
161b0753f9dSMatthew G. Knepley /*E
16287497f52SBarry Smith     PCGASMType - Type of generalized additive Schwarz method to use (differs from `PCASM` in allowing multiple processors per subdomain).
163b0753f9dSMatthew G. Knepley 
16409e53591SBarry Smith    Values:
16509e53591SBarry Smith +  `PC_GASM_BASIC`       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
16609e53591SBarry Smith                            over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
167*af27ebaaSBarry Smith                            from neighboring subdomains. Classical standard additive Schwarz {cite}`dryja1987additive`.
16809e53591SBarry Smith .  `PC_GASM_RESTRICT`    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
16909e53591SBarry Smith                            (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
17009e53591SBarry Smith                            each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
171*af27ebaaSBarry Smith                            assumption) {cite}`cs99`. Default.
17209e53591SBarry Smith .  `PC_GASM_INTERPOLATE` - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
17309e53591SBarry Smith                            applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
17409e53591SBarry Smith                            from neighboring subdomains.
175*af27ebaaSBarry Smith -  `PC_GASM_NONE`        - Residuals and corrections are zeroed out outside the local subdomains. Not very good.
176b0753f9dSMatthew G. Knepley 
177b0753f9dSMatthew G. Knepley    Level: beginner
178b0753f9dSMatthew G. Knepley 
17916a05f60SBarry Smith    Note:
18016a05f60SBarry Smith    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
18116a05f60SBarry Smith    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
18216a05f60SBarry Smith    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
18316a05f60SBarry Smith    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
18416a05f60SBarry Smith 
185*af27ebaaSBarry Smith    Developer Note:
186*af27ebaaSBarry Smith    Perhaps better to remove this since it matches `PCASMType`
187*af27ebaaSBarry Smith 
18816a05f60SBarry Smith .seealso: [](sec_pc), `PCGASM`, `PCASM`, `PC`, `PCGASMSetType()`, `PCASMType`
189b0753f9dSMatthew G. Knepley E*/
1909371c9d4SSatish Balay typedef enum {
1919371c9d4SSatish Balay   PC_GASM_BASIC       = 3,
1929371c9d4SSatish Balay   PC_GASM_RESTRICT    = 1,
1939371c9d4SSatish Balay   PC_GASM_INTERPOLATE = 2,
1949371c9d4SSatish Balay   PC_GASM_NONE        = 0
1959371c9d4SSatish Balay } PCGASMType;
196b0753f9dSMatthew G. Knepley 
197b0753f9dSMatthew G. Knepley /*E
19816a05f60SBarry Smith     PCCompositeType - Determines how two or more preconditioner are composed with the `PCType` of `PCCOMPOSITE`
199b0753f9dSMatthew G. Knepley 
20009e53591SBarry Smith   Values:
20109e53591SBarry Smith +  `PC_COMPOSITE_ADDITIVE`                 - results from application of all preconditioners are added together
20209e53591SBarry Smith .  `PC_COMPOSITE_MULTIPLICATIVE`           - preconditioners are applied sequentially to the residual freshly
20309e53591SBarry Smith                                              computed after the previous preconditioner application
20409e53591SBarry Smith .  `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly
20509e53591SBarry Smith                                              computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners)
206*af27ebaaSBarry Smith .  `PC_COMPOSITE_SPECIAL`                  - This is very special for a matrix of the form $ \alpha I + R + S$
207*af27ebaaSBarry Smith                                              where the first preconditioner is built from $\alpha I + S$ and second from $\alpha I + R$
20809e53591SBarry Smith .  `PC_COMPOSITE_SCHUR`                    - composes the Schur complement of the matrix from two blocks, see `PCFIELDSPLIT`
20909e53591SBarry Smith -  `PC_COMPOSITE_GKB`                      - the generalized Golub-Kahan bidiagonalization preconditioner, see `PCFIELDSPLIT`
210b0753f9dSMatthew G. Knepley 
211b0753f9dSMatthew G. Knepley    Level: beginner
212b0753f9dSMatthew G. Knepley 
213*af27ebaaSBarry Smith .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()`, `SNESCompositeType`, `PCCompositeSpecialSetAlpha()`
214b0753f9dSMatthew G. Knepley E*/
2159371c9d4SSatish Balay typedef enum {
2169371c9d4SSatish Balay   PC_COMPOSITE_ADDITIVE,
2179371c9d4SSatish Balay   PC_COMPOSITE_MULTIPLICATIVE,
2189371c9d4SSatish Balay   PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,
2199371c9d4SSatish Balay   PC_COMPOSITE_SPECIAL,
2209371c9d4SSatish Balay   PC_COMPOSITE_SCHUR,
2219371c9d4SSatish Balay   PC_COMPOSITE_GKB
2229371c9d4SSatish Balay } PCCompositeType;
223b0753f9dSMatthew G. Knepley 
224b0753f9dSMatthew G. Knepley /*E
22587497f52SBarry Smith     PCFieldSplitSchurPreType - Determines how to precondition a Schur complement
226b0753f9dSMatthew G. Knepley 
22716a05f60SBarry Smith     Values:
22813044ca3SPierre Jolivet +  `PC_FIELDSPLIT_SCHUR_PRE_SELF`  - the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix.
22913044ca3SPierre Jolivet                                      The only preconditioners that currently work with this symbolic representation matrix object are `PCLSC` and `PCHPDDM`
23013044ca3SPierre Jolivet .  `PC_FIELDSPLIT_SCHUR_PRE_SELFP` - the preconditioning for the Schur complement is generated from an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01.
23116a05f60SBarry Smith                                      This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be
23216a05f60SBarry Smith                                      lumped before extracting the diagonal using the additional option `-fieldsplit_1_mat_schur_complement_ainv_type lump`
23316a05f60SBarry Smith .  `PC_FIELDSPLIT_SCHUR_PRE_A11`   - the preconditioner for the Schur complement is generated from the block diagonal part of the matrix used to define the preconditioner,
23416a05f60SBarry Smith                                      associated with the Schur complement (i.e. A11), not the Schur complement matrix
23516a05f60SBarry Smith .  `PC_FIELDSPLIT_SCHUR_PRE_USER`  - the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
23616a05f60SBarry Smith                                      to this function).
23716a05f60SBarry Smith -  `PC_FIELDSPLIT_SCHUR_PRE_FULL`  - the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation
23816a05f60SBarry Smith                                      computed internally by `PCFIELDSPLIT` (this is expensive) useful mostly as a test that the Schur complement approach can work for your problem
23916a05f60SBarry Smith 
240b0753f9dSMatthew G. Knepley     Level: intermediate
241b0753f9dSMatthew G. Knepley 
24209e53591SBarry Smith .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurPre()`, `PC`
243b0753f9dSMatthew G. Knepley E*/
2449371c9d4SSatish Balay typedef enum {
2459371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_SELF,
2469371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_SELFP,
2479371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_A11,
2489371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_USER,
2499371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_FULL
2509371c9d4SSatish Balay } PCFieldSplitSchurPreType;
251b0753f9dSMatthew G. Knepley 
252b0753f9dSMatthew G. Knepley /*E
253b0753f9dSMatthew G. Knepley     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
254b0753f9dSMatthew G. Knepley 
25516a05f60SBarry Smith     Values:
25616a05f60SBarry Smith +   `PC_FIELDSPLIT_SCHUR_FACT_DIAG`  - the preconditioner is solving `D`
25716a05f60SBarry Smith .   `PC_FIELDSPLIT_SCHUR_FACT_LOWER` - the preconditioner is solving `L D`
25816a05f60SBarry Smith .   `PC_FIELDSPLIT_SCHUR_FACT_UPPER` - the preconditioner is solving `D U`
25916a05f60SBarry Smith -   `PC_FIELDSPLIT_SCHUR_FACT_FULL`  - the preconditioner is solving `L(D U)`
26016a05f60SBarry Smith 
26116a05f60SBarry Smith     where the matrix is factorized as
26216a05f60SBarry Smith .vb
26316a05f60SBarry Smith    (A   B)  = (1       0) (A   0) (1  Ainv*B)  = L D U
26416a05f60SBarry Smith    (C   E)    (C*Ainv  1) (0   S) (0       1)
26516a05f60SBarry Smith .ve
26616a05f60SBarry Smith 
267b0753f9dSMatthew G. Knepley     Level: intermediate
268b0753f9dSMatthew G. Knepley 
26909e53591SBarry Smith .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurFactType()`, `PC`
270b0753f9dSMatthew G. Knepley E*/
271b0753f9dSMatthew G. Knepley typedef enum {
272b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
273b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
274b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
275b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_FULL
276b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType;
277b0753f9dSMatthew G. Knepley 
278b0753f9dSMatthew G. Knepley /*E
27987497f52SBarry Smith     PCPARMSGlobalType - Determines the global preconditioner method in `PCPARMS`
280b0753f9dSMatthew G. Knepley 
281b0753f9dSMatthew G. Knepley     Level: intermediate
282b0753f9dSMatthew G. Knepley 
28309e53591SBarry Smith .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetGlobal()`, `PC`
284b0753f9dSMatthew G. Knepley E*/
2859371c9d4SSatish Balay typedef enum {
2869371c9d4SSatish Balay   PC_PARMS_GLOBAL_RAS,
2879371c9d4SSatish Balay   PC_PARMS_GLOBAL_SCHUR,
2889371c9d4SSatish Balay   PC_PARMS_GLOBAL_BJ
2899371c9d4SSatish Balay } PCPARMSGlobalType;
2909d8081ecSMatthew G. Knepley 
291b0753f9dSMatthew G. Knepley /*E
29287497f52SBarry Smith     PCPARMSLocalType - Determines the local preconditioner method in `PCPARMS`
293b0753f9dSMatthew G. Knepley 
294b0753f9dSMatthew G. Knepley     Level: intermediate
295b0753f9dSMatthew G. Knepley 
29609e53591SBarry Smith .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetLocal()`, `PC`
297b0753f9dSMatthew G. Knepley E*/
2989371c9d4SSatish Balay typedef enum {
2999371c9d4SSatish Balay   PC_PARMS_LOCAL_ILU0,
3009371c9d4SSatish Balay   PC_PARMS_LOCAL_ILUK,
3019371c9d4SSatish Balay   PC_PARMS_LOCAL_ILUT,
3029371c9d4SSatish Balay   PC_PARMS_LOCAL_ARMS
3039371c9d4SSatish Balay } PCPARMSLocalType;
304b0753f9dSMatthew G. Knepley 
305f0fc11ceSJed Brown /*J
30687497f52SBarry Smith     PCGAMGType - type of generalized algebraic multigrid `PCGAMG` method
307b0753f9dSMatthew G. Knepley 
30809e53591SBarry Smith    Values:
30909e53591SBarry Smith +   `PCGAMGAGG`       - (the default) smoothed aggregation algorithm, robust, very well tested
310baca6076SPierre Jolivet .   `PCGAMGGEO`       - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not supported, reference implementation (2D)
311baca6076SPierre Jolivet -   `PCGAMGCLASSICAL` - classical algebraic multigrid preconditioner, incomplete, not supported, reference implementation
31209e53591SBarry Smith 
313b0753f9dSMatthew G. Knepley      Level: intermediate
314b0753f9dSMatthew G. Knepley 
31509e53591SBarry Smith .seealso: [](sec_pc), `PCGAMG`, `PCMG`, `PC`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()`
316f0fc11ceSJed Brown J*/
317b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType;
318b0753f9dSMatthew G. Knepley #define PCGAMGAGG       "agg"
319b0753f9dSMatthew G. Knepley #define PCGAMGGEO       "geo"
320b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL "classical"
321b0753f9dSMatthew G. Knepley 
322b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType;
323b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT   "direct"
324b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard"
325b0753f9dSMatthew G. Knepley 
326b0753f9dSMatthew G. Knepley /*E
327b0753f9dSMatthew G. Knepley    PCMGType - Determines the type of multigrid method that is run.
328b0753f9dSMatthew G. Knepley 
329b0753f9dSMatthew G. Knepley    Values:
33087497f52SBarry Smith +  `PC_MG_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `PCMGSetCycleType()`
33187497f52SBarry Smith .  `PC_MG_ADDITIVE`                 - the additive multigrid preconditioner where all levels are
332b0753f9dSMatthew G. Knepley                                       smoothed before updating the residual. This only uses the
333b0753f9dSMatthew G. Knepley                                       down smoother, in the preconditioner the upper smoother is ignored
33487497f52SBarry Smith .  `PC_MG_FULL`                     - same as multiplicative except one also performs grid sequencing,
335b0753f9dSMatthew G. Knepley                                       that is starts on the coarsest grid, performs a cycle, interpolates
336b0753f9dSMatthew 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
337b0753f9dSMatthew G. Knepley                                       algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
338b0753f9dSMatthew G. Knepley                                       calls the V-cycle only on the coarser level and has a post-smoother instead.
33916a05f60SBarry Smith -  `PC_MG_KASKADE`                  - like full multigrid except one never goes back to a coarser level from a finer
34016a05f60SBarry Smith 
34116a05f60SBarry Smith    Level: beginner
342b0753f9dSMatthew G. Knepley 
34309e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()`
344b0753f9dSMatthew G. Knepley E*/
3459371c9d4SSatish Balay typedef enum {
3469371c9d4SSatish Balay   PC_MG_MULTIPLICATIVE,
3479371c9d4SSatish Balay   PC_MG_ADDITIVE,
3489371c9d4SSatish Balay   PC_MG_FULL,
3499371c9d4SSatish Balay   PC_MG_KASKADE
3509371c9d4SSatish Balay } PCMGType;
351b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE;
352b0753f9dSMatthew G. Knepley 
353b0753f9dSMatthew G. Knepley /*E
354b0753f9dSMatthew G. Knepley    PCMGCycleType - Use V-cycle or W-cycle
355b0753f9dSMatthew G. Knepley 
356b0753f9dSMatthew G. Knepley    Values:
357*af27ebaaSBarry Smith +  `PC_MG_V_CYCLE` - use the V cycle
358*af27ebaaSBarry Smith -  `PC_MG_W_CYCLE` - use the W cycle
359b0753f9dSMatthew G. Knepley 
36016a05f60SBarry Smith    Level: beginner
36116a05f60SBarry Smith 
36209e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()`
363b0753f9dSMatthew G. Knepley E*/
3649371c9d4SSatish Balay typedef enum {
3659371c9d4SSatish Balay   PC_MG_CYCLE_V = 1,
3669371c9d4SSatish Balay   PC_MG_CYCLE_W = 2
3679371c9d4SSatish Balay } PCMGCycleType;
368b0753f9dSMatthew G. Knepley 
369b0753f9dSMatthew G. Knepley /*E
3702134b1e4SBarry Smith     PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process
3712134b1e4SBarry Smith 
3722134b1e4SBarry Smith    Values:
37387497f52SBarry Smith +  `PC_MG_GALERKIN_PMAT` - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid
37487497f52SBarry Smith .  `PC_MG_GALERKIN_MAT` -  computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid
37587497f52SBarry Smith .  `PC_MG_GALERKIN_BOTH` - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once
37687497f52SBarry Smith -  `PC_MG_GALERKIN_NONE` - neither operator is computed via the Galerkin process, the user must provide the operator
3772134b1e4SBarry Smith 
37809e53591SBarry Smith    Level: beginner
37909e53591SBarry Smith 
38009e53591SBarry Smith    Note:
38187497f52SBarry Smith    Users should never set `PC_MG_GALERKIN_EXTERNAL`, it is used by `PCHYPRE` and `PCML`
3822134b1e4SBarry Smith 
38309e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()`
3842134b1e4SBarry Smith E*/
3859371c9d4SSatish Balay typedef enum {
3869371c9d4SSatish Balay   PC_MG_GALERKIN_BOTH,
3879371c9d4SSatish Balay   PC_MG_GALERKIN_PMAT,
3889371c9d4SSatish Balay   PC_MG_GALERKIN_MAT,
3899371c9d4SSatish Balay   PC_MG_GALERKIN_NONE,
3909371c9d4SSatish Balay   PC_MG_GALERKIN_EXTERNAL
3919371c9d4SSatish Balay } PCMGGalerkinType;
3922134b1e4SBarry Smith 
3932134b1e4SBarry Smith /*E
394b0753f9dSMatthew G. Knepley     PCExoticType - Face based or wirebasket based coarse grid space
395b0753f9dSMatthew G. Knepley 
396b0753f9dSMatthew G. Knepley    Level: beginner
397b0753f9dSMatthew G. Knepley 
39809e53591SBarry Smith .seealso: [](sec_pc), `PCExoticSetType()`, `PCEXOTIC`
399b0753f9dSMatthew G. Knepley E*/
4009371c9d4SSatish Balay typedef enum {
4019371c9d4SSatish Balay   PC_EXOTIC_FACE,
4029371c9d4SSatish Balay   PC_EXOTIC_WIREBASKET
4039371c9d4SSatish Balay } PCExoticType;
404b0753f9dSMatthew G. Knepley 
4058c1cd74cSHong Zhang /*E
406bc960bbfSJed Brown    PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains.
407bc960bbfSJed Brown 
408bc960bbfSJed Brown    Values:
40987497f52SBarry Smith +  `PC_BDDC_INTERFACE_EXT_DIRICHLET` - solves Dirichlet interior problem; this is the standard BDDC algorithm
41087497f52SBarry Smith -  `PC_BDDC_INTERFACE_EXT_LUMP`      - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP"
411bc960bbfSJed Brown 
41216a05f60SBarry Smith    Level: intermediate
41316a05f60SBarry Smith 
41409e53591SBarry Smith .seealso: [](sec_pc), `PCBDDC`, `PC`
415bc960bbfSJed Brown E*/
416bc960bbfSJed Brown typedef enum {
417bc960bbfSJed Brown   PC_BDDC_INTERFACE_EXT_DIRICHLET,
418bc960bbfSJed Brown   PC_BDDC_INTERFACE_EXT_LUMP
419bc960bbfSJed Brown } PCBDDCInterfaceExtType;
420bc960bbfSJed Brown 
421bc960bbfSJed Brown /*E
422f3b08a26SMatthew G. Knepley   PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation
423f3b08a26SMatthew G. Knepley 
424f3b08a26SMatthew G. Knepley   Level: beginner
425f3b08a26SMatthew G. Knepley 
42609e53591SBarry Smith .seealso: [](sec_pc), `PCMGSetAdaptCoarseSpaceType()`, `PCMG`, `PC`
427f3b08a26SMatthew G. Knepley E*/
4289371c9d4SSatish Balay typedef enum {
4299371c9d4SSatish Balay   PCMG_ADAPT_NONE,
4309371c9d4SSatish Balay   PCMG_ADAPT_POLYNOMIAL,
4319371c9d4SSatish Balay   PCMG_ADAPT_HARMONIC,
4329371c9d4SSatish Balay   PCMG_ADAPT_EIGENVECTOR,
4339371c9d4SSatish Balay   PCMG_ADAPT_GENERALIZED_EIGENVECTOR,
4349371c9d4SSatish Balay   PCMG_ADAPT_GDSW
4359371c9d4SSatish Balay } PCMGCoarseSpaceType;
436f3b08a26SMatthew G. Knepley 
437f3b08a26SMatthew G. Knepley /*E
43816a05f60SBarry Smith     PCPatchConstructType - The algorithm used to construct patches for the `PCPATCH` preconditioner
4394bbf5ea8SMatthew G. Knepley 
4404bbf5ea8SMatthew G. Knepley    Level: beginner
4414bbf5ea8SMatthew G. Knepley 
44209e53591SBarry Smith .seealso: [](sec_pc), `PCPatchSetConstructType()`, `PCPATCH`, `PC`
4434bbf5ea8SMatthew G. Knepley E*/
4449371c9d4SSatish Balay typedef enum {
4459371c9d4SSatish Balay   PC_PATCH_STAR,
4469371c9d4SSatish Balay   PC_PATCH_VANKA,
4479371c9d4SSatish Balay   PC_PATCH_PARDECOMP,
4489371c9d4SSatish Balay   PC_PATCH_USER,
4499371c9d4SSatish Balay   PC_PATCH_PYTHON
4509371c9d4SSatish Balay } PCPatchConstructType;
4514bbf5ea8SMatthew G. Knepley 
4524bbf5ea8SMatthew G. Knepley /*E
453e26ad46dSJakub Kruzik     PCDeflationSpaceType - Type of deflation
454e662bc50SJakub Kruzik 
455e662bc50SJakub Kruzik     Values:
45687497f52SBarry Smith +   `PC_DEFLATION_SPACE_HAAR`        - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off
45709e53591SBarry Smith .   `PC_DEFLATION_SPACE_DB2`         - `MATCOMPOSITE` of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet)
45887497f52SBarry Smith .   `PC_DEFLATION_SPACE_DB4`         - same as above, but with db4 (4 coefficient Daubechies)
45987497f52SBarry Smith .   `PC_DEFLATION_SPACE_DB8`         - same as above, but with db8 (8 coefficient Daubechies)
46087497f52SBarry Smith .   `PC_DEFLATION_SPACE_DB16`        - same as above, but with db16 (16 coefficient Daubechies)
46187497f52SBarry Smith .   `PC_DEFLATION_SPACE_BIORTH22`    - same as above, but with biorthogonal 2.2 (6 coefficients)
46287497f52SBarry Smith .   `PC_DEFLATION_SPACE_MEYER`       - same as above, but with Meyer/FIR (62 coefficients)
463d5b43468SJose E. Roman .   `PC_DEFLATION_SPACE_AGGREGATION` - aggregates local indices (given by operator matrix distribution) into a subdomain
46487497f52SBarry Smith -   `PC_DEFLATION_SPACE_USER`        - indicates space set by user
465e662bc50SJakub Kruzik 
4661fdb00f9SJakub Kruzik     Level: intermediate
4671fdb00f9SJakub Kruzik 
46809e53591SBarry Smith     Note:
46909e53591SBarry Smith     Wavelet-based space (except Haar) can be used in multilevel deflation.
47009e53591SBarry Smith 
47109e53591SBarry Smith .seealso: [](sec_pc), `PCDeflationSetSpaceToCompute()`, `PCDEFLATION`, `PC`
472e662bc50SJakub Kruzik E*/
473e662bc50SJakub Kruzik typedef enum {
474e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_HAAR,
475e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB2,
476e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB4,
477e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB8,
478e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB16,
479e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_BIORTH22,
480e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_MEYER,
481e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_AGGREGATION,
482e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_USER
483e662bc50SJakub Kruzik } PCDeflationSpaceType;
484e662bc50SJakub Kruzik 
485e662bc50SJakub Kruzik /*E
48687497f52SBarry Smith     PCHPDDMCoarseCorrectionType - Type of coarse correction used by `PCHPDDM`
48738cfc46eSPierre Jolivet 
48838cfc46eSPierre Jolivet     Values:
48909e53591SBarry Smith +   `PC_HPDDM_COARSE_CORRECTION_DEFLATED` (default) - eq. (1) in `PCHPDDMShellApply()`
49087497f52SBarry Smith .   `PC_HPDDM_COARSE_CORRECTION_ADDITIVE`           - eq. (2)
49187497f52SBarry Smith -   `PC_HPDDM_COARSE_CORRECTION_BALANCED`           - eq. (3)
49238cfc46eSPierre Jolivet 
49316a05f60SBarry Smith     Level: intermediate
49416a05f60SBarry Smith 
49509e53591SBarry Smith .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCSetType()`, `PCHPDDMShellApply()`
49638cfc46eSPierre Jolivet E*/
4979371c9d4SSatish Balay typedef enum {
4989371c9d4SSatish Balay   PC_HPDDM_COARSE_CORRECTION_DEFLATED,
4999371c9d4SSatish Balay   PC_HPDDM_COARSE_CORRECTION_ADDITIVE,
5009371c9d4SSatish Balay   PC_HPDDM_COARSE_CORRECTION_BALANCED
5019371c9d4SSatish Balay } PCHPDDMCoarseCorrectionType;
50238cfc46eSPierre Jolivet 
50338cfc46eSPierre Jolivet /*E
50413044ca3SPierre Jolivet     PCHPDDMSchurPreType - Type of `PCHPDDM` preconditioner for a `MATSCHURCOMPLEMENT` generated by `PCFIELDSPLIT` with `PCFieldSplitSchurPreType` set to `PC_FIELDSPLIT_SCHUR_PRE_SELF`
50513044ca3SPierre Jolivet 
50613044ca3SPierre Jolivet     Values:
50713044ca3SPierre 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
50813044ca3SPierre 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
50913044ca3SPierre Jolivet 
51013044ca3SPierre Jolivet     Level: advanced
51113044ca3SPierre Jolivet 
51213044ca3SPierre Jolivet .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCFIELDSPLIT`, `PC_FIELDSPLIT_SCHUR_PRE_SELF`, `PCFieldSplitSetSchurPre()`, `PCHPDDMSetAuxiliaryMat()`
51313044ca3SPierre Jolivet E*/
51413044ca3SPierre Jolivet typedef enum {
51513044ca3SPierre Jolivet   PC_HPDDM_SCHUR_PRE_LEAST_SQUARES,
51613044ca3SPierre Jolivet   PC_HPDDM_SCHUR_PRE_GENEO,
51713044ca3SPierre Jolivet } PCHPDDMSchurPreType;
51813044ca3SPierre Jolivet 
51913044ca3SPierre Jolivet /*E
52087497f52SBarry Smith     PCFailedReason - indicates type of `PC` failure
5218c1cd74cSHong Zhang 
5228c1cd74cSHong Zhang     Level: beginner
5238c1cd74cSHong Zhang 
52409e53591SBarry Smith .seealso: [](sec_pc), `PC`
5258c1cd74cSHong Zhang E*/
5269371c9d4SSatish Balay typedef enum {
5279371c9d4SSatish Balay   PC_SETUP_ERROR = -1,
5289371c9d4SSatish Balay   PC_NOERROR,
5299371c9d4SSatish Balay   PC_FACTOR_STRUCT_ZEROPIVOT,
5309371c9d4SSatish Balay   PC_FACTOR_NUMERIC_ZEROPIVOT,
5319371c9d4SSatish Balay   PC_FACTOR_OUTMEMORY,
5329371c9d4SSatish Balay   PC_FACTOR_OTHER,
53387b47708SBarry Smith   PC_INCONSISTENT_RHS,
5349371c9d4SSatish Balay   PC_SUBPC_ERROR
5359371c9d4SSatish Balay } PCFailedReason;
536ce7c7f2fSMark Adams 
537ce7c7f2fSMark Adams /*E
538ce7c7f2fSMark Adams     PCGAMGLayoutType - Layout for reduced grids
539ce7c7f2fSMark Adams 
540ce7c7f2fSMark Adams     Level: intermediate
541ce7c7f2fSMark Adams 
54209e53591SBarry Smith .seealso: [](sec_pc), `PCGAMG`, `PC`, `PCGAMGSetCoarseGridLayoutType()`
543ce7c7f2fSMark Adams E*/
5449371c9d4SSatish Balay typedef enum {
5459371c9d4SSatish Balay   PCGAMG_LAYOUT_COMPACT,
5469371c9d4SSatish Balay   PCGAMG_LAYOUT_SPREAD
5479371c9d4SSatish Balay } PCGAMGLayoutType;
548