xref: /petsc/include/petscpctypes.h (revision 285fb4e2b69b3de46a0633bd0adc6a7f684caa1e)
1 #if !defined(_PETSCPCTYPES_H)
2 #define _PETSCPCTYPES_H
3 
4 #include <petscdmtypes.h>
5 
6 /*S
7      PC - Abstract PETSc object that manages all preconditioners including direct solvers such as PCLU
8 
9    Level: beginner
10 
11   Concepts: preconditioners
12 
13 .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
14 S*/
15 typedef struct _p_PC* PC;
16 
17 /*J
18     PCType - String with the name of a PETSc preconditioner method.
19 
20    Level: beginner
21 
22    Notes:
23     Click on the links above to see details on a particular solver
24 
25           PCRegister() is used to register preconditioners that are then accessible via PCSetType()
26 
27 .seealso: PCSetType(), PC, PCCreate(), PCRegister(), PCSetFromOptions()
28 J*/
29 typedef const char* PCType;
30 #define PCNONE            "none"
31 #define PCJACOBI          "jacobi"
32 #define PCSOR             "sor"
33 #define PCLU              "lu"
34 #define PCSHELL           "shell"
35 #define PCBJACOBI         "bjacobi"
36 #define PCMG              "mg"
37 #define PCEISENSTAT       "eisenstat"
38 #define PCILU             "ilu"
39 #define PCICC             "icc"
40 #define PCASM             "asm"
41 #define PCGASM            "gasm"
42 #define PCKSP             "ksp"
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 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 PCSYSPFMG         "syspfmg"
63 #define PCREDISTRIBUTE    "redistribute"
64 #define PCSVD             "svd"
65 #define PCGAMG            "gamg"
66 #define PCCHOWILUVIENNACL "chowiluviennacl"
67 #define PCROWSCALINGVIENNACL "rowscalingviennacl"
68 #define PCSAVIENNACL      "saviennacl"
69 #define PCBDDC            "bddc"
70 #define PCKACZMARZ        "kaczmarz"
71 #define PCTELESCOPE       "telescope"
72 
73 /*E
74     PCSide - If the preconditioner is to be applied to the left, right
75      or symmetrically around the operator.
76 
77    Level: beginner
78 
79 .seealso:
80 E*/
81 typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
82 #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
83 PETSC_EXTERN const char *const *const PCSides;
84 
85 /*E
86     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates
87 
88    Level: advanced
89 
90    Notes:
91     this must match petsc/finclude/petscpc.h and the KSPConvergedReason values in petscksp.h
92 
93 .seealso: PCApplyRichardson()
94 E*/
95 typedef enum {
96               PCRICHARDSON_CONVERGED_RTOL               =  2,
97               PCRICHARDSON_CONVERGED_ATOL               =  3,
98               PCRICHARDSON_CONVERGED_ITS                =  4,
99               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;
100 
101 /*E
102     PCJacobiType - What elements are used to form the Jacobi preconditioner
103 
104    Level: intermediate
105 
106 .seealso:
107 E*/
108 typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType;
109 PETSC_EXTERN const char *const PCJacobiTypes[];
110 
111 /*E
112     PCASMType - Type of additive Schwarz method to use
113 
114 $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
115 $                        and computed values in ghost regions are added together.
116 $                        Classical standard additive Schwarz.
117 $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
118 $                        region are discarded.
119 $                        Default.
120 $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
121 $                        region are added back in.
122 $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
123 $                        discarded.
124 $                        Not very good.
125 
126    Level: beginner
127 
128 .seealso: PCASMSetType()
129 E*/
130 typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;
131 PETSC_EXTERN const char *const PCASMTypes[];
132 
133 /*E
134     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).
135 
136    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
137    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
138    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
139    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.
140 
141 $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
142 $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
143 $                        from neighboring subdomains.
144 $                        Classical standard additive Schwarz.
145 $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
146 $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
147 $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
148 $                        assumption).
149 $                        Default.
150 $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
151 $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
152 $                        from neighboring subdomains.
153 $
154 $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
155 $                        Not very good.
156 
157    Level: beginner
158 
159 .seealso: PCGASMSetType()
160 E*/
161 typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
162 PETSC_EXTERN const char *const PCGASMTypes[];
163 
164 /*E
165     PCCompositeType - Determines how two or more preconditioner are composed
166 
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
175 
176    Level: beginner
177 
178 .seealso: PCCompositeSetType()
179 E*/
180 typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
181 PETSC_EXTERN const char *const PCCompositeTypes[];
182 
183 /*E
184     PCFieldSplitSchurPreType - Determines how to precondition Schur complement
185 
186     Level: intermediate
187 
188 .seealso: PCFieldSplitSetSchurPre()
189 E*/
190 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 PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];
192 
193 /*E
194     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use
195 
196     Level: intermediate
197 
198 .seealso: PCFieldSplitSetSchurFactType()
199 E*/
200 typedef enum {
201   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
202   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
203   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
204   PC_FIELDSPLIT_SCHUR_FACT_FULL
205 } PCFieldSplitSchurFactType;
206 PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];
207 
208 /*E
209     PCPARMSGlobalType - Determines the global preconditioner method in PARMS
210 
211     Level: intermediate
212 
213 .seealso: PCPARMSSetGlobal()
214 E*/
215 typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
216 PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
217 /*E
218     PCPARMSLocalType - Determines the local preconditioner method in PARMS
219 
220     Level: intermediate
221 
222 .seealso: PCPARMSSetLocal()
223 E*/
224 typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
225 PETSC_EXTERN const char *const PCPARMSLocalTypes[];
226 
227 /*E
228     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method
229 
230     Level: intermediate
231 
232 .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseInterpolation()
233 E*/
234 typedef const char *PCGAMGType;
235 #define PCGAMGAGG         "agg"
236 #define PCGAMGGEO         "geo"
237 #define PCGAMGCLASSICAL   "classical"
238 
239 typedef const char *PCGAMGClassicalType;
240 #define PCGAMGCLASSICALDIRECT   "direct"
241 #define PCGAMGCLASSICALSTANDARD "standard"
242 
243 /*E
244     PCMGType - Determines the type of multigrid method that is run.
245 
246    Level: beginner
247 
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
260 
261 .seealso: PCMGSetType(), PCMGSetCycleType(), PCMGSetCycleTypeOnLevel()
262 
263 E*/
264 typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
265 PETSC_EXTERN const char *const PCMGTypes[];
266 #define PC_MG_CASCADE PC_MG_KASKADE;
267 
268 /*E
269     PCMGCycleType - Use V-cycle or W-cycle
270 
271    Level: beginner
272 
273    Values:
274 +  PC_MG_V_CYCLE
275 -  PC_MG_W_CYCLE
276 
277 .seealso: PCMGSetCycleType()
278 
279 E*/
280 typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
281 PETSC_EXTERN const char *const PCMGCycleTypes[];
282 
283 /*E
284     PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process
285 
286    Level: beginner
287 
288    Values:
289 +  PC_MG_GALERKIN_PMAT - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid
290 .  PC_MG_GALERKIN_MAT -  computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid
291 .  PC_MG_GALERKIN_BOTH - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once
292 -  PC_MG_GALERKIN_NONE - neither operator is computed via the Galerkin process, the user must provide the operator
293 
294    Users should never set PC_MG_GALERKIN_EXTERNAL, it is used by GAMG and ML
295 
296 .seealso: PCMGSetCycleType()
297 
298 E*/
299 typedef enum { PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, PC_MG_GALERKIN_EXTERNAL} PCMGGalerkinType;
300 PETSC_EXTERN const char *const PCMGGalerkinTypes[];
301 
302 /*E
303     PCExoticType - Face based or wirebasket based coarse grid space
304 
305    Level: beginner
306 
307 .seealso: PCExoticSetType(), PCEXOTIC
308 E*/
309 typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
310 PETSC_EXTERN const char *const PCExoticTypes[];
311 PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);
312 
313 /*E
314     PCFailedReason - indicates type of PC failure
315 
316     Level: beginner
317 
318     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
319 E*/
320 typedef enum {PC_NOERROR,PC_FACTOR_STRUCT_ZEROPIVOT,PC_FACTOR_NUMERIC_ZEROPIVOT,PC_FACTOR_OUTMEMORY,PC_FACTOR_OTHER,PC_SUBPC_ERROR} PCFailedReason;
321 PETSC_EXTERN const char *const PCFailedReasons[];
322 #endif
323