Actual source code: petscpctypes.h
petsc-3.9.4 2018-09-11
1: #if !defined(_PETSCPCTYPES_H)
2: #define _PETSCPCTYPES_H
4: #include <petscdmtypes.h>
6: /*S
7: PC - Abstract PETSc object that manages all preconditioners including direct solvers such as PCLU
9: Level: beginner
11: Concepts: preconditioners
13: .seealso: PCCreate(), PCSetType(), PCType (for list of available types)
14: S*/
15: typedef struct _p_PC* PC;
17: /*J
18: PCType - String with the name of a PETSc preconditioner method.
20: Level: beginner
22: Notes: Click on the links above to see details on a particular solver
24: PCRegister() is used to register preconditioners that are then accessible via PCSetType()
26: .seealso: PCSetType(), PC, PCCreate(), PCRegister(), PCSetFromOptions()
27: J*/
28: typedef const char* PCType;
29: #define PCNONE "none"
30: #define PCJACOBI "jacobi"
31: #define PCSOR "sor"
32: #define PCLU "lu"
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 PCCOMPOSITE "composite"
43: #define PCREDUNDANT "redundant"
44: #define PCSPAI "spai"
45: #define PCNN "nn"
46: #define PCCHOLESKY "cholesky"
47: #define PCPBJACOBI "pbjacobi"
48: #define PCMAT "mat"
49: #define PCHYPRE "hypre"
50: #define PCPARMS "parms"
51: #define PCFIELDSPLIT "fieldsplit"
52: #define PCTFS "tfs"
53: #define PCML "ml"
54: #define PCGALERKIN "galerkin"
55: #define PCEXOTIC "exotic"
56: #define PCCP "cp"
57: #define PCBFBT "bfbt"
58: #define PCLSC "lsc"
59: #define PCPYTHON "python"
60: #define PCPFMG "pfmg"
61: #define PCSYSPFMG "syspfmg"
62: #define PCREDISTRIBUTE "redistribute"
63: #define PCSVD "svd"
64: #define PCGAMG "gamg"
65: #define PCCHOWILUVIENNACL "chowiluviennacl"
66: #define PCROWSCALINGVIENNACL "rowscalingviennacl"
67: #define PCSAVIENNACL "saviennacl"
68: #define PCBDDC "bddc"
69: #define PCKACZMARZ "kaczmarz"
70: #define PCTELESCOPE "telescope"
72: /*E
73: PCSide - If the preconditioner is to be applied to the left, right
74: or symmetrically around the operator.
76: Level: beginner
78: .seealso:
79: E*/
80: typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
81: #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
82: PETSC_EXTERN const char *const *const PCSides;
84: /*E
85: PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
87: Level: advanced
89: Notes: this must match petsc/finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
91: .seealso: PCApplyRichardson()
92: E*/
93: typedef enum {
94: PCRICHARDSON_CONVERGED_RTOL = 2,
95: PCRICHARDSON_CONVERGED_ATOL = 3,
96: PCRICHARDSON_CONVERGED_ITS = 4,
97: PCRICHARDSON_DIVERGED_DTOL = -4} PCRichardsonConvergedReason;
99: /*E
100: PCJacobiType - What elements are used to form the Jacobi preconditioner
102: Level: intermediate
104: .seealso:
105: E*/
106: typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType;
107: PETSC_EXTERN const char *const PCJacobiTypes[];
109: /*E
110: PCASMType - Type of additive Schwarz method to use
112: $ PC_ASM_BASIC - Symmetric version where residuals from the ghost points are used
113: $ and computed values in ghost regions are added together.
114: $ Classical standard additive Schwarz.
115: $ PC_ASM_RESTRICT - Residuals from ghost points are used but computed values in ghost
116: $ region are discarded.
117: $ Default.
118: $ PC_ASM_INTERPOLATE - Residuals from ghost points are not used, computed values in ghost
119: $ region are added back in.
120: $ PC_ASM_NONE - Residuals from ghost points are not used, computed ghost values are
121: $ discarded.
122: $ Not very good.
124: Level: beginner
126: .seealso: PCASMSetType()
127: E*/
128: typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
129: PETSC_EXTERN const char *const PCASMTypes[];
131: /*E
132: PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
134: Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational
135: domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute
136: a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
137: (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
139: $ PC_GASM_BASIC - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
140: $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections
141: $ from neighboring subdomains.
142: $ Classical standard additive Schwarz.
143: $ PC_GASM_RESTRICT - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
144: $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result,
145: $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
146: $ assumption).
147: $ Default.
148: $ PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
149: $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
150: $ from neighboring subdomains.
151: $
152: $ PC_GASM_NONE - Residuals and corrections are zeroed out outside the local subdomains.
153: $ Not very good.
155: Level: beginner
157: .seealso: PCGASMSetType()
158: E*/
159: typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
160: PETSC_EXTERN const char *const PCGASMTypes[];
162: /*E
163: PCCompositeType - Determines how two or more preconditioner are composed
165: $ PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
166: $ PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
167: $ computed after the previous preconditioner application
168: $ PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
169: $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners)
170: $ PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
171: $ where first preconditioner is built from alpha I + S and second from
172: $ alpha I + R
174: Level: beginner
176: .seealso: PCCompositeSetType()
177: E*/
178: typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
179: PETSC_EXTERN const char *const PCCompositeTypes[];
181: /*E
182: PCFieldSplitSchurPreType - Determines how to precondition Schur complement
184: Level: intermediate
186: .seealso: PCFieldSplitSetSchurPre()
187: E*/
188: 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;
189: PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];
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;
204: PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];
206: /*E
207: PCPARMSGlobalType - Determines the global preconditioner method in PARMS
209: Level: intermediate
211: .seealso: PCPARMSSetGlobal()
212: E*/
213: typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
214: PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
215: /*E
216: PCPARMSLocalType - Determines the local preconditioner method in PARMS
218: Level: intermediate
220: .seealso: PCPARMSSetLocal()
221: E*/
222: typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
223: PETSC_EXTERN const char *const PCPARMSLocalTypes[];
225: /*E
226: PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method
228: Level: intermediate
230: .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseInterpolation()
231: E*/
232: typedef const char *PCGAMGType;
233: #define PCGAMGAGG "agg"
234: #define PCGAMGGEO "geo"
235: #define PCGAMGCLASSICAL "classical"
237: typedef const char *PCGAMGClassicalType;
238: #define PCGAMGCLASSICALDIRECT "direct"
239: #define PCGAMGCLASSICALSTANDARD "standard"
241: /*E
242: PCMGType - Determines the type of multigrid method that is run.
244: Level: beginner
246: Values:
247: + PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycleType()
248: . PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
249: smoothed before updating the residual. This only uses the
250: down smoother, in the preconditioner the upper smoother is ignored
251: . PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
252: that is starts on the coarsest grid, performs a cycle, interpolates
253: 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
254: algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
255: calls the V-cycle only on the coarser level and has a post-smoother instead.
256: - PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
257: from a finer
259: .seealso: PCMGSetType(), PCMGSetCycleType(), PCMGSetCycleTypeOnLevel()
261: E*/
262: typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
263: PETSC_EXTERN const char *const PCMGTypes[];
264: #define PC_MG_CASCADE PC_MG_KASKADE;
266: /*E
267: PCMGCycleType - Use V-cycle or W-cycle
269: Level: beginner
271: Values:
272: + PC_MG_V_CYCLE
273: - PC_MG_W_CYCLE
275: .seealso: PCMGSetCycleType()
277: E*/
278: typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
279: PETSC_EXTERN const char *const PCMGCycleTypes[];
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;
298: PETSC_EXTERN const char *const PCMGGalerkinTypes[];
300: /*E
301: PCExoticType - Face based or wirebasket based coarse grid space
303: Level: beginner
305: .seealso: PCExoticSetType(), PCEXOTIC
306: E*/
307: typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
308: PETSC_EXTERN const char *const PCExoticTypes[];
309: PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);
311: /*E
312: PCFailedReason - indicates type of PC failure
314: Level: beginner
316: Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
317: E*/
318: typedef enum {PC_NOERROR,PC_FACTOR_STRUCT_ZEROPIVOT,PC_FACTOR_NUMERIC_ZEROPIVOT,PC_FACTOR_OUTMEMORY,PC_FACTOR_OTHER,PC_SUBPC_ERROR} PCFailedReason;
319: PETSC_EXTERN const char *const PCFailedReasons[];
320: #endif