xref: /petsc/include/petscpctypes.h (revision 76be6f4ff3bd4e251c19fc00ebbebfd58b6e7589)
1 #if !defined(PETSCPCTYPES_H)
2 #define PETSCPCTYPES_H
3 
4 /* SUBMANSEC = PC */
5 
6 /*S
7      PC - Abstract PETSc object that manages all preconditioners including direct solvers such as PCLU
8 
9    Level: beginner
10 
11 .seealso: `PCCreate()`, `PCSetType()`, `PCType`
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 PCQR              "qr"
33 #define PCSHELL           "shell"
34 #define PCBJACOBI         "bjacobi"
35 #define PCMG              "mg"
36 #define PCEISENSTAT       "eisenstat"
37 #define PCILU             "ilu"
38 #define PCICC             "icc"
39 #define PCASM             "asm"
40 #define PCGASM            "gasm"
41 #define PCKSP             "ksp"
42 #define PCBJKOKKOS        "bjkokkos"
43 #define PCCOMPOSITE       "composite"
44 #define PCREDUNDANT       "redundant"
45 #define PCSPAI            "spai"
46 #define PCNN              "nn"
47 #define PCCHOLESKY        "cholesky"
48 #define PCPBJACOBI        "pbjacobi"
49 #define PCVPBJACOBI       "vpbjacobi"
50 #define PCMAT             "mat"
51 #define PCHYPRE           "hypre"
52 #define PCPARMS           "parms"
53 #define PCFIELDSPLIT      "fieldsplit"
54 #define PCTFS             "tfs"
55 #define PCML              "ml"
56 #define PCGALERKIN        "galerkin"
57 #define PCEXOTIC          "exotic"
58 #define PCCP              "cp"
59 #define PCBFBT            "bfbt"
60 #define PCLSC             "lsc"
61 #define PCPYTHON          "python"
62 #define PCPFMG            "pfmg"
63 #define PCSYSPFMG         "syspfmg"
64 #define PCREDISTRIBUTE    "redistribute"
65 #define PCSVD             "svd"
66 #define PCGAMG            "gamg"
67 #define PCCHOWILUVIENNACL "chowiluviennacl"
68 #define PCROWSCALINGVIENNACL "rowscalingviennacl"
69 #define PCSAVIENNACL      "saviennacl"
70 #define PCBDDC            "bddc"
71 #define PCKACZMARZ        "kaczmarz"
72 #define PCTELESCOPE       "telescope"
73 #define PCPATCH           "patch"
74 #define PCLMVM            "lmvm"
75 #define PCHMG             "hmg"
76 #define PCDEFLATION       "deflation"
77 #define PCHPDDM           "hpddm"
78 #define PCH2OPUS          "h2opus"
79 
80 /*E
81     PCSide - If the preconditioner is to be applied to the left, right
82      or symmetrically around the operator.
83 
84    Level: beginner
85 
86 .seealso:
87 E*/
88 typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
89 #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
90 
91 /*E
92     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
93 
94    Level: advanced
95 
96    Notes:
97     this must match petsc/finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
98 
99 .seealso: `PCApplyRichardson()`
100 E*/
101 typedef enum {
102               PCRICHARDSON_CONVERGED_RTOL               =  2,
103               PCRICHARDSON_CONVERGED_ATOL               =  3,
104               PCRICHARDSON_CONVERGED_ITS                =  4,
105               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
106 
107 /*E
108     PCJacobiType - What elements are used to form the Jacobi preconditioner
109 
110    Level: intermediate
111 
112 .seealso:
113 E*/
114 typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType;
115 
116 /*E
117     PCASMType - Type of additive Schwarz method to use
118 
119 $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
120 $                        and computed values in ghost regions are added together.
121 $                        Classical standard additive Schwarz.
122 $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
123 $                        region are discarded.
124 $                        Default.
125 $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
126 $                        region are added back in.
127 $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
128 $                        discarded.
129 $                        Not very good.
130 
131    Level: beginner
132 
133 .seealso: `PCASMSetType()`
134 E*/
135 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
136 
137 /*E
138     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
139 
140    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
141    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
142    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
143    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
144 
145 $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
146 $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
147 $                        from neighboring subdomains.
148 $                        Classical standard additive Schwarz.
149 $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
150 $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
151 $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
152 $                        assumption).
153 $                        Default.
154 $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
155 $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
156 $                        from neighboring subdomains.
157 $
158 $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
159 $                        Not very good.
160 
161    Level: beginner
162 
163 .seealso: `PCGASMSetType()`
164 E*/
165 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
166 
167 /*E
168     PCCompositeType - Determines how two or more preconditioner are composed
169 
170 $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
171 $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
172 $                                computed after the previous preconditioner application
173 $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
174 $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners)
175 $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
176 $                         where first preconditioner is built from alpha I + S and second from
177 $                         alpha I + R
178 
179    Level: beginner
180 
181 .seealso: `PCCompositeSetType()`
182 E*/
183 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR,PC_COMPOSITE_GKB} PCCompositeType;
184 
185 /*E
186     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
187 
188     Level: intermediate
189 
190 .seealso: `PCFieldSplitSetSchurPre()`
191 E*/
192 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;
193 
194 /*E
195     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
196 
197     Level: intermediate
198 
199 .seealso: `PCFieldSplitSetSchurFactType()`
200 E*/
201 typedef enum {
202   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
203   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
204   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
205   PC_FIELDSPLIT_SCHUR_FACT_FULL
206 } PCFieldSplitSchurFactType;
207 
208 /*E
209     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
210 
211     Level: intermediate
212 
213 .seealso: `PCPARMSSetGlobal()`
214 E*/
215 typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
216 
217 /*E
218     PCPARMSLocalType - Determines the local preconditioner method in PARMS
219 
220     Level: intermediate
221 
222 .seealso: `PCPARMSSetLocal()`
223 E*/
224 typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
225 
226 /*J
227     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method
228 
229     Level: intermediate
230 
231 $   PCGAMGAGG - (the default) smoothed aggregation algorithm, robust, very well tested
232 $   PCGAMGGEO - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not well tested
233 $   PCGAMGCLASSICAL - classical algebraic multigrid preconditioner, incomplete, poorly tested
234 
235 .seealso: `PCMG`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()`
236 J*/
237 typedef const char *PCGAMGType;
238 #define PCGAMGAGG         "agg"
239 #define PCGAMGGEO         "geo"
240 #define PCGAMGCLASSICAL   "classical"
241 
242 typedef const char *PCGAMGClassicalType;
243 #define PCGAMGCLASSICALDIRECT   "direct"
244 #define PCGAMGCLASSICALSTANDARD "standard"
245 
246 /*E
247     PCMGType - Determines the type of multigrid method that is run.
248 
249    Level: beginner
250 
251    Values:
252 +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycleType()
253 .  PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
254                 smoothed before updating the residual. This only uses the
255                 down smoother, in the preconditioner the upper smoother is ignored
256 .  PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
257             that is starts on the coarsest grid, performs a cycle, interpolates
258             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
259             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
260             calls the V-cycle only on the coarser level and has a post-smoother instead.
261 -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
262                from a finer
263 
264 .seealso: `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()`
265 
266 E*/
267 typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
268 #define PC_MG_CASCADE PC_MG_KASKADE;
269 
270 /*E
271     PCMGCycleType - Use V-cycle or W-cycle
272 
273    Level: beginner
274 
275    Values:
276 +  PC_MG_V_CYCLE - use the v cycle
277 -  PC_MG_W_CYCLE - use the w cycle
278 
279 .seealso: `PCMGSetCycleType()`
280 
281 E*/
282 typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
283 
284 /*E
285     PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process
286 
287    Level: beginner
288 
289    Values:
290 +  PC_MG_GALERKIN_PMAT - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid
291 .  PC_MG_GALERKIN_MAT -  computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid
292 .  PC_MG_GALERKIN_BOTH - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once
293 -  PC_MG_GALERKIN_NONE - neither operator is computed via the Galerkin process, the user must provide the operator
294 
295    Users should never set PC_MG_GALERKIN_EXTERNAL, it is used by GAMG and ML
296 
297 .seealso: `PCMGSetCycleType()`
298 
299 E*/
300 typedef enum { PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, PC_MG_GALERKIN_EXTERNAL} PCMGGalerkinType;
301 
302 /*E
303     PCExoticType - Face based or wirebasket based coarse grid space
304 
305    Level: beginner
306 
307 .seealso: `PCExoticSetType()`, `PCEXOTIC`
308 E*/
309 typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
310 
311 /*E
312    PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains.
313 
314    Level: intermediate
315 
316    Values:
317 +  PC_BDDC_INTERFACE_EXT_DIRICHLET - solves Dirichlet interior problem; this is the standard BDDC algorithm
318 -  PC_BDDC_INTERFACE_EXT_LUMP - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP"
319 
320 E*/
321 typedef enum {
322   PC_BDDC_INTERFACE_EXT_DIRICHLET,
323   PC_BDDC_INTERFACE_EXT_LUMP
324 } PCBDDCInterfaceExtType;
325 
326 /*E
327   PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation
328 
329   Level: beginner
330 
331 .seealso: `PCMGSetAdaptCoarseSpaceType()`, `PCMG`
332 E*/
333 typedef enum { PCMG_POLYNOMIAL, PCMG_HARMONIC, PCMG_EIGENVECTOR, PCMG_GENERALIZED_EIGENVECTOR } PCMGCoarseSpaceType;
334 
335 /*E
336     PCPatchConstructType - The algorithm used to construct patches for the preconditioner
337 
338    Level: beginner
339 
340 .seealso: `PCPatchSetConstructType()`, `PCEXOTIC`
341 E*/
342 typedef enum {PC_PATCH_STAR, PC_PATCH_VANKA, PC_PATCH_PARDECOMP, PC_PATCH_USER, PC_PATCH_PYTHON} PCPatchConstructType;
343 
344 /*E
345     PCDeflationSpaceType - Type of deflation
346 
347     Values:
348 +   PC_DEFLATION_SPACE_HAAR        - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off
349 .   PC_DEFLATION_SPACE_DB2         - MATCOMPOSITE of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet)
350 .   PC_DEFLATION_SPACE_DB4         - same as above, but with db4 (4 coefficient Daubechies)
351 .   PC_DEFLATION_SPACE_DB8         - same as above, but with db8 (8 coefficient Daubechies)
352 .   PC_DEFLATION_SPACE_DB16        - same as above, but with db16 (16 coefficient Daubechies)
353 .   PC_DEFLATION_SPACE_BIORTH22    - same as above, but with biorthogonal 2.2 (6 coefficients)
354 .   PC_DEFLATION_SPACE_MEYER       - same as above, but with Meyer/FIR (62 coefficients)
355 .   PC_DEFLATION_SPACE_AGGREGATION - aggregates local indices (given by operator matix distribution) into a subdomain
356 -   PC_DEFLATION_SPACE_USER        - indicates space set by user
357 
358     Notes:
359       Wavelet-based space (except Haar) can be used in multilevel deflation.
360 
361     Level: intermediate
362 
363 .seealso: `PCDeflationSetSpaceToCompute()`, `PCDEFLATION`
364 E*/
365 typedef enum {
366   PC_DEFLATION_SPACE_HAAR,
367   PC_DEFLATION_SPACE_DB2,
368   PC_DEFLATION_SPACE_DB4,
369   PC_DEFLATION_SPACE_DB8,
370   PC_DEFLATION_SPACE_DB16,
371   PC_DEFLATION_SPACE_BIORTH22,
372   PC_DEFLATION_SPACE_MEYER,
373   PC_DEFLATION_SPACE_AGGREGATION,
374   PC_DEFLATION_SPACE_USER
375 } PCDeflationSpaceType;
376 
377 /*E
378     PCHPDDMCoarseCorrectionType - Type of coarse correction used by PCHPDDM
379 
380     Level: intermediate
381 
382     Values:
383 +   PC_HPDDM_COARSE_CORRECTION_DEFLATED (default) - eq. (1) in PCHPDDMShellApply()
384 .   PC_HPDDM_COARSE_CORRECTION_ADDITIVE - eq. (2)
385 -   PC_HPDDM_COARSE_CORRECTION_BALANCED - eq. (3)
386 
387 .seealso: `PCHPDDM`, `PCSetType()`, `PCHPDDMShellApply()`
388 E*/
389 typedef enum { PC_HPDDM_COARSE_CORRECTION_DEFLATED, PC_HPDDM_COARSE_CORRECTION_ADDITIVE, PC_HPDDM_COARSE_CORRECTION_BALANCED } PCHPDDMCoarseCorrectionType;
390 
391 /*E
392     PCFailedReason - indicates type of PC failure
393 
394     Level: beginner
395 
396     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
397 E*/
398 typedef enum {PC_SETUP_ERROR = -1,PC_NOERROR,PC_FACTOR_STRUCT_ZEROPIVOT,PC_FACTOR_NUMERIC_ZEROPIVOT,PC_FACTOR_OUTMEMORY,PC_FACTOR_OTHER,PC_SUBPC_ERROR} PCFailedReason;
399 
400 /*E
401     PCGAMGLayoutType - Layout for reduced grids
402 
403     Level: intermediate
404 
405 .seealso: `PCGAMGSetCoarseGridLayoutType()`
406 
407     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
408 E*/
409 typedef enum {PCGAMG_LAYOUT_COMPACT,PCGAMG_LAYOUT_SPREAD} PCGAMGLayoutType;
410 
411 #endif
412