xref: /petsc/include/petscpctypes.h (revision ac09b9214d23ea9ad238aa607de9fa447fd4e91b)
126bd1501SBarry Smith #if !defined(PETSCPCTYPES_H)
226bd1501SBarry Smith #define PETSCPCTYPES_H
3b0753f9dSMatthew G. Knepley 
4*ac09b921SBarry Smith /* SUBMANSEC = PC */
5*ac09b921SBarry Smith 
6b0753f9dSMatthew G. Knepley /*S
7b0753f9dSMatthew G. Knepley      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 
2095452b02SPatrick Sanan    Notes:
2195452b02SPatrick Sanan     Click on the links above to see details on a particular solver
22b0753f9dSMatthew G. Knepley 
23b0753f9dSMatthew G. Knepley           PCRegister() is used to register preconditioners that are then accessible via PCSetType()
24b0753f9dSMatthew G. Knepley 
25db781477SPatrick Sanan .seealso: `PCSetType()`, `PC`, `PCCreate()`, `PCRegister()`, `PCSetFromOptions()`
26b0753f9dSMatthew G. Knepley J*/
27b0753f9dSMatthew G. Knepley typedef const char* PCType;
28b0753f9dSMatthew G. Knepley #define PCNONE            "none"
29b0753f9dSMatthew G. Knepley #define PCJACOBI          "jacobi"
30b0753f9dSMatthew G. Knepley #define PCSOR             "sor"
31b0753f9dSMatthew G. Knepley #define PCLU              "lu"
32a2fc1e05SToby Isaac #define PCQR              "qr"
33b0753f9dSMatthew G. Knepley #define PCSHELL           "shell"
34b0753f9dSMatthew G. Knepley #define PCBJACOBI         "bjacobi"
35b0753f9dSMatthew G. Knepley #define PCMG              "mg"
36b0753f9dSMatthew G. Knepley #define PCEISENSTAT       "eisenstat"
37b0753f9dSMatthew G. Knepley #define PCILU             "ilu"
38b0753f9dSMatthew G. Knepley #define PCICC             "icc"
39b0753f9dSMatthew G. Knepley #define PCASM             "asm"
40b0753f9dSMatthew G. Knepley #define PCGASM            "gasm"
41b0753f9dSMatthew G. Knepley #define PCKSP             "ksp"
42e607c864SMark Adams #define PCBJKOKKOS        "bjkokkos"
43b0753f9dSMatthew G. Knepley #define PCCOMPOSITE       "composite"
44b0753f9dSMatthew G. Knepley #define PCREDUNDANT       "redundant"
45b0753f9dSMatthew G. Knepley #define PCSPAI            "spai"
46b0753f9dSMatthew G. Knepley #define PCNN              "nn"
47b0753f9dSMatthew G. Knepley #define PCCHOLESKY        "cholesky"
48b0753f9dSMatthew G. Knepley #define PCPBJACOBI        "pbjacobi"
490da83c2eSBarry Smith #define PCVPBJACOBI       "vpbjacobi"
50b0753f9dSMatthew G. Knepley #define PCMAT             "mat"
51b0753f9dSMatthew G. Knepley #define PCHYPRE           "hypre"
52b0753f9dSMatthew G. Knepley #define PCPARMS           "parms"
53b0753f9dSMatthew G. Knepley #define PCFIELDSPLIT      "fieldsplit"
54b0753f9dSMatthew G. Knepley #define PCTFS             "tfs"
55b0753f9dSMatthew G. Knepley #define PCML              "ml"
56b0753f9dSMatthew G. Knepley #define PCGALERKIN        "galerkin"
57b0753f9dSMatthew G. Knepley #define PCEXOTIC          "exotic"
58b0753f9dSMatthew G. Knepley #define PCCP              "cp"
59b0753f9dSMatthew G. Knepley #define PCBFBT            "bfbt"
60b0753f9dSMatthew G. Knepley #define PCLSC             "lsc"
61b0753f9dSMatthew G. Knepley #define PCPYTHON          "python"
62b0753f9dSMatthew G. Knepley #define PCPFMG            "pfmg"
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"
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 
84b0753f9dSMatthew G. Knepley    Level: beginner
85b0753f9dSMatthew G. Knepley 
86b0753f9dSMatthew G. Knepley .seealso:
87b0753f9dSMatthew G. Knepley E*/
88b0753f9dSMatthew G. Knepley typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
89b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
90b0753f9dSMatthew G. Knepley 
91b0753f9dSMatthew G. Knepley /*E
92b0753f9dSMatthew G. Knepley     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
93b0753f9dSMatthew G. Knepley 
94b0753f9dSMatthew G. Knepley    Level: advanced
95b0753f9dSMatthew G. Knepley 
9695452b02SPatrick Sanan    Notes:
9795452b02SPatrick Sanan     this must match petsc/finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
98b0753f9dSMatthew G. Knepley 
99db781477SPatrick Sanan .seealso: `PCApplyRichardson()`
100b0753f9dSMatthew G. Knepley E*/
101b0753f9dSMatthew G. Knepley typedef enum {
102b0753f9dSMatthew G. Knepley               PCRICHARDSON_CONVERGED_RTOL               =  2,
103b0753f9dSMatthew G. Knepley               PCRICHARDSON_CONVERGED_ATOL               =  3,
104b0753f9dSMatthew G. Knepley               PCRICHARDSON_CONVERGED_ITS                =  4,
105b0753f9dSMatthew G. Knepley               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
106b0753f9dSMatthew G. Knepley 
107b0753f9dSMatthew G. Knepley /*E
108b0753f9dSMatthew G. Knepley     PCJacobiType - What elements are used to form the Jacobi preconditioner
109b0753f9dSMatthew G. Knepley 
110b0753f9dSMatthew G. Knepley    Level: intermediate
111b0753f9dSMatthew G. Knepley 
112b0753f9dSMatthew G. Knepley .seealso:
113b0753f9dSMatthew G. Knepley E*/
114b0753f9dSMatthew G. Knepley typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType;
115b0753f9dSMatthew G. Knepley 
116b0753f9dSMatthew G. Knepley /*E
117b0753f9dSMatthew G. Knepley     PCASMType - Type of additive Schwarz method to use
118b0753f9dSMatthew G. Knepley 
119b0753f9dSMatthew G. Knepley $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
120b0753f9dSMatthew G. Knepley $                        and computed values in ghost regions are added together.
121b0753f9dSMatthew G. Knepley $                        Classical standard additive Schwarz.
122b0753f9dSMatthew G. Knepley $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
123b0753f9dSMatthew G. Knepley $                        region are discarded.
124b0753f9dSMatthew G. Knepley $                        Default.
125b0753f9dSMatthew G. Knepley $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
126b0753f9dSMatthew G. Knepley $                        region are added back in.
127b0753f9dSMatthew G. Knepley $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
128b0753f9dSMatthew G. Knepley $                        discarded.
129b0753f9dSMatthew G. Knepley $                        Not very good.
130b0753f9dSMatthew G. Knepley 
131b0753f9dSMatthew G. Knepley    Level: beginner
132b0753f9dSMatthew G. Knepley 
133db781477SPatrick Sanan .seealso: `PCASMSetType()`
134b0753f9dSMatthew G. Knepley E*/
135b0753f9dSMatthew G. Knepley typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
136b0753f9dSMatthew G. Knepley 
137b0753f9dSMatthew G. Knepley /*E
138b0753f9dSMatthew G. Knepley     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
139b0753f9dSMatthew G. Knepley 
140b0753f9dSMatthew G. Knepley    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
141b0753f9dSMatthew G. Knepley    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
142b0753f9dSMatthew G. Knepley    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
143b0753f9dSMatthew G. Knepley    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
144b0753f9dSMatthew G. Knepley 
145b0753f9dSMatthew G. Knepley $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
146b0753f9dSMatthew G. Knepley $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
147b0753f9dSMatthew G. Knepley $                        from neighboring subdomains.
148b0753f9dSMatthew G. Knepley $                        Classical standard additive Schwarz.
149b0753f9dSMatthew G. Knepley $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
150b0753f9dSMatthew G. Knepley $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
151b0753f9dSMatthew G. Knepley $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
152b0753f9dSMatthew G. Knepley $                        assumption).
153b0753f9dSMatthew G. Knepley $                        Default.
154b0753f9dSMatthew G. Knepley $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
155b0753f9dSMatthew G. Knepley $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
156b0753f9dSMatthew G. Knepley $                        from neighboring subdomains.
157b0753f9dSMatthew G. Knepley $
158b0753f9dSMatthew G. Knepley $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
159b0753f9dSMatthew G. Knepley $                        Not very good.
160b0753f9dSMatthew G. Knepley 
161b0753f9dSMatthew G. Knepley    Level: beginner
162b0753f9dSMatthew G. Knepley 
163db781477SPatrick Sanan .seealso: `PCGASMSetType()`
164b0753f9dSMatthew G. Knepley E*/
165b0753f9dSMatthew G. Knepley typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
166b0753f9dSMatthew G. Knepley 
167b0753f9dSMatthew G. Knepley /*E
168b0753f9dSMatthew G. Knepley     PCCompositeType - Determines how two or more preconditioner are composed
169b0753f9dSMatthew G. Knepley 
170b0753f9dSMatthew G. Knepley $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
171b0753f9dSMatthew G. Knepley $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
172b0753f9dSMatthew G. Knepley $                                computed after the previous preconditioner application
173b0753f9dSMatthew G. Knepley $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
174d0ecd4dfSBarry Smith $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners)
175b0753f9dSMatthew G. Knepley $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
176b0753f9dSMatthew G. Knepley $                         where first preconditioner is built from alpha I + S and second from
177b0753f9dSMatthew G. Knepley $                         alpha I + R
178b0753f9dSMatthew G. Knepley 
179b0753f9dSMatthew G. Knepley    Level: beginner
180b0753f9dSMatthew G. Knepley 
181db781477SPatrick Sanan .seealso: `PCCompositeSetType()`
182b0753f9dSMatthew G. Knepley E*/
183a51937d4SCarola Kruse typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR,PC_COMPOSITE_GKB} PCCompositeType;
184b0753f9dSMatthew G. Knepley 
185b0753f9dSMatthew G. Knepley /*E
186b0753f9dSMatthew G. Knepley     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
187b0753f9dSMatthew G. Knepley 
188b0753f9dSMatthew G. Knepley     Level: intermediate
189b0753f9dSMatthew G. Knepley 
190db781477SPatrick Sanan .seealso: `PCFieldSplitSetSchurPre()`
191b0753f9dSMatthew G. Knepley E*/
192b0753f9dSMatthew G. Knepley typedef enum {PC_FIELDSPLIT_SCHUR_PRE_SELF,PC_FIELDSPLIT_SCHUR_PRE_SELFP,PC_FIELDSPLIT_SCHUR_PRE_A11,PC_FIELDSPLIT_SCHUR_PRE_USER,PC_FIELDSPLIT_SCHUR_PRE_FULL} PCFieldSplitSchurPreType;
193b0753f9dSMatthew G. Knepley 
194b0753f9dSMatthew G. Knepley /*E
195b0753f9dSMatthew G. Knepley     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
196b0753f9dSMatthew G. Knepley 
197b0753f9dSMatthew G. Knepley     Level: intermediate
198b0753f9dSMatthew G. Knepley 
199db781477SPatrick Sanan .seealso: `PCFieldSplitSetSchurFactType()`
200b0753f9dSMatthew G. Knepley E*/
201b0753f9dSMatthew G. Knepley typedef enum {
202b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
203b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
204b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
205b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_FULL
206b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType;
207b0753f9dSMatthew G. Knepley 
208b0753f9dSMatthew G. Knepley /*E
209b0753f9dSMatthew G. Knepley     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
210b0753f9dSMatthew G. Knepley 
211b0753f9dSMatthew G. Knepley     Level: intermediate
212b0753f9dSMatthew G. Knepley 
213db781477SPatrick Sanan .seealso: `PCPARMSSetGlobal()`
214b0753f9dSMatthew G. Knepley E*/
215b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
2169d8081ecSMatthew G. Knepley 
217b0753f9dSMatthew G. Knepley /*E
218b0753f9dSMatthew G. Knepley     PCPARMSLocalType - Determines the local preconditioner method in PARMS
219b0753f9dSMatthew G. Knepley 
220b0753f9dSMatthew G. Knepley     Level: intermediate
221b0753f9dSMatthew G. Knepley 
222db781477SPatrick Sanan .seealso: `PCPARMSSetLocal()`
223b0753f9dSMatthew G. Knepley E*/
224b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
225b0753f9dSMatthew G. Knepley 
226f0fc11ceSJed Brown /*J
227b0753f9dSMatthew G. Knepley     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method
228b0753f9dSMatthew G. Knepley 
229b0753f9dSMatthew G. Knepley     Level: intermediate
230b0753f9dSMatthew G. Knepley 
231f6130c64SBarry Smith $   PCGAMGAGG - (the default) smoothed aggregation algorithm, robust, very well tested
232f6130c64SBarry Smith $   PCGAMGGEO - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not well tested
233f6130c64SBarry Smith $   PCGAMGCLASSICAL - classical algebraic multigrid preconditioner, incomplete, poorly tested
234f6130c64SBarry Smith 
235db781477SPatrick Sanan .seealso: `PCMG`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()`
236f0fc11ceSJed Brown J*/
237b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType;
238b0753f9dSMatthew G. Knepley #define PCGAMGAGG         "agg"
239b0753f9dSMatthew G. Knepley #define PCGAMGGEO         "geo"
240b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL   "classical"
241b0753f9dSMatthew G. Knepley 
242b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType;
243b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT   "direct"
244b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard"
245b0753f9dSMatthew G. Knepley 
246b0753f9dSMatthew G. Knepley /*E
247b0753f9dSMatthew G. Knepley     PCMGType - Determines the type of multigrid method that is run.
248b0753f9dSMatthew G. Knepley 
249b0753f9dSMatthew G. Knepley    Level: beginner
250b0753f9dSMatthew G. Knepley 
251b0753f9dSMatthew G. Knepley    Values:
252c1cbb1deSBarry Smith +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycleType()
253b0753f9dSMatthew G. Knepley .  PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
254b0753f9dSMatthew G. Knepley                 smoothed before updating the residual. This only uses the
255b0753f9dSMatthew G. Knepley                 down smoother, in the preconditioner the upper smoother is ignored
256b0753f9dSMatthew G. Knepley .  PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
257b0753f9dSMatthew G. Knepley             that is starts on the coarsest grid, performs a cycle, interpolates
258b0753f9dSMatthew 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
259b0753f9dSMatthew G. Knepley             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
260b0753f9dSMatthew G. Knepley             calls the V-cycle only on the coarser level and has a post-smoother instead.
261b0753f9dSMatthew G. Knepley -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
262b0753f9dSMatthew G. Knepley                from a finer
263b0753f9dSMatthew G. Knepley 
264db781477SPatrick Sanan .seealso: `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()`
265b0753f9dSMatthew G. Knepley 
266b0753f9dSMatthew G. Knepley E*/
267b0753f9dSMatthew G. Knepley typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
268b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE;
269b0753f9dSMatthew G. Knepley 
270b0753f9dSMatthew G. Knepley /*E
271b0753f9dSMatthew G. Knepley     PCMGCycleType - Use V-cycle or W-cycle
272b0753f9dSMatthew G. Knepley 
273b0753f9dSMatthew G. Knepley    Level: beginner
274b0753f9dSMatthew G. Knepley 
275b0753f9dSMatthew G. Knepley    Values:
276147403d9SBarry Smith +  PC_MG_V_CYCLE - use the v cycle
277147403d9SBarry Smith -  PC_MG_W_CYCLE - use the w cycle
278b0753f9dSMatthew G. Knepley 
279db781477SPatrick Sanan .seealso: `PCMGSetCycleType()`
280b0753f9dSMatthew G. Knepley 
281b0753f9dSMatthew G. Knepley E*/
282b0753f9dSMatthew G. Knepley typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
283b0753f9dSMatthew G. Knepley 
284b0753f9dSMatthew G. Knepley /*E
2852134b1e4SBarry Smith     PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process
2862134b1e4SBarry Smith 
2872134b1e4SBarry Smith    Level: beginner
2882134b1e4SBarry Smith 
2892134b1e4SBarry Smith    Values:
2902134b1e4SBarry Smith +  PC_MG_GALERKIN_PMAT - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid
2912134b1e4SBarry Smith .  PC_MG_GALERKIN_MAT -  computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid
2922134b1e4SBarry Smith .  PC_MG_GALERKIN_BOTH - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once
2932134b1e4SBarry Smith -  PC_MG_GALERKIN_NONE - neither operator is computed via the Galerkin process, the user must provide the operator
2942134b1e4SBarry Smith 
2952134b1e4SBarry Smith    Users should never set PC_MG_GALERKIN_EXTERNAL, it is used by GAMG and ML
2962134b1e4SBarry Smith 
297db781477SPatrick Sanan .seealso: `PCMGSetCycleType()`
2982134b1e4SBarry Smith 
2992134b1e4SBarry Smith E*/
3002134b1e4SBarry Smith typedef enum { PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, PC_MG_GALERKIN_EXTERNAL} PCMGGalerkinType;
3012134b1e4SBarry Smith 
3022134b1e4SBarry Smith /*E
303b0753f9dSMatthew G. Knepley     PCExoticType - Face based or wirebasket based coarse grid space
304b0753f9dSMatthew G. Knepley 
305b0753f9dSMatthew G. Knepley    Level: beginner
306b0753f9dSMatthew G. Knepley 
307db781477SPatrick Sanan .seealso: `PCExoticSetType()`, `PCEXOTIC`
308b0753f9dSMatthew G. Knepley E*/
309b0753f9dSMatthew G. Knepley typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
310b0753f9dSMatthew G. Knepley 
3118c1cd74cSHong Zhang /*E
312bc960bbfSJed Brown    PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains.
313bc960bbfSJed Brown 
314bc960bbfSJed Brown    Level: intermediate
315bc960bbfSJed Brown 
316bc960bbfSJed Brown    Values:
317bc960bbfSJed Brown +  PC_BDDC_INTERFACE_EXT_DIRICHLET - solves Dirichlet interior problem; this is the standard BDDC algorithm
318bc960bbfSJed Brown -  PC_BDDC_INTERFACE_EXT_LUMP - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP"
319bc960bbfSJed Brown 
320bc960bbfSJed Brown E*/
321bc960bbfSJed Brown typedef enum {
322bc960bbfSJed Brown   PC_BDDC_INTERFACE_EXT_DIRICHLET,
323bc960bbfSJed Brown   PC_BDDC_INTERFACE_EXT_LUMP
324bc960bbfSJed Brown } PCBDDCInterfaceExtType;
325bc960bbfSJed Brown 
326bc960bbfSJed Brown /*E
327f3b08a26SMatthew G. Knepley   PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation
328f3b08a26SMatthew G. Knepley 
329f3b08a26SMatthew G. Knepley   Level: beginner
330f3b08a26SMatthew G. Knepley 
331db781477SPatrick Sanan .seealso: `PCMGSetAdaptCoarseSpaceType()`, `PCMG`
332f3b08a26SMatthew G. Knepley E*/
333f3b08a26SMatthew G. Knepley typedef enum { PCMG_POLYNOMIAL, PCMG_HARMONIC, PCMG_EIGENVECTOR, PCMG_GENERALIZED_EIGENVECTOR } PCMGCoarseSpaceType;
334f3b08a26SMatthew G. Knepley 
335f3b08a26SMatthew G. Knepley /*E
336ee6acaa3SMatthew G. Knepley     PCPatchConstructType - The algorithm used to construct patches for the preconditioner
3374bbf5ea8SMatthew G. Knepley 
3384bbf5ea8SMatthew G. Knepley    Level: beginner
3394bbf5ea8SMatthew G. Knepley 
340db781477SPatrick Sanan .seealso: `PCPatchSetConstructType()`, `PCEXOTIC`
3414bbf5ea8SMatthew G. Knepley E*/
342e5b9877fSPatrick Farrell typedef enum {PC_PATCH_STAR, PC_PATCH_VANKA, PC_PATCH_PARDECOMP, PC_PATCH_USER, PC_PATCH_PYTHON} PCPatchConstructType;
3434bbf5ea8SMatthew G. Knepley 
3444bbf5ea8SMatthew G. Knepley /*E
345e26ad46dSJakub Kruzik     PCDeflationSpaceType - Type of deflation
346e662bc50SJakub Kruzik 
347e662bc50SJakub Kruzik     Values:
3487e9012c5SJakub Kruzik +   PC_DEFLATION_SPACE_HAAR        - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off
349e26ad46dSJakub Kruzik .   PC_DEFLATION_SPACE_DB2         - MATCOMPOSITE of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet)
350e26ad46dSJakub Kruzik .   PC_DEFLATION_SPACE_DB4         - same as above, but with db4 (4 coefficient Daubechies)
351e26ad46dSJakub Kruzik .   PC_DEFLATION_SPACE_DB8         - same as above, but with db8 (8 coefficient Daubechies)
352e26ad46dSJakub Kruzik .   PC_DEFLATION_SPACE_DB16        - same as above, but with db16 (16 coefficient Daubechies)
353e26ad46dSJakub Kruzik .   PC_DEFLATION_SPACE_BIORTH22    - same as above, but with biorthogonal 2.2 (6 coefficients)
354e26ad46dSJakub Kruzik .   PC_DEFLATION_SPACE_MEYER       - same as above, but with Meyer/FIR (62 coefficients)
355e26ad46dSJakub Kruzik .   PC_DEFLATION_SPACE_AGGREGATION - aggregates local indices (given by operator matix distribution) into a subdomain
356e662bc50SJakub Kruzik -   PC_DEFLATION_SPACE_USER        - indicates space set by user
357e662bc50SJakub Kruzik 
358e662bc50SJakub Kruzik     Notes:
3597e9012c5SJakub Kruzik       Wavelet-based space (except Haar) can be used in multilevel deflation.
360e662bc50SJakub Kruzik 
3611fdb00f9SJakub Kruzik     Level: intermediate
3621fdb00f9SJakub Kruzik 
363db781477SPatrick Sanan .seealso: `PCDeflationSetSpaceToCompute()`, `PCDEFLATION`
364e662bc50SJakub Kruzik E*/
365e662bc50SJakub Kruzik typedef enum {
366e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_HAAR,
367e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB2,
368e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB4,
369e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB8,
370e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB16,
371e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_BIORTH22,
372e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_MEYER,
373e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_AGGREGATION,
374e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_USER
375e662bc50SJakub Kruzik } PCDeflationSpaceType;
376e662bc50SJakub Kruzik 
377e662bc50SJakub Kruzik /*E
37838cfc46eSPierre Jolivet     PCHPDDMCoarseCorrectionType - Type of coarse correction used by PCHPDDM
37938cfc46eSPierre Jolivet 
38038cfc46eSPierre Jolivet     Level: intermediate
38138cfc46eSPierre Jolivet 
38238cfc46eSPierre Jolivet     Values:
383e5cd072eSprj- +   PC_HPDDM_COARSE_CORRECTION_DEFLATED (default) - eq. (1) in PCHPDDMShellApply()
384e5cd072eSprj- .   PC_HPDDM_COARSE_CORRECTION_ADDITIVE - eq. (2)
385e5cd072eSprj- -   PC_HPDDM_COARSE_CORRECTION_BALANCED - eq. (3)
38638cfc46eSPierre Jolivet 
387db781477SPatrick Sanan .seealso: `PCHPDDM`, `PCSetType()`, `PCHPDDMShellApply()`
38838cfc46eSPierre Jolivet E*/
38938cfc46eSPierre Jolivet typedef enum { PC_HPDDM_COARSE_CORRECTION_DEFLATED, PC_HPDDM_COARSE_CORRECTION_ADDITIVE, PC_HPDDM_COARSE_CORRECTION_BALANCED } PCHPDDMCoarseCorrectionType;
39038cfc46eSPierre Jolivet 
39138cfc46eSPierre Jolivet /*E
3928c1cd74cSHong Zhang     PCFailedReason - indicates type of PC failure
3938c1cd74cSHong Zhang 
3948c1cd74cSHong Zhang     Level: beginner
3958c1cd74cSHong Zhang 
3968c1cd74cSHong Zhang     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
3978c1cd74cSHong Zhang E*/
3981b2b9847SBarry Smith typedef enum {PC_SETUP_ERROR = -1,PC_NOERROR,PC_FACTOR_STRUCT_ZEROPIVOT,PC_FACTOR_NUMERIC_ZEROPIVOT,PC_FACTOR_OUTMEMORY,PC_FACTOR_OTHER,PC_SUBPC_ERROR} PCFailedReason;
399ce7c7f2fSMark Adams 
400ce7c7f2fSMark Adams /*E
401ce7c7f2fSMark Adams     PCGAMGLayoutType - Layout for reduced grids
402ce7c7f2fSMark Adams 
403ce7c7f2fSMark Adams     Level: intermediate
404ce7c7f2fSMark Adams 
405db781477SPatrick Sanan .seealso: `PCGAMGSetCoarseGridLayoutType()`
4069969d477SPatrick Sanan 
407ce7c7f2fSMark Adams     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
408ce7c7f2fSMark Adams E*/
409a0095786SMark typedef enum {PCGAMG_LAYOUT_COMPACT,PCGAMG_LAYOUT_SPREAD} PCGAMGLayoutType;
410ce7c7f2fSMark Adams 
411b0753f9dSMatthew G. Knepley #endif
412