xref: /petsc/include/petscpctypes.h (revision e8812f7b072b3fdba59c6049ad9d46d6caf3ff6a)
126bd1501SBarry Smith #if !defined(PETSCPCTYPES_H)
226bd1501SBarry Smith #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 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
10b0753f9dSMatthew G. Knepley S*/
11b0753f9dSMatthew G. Knepley typedef struct _p_PC* PC;
12b0753f9dSMatthew G. Knepley 
13b0753f9dSMatthew G. Knepley /*J
14b0753f9dSMatthew G. Knepley     PCType - String with the name of a PETSc preconditioner method.
15b0753f9dSMatthew G. Knepley 
16b0753f9dSMatthew G. Knepley    Level: beginner
17b0753f9dSMatthew G. Knepley 
1895452b02SPatrick Sanan    Notes:
1995452b02SPatrick Sanan     Click on the links above to see details on a particular solver
20b0753f9dSMatthew G. Knepley 
21b0753f9dSMatthew G. Knepley           PCRegister() is used to register preconditioners that are then accessible via PCSetType()
22b0753f9dSMatthew G. Knepley 
23b0753f9dSMatthew G. Knepley .seealso: PCSetType(), PC, PCCreate(), PCRegister(), PCSetFromOptions()
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"
30b0753f9dSMatthew G. Knepley #define PCSHELL           "shell"
31b0753f9dSMatthew G. Knepley #define PCBJACOBI         "bjacobi"
32b0753f9dSMatthew G. Knepley #define PCMG              "mg"
33b0753f9dSMatthew G. Knepley #define PCEISENSTAT       "eisenstat"
34b0753f9dSMatthew G. Knepley #define PCILU             "ilu"
35b0753f9dSMatthew G. Knepley #define PCICC             "icc"
36b0753f9dSMatthew G. Knepley #define PCASM             "asm"
37b0753f9dSMatthew G. Knepley #define PCGASM            "gasm"
38b0753f9dSMatthew G. Knepley #define PCKSP             "ksp"
39b0753f9dSMatthew G. Knepley #define PCCOMPOSITE       "composite"
40b0753f9dSMatthew G. Knepley #define PCREDUNDANT       "redundant"
41b0753f9dSMatthew G. Knepley #define PCSPAI            "spai"
42b0753f9dSMatthew G. Knepley #define PCNN              "nn"
43b0753f9dSMatthew G. Knepley #define PCCHOLESKY        "cholesky"
44b0753f9dSMatthew G. Knepley #define PCPBJACOBI        "pbjacobi"
450da83c2eSBarry Smith #define PCVPBJACOBI       "vpbjacobi"
46b0753f9dSMatthew G. Knepley #define PCMAT             "mat"
47b0753f9dSMatthew G. Knepley #define PCHYPRE           "hypre"
48b0753f9dSMatthew G. Knepley #define PCPARMS           "parms"
49b0753f9dSMatthew G. Knepley #define PCFIELDSPLIT      "fieldsplit"
50b0753f9dSMatthew G. Knepley #define PCTFS             "tfs"
51b0753f9dSMatthew G. Knepley #define PCML              "ml"
52b0753f9dSMatthew G. Knepley #define PCGALERKIN        "galerkin"
53b0753f9dSMatthew G. Knepley #define PCEXOTIC          "exotic"
54b0753f9dSMatthew G. Knepley #define PCCP              "cp"
55b0753f9dSMatthew G. Knepley #define PCBFBT            "bfbt"
56b0753f9dSMatthew G. Knepley #define PCLSC             "lsc"
57b0753f9dSMatthew G. Knepley #define PCPYTHON          "python"
58b0753f9dSMatthew G. Knepley #define PCPFMG            "pfmg"
59b0753f9dSMatthew G. Knepley #define PCSYSPFMG         "syspfmg"
60b0753f9dSMatthew G. Knepley #define PCREDISTRIBUTE    "redistribute"
61b0753f9dSMatthew G. Knepley #define PCSVD             "svd"
62b0753f9dSMatthew G. Knepley #define PCGAMG            "gamg"
634b3f184cSKarl Rupp #define PCCHOWILUVIENNACL "chowiluviennacl"
6470baa948SKarl Rupp #define PCROWSCALINGVIENNACL "rowscalingviennacl"
6507598726SKarl Rupp #define PCSAVIENNACL      "saviennacl"
66b0753f9dSMatthew G. Knepley #define PCBDDC            "bddc"
67b0753f9dSMatthew G. Knepley #define PCKACZMARZ        "kaczmarz"
6868ddcbeaSBarry Smith #define PCTELESCOPE       "telescope"
694bbf5ea8SMatthew G. Knepley #define PCPATCH           "patch"
70b9ac7092SAlp Dener #define PCLMVM            "lmvm"
71360ee056SFande Kong #define PCHMG             "hmg"
7237eeb815SJakub Kruzik #define PCDEFLATION       "deflation"
7338cfc46eSPierre Jolivet #define PCHPDDM           "hpddm"
74*e8812f7bSStefano Zampini #define PCHARA            "hara"
75b0753f9dSMatthew G. Knepley 
76b0753f9dSMatthew G. Knepley /*E
77b0753f9dSMatthew G. Knepley     PCSide - If the preconditioner is to be applied to the left, right
78b0753f9dSMatthew G. Knepley      or symmetrically around the operator.
79b0753f9dSMatthew G. Knepley 
80b0753f9dSMatthew G. Knepley    Level: beginner
81b0753f9dSMatthew G. Knepley 
82b0753f9dSMatthew G. Knepley .seealso:
83b0753f9dSMatthew G. Knepley E*/
84b0753f9dSMatthew G. Knepley typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
85b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
86b0753f9dSMatthew G. Knepley 
87b0753f9dSMatthew G. Knepley /*E
88b0753f9dSMatthew G. Knepley     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
89b0753f9dSMatthew G. Knepley 
90b0753f9dSMatthew G. Knepley    Level: advanced
91b0753f9dSMatthew G. Knepley 
9295452b02SPatrick Sanan    Notes:
9395452b02SPatrick Sanan     this must match petsc/finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
94b0753f9dSMatthew G. Knepley 
95b0753f9dSMatthew G. Knepley .seealso: PCApplyRichardson()
96b0753f9dSMatthew G. Knepley E*/
97b0753f9dSMatthew G. Knepley typedef enum {
98b0753f9dSMatthew G. Knepley               PCRICHARDSON_CONVERGED_RTOL               =  2,
99b0753f9dSMatthew G. Knepley               PCRICHARDSON_CONVERGED_ATOL               =  3,
100b0753f9dSMatthew G. Knepley               PCRICHARDSON_CONVERGED_ITS                =  4,
101b0753f9dSMatthew G. Knepley               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
102b0753f9dSMatthew G. Knepley 
103b0753f9dSMatthew G. Knepley /*E
104b0753f9dSMatthew G. Knepley     PCJacobiType - What elements are used to form the Jacobi preconditioner
105b0753f9dSMatthew G. Knepley 
106b0753f9dSMatthew G. Knepley    Level: intermediate
107b0753f9dSMatthew G. Knepley 
108b0753f9dSMatthew G. Knepley .seealso:
109b0753f9dSMatthew G. Knepley E*/
110b0753f9dSMatthew G. Knepley typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType;
111b0753f9dSMatthew G. Knepley 
112b0753f9dSMatthew G. Knepley /*E
113b0753f9dSMatthew G. Knepley     PCASMType - Type of additive Schwarz method to use
114b0753f9dSMatthew G. Knepley 
115b0753f9dSMatthew G. Knepley $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
116b0753f9dSMatthew G. Knepley $                        and computed values in ghost regions are added together.
117b0753f9dSMatthew G. Knepley $                        Classical standard additive Schwarz.
118b0753f9dSMatthew G. Knepley $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
119b0753f9dSMatthew G. Knepley $                        region are discarded.
120b0753f9dSMatthew G. Knepley $                        Default.
121b0753f9dSMatthew G. Knepley $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
122b0753f9dSMatthew G. Knepley $                        region are added back in.
123b0753f9dSMatthew G. Knepley $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
124b0753f9dSMatthew G. Knepley $                        discarded.
125b0753f9dSMatthew G. Knepley $                        Not very good.
126b0753f9dSMatthew G. Knepley 
127b0753f9dSMatthew G. Knepley    Level: beginner
128b0753f9dSMatthew G. Knepley 
129b0753f9dSMatthew G. Knepley .seealso: PCASMSetType()
130b0753f9dSMatthew G. Knepley E*/
131b0753f9dSMatthew G. Knepley typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
132b0753f9dSMatthew G. Knepley 
133b0753f9dSMatthew G. Knepley /*E
134b0753f9dSMatthew G. Knepley     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
135b0753f9dSMatthew G. Knepley 
136b0753f9dSMatthew G. Knepley    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
137b0753f9dSMatthew G. Knepley    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
138b0753f9dSMatthew G. Knepley    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
139b0753f9dSMatthew G. Knepley    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
140b0753f9dSMatthew G. Knepley 
141b0753f9dSMatthew G. Knepley $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
142b0753f9dSMatthew G. Knepley $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
143b0753f9dSMatthew G. Knepley $                        from neighboring subdomains.
144b0753f9dSMatthew G. Knepley $                        Classical standard additive Schwarz.
145b0753f9dSMatthew G. Knepley $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
146b0753f9dSMatthew G. Knepley $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
147b0753f9dSMatthew G. Knepley $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
148b0753f9dSMatthew G. Knepley $                        assumption).
149b0753f9dSMatthew G. Knepley $                        Default.
150b0753f9dSMatthew G. Knepley $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
151b0753f9dSMatthew G. Knepley $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
152b0753f9dSMatthew G. Knepley $                        from neighboring subdomains.
153b0753f9dSMatthew G. Knepley $
154b0753f9dSMatthew G. Knepley $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
155b0753f9dSMatthew G. Knepley $                        Not very good.
156b0753f9dSMatthew G. Knepley 
157b0753f9dSMatthew G. Knepley    Level: beginner
158b0753f9dSMatthew G. Knepley 
159b0753f9dSMatthew G. Knepley .seealso: PCGASMSetType()
160b0753f9dSMatthew G. Knepley E*/
161b0753f9dSMatthew G. Knepley typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
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
170d0ecd4dfSBarry Smith $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners)
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*/
179a51937d4SCarola Kruse typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR,PC_COMPOSITE_GKB} PCCompositeType;
180b0753f9dSMatthew G. Knepley 
181b0753f9dSMatthew G. Knepley /*E
182b0753f9dSMatthew G. Knepley     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
183b0753f9dSMatthew G. Knepley 
184b0753f9dSMatthew G. Knepley     Level: intermediate
185b0753f9dSMatthew G. Knepley 
186b0753f9dSMatthew G. Knepley .seealso: PCFieldSplitSetSchurPre()
187b0753f9dSMatthew G. Knepley E*/
188b0753f9dSMatthew 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;
189b0753f9dSMatthew G. Knepley 
190b0753f9dSMatthew G. Knepley /*E
191b0753f9dSMatthew G. Knepley     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
192b0753f9dSMatthew G. Knepley 
193b0753f9dSMatthew G. Knepley     Level: intermediate
194b0753f9dSMatthew G. Knepley 
195b0753f9dSMatthew G. Knepley .seealso: PCFieldSplitSetSchurFactType()
196b0753f9dSMatthew G. Knepley E*/
197b0753f9dSMatthew G. Knepley typedef enum {
198b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
199b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
200b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
201b0753f9dSMatthew G. Knepley   PC_FIELDSPLIT_SCHUR_FACT_FULL
202b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType;
203b0753f9dSMatthew G. Knepley 
204b0753f9dSMatthew G. Knepley /*E
205b0753f9dSMatthew G. Knepley     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
206b0753f9dSMatthew G. Knepley 
207b0753f9dSMatthew G. Knepley     Level: intermediate
208b0753f9dSMatthew G. Knepley 
209b0753f9dSMatthew G. Knepley .seealso: PCPARMSSetGlobal()
210b0753f9dSMatthew G. Knepley E*/
211b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
2129d8081ecSMatthew G. Knepley 
213b0753f9dSMatthew G. Knepley /*E
214b0753f9dSMatthew G. Knepley     PCPARMSLocalType - Determines the local preconditioner method in PARMS
215b0753f9dSMatthew G. Knepley 
216b0753f9dSMatthew G. Knepley     Level: intermediate
217b0753f9dSMatthew G. Knepley 
218b0753f9dSMatthew G. Knepley .seealso: PCPARMSSetLocal()
219b0753f9dSMatthew G. Knepley E*/
220b0753f9dSMatthew G. Knepley typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
221b0753f9dSMatthew G. Knepley 
222f0fc11ceSJed Brown /*J
223b0753f9dSMatthew G. Knepley     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method
224b0753f9dSMatthew G. Knepley 
225b0753f9dSMatthew G. Knepley     Level: intermediate
226b0753f9dSMatthew G. Knepley 
227f6130c64SBarry Smith $   PCGAMGAGG - (the default) smoothed aggregation algorithm, robust, very well tested
228f6130c64SBarry Smith $   PCGAMGGEO - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not well tested
229f6130c64SBarry Smith $   PCGAMGCLASSICAL - classical algebraic multigrid preconditioner, incomplete, poorly tested
230f6130c64SBarry Smith 
231b0753f9dSMatthew G. Knepley .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseInterpolation()
232f0fc11ceSJed Brown J*/
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:
248c1cbb1deSBarry Smith +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycleType()
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 
260c1cbb1deSBarry Smith .seealso: PCMGSetType(), PCMGSetCycleType(), PCMGSetCycleTypeOnLevel()
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 #define PC_MG_CASCADE PC_MG_KASKADE;
265b0753f9dSMatthew G. Knepley 
266b0753f9dSMatthew G. Knepley /*E
267b0753f9dSMatthew G. Knepley     PCMGCycleType - Use V-cycle or W-cycle
268b0753f9dSMatthew G. Knepley 
269b0753f9dSMatthew G. Knepley    Level: beginner
270b0753f9dSMatthew G. Knepley 
271b0753f9dSMatthew G. Knepley    Values:
272b0753f9dSMatthew G. Knepley +  PC_MG_V_CYCLE
273b0753f9dSMatthew G. Knepley -  PC_MG_W_CYCLE
274b0753f9dSMatthew G. Knepley 
275b0753f9dSMatthew G. Knepley .seealso: PCMGSetCycleType()
276b0753f9dSMatthew G. Knepley 
277b0753f9dSMatthew G. Knepley E*/
278b0753f9dSMatthew G. Knepley typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
279b0753f9dSMatthew G. Knepley 
280b0753f9dSMatthew G. Knepley /*E
2812134b1e4SBarry Smith     PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process
2822134b1e4SBarry Smith 
2832134b1e4SBarry Smith    Level: beginner
2842134b1e4SBarry Smith 
2852134b1e4SBarry Smith    Values:
2862134b1e4SBarry Smith +  PC_MG_GALERKIN_PMAT - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid
2872134b1e4SBarry Smith .  PC_MG_GALERKIN_MAT -  computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid
2882134b1e4SBarry Smith .  PC_MG_GALERKIN_BOTH - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once
2892134b1e4SBarry Smith -  PC_MG_GALERKIN_NONE - neither operator is computed via the Galerkin process, the user must provide the operator
2902134b1e4SBarry Smith 
2912134b1e4SBarry Smith    Users should never set PC_MG_GALERKIN_EXTERNAL, it is used by GAMG and ML
2922134b1e4SBarry Smith 
2932134b1e4SBarry Smith .seealso: PCMGSetCycleType()
2942134b1e4SBarry Smith 
2952134b1e4SBarry Smith E*/
2962134b1e4SBarry Smith typedef enum { PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, PC_MG_GALERKIN_EXTERNAL} PCMGGalerkinType;
2972134b1e4SBarry Smith 
2982134b1e4SBarry Smith /*E
299b0753f9dSMatthew G. Knepley     PCExoticType - Face based or wirebasket based coarse grid space
300b0753f9dSMatthew G. Knepley 
301b0753f9dSMatthew G. Knepley    Level: beginner
302b0753f9dSMatthew G. Knepley 
303b0753f9dSMatthew G. Knepley .seealso: PCExoticSetType(), PCEXOTIC
304b0753f9dSMatthew G. Knepley E*/
305b0753f9dSMatthew G. Knepley typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
306b0753f9dSMatthew G. Knepley 
3078c1cd74cSHong Zhang /*E
308bc960bbfSJed Brown    PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains.
309bc960bbfSJed Brown 
310bc960bbfSJed Brown    Level: intermediate
311bc960bbfSJed Brown 
312bc960bbfSJed Brown    Values:
313bc960bbfSJed Brown +  PC_BDDC_INTERFACE_EXT_DIRICHLET - solves Dirichlet interior problem; this is the standard BDDC algorithm
314bc960bbfSJed Brown -  PC_BDDC_INTERFACE_EXT_LUMP - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP"
315bc960bbfSJed Brown 
316bc960bbfSJed Brown E*/
317bc960bbfSJed Brown typedef enum {
318bc960bbfSJed Brown   PC_BDDC_INTERFACE_EXT_DIRICHLET,
319bc960bbfSJed Brown   PC_BDDC_INTERFACE_EXT_LUMP
320bc960bbfSJed Brown } PCBDDCInterfaceExtType;
321bc960bbfSJed Brown 
322bc960bbfSJed Brown /*E
323ee6acaa3SMatthew G. Knepley     PCPatchConstructType - The algorithm used to construct patches for the preconditioner
3244bbf5ea8SMatthew G. Knepley 
3254bbf5ea8SMatthew G. Knepley    Level: beginner
3264bbf5ea8SMatthew G. Knepley 
3274bbf5ea8SMatthew G. Knepley .seealso: PCPatchSetConstructType(), PCEXOTIC
3284bbf5ea8SMatthew G. Knepley E*/
329e5b9877fSPatrick Farrell typedef enum {PC_PATCH_STAR, PC_PATCH_VANKA, PC_PATCH_PARDECOMP, PC_PATCH_USER, PC_PATCH_PYTHON} PCPatchConstructType;
3304bbf5ea8SMatthew G. Knepley 
3314bbf5ea8SMatthew G. Knepley /*E
332e26ad46dSJakub Kruzik     PCDeflationSpaceType - Type of deflation
333e662bc50SJakub Kruzik 
334e662bc50SJakub Kruzik     Values:
3357e9012c5SJakub Kruzik +   PC_DEFLATION_SPACE_HAAR        - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off
336e26ad46dSJakub Kruzik .   PC_DEFLATION_SPACE_DB2         - MATCOMPOSITE of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet)
337e26ad46dSJakub Kruzik .   PC_DEFLATION_SPACE_DB4         - same as above, but with db4 (4 coefficient Daubechies)
338e26ad46dSJakub Kruzik .   PC_DEFLATION_SPACE_DB8         - same as above, but with db8 (8 coefficient Daubechies)
339e26ad46dSJakub Kruzik .   PC_DEFLATION_SPACE_DB16        - same as above, but with db16 (16 coefficient Daubechies)
340e26ad46dSJakub Kruzik .   PC_DEFLATION_SPACE_BIORTH22    - same as above, but with biorthogonal 2.2 (6 coefficients)
341e26ad46dSJakub Kruzik .   PC_DEFLATION_SPACE_MEYER       - same as above, but with Meyer/FIR (62 coefficients)
342e26ad46dSJakub Kruzik .   PC_DEFLATION_SPACE_AGGREGATION - aggregates local indices (given by operator matix distribution) into a subdomain
343e662bc50SJakub Kruzik -   PC_DEFLATION_SPACE_USER        - indicates space set by user
344e662bc50SJakub Kruzik 
345e662bc50SJakub Kruzik     Notes:
3467e9012c5SJakub Kruzik       Wavelet-based space (except Haar) can be used in multilevel deflation.
347e662bc50SJakub Kruzik 
3481fdb00f9SJakub Kruzik     Level: intermediate
3491fdb00f9SJakub Kruzik 
3501fdb00f9SJakub Kruzik .seealso: PCDeflationSetSpaceToCompute(), PCDEFLATION
351e662bc50SJakub Kruzik E*/
352e662bc50SJakub Kruzik typedef enum {
353e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_HAAR,
354e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB2,
355e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB4,
356e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB8,
357e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_DB16,
358e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_BIORTH22,
359e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_MEYER,
360e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_AGGREGATION,
361e662bc50SJakub Kruzik   PC_DEFLATION_SPACE_USER
362e662bc50SJakub Kruzik } PCDeflationSpaceType;
363e662bc50SJakub Kruzik 
364e662bc50SJakub Kruzik /*E
36538cfc46eSPierre Jolivet     PCHPDDMCoarseCorrectionType - Type of coarse correction used by PCHPDDM
36638cfc46eSPierre Jolivet 
36738cfc46eSPierre Jolivet     Level: intermediate
36838cfc46eSPierre Jolivet 
36938cfc46eSPierre Jolivet     Values:
370e5cd072eSprj- +   PC_HPDDM_COARSE_CORRECTION_DEFLATED (default) - eq. (1) in PCHPDDMShellApply()
371e5cd072eSprj- .   PC_HPDDM_COARSE_CORRECTION_ADDITIVE - eq. (2)
372e5cd072eSprj- -   PC_HPDDM_COARSE_CORRECTION_BALANCED - eq. (3)
37338cfc46eSPierre Jolivet 
37438cfc46eSPierre Jolivet .seealso: PCHPDDM, PCSetType(), PCHPDDMShellApply()
37538cfc46eSPierre Jolivet E*/
37638cfc46eSPierre Jolivet typedef enum { PC_HPDDM_COARSE_CORRECTION_DEFLATED, PC_HPDDM_COARSE_CORRECTION_ADDITIVE, PC_HPDDM_COARSE_CORRECTION_BALANCED } PCHPDDMCoarseCorrectionType;
37738cfc46eSPierre Jolivet 
37838cfc46eSPierre Jolivet /*E
3798c1cd74cSHong Zhang     PCFailedReason - indicates type of PC failure
3808c1cd74cSHong Zhang 
3818c1cd74cSHong Zhang     Level: beginner
3828c1cd74cSHong Zhang 
3838c1cd74cSHong Zhang     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
3848c1cd74cSHong Zhang E*/
385284f0211SHong Zhang typedef enum {PC_NOERROR,PC_FACTOR_STRUCT_ZEROPIVOT,PC_FACTOR_NUMERIC_ZEROPIVOT,PC_FACTOR_OUTMEMORY,PC_FACTOR_OTHER,PC_SUBPC_ERROR} PCFailedReason;
386ce7c7f2fSMark Adams 
387ce7c7f2fSMark Adams /*E
388ce7c7f2fSMark Adams     PCGAMGLayoutType - Layout for reduced grids
389ce7c7f2fSMark Adams 
390ce7c7f2fSMark Adams     Level: intermediate
391ce7c7f2fSMark Adams 
392f73ac439SMark Adams .seealso: PCGAMGSetCoarseGridLayoutType()
393ce7c7f2fSMark Adams     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
394ce7c7f2fSMark Adams E*/
395a0095786SMark typedef enum {PCGAMG_LAYOUT_COMPACT,PCGAMG_LAYOUT_SPREAD} PCGAMGLayoutType;
396ce7c7f2fSMark Adams 
397b0753f9dSMatthew G. Knepley #endif
398