xref: /petsc/include/petscpctypes.h (revision 284f0211dbdbe487df0cc78fa00fd614e5cc2783)
1b0753f9dSMatthew G. Knepley #if !defined(_PETSCPCTYPES_H)
2b0753f9dSMatthew G. Knepley #define _PETSCPCTYPES_H
3b0753f9dSMatthew G. Knepley 
4b0753f9dSMatthew G. Knepley #include <petscdmtypes.h>
5b0753f9dSMatthew G. Knepley 
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 
11b0753f9dSMatthew G. Knepley   Concepts: preconditioners
12b0753f9dSMatthew G. Knepley 
13b0753f9dSMatthew G. Knepley .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
14b0753f9dSMatthew G. Knepley S*/
15b0753f9dSMatthew G. Knepley typedef struct _p_PC* PC;
16b0753f9dSMatthew G. Knepley 
17b0753f9dSMatthew G. Knepley /*J
18b0753f9dSMatthew G. Knepley     PCType - String with the name of a PETSc preconditioner method.
19b0753f9dSMatthew G. Knepley 
20b0753f9dSMatthew G. Knepley    Level: beginner
21b0753f9dSMatthew G. Knepley 
22b0753f9dSMatthew G. Knepley    Notes: Click on the links below to see details on a particular solver
23b0753f9dSMatthew G. Knepley 
24b0753f9dSMatthew G. Knepley           PCRegister() is used to register preconditioners that are then accessible via PCSetType()
25b0753f9dSMatthew G. Knepley 
26b0753f9dSMatthew G. Knepley .seealso: PCSetType(), PC, PCCreate(), PCRegister(), PCSetFromOptions()
27b0753f9dSMatthew G. Knepley J*/
28b0753f9dSMatthew G. Knepley typedef const char* PCType;
29b0753f9dSMatthew G. Knepley #define PCNONE            "none"
30b0753f9dSMatthew G. Knepley #define PCJACOBI          "jacobi"
31b0753f9dSMatthew G. Knepley #define PCSOR             "sor"
32b0753f9dSMatthew G. Knepley #define PCLU              "lu"
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"
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"
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"
65b0753f9dSMatthew G. Knepley #define PCSACUSP          "sacusp"        /* these four run on NVIDIA GPUs using CUSP */
66b0753f9dSMatthew G. Knepley #define PCSACUSPPOLY      "sacusppoly"
67b0753f9dSMatthew G. Knepley #define PCBICGSTABCUSP    "bicgstabcusp"
68b0753f9dSMatthew G. Knepley #define PCAINVCUSP        "ainvcusp"
69b0753f9dSMatthew G. Knepley #define PCBDDC            "bddc"
70b0753f9dSMatthew G. Knepley #define PCKACZMARZ        "kaczmarz"
7168ddcbeaSBarry Smith #define PCTELESCOPE       "telescope"
72b0753f9dSMatthew G. Knepley 
73b0753f9dSMatthew G. Knepley /*E
74b0753f9dSMatthew G. Knepley     PCSide - If the preconditioner is to be applied to the left, right
75b0753f9dSMatthew G. Knepley      or symmetrically around the operator.
76b0753f9dSMatthew G. Knepley 
77b0753f9dSMatthew G. Knepley    Level: beginner
78b0753f9dSMatthew G. Knepley 
79b0753f9dSMatthew G. Knepley .seealso:
80b0753f9dSMatthew G. Knepley E*/
81b0753f9dSMatthew G. Knepley typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
82b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
83b0753f9dSMatthew G. Knepley PETSC_EXTERN const char *const *const PCSides;
84b0753f9dSMatthew G. Knepley 
85b0753f9dSMatthew G. Knepley /*E
86b0753f9dSMatthew G. Knepley     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
87b0753f9dSMatthew G. Knepley 
88b0753f9dSMatthew G. Knepley    Level: advanced
89b0753f9dSMatthew G. Knepley 
90b0753f9dSMatthew G. Knepley    Notes: this must match petsc/finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
91b0753f9dSMatthew G. Knepley 
92b0753f9dSMatthew G. Knepley .seealso: PCApplyRichardson()
93b0753f9dSMatthew G. Knepley E*/
94b0753f9dSMatthew G. Knepley typedef enum {
95b0753f9dSMatthew G. Knepley               PCRICHARDSON_CONVERGED_RTOL               =  2,
96b0753f9dSMatthew G. Knepley               PCRICHARDSON_CONVERGED_ATOL               =  3,
97b0753f9dSMatthew G. Knepley               PCRICHARDSON_CONVERGED_ITS                =  4,
98b0753f9dSMatthew G. Knepley               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
99b0753f9dSMatthew G. Knepley 
100b0753f9dSMatthew G. Knepley /*E
101b0753f9dSMatthew G. Knepley     PCJacobiType - What elements are used to form the Jacobi preconditioner
102b0753f9dSMatthew G. Knepley 
103b0753f9dSMatthew G. Knepley    Level: intermediate
104b0753f9dSMatthew G. Knepley 
105b0753f9dSMatthew G. Knepley .seealso:
106b0753f9dSMatthew G. Knepley E*/
107b0753f9dSMatthew G. Knepley typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType;
108b0753f9dSMatthew G. Knepley PETSC_EXTERN const char *const PCJacobiTypes[];
109b0753f9dSMatthew G. Knepley 
110b0753f9dSMatthew G. Knepley /*E
111b0753f9dSMatthew G. Knepley     PCASMType - Type of additive Schwarz method to use
112b0753f9dSMatthew G. Knepley 
113b0753f9dSMatthew G. Knepley $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
114b0753f9dSMatthew G. Knepley $                        and computed values in ghost regions are added together.
115b0753f9dSMatthew G. Knepley $                        Classical standard additive Schwarz.
116b0753f9dSMatthew G. Knepley $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
117b0753f9dSMatthew G. Knepley $                        region are discarded.
118b0753f9dSMatthew G. Knepley $                        Default.
119b0753f9dSMatthew G. Knepley $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
120b0753f9dSMatthew G. Knepley $                        region are added back in.
121b0753f9dSMatthew G. Knepley $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
122b0753f9dSMatthew G. Knepley $                        discarded.
123b0753f9dSMatthew G. Knepley $                        Not very good.
124b0753f9dSMatthew G. Knepley 
125b0753f9dSMatthew G. Knepley    Level: beginner
126b0753f9dSMatthew G. Knepley 
127b0753f9dSMatthew G. Knepley .seealso: PCASMSetType()
128b0753f9dSMatthew G. Knepley E*/
129b0753f9dSMatthew G. Knepley typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
130b0753f9dSMatthew G. Knepley PETSC_EXTERN const char *const PCASMTypes[];
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 PETSC_EXTERN const char *const PCGASMTypes[];
162b0753f9dSMatthew G. Knepley 
163b0753f9dSMatthew G. Knepley /*E
164b0753f9dSMatthew G. Knepley     PCCompositeType - Determines how two or more preconditioner are composed
165b0753f9dSMatthew G. Knepley 
166b0753f9dSMatthew G. Knepley $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
167b0753f9dSMatthew G. Knepley $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
168b0753f9dSMatthew G. Knepley $                                computed after the previous preconditioner application
169b0753f9dSMatthew G. Knepley $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
170b0753f9dSMatthew G. Knepley $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditions)
171b0753f9dSMatthew G. Knepley $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
172b0753f9dSMatthew G. Knepley $                         where first preconditioner is built from alpha I + S and second from
173b0753f9dSMatthew G. Knepley $                         alpha I + R
174b0753f9dSMatthew G. Knepley 
175b0753f9dSMatthew G. Knepley    Level: beginner
176b0753f9dSMatthew G. Knepley 
177b0753f9dSMatthew G. Knepley .seealso: PCCompositeSetType()
178b0753f9dSMatthew G. Knepley E*/
179b0753f9dSMatthew G. Knepley typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
180b0753f9dSMatthew G. Knepley PETSC_EXTERN const char *const PCCompositeTypes[];
181b0753f9dSMatthew G. Knepley 
182b0753f9dSMatthew G. Knepley /*E
183b0753f9dSMatthew G. Knepley     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
184b0753f9dSMatthew G. Knepley 
185b0753f9dSMatthew G. Knepley     Level: intermediate
186b0753f9dSMatthew G. Knepley 
187b0753f9dSMatthew G. Knepley .seealso: PCFieldSplitSetSchurPre()
188b0753f9dSMatthew G. Knepley E*/
189b0753f9dSMatthew 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;
190b0753f9dSMatthew G. Knepley PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];
191b0753f9dSMatthew G. Knepley 
192b0753f9dSMatthew G. Knepley /*E
193b0753f9dSMatthew G. Knepley     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
194b0753f9dSMatthew G. Knepley 
195b0753f9dSMatthew G. Knepley     Level: intermediate
196b0753f9dSMatthew G. Knepley 
197b0753f9dSMatthew G. Knepley .seealso: PCFieldSplitSetSchurFactType()
198b0753f9dSMatthew G. Knepley E*/
199b0753f9dSMatthew G. Knepley typedef enum {
200b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
201b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
202b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
203b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_FULL
204b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType;
205b0753f9dSMatthew G. Knepley PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];
206b0753f9dSMatthew G. Knepley 
207b0753f9dSMatthew G. Knepley /*E
208b0753f9dSMatthew G. Knepley     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
209b0753f9dSMatthew G. Knepley 
210b0753f9dSMatthew G. Knepley     Level: intermediate
211b0753f9dSMatthew G. Knepley 
212b0753f9dSMatthew G. Knepley .seealso: PCPARMSSetGlobal()
213b0753f9dSMatthew G. Knepley E*/
214b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
215b0753f9dSMatthew G. Knepley PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
216b0753f9dSMatthew G. Knepley /*E
217b0753f9dSMatthew G. Knepley     PCPARMSLocalType - Determines the local preconditioner method in PARMS
218b0753f9dSMatthew G. Knepley 
219b0753f9dSMatthew G. Knepley     Level: intermediate
220b0753f9dSMatthew G. Knepley 
221b0753f9dSMatthew G. Knepley .seealso: PCPARMSSetLocal()
222b0753f9dSMatthew G. Knepley E*/
223b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
224b0753f9dSMatthew G. Knepley PETSC_EXTERN const char *const PCPARMSLocalTypes[];
225b0753f9dSMatthew G. Knepley 
226b0753f9dSMatthew G. Knepley /*E
227b0753f9dSMatthew G. Knepley     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method
228b0753f9dSMatthew G. Knepley 
229b0753f9dSMatthew G. Knepley     Level: intermediate
230b0753f9dSMatthew G. Knepley 
231b0753f9dSMatthew G. Knepley .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseInterpolation()
232b0753f9dSMatthew G. Knepley E*/
233b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType;
234b0753f9dSMatthew G. Knepley #define PCGAMGAGG         "agg"
235b0753f9dSMatthew G. Knepley #define PCGAMGGEO         "geo"
236b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL   "classical"
237b0753f9dSMatthew G. Knepley 
238b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType;
239b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT   "direct"
240b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard"
241b0753f9dSMatthew G. Knepley 
242b0753f9dSMatthew G. Knepley /*E
243b0753f9dSMatthew G. Knepley     PCMGType - Determines the type of multigrid method that is run.
244b0753f9dSMatthew G. Knepley 
245b0753f9dSMatthew G. Knepley    Level: beginner
246b0753f9dSMatthew G. Knepley 
247b0753f9dSMatthew G. Knepley    Values:
248b0753f9dSMatthew G. Knepley +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycles()
249b0753f9dSMatthew G. Knepley .  PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
250b0753f9dSMatthew G. Knepley                 smoothed before updating the residual. This only uses the
251b0753f9dSMatthew G. Knepley                 down smoother, in the preconditioner the upper smoother is ignored
252b0753f9dSMatthew G. Knepley .  PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
253b0753f9dSMatthew G. Knepley             that is starts on the coarsest grid, performs a cycle, interpolates
254b0753f9dSMatthew 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
255b0753f9dSMatthew G. Knepley             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
256b0753f9dSMatthew G. Knepley             calls the V-cycle only on the coarser level and has a post-smoother instead.
257b0753f9dSMatthew G. Knepley -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
258b0753f9dSMatthew G. Knepley                from a finer
259b0753f9dSMatthew G. Knepley 
260b0753f9dSMatthew G. Knepley .seealso: PCMGSetType()
261b0753f9dSMatthew G. Knepley 
262b0753f9dSMatthew G. Knepley E*/
263b0753f9dSMatthew G. Knepley typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
264b0753f9dSMatthew G. Knepley PETSC_EXTERN const char *const PCMGTypes[];
265b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE;
266b0753f9dSMatthew G. Knepley 
267b0753f9dSMatthew G. Knepley /*E
268b0753f9dSMatthew G. Knepley     PCMGCycleType - Use V-cycle or W-cycle
269b0753f9dSMatthew G. Knepley 
270b0753f9dSMatthew G. Knepley    Level: beginner
271b0753f9dSMatthew G. Knepley 
272b0753f9dSMatthew G. Knepley    Values:
273b0753f9dSMatthew G. Knepley +  PC_MG_V_CYCLE
274b0753f9dSMatthew G. Knepley -  PC_MG_W_CYCLE
275b0753f9dSMatthew G. Knepley 
276b0753f9dSMatthew G. Knepley .seealso: PCMGSetCycleType()
277b0753f9dSMatthew G. Knepley 
278b0753f9dSMatthew G. Knepley E*/
279b0753f9dSMatthew G. Knepley typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
280b0753f9dSMatthew G. Knepley PETSC_EXTERN const char *const PCMGCycleTypes[];
281b0753f9dSMatthew G. Knepley 
282b0753f9dSMatthew G. Knepley /*E
283b0753f9dSMatthew G. Knepley     PCExoticType - Face based or wirebasket based coarse grid space
284b0753f9dSMatthew G. Knepley 
285b0753f9dSMatthew G. Knepley    Level: beginner
286b0753f9dSMatthew G. Knepley 
287b0753f9dSMatthew G. Knepley .seealso: PCExoticSetType(), PCEXOTIC
288b0753f9dSMatthew G. Knepley E*/
289b0753f9dSMatthew G. Knepley typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
290b0753f9dSMatthew G. Knepley PETSC_EXTERN const char *const PCExoticTypes[];
291b0753f9dSMatthew G. Knepley PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);
292b0753f9dSMatthew G. Knepley 
2938c1cd74cSHong Zhang /*E
2948c1cd74cSHong Zhang     PCFailedReason - indicates type of PC failure
2958c1cd74cSHong Zhang 
2968c1cd74cSHong Zhang     Level: beginner
2978c1cd74cSHong Zhang 
2988c1cd74cSHong Zhang     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
2998c1cd74cSHong Zhang E*/
300*284f0211SHong Zhang typedef enum {PC_NOERROR,PC_FACTOR_STRUCT_ZEROPIVOT,PC_FACTOR_NUMERIC_ZEROPIVOT,PC_FACTOR_OUTMEMORY,PC_FACTOR_OTHER,PC_SUBPC_ERROR} PCFailedReason;
3011c852feaSHong Zhang PETSC_EXTERN const char *const PCFailedReasons[];
302b0753f9dSMatthew G. Knepley #endif
303