xref: /petsc/include/petscpctypes.h (revision 58c0e5077dcf40d6a880c19c87f1075ae1d22c8e)
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 #define PCHMG             "hmg"
72 #define PCDEFLATION       "deflation"
73 
74 /*E
75     PCSide - If the preconditioner is to be applied to the left, right
76      or symmetrically around the operator.
77 
78    Level: beginner
79 
80 .seealso:
81 E*/
82 typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
83 #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
84 
85 /*E
86     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
87 
88    Level: advanced
89 
90    Notes:
91     this must match petsc/finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
92 
93 .seealso: PCApplyRichardson()
94 E*/
95 typedef enum {
96               PCRICHARDSON_CONVERGED_RTOL               =  2,
97               PCRICHARDSON_CONVERGED_ATOL               =  3,
98               PCRICHARDSON_CONVERGED_ITS                =  4,
99               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
100 
101 /*E
102     PCJacobiType - What elements are used to form the Jacobi preconditioner
103 
104    Level: intermediate
105 
106 .seealso:
107 E*/
108 typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType;
109 
110 /*E
111     PCASMType - Type of additive Schwarz method to use
112 
113 $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
114 $                        and computed values in ghost regions are added together.
115 $                        Classical standard additive Schwarz.
116 $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
117 $                        region are discarded.
118 $                        Default.
119 $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
120 $                        region are added back in.
121 $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
122 $                        discarded.
123 $                        Not very good.
124 
125    Level: beginner
126 
127 .seealso: PCASMSetType()
128 E*/
129 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
130 
131 /*E
132     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
133 
134    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
135    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
136    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
137    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
138 
139 $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
140 $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
141 $                        from neighboring subdomains.
142 $                        Classical standard additive Schwarz.
143 $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
144 $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
145 $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
146 $                        assumption).
147 $                        Default.
148 $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
149 $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
150 $                        from neighboring subdomains.
151 $
152 $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
153 $                        Not very good.
154 
155    Level: beginner
156 
157 .seealso: PCGASMSetType()
158 E*/
159 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
160 
161 /*E
162     PCCompositeType - Determines how two or more preconditioner are composed
163 
164 $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
165 $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
166 $                                computed after the previous preconditioner application
167 $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
168 $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners)
169 $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
170 $                         where first preconditioner is built from alpha I + S and second from
171 $                         alpha I + R
172 
173    Level: beginner
174 
175 .seealso: PCCompositeSetType()
176 E*/
177 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR,PC_COMPOSITE_GKB} PCCompositeType;
178 
179 /*E
180     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
181 
182     Level: intermediate
183 
184 .seealso: PCFieldSplitSetSchurPre()
185 E*/
186 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;
187 
188 /*E
189     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
190 
191     Level: intermediate
192 
193 .seealso: PCFieldSplitSetSchurFactType()
194 E*/
195 typedef enum {
196   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
197   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
198   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
199   PC_FIELDSPLIT_SCHUR_FACT_FULL
200 } PCFieldSplitSchurFactType;
201 
202 /*E
203     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
204 
205     Level: intermediate
206 
207 .seealso: PCPARMSSetGlobal()
208 E*/
209 typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
210 
211 /*E
212     PCPARMSLocalType - Determines the local preconditioner method in PARMS
213 
214     Level: intermediate
215 
216 .seealso: PCPARMSSetLocal()
217 E*/
218 typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
219 
220 /*E
221     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method
222 
223     Level: intermediate
224 
225 $   PCGAMGAGG - (the default) smoothed aggregation algorithm, robust, very well tested
226 $   PCGAMGGEO - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not well tested
227 $   PCGAMGCLASSICAL - classical algebraic multigrid preconditioner, incomplete, poorly tested
228 
229 .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseInterpolation()
230 E*/
231 typedef const char *PCGAMGType;
232 #define PCGAMGAGG         "agg"
233 #define PCGAMGGEO         "geo"
234 #define PCGAMGCLASSICAL   "classical"
235 
236 typedef const char *PCGAMGClassicalType;
237 #define PCGAMGCLASSICALDIRECT   "direct"
238 #define PCGAMGCLASSICALSTANDARD "standard"
239 
240 /*E
241     PCMGType - Determines the type of multigrid method that is run.
242 
243    Level: beginner
244 
245    Values:
246 +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycleType()
247 .  PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
248                 smoothed before updating the residual. This only uses the
249                 down smoother, in the preconditioner the upper smoother is ignored
250 .  PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
251             that is starts on the coarsest grid, performs a cycle, interpolates
252             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
253             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
254             calls the V-cycle only on the coarser level and has a post-smoother instead.
255 -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
256                from a finer
257 
258 .seealso: PCMGSetType(), PCMGSetCycleType(), PCMGSetCycleTypeOnLevel()
259 
260 E*/
261 typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
262 #define PC_MG_CASCADE PC_MG_KASKADE;
263 
264 /*E
265     PCMGCycleType - Use V-cycle or W-cycle
266 
267    Level: beginner
268 
269    Values:
270 +  PC_MG_V_CYCLE
271 -  PC_MG_W_CYCLE
272 
273 .seealso: PCMGSetCycleType()
274 
275 E*/
276 typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
277 
278 /*E
279     PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process
280 
281    Level: beginner
282 
283    Values:
284 +  PC_MG_GALERKIN_PMAT - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid
285 .  PC_MG_GALERKIN_MAT -  computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid
286 .  PC_MG_GALERKIN_BOTH - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once
287 -  PC_MG_GALERKIN_NONE - neither operator is computed via the Galerkin process, the user must provide the operator
288 
289    Users should never set PC_MG_GALERKIN_EXTERNAL, it is used by GAMG and ML
290 
291 .seealso: PCMGSetCycleType()
292 
293 E*/
294 typedef enum { PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, PC_MG_GALERKIN_EXTERNAL} PCMGGalerkinType;
295 
296 /*E
297     PCExoticType - Face based or wirebasket based coarse grid space
298 
299    Level: beginner
300 
301 .seealso: PCExoticSetType(), PCEXOTIC
302 E*/
303 typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
304 
305 /*E
306     PCPatchConstructType - The algorithm used to construct patches for the preconditioner
307 
308    Level: beginner
309 
310 .seealso: PCPatchSetConstructType(), PCEXOTIC
311 E*/
312 typedef enum {PC_PATCH_STAR, PC_PATCH_VANKA, PC_PATCH_PARDECOMP, PC_PATCH_USER, PC_PATCH_PYTHON} PCPatchConstructType;
313 
314 /*E
315     PCDeflationSpaceType - Type of deflation
316 
317     Values:
318 +   PC_DEFLATION_SPACE_HAAR        - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off
319 .   PC_DEFLATION_SPACE_DB2         - MATCOMPOSITE of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet)
320 .   PC_DEFLATION_SPACE_DB4         - same as above, but with db4 (4 coefficient Daubechies)
321 .   PC_DEFLATION_SPACE_DB8         - same as above, but with db8 (8 coefficient Daubechies)
322 .   PC_DEFLATION_SPACE_DB16        - same as above, but with db16 (16 coefficient Daubechies)
323 .   PC_DEFLATION_SPACE_BIORTH22    - same as above, but with biorthogonal 2.2 (6 coefficients)
324 .   PC_DEFLATION_SPACE_MEYER       - same as above, but with Meyer/FIR (62 coefficients)
325 .   PC_DEFLATION_SPACE_AGGREGATION - aggregates local indices (given by operator matix distribution) into a subdomain
326 -   PC_DEFLATION_SPACE_USER        - indicates space set by user
327 
328     Notes:
329       Wavelet-based space (except Haar) can be used in multilevel deflation.
330 
331     Level: intermediate
332 
333 .seealso: PCDeflationSetSpaceToCompute(), PCDEFLATION
334 E*/
335 typedef enum {
336   PC_DEFLATION_SPACE_HAAR,
337   PC_DEFLATION_SPACE_DB2,
338   PC_DEFLATION_SPACE_DB4,
339   PC_DEFLATION_SPACE_DB8,
340   PC_DEFLATION_SPACE_DB16,
341   PC_DEFLATION_SPACE_BIORTH22,
342   PC_DEFLATION_SPACE_MEYER,
343   PC_DEFLATION_SPACE_AGGREGATION,
344   PC_DEFLATION_SPACE_USER
345 } PCDeflationSpaceType;
346 
347 /*E
348     PCFailedReason - indicates type of PC failure
349 
350     Level: beginner
351 
352     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
353 E*/
354 typedef enum {PC_NOERROR,PC_FACTOR_STRUCT_ZEROPIVOT,PC_FACTOR_NUMERIC_ZEROPIVOT,PC_FACTOR_OUTMEMORY,PC_FACTOR_OTHER,PC_SUBPC_ERROR} PCFailedReason;
355 #endif
356