xref: /petsc/include/petscpctypes.h (revision d5b43468fb8780a8feea140ccd6fa3e6a50411cc)
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 
11db781477SPatrick Sanan .seealso: `PCCreate()`, `PCSetType()`, `PCType`
12b0753f9dSMatthew G. Knepley S*/
13b0753f9dSMatthew G. Knepley typedef struct _p_PC *PC;
14b0753f9dSMatthew G. Knepley 
15b0753f9dSMatthew G. Knepley /*J
16b0753f9dSMatthew G. Knepley     PCType - String with the name of a PETSc preconditioner method.
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 
2387497f52SBarry Smith .seealso: `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 
85b0753f9dSMatthew G. Knepley    Level: beginner
86b0753f9dSMatthew G. Knepley 
87b0753f9dSMatthew G. Knepley .seealso:
88b0753f9dSMatthew G. Knepley E*/
899371c9d4SSatish Balay typedef enum {
909371c9d4SSatish Balay   PC_SIDE_DEFAULT = -1,
919371c9d4SSatish Balay   PC_LEFT,
929371c9d4SSatish Balay   PC_RIGHT,
939371c9d4SSatish Balay   PC_SYMMETRIC
949371c9d4SSatish Balay } PCSide;
95b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
96b0753f9dSMatthew G. Knepley 
97b0753f9dSMatthew G. Knepley /*E
9887497f52SBarry Smith     PCRichardsonConvergedReason - reason a `PCApplyRichardson() method terminated
99b0753f9dSMatthew G. Knepley 
100b0753f9dSMatthew G. Knepley    Level: advanced
101b0753f9dSMatthew G. Knepley 
10287497f52SBarry Smith    Developer Note:
10387497f52SBarry Smith     this must match petsc/finclude/petscpc.h and the `KSPConvergedReason` values in petscksp.h
104b0753f9dSMatthew G. Knepley 
105db781477SPatrick Sanan .seealso: `PCApplyRichardson()`
106b0753f9dSMatthew G. Knepley E*/
107b0753f9dSMatthew G. Knepley typedef enum {
108b0753f9dSMatthew G. Knepley   PCRICHARDSON_CONVERGED_RTOL = 2,
109b0753f9dSMatthew G. Knepley   PCRICHARDSON_CONVERGED_ATOL = 3,
110b0753f9dSMatthew G. Knepley   PCRICHARDSON_CONVERGED_ITS  = 4,
1119371c9d4SSatish Balay   PCRICHARDSON_DIVERGED_DTOL  = -4
1129371c9d4SSatish Balay } PCRichardsonConvergedReason;
113b0753f9dSMatthew G. Knepley 
114b0753f9dSMatthew G. Knepley /*E
115b0753f9dSMatthew G. Knepley     PCJacobiType - What elements are used to form the Jacobi preconditioner
116b0753f9dSMatthew G. Knepley 
117b0753f9dSMatthew G. Knepley    Level: intermediate
118b0753f9dSMatthew G. Knepley 
119b0753f9dSMatthew G. Knepley .seealso:
120b0753f9dSMatthew G. Knepley E*/
1219371c9d4SSatish Balay typedef enum {
1229371c9d4SSatish Balay   PC_JACOBI_DIAGONAL,
1239371c9d4SSatish Balay   PC_JACOBI_ROWMAX,
1249371c9d4SSatish Balay   PC_JACOBI_ROWSUM
1259371c9d4SSatish Balay } PCJacobiType;
126b0753f9dSMatthew G. Knepley 
127b0753f9dSMatthew G. Knepley /*E
128b0753f9dSMatthew G. Knepley     PCASMType - Type of additive Schwarz method to use
129b0753f9dSMatthew G. Knepley 
13087497f52SBarry Smith $  `PC_ASM_BASIC`        - Symmetric version where residuals from the ghost points are used
131b0753f9dSMatthew G. Knepley $                        and computed values in ghost regions are added together.
132b0753f9dSMatthew G. Knepley $                        Classical standard additive Schwarz.
13387497f52SBarry Smith $  `PC_ASM_RESTRICT`     - Residuals from ghost points are used but computed values in ghost
134b0753f9dSMatthew G. Knepley $                        region are discarded.
135b0753f9dSMatthew G. Knepley $                        Default.
13687497f52SBarry Smith $  `PC_ASM_INTERPOLATE`  - Residuals from ghost points are not used, computed values in ghost
137b0753f9dSMatthew G. Knepley $                        region are added back in.
13887497f52SBarry Smith $  `PC_ASM_NONE`         - Residuals from ghost points are not used, computed ghost values are
139b0753f9dSMatthew G. Knepley $                        discarded.
140b0753f9dSMatthew G. Knepley $                        Not very good.
141b0753f9dSMatthew G. Knepley 
142b0753f9dSMatthew G. Knepley    Level: beginner
143b0753f9dSMatthew G. Knepley 
144db781477SPatrick Sanan .seealso: `PCASMSetType()`
145b0753f9dSMatthew G. Knepley E*/
1469371c9d4SSatish Balay typedef enum {
1479371c9d4SSatish Balay   PC_ASM_BASIC       = 3,
1489371c9d4SSatish Balay   PC_ASM_RESTRICT    = 1,
1499371c9d4SSatish Balay   PC_ASM_INTERPOLATE = 2,
1509371c9d4SSatish Balay   PC_ASM_NONE        = 0
1519371c9d4SSatish Balay } PCASMType;
152b0753f9dSMatthew G. Knepley 
153b0753f9dSMatthew G. Knepley /*E
15487497f52SBarry Smith     PCGASMType - Type of generalized additive Schwarz method to use (differs from `PCASM` in allowing multiple processors per subdomain).
155b0753f9dSMatthew G. Knepley 
156b0753f9dSMatthew G. Knepley    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
157b0753f9dSMatthew G. Knepley    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
158b0753f9dSMatthew G. Knepley    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
159b0753f9dSMatthew G. Knepley    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
160b0753f9dSMatthew G. Knepley 
16187497f52SBarry Smith $  `PC_GASM_BASIC`       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
162b0753f9dSMatthew G. Knepley $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
163b0753f9dSMatthew G. Knepley $                        from neighboring subdomains.
164b0753f9dSMatthew G. Knepley $                        Classical standard additive Schwarz.
16587497f52SBarry Smith $  `PC_GASM_RESTRICT`    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
166b0753f9dSMatthew G. Knepley $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
167b0753f9dSMatthew G. Knepley $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
168b0753f9dSMatthew G. Knepley $                        assumption).
169b0753f9dSMatthew G. Knepley $                        Default.
17087497f52SBarry Smith $  `PC_GASM_INTERPOLATE` - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
171b0753f9dSMatthew G. Knepley $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
172b0753f9dSMatthew G. Knepley $                        from neighboring subdomains.
173b0753f9dSMatthew G. Knepley $
17487497f52SBarry Smith $  `PC_GASM_NONE`       - Residuals and corrections are zeroed out outside the local subdomains.
175b0753f9dSMatthew G. Knepley $                        Not very good.
176b0753f9dSMatthew G. Knepley 
177b0753f9dSMatthew G. Knepley    Level: beginner
178b0753f9dSMatthew G. Knepley 
179db781477SPatrick Sanan .seealso: `PCGASMSetType()`
180b0753f9dSMatthew G. Knepley E*/
1819371c9d4SSatish Balay typedef enum {
1829371c9d4SSatish Balay   PC_GASM_BASIC       = 3,
1839371c9d4SSatish Balay   PC_GASM_RESTRICT    = 1,
1849371c9d4SSatish Balay   PC_GASM_INTERPOLATE = 2,
1859371c9d4SSatish Balay   PC_GASM_NONE        = 0
1869371c9d4SSatish Balay } PCGASMType;
187b0753f9dSMatthew G. Knepley 
188b0753f9dSMatthew G. Knepley /*E
189b0753f9dSMatthew G. Knepley     PCCompositeType - Determines how two or more preconditioner are composed
190b0753f9dSMatthew G. Knepley 
19187497f52SBarry Smith $  `PC_COMPOSITE_ADDITIVE` - results from application of all preconditioners are added together
19287497f52SBarry Smith $  `PC_COMPOSITE_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly
193b0753f9dSMatthew G. Knepley $                                computed after the previous preconditioner application
19487497f52SBarry Smith $  `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly
195d0ecd4dfSBarry Smith $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners)
19687497f52SBarry Smith $  `PC_COMPOSITE_SPECIAL` - This is very special for a matrix of the form alpha I + R + S
197b0753f9dSMatthew G. Knepley $                         where first preconditioner is built from alpha I + S and second from
198b0753f9dSMatthew G. Knepley $                         alpha I + R
199b0753f9dSMatthew G. Knepley 
200b0753f9dSMatthew G. Knepley    Level: beginner
201b0753f9dSMatthew G. Knepley 
202db781477SPatrick Sanan .seealso: `PCCompositeSetType()`
203b0753f9dSMatthew G. Knepley E*/
2049371c9d4SSatish Balay typedef enum {
2059371c9d4SSatish Balay   PC_COMPOSITE_ADDITIVE,
2069371c9d4SSatish Balay   PC_COMPOSITE_MULTIPLICATIVE,
2079371c9d4SSatish Balay   PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,
2089371c9d4SSatish Balay   PC_COMPOSITE_SPECIAL,
2099371c9d4SSatish Balay   PC_COMPOSITE_SCHUR,
2109371c9d4SSatish Balay   PC_COMPOSITE_GKB
2119371c9d4SSatish Balay } PCCompositeType;
212b0753f9dSMatthew G. Knepley 
213b0753f9dSMatthew G. Knepley /*E
21487497f52SBarry Smith     PCFieldSplitSchurPreType - Determines how to precondition a Schur complement
215b0753f9dSMatthew G. Knepley 
216b0753f9dSMatthew G. Knepley     Level: intermediate
217b0753f9dSMatthew G. Knepley 
218db781477SPatrick Sanan .seealso: `PCFieldSplitSetSchurPre()`
219b0753f9dSMatthew G. Knepley E*/
2209371c9d4SSatish Balay typedef enum {
2219371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_SELF,
2229371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_SELFP,
2239371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_A11,
2249371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_USER,
2259371c9d4SSatish Balay   PC_FIELDSPLIT_SCHUR_PRE_FULL
2269371c9d4SSatish Balay } PCFieldSplitSchurPreType;
227b0753f9dSMatthew G. Knepley 
228b0753f9dSMatthew G. Knepley /*E
229b0753f9dSMatthew G. Knepley     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
230b0753f9dSMatthew G. Knepley 
231b0753f9dSMatthew G. Knepley     Level: intermediate
232b0753f9dSMatthew G. Knepley 
233db781477SPatrick Sanan .seealso: `PCFieldSplitSetSchurFactType()`
234b0753f9dSMatthew G. Knepley E*/
235b0753f9dSMatthew G. Knepley typedef enum {
236b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
237b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
238b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
239b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_FULL
240b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType;
241b0753f9dSMatthew G. Knepley 
242b0753f9dSMatthew G. Knepley /*E
24387497f52SBarry Smith     PCPARMSGlobalType - Determines the global preconditioner method in `PCPARMS`
244b0753f9dSMatthew G. Knepley 
245b0753f9dSMatthew G. Knepley     Level: intermediate
246b0753f9dSMatthew G. Knepley 
247db781477SPatrick Sanan .seealso: `PCPARMSSetGlobal()`
248b0753f9dSMatthew G. Knepley E*/
2499371c9d4SSatish Balay typedef enum {
2509371c9d4SSatish Balay   PC_PARMS_GLOBAL_RAS,
2519371c9d4SSatish Balay   PC_PARMS_GLOBAL_SCHUR,
2529371c9d4SSatish Balay   PC_PARMS_GLOBAL_BJ
2539371c9d4SSatish Balay } PCPARMSGlobalType;
2549d8081ecSMatthew G. Knepley 
255b0753f9dSMatthew G. Knepley /*E
25687497f52SBarry Smith     PCPARMSLocalType - Determines the local preconditioner method in `PCPARMS`
257b0753f9dSMatthew G. Knepley 
258b0753f9dSMatthew G. Knepley     Level: intermediate
259b0753f9dSMatthew G. Knepley 
260db781477SPatrick Sanan .seealso: `PCPARMSSetLocal()`
261b0753f9dSMatthew G. Knepley E*/
2629371c9d4SSatish Balay typedef enum {
2639371c9d4SSatish Balay   PC_PARMS_LOCAL_ILU0,
2649371c9d4SSatish Balay   PC_PARMS_LOCAL_ILUK,
2659371c9d4SSatish Balay   PC_PARMS_LOCAL_ILUT,
2669371c9d4SSatish Balay   PC_PARMS_LOCAL_ARMS
2679371c9d4SSatish Balay } PCPARMSLocalType;
268b0753f9dSMatthew G. Knepley 
269f0fc11ceSJed Brown /*J
27087497f52SBarry Smith     PCGAMGType - type of generalized algebraic multigrid `PCGAMG` method
271b0753f9dSMatthew G. Knepley 
272b0753f9dSMatthew G. Knepley     Level: intermediate
273b0753f9dSMatthew G. Knepley 
27487497f52SBarry Smith $   `PCGAMGAGG` - (the default) smoothed aggregation algorithm, robust, very well tested
27587497f52SBarry Smith $   `PCGAMGGEO` - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not well tested
27687497f52SBarry Smith $   `PCGAMGCLASSICAL` - classical algebraic multigrid preconditioner, incomplete, poorly tested
277f6130c64SBarry Smith 
278db781477SPatrick Sanan .seealso: `PCMG`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()`
279f0fc11ceSJed Brown J*/
280b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType;
281b0753f9dSMatthew G. Knepley #define PCGAMGAGG       "agg"
282b0753f9dSMatthew G. Knepley #define PCGAMGGEO       "geo"
283b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL "classical"
284b0753f9dSMatthew G. Knepley 
285b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType;
286b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT   "direct"
287b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard"
288b0753f9dSMatthew G. Knepley 
289b0753f9dSMatthew G. Knepley /*E
290b0753f9dSMatthew G. Knepley     PCMGType - Determines the type of multigrid method that is run.
291b0753f9dSMatthew G. Knepley 
292b0753f9dSMatthew G. Knepley    Level: beginner
293b0753f9dSMatthew G. Knepley 
294b0753f9dSMatthew G. Knepley    Values:
29587497f52SBarry Smith +  `PC_MG_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `PCMGSetCycleType()`
29687497f52SBarry Smith .  `PC_MG_ADDITIVE` - the additive multigrid preconditioner where all levels are
297b0753f9dSMatthew G. Knepley                 smoothed before updating the residual. This only uses the
298b0753f9dSMatthew G. Knepley                 down smoother, in the preconditioner the upper smoother is ignored
29987497f52SBarry Smith .  `PC_MG_FULL` - same as multiplicative except one also performs grid sequencing,
300b0753f9dSMatthew G. Knepley             that is starts on the coarsest grid, performs a cycle, interpolates
301b0753f9dSMatthew 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
302b0753f9dSMatthew G. Knepley             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
303b0753f9dSMatthew G. Knepley             calls the V-cycle only on the coarser level and has a post-smoother instead.
30487497f52SBarry Smith -  `PC_MG_KASKADE` - like full multigrid except one never goes back to a coarser level
305b0753f9dSMatthew G. Knepley                from a finer
306b0753f9dSMatthew G. Knepley 
307db781477SPatrick Sanan .seealso: `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()`
308b0753f9dSMatthew G. Knepley 
309b0753f9dSMatthew G. Knepley E*/
3109371c9d4SSatish Balay typedef enum {
3119371c9d4SSatish Balay   PC_MG_MULTIPLICATIVE,
3129371c9d4SSatish Balay   PC_MG_ADDITIVE,
3139371c9d4SSatish Balay   PC_MG_FULL,
3149371c9d4SSatish Balay   PC_MG_KASKADE
3159371c9d4SSatish Balay } PCMGType;
316b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE;
317b0753f9dSMatthew G. Knepley 
318b0753f9dSMatthew G. Knepley /*E
319b0753f9dSMatthew G. Knepley     PCMGCycleType - Use V-cycle or W-cycle
320b0753f9dSMatthew G. Knepley 
321b0753f9dSMatthew G. Knepley    Level: beginner
322b0753f9dSMatthew G. Knepley 
323b0753f9dSMatthew G. Knepley    Values:
32487497f52SBarry Smith +  `PC_MG_V_CYCLE` - use the v cycle
32587497f52SBarry Smith -  `PC_MG_W_CYCLE` - use the w cycle
326b0753f9dSMatthew G. Knepley 
327db781477SPatrick Sanan .seealso: `PCMGSetCycleType()`
328b0753f9dSMatthew G. Knepley 
329b0753f9dSMatthew G. Knepley E*/
3309371c9d4SSatish Balay typedef enum {
3319371c9d4SSatish Balay   PC_MG_CYCLE_V = 1,
3329371c9d4SSatish Balay   PC_MG_CYCLE_W = 2
3339371c9d4SSatish Balay } PCMGCycleType;
334b0753f9dSMatthew G. Knepley 
335b0753f9dSMatthew G. Knepley /*E
3362134b1e4SBarry Smith     PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process
3372134b1e4SBarry Smith 
3382134b1e4SBarry Smith    Level: beginner
3392134b1e4SBarry Smith 
3402134b1e4SBarry Smith    Values:
34187497f52SBarry Smith +  `PC_MG_GALERKIN_PMAT` - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid
34287497f52SBarry Smith .  `PC_MG_GALERKIN_MAT` -  computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid
34387497f52SBarry Smith .  `PC_MG_GALERKIN_BOTH` - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once
34487497f52SBarry Smith -  `PC_MG_GALERKIN_NONE` - neither operator is computed via the Galerkin process, the user must provide the operator
3452134b1e4SBarry Smith 
34687497f52SBarry Smith    Users should never set `PC_MG_GALERKIN_EXTERNAL`, it is used by `PCHYPRE` and `PCML`
3472134b1e4SBarry Smith 
348db781477SPatrick Sanan .seealso: `PCMGSetCycleType()`
3492134b1e4SBarry Smith 
3502134b1e4SBarry Smith E*/
3519371c9d4SSatish Balay typedef enum {
3529371c9d4SSatish Balay   PC_MG_GALERKIN_BOTH,
3539371c9d4SSatish Balay   PC_MG_GALERKIN_PMAT,
3549371c9d4SSatish Balay   PC_MG_GALERKIN_MAT,
3559371c9d4SSatish Balay   PC_MG_GALERKIN_NONE,
3569371c9d4SSatish Balay   PC_MG_GALERKIN_EXTERNAL
3579371c9d4SSatish Balay } PCMGGalerkinType;
3582134b1e4SBarry Smith 
3592134b1e4SBarry Smith /*E
360b0753f9dSMatthew G. Knepley     PCExoticType - Face based or wirebasket based coarse grid space
361b0753f9dSMatthew G. Knepley 
362b0753f9dSMatthew G. Knepley    Level: beginner
363b0753f9dSMatthew G. Knepley 
364db781477SPatrick Sanan .seealso: `PCExoticSetType()`, `PCEXOTIC`
365b0753f9dSMatthew G. Knepley E*/
3669371c9d4SSatish Balay typedef enum {
3679371c9d4SSatish Balay   PC_EXOTIC_FACE,
3689371c9d4SSatish Balay   PC_EXOTIC_WIREBASKET
3699371c9d4SSatish Balay } PCExoticType;
370b0753f9dSMatthew G. Knepley 
3718c1cd74cSHong Zhang /*E
372bc960bbfSJed Brown    PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains.
373bc960bbfSJed Brown 
374bc960bbfSJed Brown    Level: intermediate
375bc960bbfSJed Brown 
376bc960bbfSJed Brown    Values:
37787497f52SBarry Smith +  `PC_BDDC_INTERFACE_EXT_DIRICHLET` - solves Dirichlet interior problem; this is the standard BDDC algorithm
37887497f52SBarry Smith -  `PC_BDDC_INTERFACE_EXT_LUMP` - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP"
379bc960bbfSJed Brown 
380bc960bbfSJed Brown E*/
381bc960bbfSJed Brown typedef enum {
382bc960bbfSJed Brown   PC_BDDC_INTERFACE_EXT_DIRICHLET,
383bc960bbfSJed Brown   PC_BDDC_INTERFACE_EXT_LUMP
384bc960bbfSJed Brown } PCBDDCInterfaceExtType;
385bc960bbfSJed Brown 
386bc960bbfSJed Brown /*E
387f3b08a26SMatthew G. Knepley   PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation
388f3b08a26SMatthew G. Knepley 
389f3b08a26SMatthew G. Knepley   Level: beginner
390f3b08a26SMatthew G. Knepley 
391db781477SPatrick Sanan .seealso: `PCMGSetAdaptCoarseSpaceType()`, `PCMG`
392f3b08a26SMatthew G. Knepley E*/
3939371c9d4SSatish Balay typedef enum {
3949371c9d4SSatish Balay   PCMG_ADAPT_NONE,
3959371c9d4SSatish Balay   PCMG_ADAPT_POLYNOMIAL,
3969371c9d4SSatish Balay   PCMG_ADAPT_HARMONIC,
3979371c9d4SSatish Balay   PCMG_ADAPT_EIGENVECTOR,
3989371c9d4SSatish Balay   PCMG_ADAPT_GENERALIZED_EIGENVECTOR,
3999371c9d4SSatish Balay   PCMG_ADAPT_GDSW
4009371c9d4SSatish Balay } PCMGCoarseSpaceType;
401f3b08a26SMatthew G. Knepley 
402f3b08a26SMatthew G. Knepley /*E
403ee6acaa3SMatthew G. Knepley     PCPatchConstructType - The algorithm used to construct patches for the preconditioner
4044bbf5ea8SMatthew G. Knepley 
4054bbf5ea8SMatthew G. Knepley    Level: beginner
4064bbf5ea8SMatthew G. Knepley 
407db781477SPatrick Sanan .seealso: `PCPatchSetConstructType()`, `PCEXOTIC`
4084bbf5ea8SMatthew G. Knepley E*/
4099371c9d4SSatish Balay typedef enum {
4109371c9d4SSatish Balay   PC_PATCH_STAR,
4119371c9d4SSatish Balay   PC_PATCH_VANKA,
4129371c9d4SSatish Balay   PC_PATCH_PARDECOMP,
4139371c9d4SSatish Balay   PC_PATCH_USER,
4149371c9d4SSatish Balay   PC_PATCH_PYTHON
4159371c9d4SSatish Balay } PCPatchConstructType;
4164bbf5ea8SMatthew G. Knepley 
4174bbf5ea8SMatthew G. Knepley /*E
418e26ad46dSJakub Kruzik     PCDeflationSpaceType - Type of deflation
419e662bc50SJakub Kruzik 
420e662bc50SJakub Kruzik     Values:
42187497f52SBarry Smith +   `PC_DEFLATION_SPACE_HAAR`        - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off
42287497f52SBarry Smith .   `PC_DEFLATION_SPACE_DB2`         - MATCOMPOSITE of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet)
42387497f52SBarry Smith .   `PC_DEFLATION_SPACE_DB4`         - same as above, but with db4 (4 coefficient Daubechies)
42487497f52SBarry Smith .   `PC_DEFLATION_SPACE_DB8`         - same as above, but with db8 (8 coefficient Daubechies)
42587497f52SBarry Smith .   `PC_DEFLATION_SPACE_DB16`        - same as above, but with db16 (16 coefficient Daubechies)
42687497f52SBarry Smith .   `PC_DEFLATION_SPACE_BIORTH22`    - same as above, but with biorthogonal 2.2 (6 coefficients)
42787497f52SBarry Smith .   `PC_DEFLATION_SPACE_MEYER`       - same as above, but with Meyer/FIR (62 coefficients)
428*d5b43468SJose E. Roman .   `PC_DEFLATION_SPACE_AGGREGATION` - aggregates local indices (given by operator matrix distribution) into a subdomain
42987497f52SBarry Smith -   `PC_DEFLATION_SPACE_USER`        - indicates space set by user
430e662bc50SJakub Kruzik 
431e662bc50SJakub Kruzik     Notes:
4327e9012c5SJakub Kruzik       Wavelet-based space (except Haar) can be used in multilevel deflation.
433e662bc50SJakub Kruzik 
4341fdb00f9SJakub Kruzik     Level: intermediate
4351fdb00f9SJakub Kruzik 
436db781477SPatrick Sanan .seealso: `PCDeflationSetSpaceToCompute()`, `PCDEFLATION`
437e662bc50SJakub Kruzik E*/
438e662bc50SJakub Kruzik typedef enum {
439e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_HAAR,
440e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB2,
441e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB4,
442e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB8,
443e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB16,
444e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_BIORTH22,
445e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_MEYER,
446e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_AGGREGATION,
447e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_USER
448e662bc50SJakub Kruzik } PCDeflationSpaceType;
449e662bc50SJakub Kruzik 
450e662bc50SJakub Kruzik /*E
45187497f52SBarry Smith     PCHPDDMCoarseCorrectionType - Type of coarse correction used by `PCHPDDM`
45238cfc46eSPierre Jolivet 
45338cfc46eSPierre Jolivet     Level: intermediate
45438cfc46eSPierre Jolivet 
45538cfc46eSPierre Jolivet     Values:
45687497f52SBarry Smith +   `PC_HPDDM_COARSE_CORRECTION_DEFLATED` (default) - eq. (1) in PCHPDDMShellApply()
45787497f52SBarry Smith .   `PC_HPDDM_COARSE_CORRECTION_ADDITIVE` - eq. (2)
45887497f52SBarry Smith -   `PC_HPDDM_COARSE_CORRECTION_BALANCED` - eq. (3)
45938cfc46eSPierre Jolivet 
460db781477SPatrick Sanan .seealso: `PCHPDDM`, `PCSetType()`, `PCHPDDMShellApply()`
46138cfc46eSPierre Jolivet E*/
4629371c9d4SSatish Balay typedef enum {
4639371c9d4SSatish Balay   PC_HPDDM_COARSE_CORRECTION_DEFLATED,
4649371c9d4SSatish Balay   PC_HPDDM_COARSE_CORRECTION_ADDITIVE,
4659371c9d4SSatish Balay   PC_HPDDM_COARSE_CORRECTION_BALANCED
4669371c9d4SSatish Balay } PCHPDDMCoarseCorrectionType;
46738cfc46eSPierre Jolivet 
46838cfc46eSPierre Jolivet /*E
46987497f52SBarry Smith     PCFailedReason - indicates type of `PC` failure
4708c1cd74cSHong Zhang 
4718c1cd74cSHong Zhang     Level: beginner
4728c1cd74cSHong Zhang 
47387497f52SBarry Smith     Developer Note:
4748c1cd74cSHong Zhang     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
4758c1cd74cSHong Zhang E*/
4769371c9d4SSatish Balay typedef enum {
4779371c9d4SSatish Balay   PC_SETUP_ERROR = -1,
4789371c9d4SSatish Balay   PC_NOERROR,
4799371c9d4SSatish Balay   PC_FACTOR_STRUCT_ZEROPIVOT,
4809371c9d4SSatish Balay   PC_FACTOR_NUMERIC_ZEROPIVOT,
4819371c9d4SSatish Balay   PC_FACTOR_OUTMEMORY,
4829371c9d4SSatish Balay   PC_FACTOR_OTHER,
48387b47708SBarry Smith   PC_INCONSISTENT_RHS,
4849371c9d4SSatish Balay   PC_SUBPC_ERROR
4859371c9d4SSatish Balay } PCFailedReason;
486ce7c7f2fSMark Adams 
487ce7c7f2fSMark Adams /*E
488ce7c7f2fSMark Adams     PCGAMGLayoutType - Layout for reduced grids
489ce7c7f2fSMark Adams 
490ce7c7f2fSMark Adams     Level: intermediate
491ce7c7f2fSMark Adams 
492db781477SPatrick Sanan .seealso: `PCGAMGSetCoarseGridLayoutType()`
4939969d477SPatrick Sanan 
494ce7c7f2fSMark Adams     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
495ce7c7f2fSMark Adams E*/
4969371c9d4SSatish Balay typedef enum {
4979371c9d4SSatish Balay   PCGAMG_LAYOUT_COMPACT,
4989371c9d4SSatish Balay   PCGAMG_LAYOUT_SPREAD
4999371c9d4SSatish Balay } PCGAMGLayoutType;
500ce7c7f2fSMark Adams 
501b0753f9dSMatthew G. Knepley #endif
502