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