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