xref: /petsc/include/petscpctypes.h (revision 360ee056ff647363bf77f2811d7ea75daaaacc81)
1b0753f9dSMatthew G. Knepley #if !defined(_PETSCPCTYPES_H)
2b0753f9dSMatthew G. Knepley #define _PETSCPCTYPES_H
3b0753f9dSMatthew G. Knepley 
4b0753f9dSMatthew G. Knepley /*S
5b0753f9dSMatthew G. Knepley      PC - Abstract PETSc object that manages all preconditioners including direct solvers such as PCLU
6b0753f9dSMatthew G. Knepley 
7b0753f9dSMatthew G. Knepley    Level: beginner
8b0753f9dSMatthew G. Knepley 
9b0753f9dSMatthew G. Knepley   Concepts: preconditioners
10b0753f9dSMatthew G. Knepley 
11b0753f9dSMatthew G. Knepley .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
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 
25b0753f9dSMatthew G. Knepley .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"
32b0753f9dSMatthew G. Knepley #define PCSHELL           "shell"
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"
41b0753f9dSMatthew G. Knepley #define PCCOMPOSITE       "composite"
42b0753f9dSMatthew G. Knepley #define PCREDUNDANT       "redundant"
43b0753f9dSMatthew G. Knepley #define PCSPAI            "spai"
44b0753f9dSMatthew G. Knepley #define PCNN              "nn"
45b0753f9dSMatthew G. Knepley #define PCCHOLESKY        "cholesky"
46b0753f9dSMatthew G. Knepley #define PCPBJACOBI        "pbjacobi"
470da83c2eSBarry Smith #define PCVPBJACOBI       "vpbjacobi"
48b0753f9dSMatthew G. Knepley #define PCMAT             "mat"
49b0753f9dSMatthew G. Knepley #define PCHYPRE           "hypre"
50b0753f9dSMatthew G. Knepley #define PCPARMS           "parms"
51b0753f9dSMatthew G. Knepley #define PCFIELDSPLIT      "fieldsplit"
52b0753f9dSMatthew G. Knepley #define PCTFS             "tfs"
53b0753f9dSMatthew G. Knepley #define PCML              "ml"
54b0753f9dSMatthew G. Knepley #define PCGALERKIN        "galerkin"
55b0753f9dSMatthew G. Knepley #define PCEXOTIC          "exotic"
56b0753f9dSMatthew G. Knepley #define PCCP              "cp"
57b0753f9dSMatthew G. Knepley #define PCBFBT            "bfbt"
58b0753f9dSMatthew G. Knepley #define PCLSC             "lsc"
59b0753f9dSMatthew G. Knepley #define PCPYTHON          "python"
60b0753f9dSMatthew G. Knepley #define PCPFMG            "pfmg"
61b0753f9dSMatthew G. Knepley #define PCSYSPFMG         "syspfmg"
62b0753f9dSMatthew G. Knepley #define PCREDISTRIBUTE    "redistribute"
63b0753f9dSMatthew G. Knepley #define PCSVD             "svd"
64b0753f9dSMatthew G. Knepley #define PCGAMG            "gamg"
654b3f184cSKarl Rupp #define PCCHOWILUVIENNACL "chowiluviennacl"
6670baa948SKarl Rupp #define PCROWSCALINGVIENNACL "rowscalingviennacl"
6707598726SKarl Rupp #define PCSAVIENNACL      "saviennacl"
68b0753f9dSMatthew G. Knepley #define PCBDDC            "bddc"
69b0753f9dSMatthew G. Knepley #define PCKACZMARZ        "kaczmarz"
7068ddcbeaSBarry Smith #define PCTELESCOPE       "telescope"
714bbf5ea8SMatthew G. Knepley #define PCPATCH           "patch"
72b9ac7092SAlp Dener #define PCLMVM            "lmvm"
73*360ee056SFande Kong #define PCHMG             "hmg"
74b0753f9dSMatthew G. Knepley 
75b0753f9dSMatthew G. Knepley /*E
76b0753f9dSMatthew G. Knepley     PCSide - If the preconditioner is to be applied to the left, right
77b0753f9dSMatthew G. Knepley      or symmetrically around the operator.
78b0753f9dSMatthew G. Knepley 
79b0753f9dSMatthew G. Knepley    Level: beginner
80b0753f9dSMatthew G. Knepley 
81b0753f9dSMatthew G. Knepley .seealso:
82b0753f9dSMatthew G. Knepley E*/
83b0753f9dSMatthew G. Knepley typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
84b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
85b0753f9dSMatthew G. Knepley 
86b0753f9dSMatthew G. Knepley /*E
87b0753f9dSMatthew G. Knepley     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
88b0753f9dSMatthew G. Knepley 
89b0753f9dSMatthew G. Knepley    Level: advanced
90b0753f9dSMatthew G. Knepley 
9195452b02SPatrick Sanan    Notes:
9295452b02SPatrick Sanan     this must match petsc/finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
93b0753f9dSMatthew G. Knepley 
94b0753f9dSMatthew G. Knepley .seealso: PCApplyRichardson()
95b0753f9dSMatthew G. Knepley E*/
96b0753f9dSMatthew G. Knepley typedef enum {
97b0753f9dSMatthew G. Knepley               PCRICHARDSON_CONVERGED_RTOL               =  2,
98b0753f9dSMatthew G. Knepley               PCRICHARDSON_CONVERGED_ATOL               =  3,
99b0753f9dSMatthew G. Knepley               PCRICHARDSON_CONVERGED_ITS                =  4,
100b0753f9dSMatthew G. Knepley               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
101b0753f9dSMatthew G. Knepley 
102b0753f9dSMatthew G. Knepley /*E
103b0753f9dSMatthew G. Knepley     PCJacobiType - What elements are used to form the Jacobi preconditioner
104b0753f9dSMatthew G. Knepley 
105b0753f9dSMatthew G. Knepley    Level: intermediate
106b0753f9dSMatthew G. Knepley 
107b0753f9dSMatthew G. Knepley .seealso:
108b0753f9dSMatthew G. Knepley E*/
109b0753f9dSMatthew G. Knepley typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType;
110b0753f9dSMatthew G. Knepley 
111b0753f9dSMatthew G. Knepley /*E
112b0753f9dSMatthew G. Knepley     PCASMType - Type of additive Schwarz method to use
113b0753f9dSMatthew G. Knepley 
114b0753f9dSMatthew G. Knepley $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
115b0753f9dSMatthew G. Knepley $                        and computed values in ghost regions are added together.
116b0753f9dSMatthew G. Knepley $                        Classical standard additive Schwarz.
117b0753f9dSMatthew G. Knepley $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
118b0753f9dSMatthew G. Knepley $                        region are discarded.
119b0753f9dSMatthew G. Knepley $                        Default.
120b0753f9dSMatthew G. Knepley $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
121b0753f9dSMatthew G. Knepley $                        region are added back in.
122b0753f9dSMatthew G. Knepley $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
123b0753f9dSMatthew G. Knepley $                        discarded.
124b0753f9dSMatthew G. Knepley $                        Not very good.
125b0753f9dSMatthew G. Knepley 
126b0753f9dSMatthew G. Knepley    Level: beginner
127b0753f9dSMatthew G. Knepley 
128b0753f9dSMatthew G. Knepley .seealso: PCASMSetType()
129b0753f9dSMatthew G. Knepley E*/
130b0753f9dSMatthew G. Knepley typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
131b0753f9dSMatthew G. Knepley 
132b0753f9dSMatthew G. Knepley /*E
133b0753f9dSMatthew G. Knepley     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
134b0753f9dSMatthew G. Knepley 
135b0753f9dSMatthew G. Knepley    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
136b0753f9dSMatthew G. Knepley    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
137b0753f9dSMatthew G. Knepley    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
138b0753f9dSMatthew G. Knepley    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
139b0753f9dSMatthew G. Knepley 
140b0753f9dSMatthew G. Knepley $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
141b0753f9dSMatthew G. Knepley $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
142b0753f9dSMatthew G. Knepley $                        from neighboring subdomains.
143b0753f9dSMatthew G. Knepley $                        Classical standard additive Schwarz.
144b0753f9dSMatthew G. Knepley $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
145b0753f9dSMatthew G. Knepley $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
146b0753f9dSMatthew G. Knepley $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
147b0753f9dSMatthew G. Knepley $                        assumption).
148b0753f9dSMatthew G. Knepley $                        Default.
149b0753f9dSMatthew G. Knepley $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
150b0753f9dSMatthew G. Knepley $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
151b0753f9dSMatthew G. Knepley $                        from neighboring subdomains.
152b0753f9dSMatthew G. Knepley $
153b0753f9dSMatthew G. Knepley $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
154b0753f9dSMatthew G. Knepley $                        Not very good.
155b0753f9dSMatthew G. Knepley 
156b0753f9dSMatthew G. Knepley    Level: beginner
157b0753f9dSMatthew G. Knepley 
158b0753f9dSMatthew G. Knepley .seealso: PCGASMSetType()
159b0753f9dSMatthew G. Knepley E*/
160b0753f9dSMatthew G. Knepley typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
161b0753f9dSMatthew G. Knepley 
162b0753f9dSMatthew G. Knepley /*E
163b0753f9dSMatthew G. Knepley     PCCompositeType - Determines how two or more preconditioner are composed
164b0753f9dSMatthew G. Knepley 
165b0753f9dSMatthew G. Knepley $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
166b0753f9dSMatthew G. Knepley $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
167b0753f9dSMatthew G. Knepley $                                computed after the previous preconditioner application
168b0753f9dSMatthew G. Knepley $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
169d0ecd4dfSBarry Smith $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners)
170b0753f9dSMatthew G. Knepley $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
171b0753f9dSMatthew G. Knepley $                         where first preconditioner is built from alpha I + S and second from
172b0753f9dSMatthew G. Knepley $                         alpha I + R
173b0753f9dSMatthew G. Knepley 
174b0753f9dSMatthew G. Knepley    Level: beginner
175b0753f9dSMatthew G. Knepley 
176b0753f9dSMatthew G. Knepley .seealso: PCCompositeSetType()
177b0753f9dSMatthew G. Knepley E*/
178a51937d4SCarola Kruse typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR,PC_COMPOSITE_GKB} PCCompositeType;
179b0753f9dSMatthew G. Knepley 
180b0753f9dSMatthew G. Knepley /*E
181b0753f9dSMatthew G. Knepley     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
182b0753f9dSMatthew G. Knepley 
183b0753f9dSMatthew G. Knepley     Level: intermediate
184b0753f9dSMatthew G. Knepley 
185b0753f9dSMatthew G. Knepley .seealso: PCFieldSplitSetSchurPre()
186b0753f9dSMatthew G. Knepley E*/
187b0753f9dSMatthew 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;
188b0753f9dSMatthew G. Knepley 
189b0753f9dSMatthew G. Knepley /*E
190b0753f9dSMatthew G. Knepley     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
191b0753f9dSMatthew G. Knepley 
192b0753f9dSMatthew G. Knepley     Level: intermediate
193b0753f9dSMatthew G. Knepley 
194b0753f9dSMatthew G. Knepley .seealso: PCFieldSplitSetSchurFactType()
195b0753f9dSMatthew G. Knepley E*/
196b0753f9dSMatthew G. Knepley typedef enum {
197b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
198b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
199b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
200b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_FULL
201b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType;
202b0753f9dSMatthew G. Knepley 
203b0753f9dSMatthew G. Knepley /*E
204b0753f9dSMatthew G. Knepley     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
205b0753f9dSMatthew G. Knepley 
206b0753f9dSMatthew G. Knepley     Level: intermediate
207b0753f9dSMatthew G. Knepley 
208b0753f9dSMatthew G. Knepley .seealso: PCPARMSSetGlobal()
209b0753f9dSMatthew G. Knepley E*/
210b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
2119d8081ecSMatthew G. Knepley 
212b0753f9dSMatthew G. Knepley /*E
213b0753f9dSMatthew G. Knepley     PCPARMSLocalType - Determines the local preconditioner method in PARMS
214b0753f9dSMatthew G. Knepley 
215b0753f9dSMatthew G. Knepley     Level: intermediate
216b0753f9dSMatthew G. Knepley 
217b0753f9dSMatthew G. Knepley .seealso: PCPARMSSetLocal()
218b0753f9dSMatthew G. Knepley E*/
219b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
220b0753f9dSMatthew G. Knepley 
221b0753f9dSMatthew G. Knepley /*E
222b0753f9dSMatthew G. Knepley     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method
223b0753f9dSMatthew G. Knepley 
224b0753f9dSMatthew G. Knepley     Level: intermediate
225b0753f9dSMatthew G. Knepley 
226b0753f9dSMatthew G. Knepley .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseInterpolation()
227b0753f9dSMatthew G. Knepley E*/
228b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType;
229b0753f9dSMatthew G. Knepley #define PCGAMGAGG         "agg"
230b0753f9dSMatthew G. Knepley #define PCGAMGGEO         "geo"
231b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL   "classical"
232b0753f9dSMatthew G. Knepley 
233b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType;
234b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT   "direct"
235b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard"
236b0753f9dSMatthew G. Knepley 
237b0753f9dSMatthew G. Knepley /*E
238b0753f9dSMatthew G. Knepley     PCMGType - Determines the type of multigrid method that is run.
239b0753f9dSMatthew G. Knepley 
240b0753f9dSMatthew G. Knepley    Level: beginner
241b0753f9dSMatthew G. Knepley 
242b0753f9dSMatthew G. Knepley    Values:
243c1cbb1deSBarry Smith +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycleType()
244b0753f9dSMatthew G. Knepley .  PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
245b0753f9dSMatthew G. Knepley                 smoothed before updating the residual. This only uses the
246b0753f9dSMatthew G. Knepley                 down smoother, in the preconditioner the upper smoother is ignored
247b0753f9dSMatthew G. Knepley .  PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
248b0753f9dSMatthew G. Knepley             that is starts on the coarsest grid, performs a cycle, interpolates
249b0753f9dSMatthew 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
250b0753f9dSMatthew G. Knepley             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
251b0753f9dSMatthew G. Knepley             calls the V-cycle only on the coarser level and has a post-smoother instead.
252b0753f9dSMatthew G. Knepley -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
253b0753f9dSMatthew G. Knepley                from a finer
254b0753f9dSMatthew G. Knepley 
255c1cbb1deSBarry Smith .seealso: PCMGSetType(), PCMGSetCycleType(), PCMGSetCycleTypeOnLevel()
256b0753f9dSMatthew G. Knepley 
257b0753f9dSMatthew G. Knepley E*/
258b0753f9dSMatthew G. Knepley typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
259b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE;
260b0753f9dSMatthew G. Knepley 
261b0753f9dSMatthew G. Knepley /*E
262b0753f9dSMatthew G. Knepley     PCMGCycleType - Use V-cycle or W-cycle
263b0753f9dSMatthew G. Knepley 
264b0753f9dSMatthew G. Knepley    Level: beginner
265b0753f9dSMatthew G. Knepley 
266b0753f9dSMatthew G. Knepley    Values:
267b0753f9dSMatthew G. Knepley +  PC_MG_V_CYCLE
268b0753f9dSMatthew G. Knepley -  PC_MG_W_CYCLE
269b0753f9dSMatthew G. Knepley 
270b0753f9dSMatthew G. Knepley .seealso: PCMGSetCycleType()
271b0753f9dSMatthew G. Knepley 
272b0753f9dSMatthew G. Knepley E*/
273b0753f9dSMatthew G. Knepley typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
274b0753f9dSMatthew G. Knepley 
275b0753f9dSMatthew G. Knepley /*E
2762134b1e4SBarry Smith     PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process
2772134b1e4SBarry Smith 
2782134b1e4SBarry Smith    Level: beginner
2792134b1e4SBarry Smith 
2802134b1e4SBarry Smith    Values:
2812134b1e4SBarry Smith +  PC_MG_GALERKIN_PMAT - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid
2822134b1e4SBarry Smith .  PC_MG_GALERKIN_MAT -  computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid
2832134b1e4SBarry Smith .  PC_MG_GALERKIN_BOTH - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once
2842134b1e4SBarry Smith -  PC_MG_GALERKIN_NONE - neither operator is computed via the Galerkin process, the user must provide the operator
2852134b1e4SBarry Smith 
2862134b1e4SBarry Smith    Users should never set PC_MG_GALERKIN_EXTERNAL, it is used by GAMG and ML
2872134b1e4SBarry Smith 
2882134b1e4SBarry Smith .seealso: PCMGSetCycleType()
2892134b1e4SBarry Smith 
2902134b1e4SBarry Smith E*/
2912134b1e4SBarry Smith typedef enum { PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, PC_MG_GALERKIN_EXTERNAL} PCMGGalerkinType;
2922134b1e4SBarry Smith 
2932134b1e4SBarry Smith /*E
294b0753f9dSMatthew G. Knepley     PCExoticType - Face based or wirebasket based coarse grid space
295b0753f9dSMatthew G. Knepley 
296b0753f9dSMatthew G. Knepley    Level: beginner
297b0753f9dSMatthew G. Knepley 
298b0753f9dSMatthew G. Knepley .seealso: PCExoticSetType(), PCEXOTIC
299b0753f9dSMatthew G. Knepley E*/
300b0753f9dSMatthew G. Knepley typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
301b0753f9dSMatthew G. Knepley 
3028c1cd74cSHong Zhang /*E
303ee6acaa3SMatthew G. Knepley     PCPatchConstructType - The algorithm used to construct patches for the preconditioner
3044bbf5ea8SMatthew G. Knepley 
3054bbf5ea8SMatthew G. Knepley    Level: beginner
3064bbf5ea8SMatthew G. Knepley 
3074bbf5ea8SMatthew G. Knepley .seealso: PCPatchSetConstructType(), PCEXOTIC
3084bbf5ea8SMatthew G. Knepley E*/
309e5b9877fSPatrick Farrell typedef enum {PC_PATCH_STAR, PC_PATCH_VANKA, PC_PATCH_PARDECOMP, PC_PATCH_USER, PC_PATCH_PYTHON} PCPatchConstructType;
3104bbf5ea8SMatthew G. Knepley 
3114bbf5ea8SMatthew G. Knepley /*E
3128c1cd74cSHong Zhang     PCFailedReason - indicates type of PC failure
3138c1cd74cSHong Zhang 
3148c1cd74cSHong Zhang     Level: beginner
3158c1cd74cSHong Zhang 
3168c1cd74cSHong Zhang     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
3178c1cd74cSHong Zhang E*/
318284f0211SHong Zhang typedef enum {PC_NOERROR,PC_FACTOR_STRUCT_ZEROPIVOT,PC_FACTOR_NUMERIC_ZEROPIVOT,PC_FACTOR_OUTMEMORY,PC_FACTOR_OTHER,PC_SUBPC_ERROR} PCFailedReason;
319b0753f9dSMatthew G. Knepley #endif
320