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