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