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