xref: /petsc/include/petscpctypes.h (revision 5e4ac8c8fe9feb46e07e9b4aa07740c81730fddd)
16524c165SJacob Faibussowitsch #ifndef PETSCPCTYPES_H
226bd1501SBarry Smith #define PETSCPCTYPES_H
3b0753f9dSMatthew G. Knepley 
4ac09b921SBarry Smith /* SUBMANSEC = PC */
5ac09b921SBarry Smith 
6b0753f9dSMatthew G. Knepley /*S
787497f52SBarry Smith      PC - Abstract PETSc object that manages all preconditioners including direct solvers such as `PCLU`
8b0753f9dSMatthew G. Knepley 
9b0753f9dSMatthew G. Knepley    Level: beginner
10b0753f9dSMatthew G. Knepley 
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
1616a05f60SBarry Smith     PCType - String with the name of a PETSc preconditioner
17b0753f9dSMatthew G. Knepley 
18b0753f9dSMatthew G. Knepley    Level: beginner
19b0753f9dSMatthew G. Knepley 
2087497f52SBarry Smith    Note:
2187497f52SBarry Smith    `PCRegister()` is used to register preconditioners that are then accessible via `PCSetType()`
22b0753f9dSMatthew G. Knepley 
2316a05f60SBarry Smith .seealso: [](doc_linsolve), [](sec_pc), `PCSetType()`, `PC`, `PCCreate()`, `PCRegister()`, `PCSetFromOptions()`, `PCLU`, `PCJACOBI`, `PCBJACOBI`
24b0753f9dSMatthew G. Knepley J*/
25b0753f9dSMatthew G. Knepley typedef const char *PCType;
26b0753f9dSMatthew G. Knepley #define PCNONE               "none"
27b0753f9dSMatthew G. Knepley #define PCJACOBI             "jacobi"
28b0753f9dSMatthew G. Knepley #define PCSOR                "sor"
29b0753f9dSMatthew G. Knepley #define PCLU                 "lu"
30a2fc1e05SToby Isaac #define PCQR                 "qr"
31b0753f9dSMatthew G. Knepley #define PCSHELL              "shell"
32e6f8f311SMark Adams #define PCAMGX               "amgx"
33b0753f9dSMatthew G. Knepley #define PCBJACOBI            "bjacobi"
34b0753f9dSMatthew G. Knepley #define PCMG                 "mg"
35b0753f9dSMatthew G. Knepley #define PCEISENSTAT          "eisenstat"
36b0753f9dSMatthew G. Knepley #define PCILU                "ilu"
37b0753f9dSMatthew G. Knepley #define PCICC                "icc"
38b0753f9dSMatthew G. Knepley #define PCASM                "asm"
39b0753f9dSMatthew G. Knepley #define PCGASM               "gasm"
40b0753f9dSMatthew G. Knepley #define PCKSP                "ksp"
41e607c864SMark Adams #define PCBJKOKKOS           "bjkokkos"
42b0753f9dSMatthew G. Knepley #define PCCOMPOSITE          "composite"
43b0753f9dSMatthew G. Knepley #define PCREDUNDANT          "redundant"
44b0753f9dSMatthew G. Knepley #define PCSPAI               "spai"
45b0753f9dSMatthew G. Knepley #define PCNN                 "nn"
46b0753f9dSMatthew G. Knepley #define PCCHOLESKY           "cholesky"
47b0753f9dSMatthew G. Knepley #define PCPBJACOBI           "pbjacobi"
480da83c2eSBarry Smith #define PCVPBJACOBI          "vpbjacobi"
49b0753f9dSMatthew G. Knepley #define PCMAT                "mat"
50b0753f9dSMatthew G. Knepley #define PCHYPRE              "hypre"
51b0753f9dSMatthew G. Knepley #define PCPARMS              "parms"
52b0753f9dSMatthew G. Knepley #define PCFIELDSPLIT         "fieldsplit"
53b0753f9dSMatthew G. Knepley #define PCTFS                "tfs"
54b0753f9dSMatthew G. Knepley #define PCML                 "ml"
55b0753f9dSMatthew G. Knepley #define PCGALERKIN           "galerkin"
56b0753f9dSMatthew G. Knepley #define PCEXOTIC             "exotic"
57b0753f9dSMatthew G. Knepley #define PCCP                 "cp"
58b0753f9dSMatthew G. Knepley #define PCBFBT               "bfbt"
59b0753f9dSMatthew G. Knepley #define PCLSC                "lsc"
60b0753f9dSMatthew G. Knepley #define PCPYTHON             "python"
61b0753f9dSMatthew G. Knepley #define PCPFMG               "pfmg"
621c188c59Sftrigaux #define PCSMG                "smg"
63b0753f9dSMatthew G. Knepley #define PCSYSPFMG            "syspfmg"
64b0753f9dSMatthew G. Knepley #define PCREDISTRIBUTE       "redistribute"
65b0753f9dSMatthew G. Knepley #define PCSVD                "svd"
66b0753f9dSMatthew G. Knepley #define PCGAMG               "gamg"
674b3f184cSKarl Rupp #define PCCHOWILUVIENNACL    "chowiluviennacl"
6870baa948SKarl Rupp #define PCROWSCALINGVIENNACL "rowscalingviennacl"
6907598726SKarl Rupp #define PCSAVIENNACL         "saviennacl"
70b0753f9dSMatthew G. Knepley #define PCBDDC               "bddc"
71b0753f9dSMatthew G. Knepley #define PCKACZMARZ           "kaczmarz"
7268ddcbeaSBarry Smith #define PCTELESCOPE          "telescope"
734bbf5ea8SMatthew G. Knepley #define PCPATCH              "patch"
74b9ac7092SAlp Dener #define PCLMVM               "lmvm"
75360ee056SFande Kong #define PCHMG                "hmg"
7637eeb815SJakub Kruzik #define PCDEFLATION          "deflation"
7738cfc46eSPierre Jolivet #define PCHPDDM              "hpddm"
7853022affSStefano Zampini #define PCH2OPUS             "h2opus"
79f1f2ae84SBarry Smith #define PCMPI                "mpi"
80b0753f9dSMatthew G. Knepley 
81b0753f9dSMatthew G. Knepley /*E
82b0753f9dSMatthew G. Knepley     PCSide - If the preconditioner is to be applied to the left, right
83b0753f9dSMatthew G. Knepley      or symmetrically around the operator.
84b0753f9dSMatthew G. Knepley 
8516a05f60SBarry Smith    Values:
8616a05f60SBarry Smith +  `PC_LEFT` - applied after the operator is applied
8716a05f60SBarry Smith .  `PC_RIGHT` - applied before the operator is applied
8816a05f60SBarry 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.
8916a05f60SBarry Smith 
90b0753f9dSMatthew G. Knepley    Level: beginner
91b0753f9dSMatthew G. Knepley 
9216a05f60SBarry Smith    Note:
9316a05f60SBarry Smith    Certain `KSPType` support only a subset of `PCSide` values
9416a05f60SBarry Smith 
9516a05f60SBarry Smith .seealso: [](sec_pc), `PC`, `KSPSetPCSide()`
96b0753f9dSMatthew G. Knepley E*/
979371c9d4SSatish Balay typedef enum {
989371c9d4SSatish Balay   PC_SIDE_DEFAULT = -1,
999371c9d4SSatish Balay   PC_LEFT,
1009371c9d4SSatish Balay   PC_RIGHT,
1019371c9d4SSatish Balay   PC_SYMMETRIC
1029371c9d4SSatish Balay } PCSide;
103b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
104b0753f9dSMatthew G. Knepley 
105b0753f9dSMatthew G. Knepley /*E
10609e53591SBarry Smith     PCRichardsonConvergedReason - reason a `PCRICHARDSON` `PCApplyRichardson()` method terminated
107b0753f9dSMatthew G. Knepley 
108b0753f9dSMatthew G. Knepley    Level: advanced
109b0753f9dSMatthew G. Knepley 
11087497f52SBarry Smith    Developer Note:
11116a05f60SBarry Smith   This must match `include/petsc/finclude/petscpc.h` and the `KSPConvergedReason` values in `include/petscksp.h
112b0753f9dSMatthew G. Knepley 
11309e53591SBarry Smith .seealso: [](sec_pc), `PCRICHARDSON`, `PC`, `PCApplyRichardson()`
114b0753f9dSMatthew G. Knepley E*/
115b0753f9dSMatthew G. Knepley typedef enum {
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
12316a05f60SBarry Smith     PCJacobiType - What elements of the matrix are used to form the Jacobi preconditioner
124b0753f9dSMatthew G. Knepley 
12509e53591SBarry Smith    Values:
12609e53591SBarry Smith +  `PC_JACOBI_DIAGONAL` - use the diagonal entry, if it is zero use one
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,
1369371c9d4SSatish Balay   PC_JACOBI_ROWMAX,
1379371c9d4SSatish Balay   PC_JACOBI_ROWSUM
1389371c9d4SSatish Balay } PCJacobiType;
139b0753f9dSMatthew G. Knepley 
140b0753f9dSMatthew G. Knepley /*E
141b0753f9dSMatthew G. Knepley     PCASMType - Type of additive Schwarz method to use
142b0753f9dSMatthew G. Knepley 
14309e53591SBarry Smith    Values:
14409e53591SBarry Smith +  `PC_ASM_BASIC`        - Symmetric version where residuals from the ghost points are used
14509e53591SBarry Smith                         and computed values in ghost regions are added together.
14609e53591SBarry Smith                         Classical standard additive Schwarz.
14709e53591SBarry Smith .  `PC_ASM_RESTRICT`     - Residuals from ghost points are used but computed values in ghost
14809e53591SBarry Smith                         region are discarded.
14909e53591SBarry Smith                         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
15309e53591SBarry Smith                         discarded.
15409e53591SBarry Smith                         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
16887497f52SBarry Smith     PCGASMType - Type of generalized additive Schwarz method to use (differs from `PCASM` in allowing multiple processors per subdomain).
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
17309e53591SBarry Smith                         from neighboring subdomains.
17409e53591SBarry Smith                         Classical standard additive Schwarz.
17509e53591SBarry Smith .  `PC_GASM_RESTRICT`    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
17609e53591SBarry Smith                         (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
17709e53591SBarry Smith                         each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
17809e53591SBarry Smith                         assumption).
17909e53591SBarry Smith                         Default.
18009e53591SBarry Smith .  `PC_GASM_INTERPOLATE` - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
18109e53591SBarry Smith                         applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
18209e53591SBarry Smith                         from neighboring subdomains.
18309e53591SBarry Smith -  `PC_GASM_NONE`       - Residuals and corrections are zeroed out outside the local subdomains.
18409e53591SBarry Smith                         Not very good.
185b0753f9dSMatthew G. Knepley 
186b0753f9dSMatthew G. Knepley    Level: beginner
187b0753f9dSMatthew G. Knepley 
18816a05f60SBarry Smith    Note:
18916a05f60SBarry Smith      Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
19016a05f60SBarry Smith    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
19116a05f60SBarry Smith    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
19216a05f60SBarry Smith    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
19316a05f60SBarry 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)
21209e53591SBarry Smith .  `PC_COMPOSITE_SPECIAL` - This is very special for a matrix of the form alpha I + R + S
21309e53591SBarry Smith                          where first preconditioner is built from alpha I + S and second from
21409e53591SBarry Smith                          alpha I + R
21509e53591SBarry Smith .  `PC_COMPOSITE_SCHUR` -  composes the Schur complement of the matrix from two blocks, see `PCFIELDSPLIT`
21609e53591SBarry Smith -  `PC_COMPOSITE_GKB` - the generalized Golub-Kahan bidiagonalization preconditioner, see `PCFIELDSPLIT`
217b0753f9dSMatthew G. Knepley 
218b0753f9dSMatthew G. Knepley    Level: beginner
219b0753f9dSMatthew G. Knepley 
22009e53591SBarry Smith .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()`
221b0753f9dSMatthew G. Knepley E*/
2229371c9d4SSatish Balay typedef enum {
2239371c9d4SSatish Balay   PC_COMPOSITE_ADDITIVE,
2249371c9d4SSatish Balay   PC_COMPOSITE_MULTIPLICATIVE,
2259371c9d4SSatish Balay   PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,
2269371c9d4SSatish Balay   PC_COMPOSITE_SPECIAL,
2279371c9d4SSatish Balay   PC_COMPOSITE_SCHUR,
2289371c9d4SSatish Balay   PC_COMPOSITE_GKB
2299371c9d4SSatish Balay } PCCompositeType;
230b0753f9dSMatthew G. Knepley 
231b0753f9dSMatthew G. Knepley /*E
23287497f52SBarry Smith     PCFieldSplitSchurPreType - Determines how to precondition a Schur complement
233b0753f9dSMatthew G. Knepley 
23416a05f60SBarry Smith     Values:
23516a05f60SBarry Smith +  `PC_FIELDSPLIT_SCHUR_PRE_SELF` - the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix:
23616a05f60SBarry Smith           The only preconditioner that currently works with this symbolic representation matrix object is `PCLSC`
23716a05f60SBarry Smith .  `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
23816a05f60SBarry Smith           This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be
23916a05f60SBarry Smith           lumped before extracting the diagonal using the additional option `-fieldsplit_1_mat_schur_complement_ainv_type lump`
24016a05f60SBarry 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,
24116a05f60SBarry Smith                                  associated with the Schur complement (i.e. A11), not the Schur complement matrix
24216a05f60SBarry Smith .  `PC_FIELDSPLIT_SCHUR_PRE_USER` - the preconditioner for the Schur complement is generated from the user provided matrix (pre argument
24316a05f60SBarry Smith           to this function).
24416a05f60SBarry Smith -  `PC_FIELDSPLIT_SCHUR_PRE_FULL` -  the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation
24516a05f60SBarry Smith       computed internally by `PCFIELDSPLIT` (this is expensive) useful mostly as a test that the Schur complement approach can work for your problem
24616a05f60SBarry Smith 
247b0753f9dSMatthew G. Knepley     Level: intermediate
248b0753f9dSMatthew G. Knepley 
24909e53591SBarry Smith .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurPre()`, `PC`
250b0753f9dSMatthew G. Knepley E*/
2519371c9d4SSatish Balay typedef enum {
2529371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_SELF,
2539371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_SELFP,
2549371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_A11,
2559371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_USER,
2569371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_FULL
2579371c9d4SSatish Balay } PCFieldSplitSchurPreType;
258b0753f9dSMatthew G. Knepley 
259b0753f9dSMatthew G. Knepley /*E
260b0753f9dSMatthew G. Knepley     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
261b0753f9dSMatthew G. Knepley 
26216a05f60SBarry Smith     Values:
26316a05f60SBarry Smith +   `PC_FIELDSPLIT_SCHUR_FACT_DIAG` - the preconditioner is solving `D`
26416a05f60SBarry Smith .   `PC_FIELDSPLIT_SCHUR_FACT_LOWER` - the preconditioner is solving `L D`
26516a05f60SBarry Smith .   `PC_FIELDSPLIT_SCHUR_FACT_UPPER` - the preconditioner is solving `D U`
26616a05f60SBarry Smith -   `PC_FIELDSPLIT_SCHUR_FACT_FULL` - the preconditioner is solving `L(D U)`
26716a05f60SBarry Smith 
26816a05f60SBarry Smith     where the matrix is factorized as
26916a05f60SBarry Smith .vb
27016a05f60SBarry Smith    (A   B)  = (1       0) (A   0) (1  Ainv*B)  = L D U
27116a05f60SBarry Smith    (C   E)    (C*Ainv  1) (0   S) (0       1)
27216a05f60SBarry Smith .ve
27316a05f60SBarry Smith 
274b0753f9dSMatthew G. Knepley     Level: intermediate
275b0753f9dSMatthew G. Knepley 
27609e53591SBarry Smith .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurFactType()`, `PC`
277b0753f9dSMatthew G. Knepley E*/
278b0753f9dSMatthew G. Knepley typedef enum {
279b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
280b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
281b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
282b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_FULL
283b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType;
284b0753f9dSMatthew G. Knepley 
285b0753f9dSMatthew G. Knepley /*E
28687497f52SBarry Smith     PCPARMSGlobalType - Determines the global preconditioner method in `PCPARMS`
287b0753f9dSMatthew G. Knepley 
288b0753f9dSMatthew G. Knepley     Level: intermediate
289b0753f9dSMatthew G. Knepley 
29009e53591SBarry Smith .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetGlobal()`, `PC`
291b0753f9dSMatthew G. Knepley E*/
2929371c9d4SSatish Balay typedef enum {
2939371c9d4SSatish Balay   PC_PARMS_GLOBAL_RAS,
2949371c9d4SSatish Balay   PC_PARMS_GLOBAL_SCHUR,
2959371c9d4SSatish Balay   PC_PARMS_GLOBAL_BJ
2969371c9d4SSatish Balay } PCPARMSGlobalType;
2979d8081ecSMatthew G. Knepley 
298b0753f9dSMatthew G. Knepley /*E
29987497f52SBarry Smith     PCPARMSLocalType - Determines the local preconditioner method in `PCPARMS`
300b0753f9dSMatthew G. Knepley 
301b0753f9dSMatthew G. Knepley     Level: intermediate
302b0753f9dSMatthew G. Knepley 
30309e53591SBarry Smith .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetLocal()`, `PC`
304b0753f9dSMatthew G. Knepley E*/
3059371c9d4SSatish Balay typedef enum {
3069371c9d4SSatish Balay   PC_PARMS_LOCAL_ILU0,
3079371c9d4SSatish Balay   PC_PARMS_LOCAL_ILUK,
3089371c9d4SSatish Balay   PC_PARMS_LOCAL_ILUT,
3099371c9d4SSatish Balay   PC_PARMS_LOCAL_ARMS
3109371c9d4SSatish Balay } PCPARMSLocalType;
311b0753f9dSMatthew G. Knepley 
312f0fc11ceSJed Brown /*J
31387497f52SBarry Smith     PCGAMGType - type of generalized algebraic multigrid `PCGAMG` method
314b0753f9dSMatthew G. Knepley 
31509e53591SBarry Smith    Values:
31609e53591SBarry Smith +   `PCGAMGAGG` - (the default) smoothed aggregation algorithm, robust, very well tested
317*5e4ac8c8Smarkadams4 .   `PCGAMGGEO` - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not supported, referance implementation (2D)
318*5e4ac8c8Smarkadams4 -   `PCGAMGCLASSICAL` - classical algebraic multigrid preconditioner, incomplete, not supported, referance implementation
31909e53591SBarry Smith 
320b0753f9dSMatthew G. Knepley      Level: intermediate
321b0753f9dSMatthew G. Knepley 
32209e53591SBarry Smith .seealso: [](sec_pc), `PCGAMG`, `PCMG`, `PC`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()`
323f0fc11ceSJed Brown J*/
324b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType;
325b0753f9dSMatthew G. Knepley #define PCGAMGAGG       "agg"
326b0753f9dSMatthew G. Knepley #define PCGAMGGEO       "geo"
327b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL "classical"
328b0753f9dSMatthew G. Knepley 
329b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType;
330b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT   "direct"
331b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard"
332b0753f9dSMatthew G. Knepley 
333b0753f9dSMatthew G. Knepley /*E
334b0753f9dSMatthew G. Knepley     PCMGType - Determines the type of multigrid method that is run.
335b0753f9dSMatthew G. Knepley 
336b0753f9dSMatthew G. Knepley    Values:
33787497f52SBarry Smith +  `PC_MG_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `PCMGSetCycleType()`
33887497f52SBarry Smith .  `PC_MG_ADDITIVE` - the additive multigrid preconditioner where all levels are
339b0753f9dSMatthew G. Knepley                 smoothed before updating the residual. This only uses the
340b0753f9dSMatthew G. Knepley                 down smoother, in the preconditioner the upper smoother is ignored
34187497f52SBarry Smith .  `PC_MG_FULL` - same as multiplicative except one also performs grid sequencing,
342b0753f9dSMatthew G. Knepley             that is starts on the coarsest grid, performs a cycle, interpolates
343b0753f9dSMatthew 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
344b0753f9dSMatthew G. Knepley             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
345b0753f9dSMatthew G. Knepley             calls the V-cycle only on the coarser level and has a post-smoother instead.
34616a05f60SBarry Smith -  `PC_MG_KASKADE` - like full multigrid except one never goes back to a coarser level from a finer
34716a05f60SBarry Smith 
34816a05f60SBarry Smith    Level: beginner
349b0753f9dSMatthew G. Knepley 
35009e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()`
351b0753f9dSMatthew G. Knepley E*/
3529371c9d4SSatish Balay typedef enum {
3539371c9d4SSatish Balay   PC_MG_MULTIPLICATIVE,
3549371c9d4SSatish Balay   PC_MG_ADDITIVE,
3559371c9d4SSatish Balay   PC_MG_FULL,
3569371c9d4SSatish Balay   PC_MG_KASKADE
3579371c9d4SSatish Balay } PCMGType;
358b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE;
359b0753f9dSMatthew G. Knepley 
360b0753f9dSMatthew G. Knepley /*E
361b0753f9dSMatthew G. Knepley     PCMGCycleType - Use V-cycle or W-cycle
362b0753f9dSMatthew G. Knepley 
363b0753f9dSMatthew G. Knepley    Values:
36487497f52SBarry Smith +  `PC_MG_V_CYCLE` - use the v cycle
36587497f52SBarry Smith -  `PC_MG_W_CYCLE` - use the w cycle
366b0753f9dSMatthew G. Knepley 
36716a05f60SBarry Smith    Level: beginner
36816a05f60SBarry Smith 
36909e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()`
370b0753f9dSMatthew G. Knepley E*/
3719371c9d4SSatish Balay typedef enum {
3729371c9d4SSatish Balay   PC_MG_CYCLE_V = 1,
3739371c9d4SSatish Balay   PC_MG_CYCLE_W = 2
3749371c9d4SSatish Balay } PCMGCycleType;
375b0753f9dSMatthew G. Knepley 
376b0753f9dSMatthew G. Knepley /*E
3772134b1e4SBarry Smith     PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process
3782134b1e4SBarry Smith 
3792134b1e4SBarry Smith    Values:
38087497f52SBarry Smith +  `PC_MG_GALERKIN_PMAT` - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid
38187497f52SBarry Smith .  `PC_MG_GALERKIN_MAT` -  computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid
38287497f52SBarry Smith .  `PC_MG_GALERKIN_BOTH` - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once
38387497f52SBarry Smith -  `PC_MG_GALERKIN_NONE` - neither operator is computed via the Galerkin process, the user must provide the operator
3842134b1e4SBarry Smith 
38509e53591SBarry Smith    Level: beginner
38609e53591SBarry Smith 
38709e53591SBarry Smith    Note:
38887497f52SBarry Smith    Users should never set `PC_MG_GALERKIN_EXTERNAL`, it is used by `PCHYPRE` and `PCML`
3892134b1e4SBarry Smith 
39009e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()`
3912134b1e4SBarry Smith E*/
3929371c9d4SSatish Balay typedef enum {
3939371c9d4SSatish Balay   PC_MG_GALERKIN_BOTH,
3949371c9d4SSatish Balay   PC_MG_GALERKIN_PMAT,
3959371c9d4SSatish Balay   PC_MG_GALERKIN_MAT,
3969371c9d4SSatish Balay   PC_MG_GALERKIN_NONE,
3979371c9d4SSatish Balay   PC_MG_GALERKIN_EXTERNAL
3989371c9d4SSatish Balay } PCMGGalerkinType;
3992134b1e4SBarry Smith 
4002134b1e4SBarry Smith /*E
401b0753f9dSMatthew G. Knepley     PCExoticType - Face based or wirebasket based coarse grid space
402b0753f9dSMatthew G. Knepley 
403b0753f9dSMatthew G. Knepley    Level: beginner
404b0753f9dSMatthew G. Knepley 
40509e53591SBarry Smith .seealso: [](sec_pc), `PCExoticSetType()`, `PCEXOTIC`
406b0753f9dSMatthew G. Knepley E*/
4079371c9d4SSatish Balay typedef enum {
4089371c9d4SSatish Balay   PC_EXOTIC_FACE,
4099371c9d4SSatish Balay   PC_EXOTIC_WIREBASKET
4109371c9d4SSatish Balay } PCExoticType;
411b0753f9dSMatthew G. Knepley 
4128c1cd74cSHong Zhang /*E
413bc960bbfSJed Brown    PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains.
414bc960bbfSJed Brown 
415bc960bbfSJed Brown    Values:
41687497f52SBarry Smith +  `PC_BDDC_INTERFACE_EXT_DIRICHLET` - solves Dirichlet interior problem; this is the standard BDDC algorithm
41787497f52SBarry Smith -  `PC_BDDC_INTERFACE_EXT_LUMP` - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP"
418bc960bbfSJed Brown 
41916a05f60SBarry Smith    Level: intermediate
42016a05f60SBarry Smith 
42109e53591SBarry Smith .seealso: [](sec_pc), `PCBDDC`, `PC`
422bc960bbfSJed Brown E*/
423bc960bbfSJed Brown typedef enum {
424bc960bbfSJed Brown   PC_BDDC_INTERFACE_EXT_DIRICHLET,
425bc960bbfSJed Brown   PC_BDDC_INTERFACE_EXT_LUMP
426bc960bbfSJed Brown } PCBDDCInterfaceExtType;
427bc960bbfSJed Brown 
428bc960bbfSJed Brown /*E
429f3b08a26SMatthew G. Knepley   PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation
430f3b08a26SMatthew G. Knepley 
431f3b08a26SMatthew G. Knepley   Level: beginner
432f3b08a26SMatthew G. Knepley 
43309e53591SBarry Smith .seealso: [](sec_pc), `PCMGSetAdaptCoarseSpaceType()`, `PCMG`, `PC`
434f3b08a26SMatthew G. Knepley E*/
4359371c9d4SSatish Balay typedef enum {
4369371c9d4SSatish Balay   PCMG_ADAPT_NONE,
4379371c9d4SSatish Balay   PCMG_ADAPT_POLYNOMIAL,
4389371c9d4SSatish Balay   PCMG_ADAPT_HARMONIC,
4399371c9d4SSatish Balay   PCMG_ADAPT_EIGENVECTOR,
4409371c9d4SSatish Balay   PCMG_ADAPT_GENERALIZED_EIGENVECTOR,
4419371c9d4SSatish Balay   PCMG_ADAPT_GDSW
4429371c9d4SSatish Balay } PCMGCoarseSpaceType;
443f3b08a26SMatthew G. Knepley 
444f3b08a26SMatthew G. Knepley /*E
44516a05f60SBarry Smith     PCPatchConstructType - The algorithm used to construct patches for the `PCPATCH` preconditioner
4464bbf5ea8SMatthew G. Knepley 
4474bbf5ea8SMatthew G. Knepley    Level: beginner
4484bbf5ea8SMatthew G. Knepley 
44909e53591SBarry Smith .seealso: [](sec_pc), `PCPatchSetConstructType()`, `PCPATCH`, `PC`
4504bbf5ea8SMatthew G. Knepley E*/
4519371c9d4SSatish Balay typedef enum {
4529371c9d4SSatish Balay   PC_PATCH_STAR,
4539371c9d4SSatish Balay   PC_PATCH_VANKA,
4549371c9d4SSatish Balay   PC_PATCH_PARDECOMP,
4559371c9d4SSatish Balay   PC_PATCH_USER,
4569371c9d4SSatish Balay   PC_PATCH_PYTHON
4579371c9d4SSatish Balay } PCPatchConstructType;
4584bbf5ea8SMatthew G. Knepley 
4594bbf5ea8SMatthew G. Knepley /*E
460e26ad46dSJakub Kruzik     PCDeflationSpaceType - Type of deflation
461e662bc50SJakub Kruzik 
462e662bc50SJakub Kruzik     Values:
46387497f52SBarry Smith +   `PC_DEFLATION_SPACE_HAAR`        - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off
46409e53591SBarry Smith .   `PC_DEFLATION_SPACE_DB2`         - `MATCOMPOSITE` of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet)
46587497f52SBarry Smith .   `PC_DEFLATION_SPACE_DB4`         - same as above, but with db4 (4 coefficient Daubechies)
46687497f52SBarry Smith .   `PC_DEFLATION_SPACE_DB8`         - same as above, but with db8 (8 coefficient Daubechies)
46787497f52SBarry Smith .   `PC_DEFLATION_SPACE_DB16`        - same as above, but with db16 (16 coefficient Daubechies)
46887497f52SBarry Smith .   `PC_DEFLATION_SPACE_BIORTH22`    - same as above, but with biorthogonal 2.2 (6 coefficients)
46987497f52SBarry Smith .   `PC_DEFLATION_SPACE_MEYER`       - same as above, but with Meyer/FIR (62 coefficients)
470d5b43468SJose E. Roman .   `PC_DEFLATION_SPACE_AGGREGATION` - aggregates local indices (given by operator matrix distribution) into a subdomain
47187497f52SBarry Smith -   `PC_DEFLATION_SPACE_USER`        - indicates space set by user
472e662bc50SJakub Kruzik 
4731fdb00f9SJakub Kruzik     Level: intermediate
4741fdb00f9SJakub Kruzik 
47509e53591SBarry Smith     Note:
47609e53591SBarry Smith     Wavelet-based space (except Haar) can be used in multilevel deflation.
47709e53591SBarry Smith 
47809e53591SBarry Smith .seealso: [](sec_pc), `PCDeflationSetSpaceToCompute()`, `PCDEFLATION`, `PC`
479e662bc50SJakub Kruzik E*/
480e662bc50SJakub Kruzik typedef enum {
481e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_HAAR,
482e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB2,
483e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB4,
484e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB8,
485e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB16,
486e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_BIORTH22,
487e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_MEYER,
488e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_AGGREGATION,
489e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_USER
490e662bc50SJakub Kruzik } PCDeflationSpaceType;
491e662bc50SJakub Kruzik 
492e662bc50SJakub Kruzik /*E
49387497f52SBarry Smith     PCHPDDMCoarseCorrectionType - Type of coarse correction used by `PCHPDDM`
49438cfc46eSPierre Jolivet 
49538cfc46eSPierre Jolivet     Values:
49609e53591SBarry Smith +   `PC_HPDDM_COARSE_CORRECTION_DEFLATED` (default) - eq. (1) in `PCHPDDMShellApply()`
49787497f52SBarry Smith .   `PC_HPDDM_COARSE_CORRECTION_ADDITIVE` - eq. (2)
49887497f52SBarry Smith -   `PC_HPDDM_COARSE_CORRECTION_BALANCED` - eq. (3)
49938cfc46eSPierre Jolivet 
50016a05f60SBarry Smith     Level: intermediate
50116a05f60SBarry Smith 
50209e53591SBarry Smith .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCSetType()`, `PCHPDDMShellApply()`
50338cfc46eSPierre Jolivet E*/
5049371c9d4SSatish Balay typedef enum {
5059371c9d4SSatish Balay   PC_HPDDM_COARSE_CORRECTION_DEFLATED,
5069371c9d4SSatish Balay   PC_HPDDM_COARSE_CORRECTION_ADDITIVE,
5079371c9d4SSatish Balay   PC_HPDDM_COARSE_CORRECTION_BALANCED
5089371c9d4SSatish Balay } PCHPDDMCoarseCorrectionType;
50938cfc46eSPierre Jolivet 
51038cfc46eSPierre Jolivet /*E
51187497f52SBarry Smith     PCFailedReason - indicates type of `PC` failure
5128c1cd74cSHong Zhang 
5138c1cd74cSHong Zhang     Level: beginner
5148c1cd74cSHong Zhang 
51587497f52SBarry Smith     Developer Note:
51616a05f60SBarry Smith     Any additions/changes here MUST also be made in `include/petsc/finclude/petscpc.h`
51709e53591SBarry Smith 
51809e53591SBarry Smith .seealso: [](sec_pc), `PC`
5198c1cd74cSHong Zhang E*/
5209371c9d4SSatish Balay typedef enum {
5219371c9d4SSatish Balay   PC_SETUP_ERROR = -1,
5229371c9d4SSatish Balay   PC_NOERROR,
5239371c9d4SSatish Balay   PC_FACTOR_STRUCT_ZEROPIVOT,
5249371c9d4SSatish Balay   PC_FACTOR_NUMERIC_ZEROPIVOT,
5259371c9d4SSatish Balay   PC_FACTOR_OUTMEMORY,
5269371c9d4SSatish Balay   PC_FACTOR_OTHER,
52787b47708SBarry Smith   PC_INCONSISTENT_RHS,
5289371c9d4SSatish Balay   PC_SUBPC_ERROR
5299371c9d4SSatish Balay } PCFailedReason;
530ce7c7f2fSMark Adams 
531ce7c7f2fSMark Adams /*E
532ce7c7f2fSMark Adams     PCGAMGLayoutType - Layout for reduced grids
533ce7c7f2fSMark Adams 
534ce7c7f2fSMark Adams     Level: intermediate
535ce7c7f2fSMark Adams 
53609e53591SBarry Smith     Developer Note:
53716a05f60SBarry Smith     Any additions/changes here MUST also be made in `include/petsc/finclude/petscpc.h`
53809e53591SBarry Smith 
53909e53591SBarry Smith .seealso: [](sec_pc), `PCGAMG`, `PC`, `PCGAMGSetCoarseGridLayoutType()`
540ce7c7f2fSMark Adams E*/
5419371c9d4SSatish Balay typedef enum {
5429371c9d4SSatish Balay   PCGAMG_LAYOUT_COMPACT,
5439371c9d4SSatish Balay   PCGAMG_LAYOUT_SPREAD
5449371c9d4SSatish Balay } PCGAMGLayoutType;
545ce7c7f2fSMark Adams 
546b0753f9dSMatthew G. Knepley #endif
547