xref: /petsc/include/petscpctypes.h (revision 37d05b0256c1e9ba4bc423c4eccb3df226931ef0)
1 #ifndef 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: [](sec_pc), `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    Note:
21    `PCRegister()` is used to register preconditioners that are then accessible via `PCSetType()`
22 
23 .seealso: [](sec_pc), `PCSetType()`, `PC`, `PCCreate()`, `PCRegister()`, `PCSetFromOptions()`, `PCLU`, `PCJACOBI`, `PCBJACOBI`
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 PCQR                 "qr"
31 #define PCSHELL              "shell"
32 #define PCAMGX               "amgx"
33 #define PCBJACOBI            "bjacobi"
34 #define PCMG                 "mg"
35 #define PCEISENSTAT          "eisenstat"
36 #define PCILU                "ilu"
37 #define PCICC                "icc"
38 #define PCASM                "asm"
39 #define PCGASM               "gasm"
40 #define PCKSP                "ksp"
41 #define PCBJKOKKOS           "bjkokkos"
42 #define PCCOMPOSITE          "composite"
43 #define PCREDUNDANT          "redundant"
44 #define PCSPAI               "spai"
45 #define PCNN                 "nn"
46 #define PCCHOLESKY           "cholesky"
47 #define PCPBJACOBI           "pbjacobi"
48 #define PCVPBJACOBI          "vpbjacobi"
49 #define PCMAT                "mat"
50 #define PCHYPRE              "hypre"
51 #define PCPARMS              "parms"
52 #define PCFIELDSPLIT         "fieldsplit"
53 #define PCTFS                "tfs"
54 #define PCML                 "ml"
55 #define PCGALERKIN           "galerkin"
56 #define PCEXOTIC             "exotic"
57 #define PCCP                 "cp"
58 #define PCBFBT               "bfbt"
59 #define PCLSC                "lsc"
60 #define PCPYTHON             "python"
61 #define PCPFMG               "pfmg"
62 #define PCSMG                "smg"
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: [](sec_pc), `PC`
88 E*/
89 typedef enum {
90   PC_SIDE_DEFAULT = -1,
91   PC_LEFT,
92   PC_RIGHT,
93   PC_SYMMETRIC
94 } PCSide;
95 #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
96 
97 /*E
98     PCRichardsonConvergedReason - reason a `PCRICHARDSON` `PCApplyRichardson()` method terminated
99 
100    Level: advanced
101 
102    Developer Note:
103     this must match petsc/finclude/petscpc.h and the `KSPConvergedReason` values in petscksp.h
104 
105 .seealso: [](sec_pc), `PCRICHARDSON`, `PC`, `PCApplyRichardson()`
106 E*/
107 typedef enum {
108   PCRICHARDSON_CONVERGED_RTOL = 2,
109   PCRICHARDSON_CONVERGED_ATOL = 3,
110   PCRICHARDSON_CONVERGED_ITS  = 4,
111   PCRICHARDSON_DIVERGED_DTOL  = -4
112 } PCRichardsonConvergedReason;
113 
114 /*E
115     PCJacobiType - What elements are used to form the Jacobi preconditioner
116 
117    Values:
118 +  `PC_JACOBI_DIAGONAL` - use the diagonal entry, if it is zero use one
119 .  `PC_JACOBI_ROWMAX` - use the maximum absolute value in the row
120 -  `PC_JACOBI_ROWSUM` - use the sum of the values in the row (not the absolute values)
121 
122    Level: intermediate
123 
124 .seealso: [](sec_pc), `PCJACOBI`, `PC`
125 E*/
126 typedef enum {
127   PC_JACOBI_DIAGONAL,
128   PC_JACOBI_ROWMAX,
129   PC_JACOBI_ROWSUM
130 } PCJacobiType;
131 
132 /*E
133     PCASMType - Type of additive Schwarz method to use
134 
135    Values:
136 +  `PC_ASM_BASIC`        - Symmetric version where residuals from the ghost points are used
137                         and computed values in ghost regions are added together.
138                         Classical standard additive Schwarz.
139 .  `PC_ASM_RESTRICT`     - Residuals from ghost points are used but computed values in ghost
140                         region are discarded.
141                         Default.
142 .  `PC_ASM_INTERPOLATE`  - Residuals from ghost points are not used, computed values in ghost
143                         region are added back in.
144 -  `PC_ASM_NONE`         - Residuals from ghost points are not used, computed ghost values are
145                         discarded.
146                         Not very good.
147 
148    Level: beginner
149 
150 .seealso: [](sec_pc), `PC`, `PCASM`, `PCASMSetType()`
151 E*/
152 typedef enum {
153   PC_ASM_BASIC       = 3,
154   PC_ASM_RESTRICT    = 1,
155   PC_ASM_INTERPOLATE = 2,
156   PC_ASM_NONE        = 0
157 } PCASMType;
158 
159 /*E
160     PCGASMType - Type of generalized additive Schwarz method to use (differs from `PCASM` in allowing multiple processors per subdomain).
161 
162    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
163    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
164    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
165    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
166 
167    Values:
168 +  `PC_GASM_BASIC`      - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
169                         over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
170                         from neighboring subdomains.
171                         Classical standard additive Schwarz.
172 .  `PC_GASM_RESTRICT`    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
173                         (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
174                         each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
175                         assumption).
176                         Default.
177 .  `PC_GASM_INTERPOLATE` - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
178                         applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
179                         from neighboring subdomains.
180 
181 -  `PC_GASM_NONE`       - Residuals and corrections are zeroed out outside the local subdomains.
182                         Not very good.
183 
184    Level: beginner
185 
186 .seealso: [](sec_pc), `PCGASM`, `PCASM`, `PC`, `PCGASMSetType()`
187 E*/
188 typedef enum {
189   PC_GASM_BASIC       = 3,
190   PC_GASM_RESTRICT    = 1,
191   PC_GASM_INTERPOLATE = 2,
192   PC_GASM_NONE        = 0
193 } PCGASMType;
194 
195 /*E
196     PCCompositeType - Determines how two or more preconditioner are composed
197 
198   Values:
199 +  `PC_COMPOSITE_ADDITIVE` - results from application of all preconditioners are added together
200 .  `PC_COMPOSITE_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly
201                                 computed after the previous preconditioner application
202 .  `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly
203                                 computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners)
204 .  `PC_COMPOSITE_SPECIAL` - This is very special for a matrix of the form alpha I + R + S
205                          where first preconditioner is built from alpha I + S and second from
206                          alpha I + R
207 .  `PC_COMPOSITE_SCHUR` -  composes the Schur complement of the matrix from two blocks, see `PCFIELDSPLIT`
208 -  `PC_COMPOSITE_GKB` - the generalized Golub-Kahan bidiagonalization preconditioner, see `PCFIELDSPLIT`
209 
210    Level: beginner
211 
212 .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()`
213 E*/
214 typedef enum {
215   PC_COMPOSITE_ADDITIVE,
216   PC_COMPOSITE_MULTIPLICATIVE,
217   PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,
218   PC_COMPOSITE_SPECIAL,
219   PC_COMPOSITE_SCHUR,
220   PC_COMPOSITE_GKB
221 } PCCompositeType;
222 
223 /*E
224     PCFieldSplitSchurPreType - Determines how to precondition a Schur complement
225 
226     Level: intermediate
227 
228 .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurPre()`, `PC`
229 E*/
230 typedef enum {
231   PC_FIELDSPLIT_SCHUR_PRE_SELF,
232   PC_FIELDSPLIT_SCHUR_PRE_SELFP,
233   PC_FIELDSPLIT_SCHUR_PRE_A11,
234   PC_FIELDSPLIT_SCHUR_PRE_USER,
235   PC_FIELDSPLIT_SCHUR_PRE_FULL
236 } PCFieldSplitSchurPreType;
237 
238 /*E
239     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
240 
241     Level: intermediate
242 
243  .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurFactType()`, `PC`
244 E*/
245 typedef enum {
246   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
247   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
248   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
249   PC_FIELDSPLIT_SCHUR_FACT_FULL
250 } PCFieldSplitSchurFactType;
251 
252 /*E
253     PCPARMSGlobalType - Determines the global preconditioner method in `PCPARMS`
254 
255     Level: intermediate
256 
257  .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetGlobal()`, `PC`
258 E*/
259 typedef enum {
260   PC_PARMS_GLOBAL_RAS,
261   PC_PARMS_GLOBAL_SCHUR,
262   PC_PARMS_GLOBAL_BJ
263 } PCPARMSGlobalType;
264 
265 /*E
266     PCPARMSLocalType - Determines the local preconditioner method in `PCPARMS`
267 
268     Level: intermediate
269 
270 .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetLocal()`, `PC`
271 E*/
272 typedef enum {
273   PC_PARMS_LOCAL_ILU0,
274   PC_PARMS_LOCAL_ILUK,
275   PC_PARMS_LOCAL_ILUT,
276   PC_PARMS_LOCAL_ARMS
277 } PCPARMSLocalType;
278 
279 /*J
280     PCGAMGType - type of generalized algebraic multigrid `PCGAMG` method
281 
282    Values:
283 +   `PCGAMGAGG` - (the default) smoothed aggregation algorithm, robust, very well tested
284 .   `PCGAMGGEO` - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not well tested
285 -   `PCGAMGCLASSICAL` - classical algebraic multigrid preconditioner, incomplete, poorly tested
286 
287      Level: intermediate
288 
289 .seealso: [](sec_pc), `PCGAMG`, `PCMG`, `PC`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()`
290 J*/
291 typedef const char *PCGAMGType;
292 #define PCGAMGAGG       "agg"
293 #define PCGAMGGEO       "geo"
294 #define PCGAMGCLASSICAL "classical"
295 
296 typedef const char *PCGAMGClassicalType;
297 #define PCGAMGCLASSICALDIRECT   "direct"
298 #define PCGAMGCLASSICALSTANDARD "standard"
299 
300 /*E
301     PCMGType - Determines the type of multigrid method that is run.
302 
303    Level: beginner
304 
305    Values:
306 +  `PC_MG_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `PCMGSetCycleType()`
307 .  `PC_MG_ADDITIVE` - the additive multigrid preconditioner where all levels are
308                 smoothed before updating the residual. This only uses the
309                 down smoother, in the preconditioner the upper smoother is ignored
310 .  `PC_MG_FULL` - same as multiplicative except one also performs grid sequencing,
311             that is starts on the coarsest grid, performs a cycle, interpolates
312             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
313             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
314             calls the V-cycle only on the coarser level and has a post-smoother instead.
315 -  `PC_MG_KASKADE` - like full multigrid except one never goes back to a coarser level
316                from a finer
317 
318 .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()`
319 E*/
320 typedef enum {
321   PC_MG_MULTIPLICATIVE,
322   PC_MG_ADDITIVE,
323   PC_MG_FULL,
324   PC_MG_KASKADE
325 } PCMGType;
326 #define PC_MG_CASCADE PC_MG_KASKADE;
327 
328 /*E
329     PCMGCycleType - Use V-cycle or W-cycle
330 
331    Level: beginner
332 
333    Values:
334 +  `PC_MG_V_CYCLE` - use the v cycle
335 -  `PC_MG_W_CYCLE` - use the w cycle
336 
337 .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()`
338 E*/
339 typedef enum {
340   PC_MG_CYCLE_V = 1,
341   PC_MG_CYCLE_W = 2
342 } PCMGCycleType;
343 
344 /*E
345     PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process
346 
347    Values:
348 +  `PC_MG_GALERKIN_PMAT` - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid
349 .  `PC_MG_GALERKIN_MAT` -  computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid
350 .  `PC_MG_GALERKIN_BOTH` - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once
351 -  `PC_MG_GALERKIN_NONE` - neither operator is computed via the Galerkin process, the user must provide the operator
352 
353    Level: beginner
354 
355    Note:
356    Users should never set `PC_MG_GALERKIN_EXTERNAL`, it is used by `PCHYPRE` and `PCML`
357 
358 .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()`
359 E*/
360 typedef enum {
361   PC_MG_GALERKIN_BOTH,
362   PC_MG_GALERKIN_PMAT,
363   PC_MG_GALERKIN_MAT,
364   PC_MG_GALERKIN_NONE,
365   PC_MG_GALERKIN_EXTERNAL
366 } PCMGGalerkinType;
367 
368 /*E
369     PCExoticType - Face based or wirebasket based coarse grid space
370 
371    Level: beginner
372 
373 .seealso: [](sec_pc), `PCExoticSetType()`, `PCEXOTIC`
374 E*/
375 typedef enum {
376   PC_EXOTIC_FACE,
377   PC_EXOTIC_WIREBASKET
378 } PCExoticType;
379 
380 /*E
381    PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains.
382 
383    Level: intermediate
384 
385    Values:
386 +  `PC_BDDC_INTERFACE_EXT_DIRICHLET` - solves Dirichlet interior problem; this is the standard BDDC algorithm
387 -  `PC_BDDC_INTERFACE_EXT_LUMP` - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP"
388 
389 .seealso: [](sec_pc), `PCBDDC`, `PC`
390 E*/
391 typedef enum {
392   PC_BDDC_INTERFACE_EXT_DIRICHLET,
393   PC_BDDC_INTERFACE_EXT_LUMP
394 } PCBDDCInterfaceExtType;
395 
396 /*E
397   PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation
398 
399   Level: beginner
400 
401 .seealso: [](sec_pc), `PCMGSetAdaptCoarseSpaceType()`, `PCMG`, `PC`
402 E*/
403 typedef enum {
404   PCMG_ADAPT_NONE,
405   PCMG_ADAPT_POLYNOMIAL,
406   PCMG_ADAPT_HARMONIC,
407   PCMG_ADAPT_EIGENVECTOR,
408   PCMG_ADAPT_GENERALIZED_EIGENVECTOR,
409   PCMG_ADAPT_GDSW
410 } PCMGCoarseSpaceType;
411 
412 /*E
413     PCPatchConstructType - The algorithm used to construct patches for the preconditioner
414 
415    Level: beginner
416 
417 .seealso: [](sec_pc), `PCPatchSetConstructType()`, `PCPATCH`, `PC`
418 E*/
419 typedef enum {
420   PC_PATCH_STAR,
421   PC_PATCH_VANKA,
422   PC_PATCH_PARDECOMP,
423   PC_PATCH_USER,
424   PC_PATCH_PYTHON
425 } PCPatchConstructType;
426 
427 /*E
428     PCDeflationSpaceType - Type of deflation
429 
430     Values:
431 +   `PC_DEFLATION_SPACE_HAAR`        - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off
432 .   `PC_DEFLATION_SPACE_DB2`         - `MATCOMPOSITE` of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet)
433 .   `PC_DEFLATION_SPACE_DB4`         - same as above, but with db4 (4 coefficient Daubechies)
434 .   `PC_DEFLATION_SPACE_DB8`         - same as above, but with db8 (8 coefficient Daubechies)
435 .   `PC_DEFLATION_SPACE_DB16`        - same as above, but with db16 (16 coefficient Daubechies)
436 .   `PC_DEFLATION_SPACE_BIORTH22`    - same as above, but with biorthogonal 2.2 (6 coefficients)
437 .   `PC_DEFLATION_SPACE_MEYER`       - same as above, but with Meyer/FIR (62 coefficients)
438 .   `PC_DEFLATION_SPACE_AGGREGATION` - aggregates local indices (given by operator matrix distribution) into a subdomain
439 -   `PC_DEFLATION_SPACE_USER`        - indicates space set by user
440 
441     Level: intermediate
442 
443     Note:
444     Wavelet-based space (except Haar) can be used in multilevel deflation.
445 
446 .seealso: [](sec_pc), `PCDeflationSetSpaceToCompute()`, `PCDEFLATION`, `PC`
447 E*/
448 typedef enum {
449   PC_DEFLATION_SPACE_HAAR,
450   PC_DEFLATION_SPACE_DB2,
451   PC_DEFLATION_SPACE_DB4,
452   PC_DEFLATION_SPACE_DB8,
453   PC_DEFLATION_SPACE_DB16,
454   PC_DEFLATION_SPACE_BIORTH22,
455   PC_DEFLATION_SPACE_MEYER,
456   PC_DEFLATION_SPACE_AGGREGATION,
457   PC_DEFLATION_SPACE_USER
458 } PCDeflationSpaceType;
459 
460 /*E
461     PCHPDDMCoarseCorrectionType - Type of coarse correction used by `PCHPDDM`
462 
463     Level: intermediate
464 
465     Values:
466 +   `PC_HPDDM_COARSE_CORRECTION_DEFLATED` (default) - eq. (1) in `PCHPDDMShellApply()`
467 .   `PC_HPDDM_COARSE_CORRECTION_ADDITIVE` - eq. (2)
468 -   `PC_HPDDM_COARSE_CORRECTION_BALANCED` - eq. (3)
469 
470 .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCSetType()`, `PCHPDDMShellApply()`
471 E*/
472 typedef enum {
473   PC_HPDDM_COARSE_CORRECTION_DEFLATED,
474   PC_HPDDM_COARSE_CORRECTION_ADDITIVE,
475   PC_HPDDM_COARSE_CORRECTION_BALANCED
476 } PCHPDDMCoarseCorrectionType;
477 
478 /*E
479     PCFailedReason - indicates type of `PC` failure
480 
481     Level: beginner
482 
483     Developer Note:
484     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
485 
486 .seealso: [](sec_pc), `PC`
487 E*/
488 typedef enum {
489   PC_SETUP_ERROR = -1,
490   PC_NOERROR,
491   PC_FACTOR_STRUCT_ZEROPIVOT,
492   PC_FACTOR_NUMERIC_ZEROPIVOT,
493   PC_FACTOR_OUTMEMORY,
494   PC_FACTOR_OTHER,
495   PC_INCONSISTENT_RHS,
496   PC_SUBPC_ERROR
497 } PCFailedReason;
498 
499 /*E
500     PCGAMGLayoutType - Layout for reduced grids
501 
502     Level: intermediate
503 
504     Developer Note:
505     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
506 
507 .seealso: [](sec_pc), `PCGAMG`, `PC`, `PCGAMGSetCoarseGridLayoutType()`
508 E*/
509 typedef enum {
510   PCGAMG_LAYOUT_COMPACT,
511   PCGAMG_LAYOUT_SPREAD
512 } PCGAMGLayoutType;
513 
514 #endif
515