Actual source code: petscpctypes.h

petsc-3.11.4 2019-09-28
Report Typos and Errors
  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:   Concepts: preconditioners

 11: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types)
 12: S*/
 13: typedef struct _p_PC* PC;

 15: /*J
 16:     PCType - String with the name of a PETSc preconditioner method.

 18:    Level: beginner

 20:    Notes:
 21:     Click on the links above to see details on a particular solver

 23:           PCRegister() is used to register preconditioners that are then accessible via PCSetType()

 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 PCSHELL           "shell"
 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 PCCOMPOSITE       "composite"
 42: #define PCREDUNDANT       "redundant"
 43: #define PCSPAI            "spai"
 44: #define PCNN              "nn"
 45: #define PCCHOLESKY        "cholesky"
 46: #define PCPBJACOBI        "pbjacobi"
 47: #define PCVPBJACOBI       "vpbjacobi"
 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"
 71: #define PCPATCH           "patch"
 72: #define PCLMVM            "lmvm"

 74: /*E
 75:     PCSide - If the preconditioner is to be applied to the left, right
 76:      or symmetrically around the operator.

 78:    Level: beginner

 80: .seealso:
 81: E*/
 82: typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
 83: #define PC_SIDE_MAX (PC_SYMMETRIC + 1)

 85: /*E
 86:     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates

 88:    Level: advanced

 90:    Notes:
 91:     this must match petsc/finclude/petscpc.h and the KSPConvergedReason values in petscksp.h

 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;

101: /*E
102:     PCJacobiType - What elements are used to form the Jacobi preconditioner

104:    Level: intermediate

106: .seealso:
107: E*/
108: typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType;

110: /*E
111:     PCASMType - Type of additive Schwarz method to use

113: $  PC_ASM_BASIC        - Symmetric version where residuals from the ghost points are used
114: $                        and computed values in ghost regions are added together.
115: $                        Classical standard additive Schwarz.
116: $  PC_ASM_RESTRICT     - Residuals from ghost points are used but computed values in ghost
117: $                        region are discarded.
118: $                        Default.
119: $  PC_ASM_INTERPOLATE  - Residuals from ghost points are not used, computed values in ghost
120: $                        region are added back in.
121: $  PC_ASM_NONE         - Residuals from ghost points are not used, computed ghost values are
122: $                        discarded.
123: $                        Not very good.

125:    Level: beginner

127: .seealso: PCASMSetType()
128: E*/
129: typedef enum {PC_ASM_BASIC = 3,PC_ASM_RESTRICT = 1,PC_ASM_INTERPOLATE = 2,PC_ASM_NONE = 0} PCASMType;

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;

161: /*E
162:     PCCompositeType - Determines how two or more preconditioner are composed

164: $  PC_COMPOSITE_ADDITIVE - results from application of all preconditioners are added together
165: $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
166: $                                computed after the previous preconditioner application
167: $  PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
168: $                                computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners)
169: $  PC_COMPOSITE_SPECIAL - This is very special for a matrix of the form alpha I + R + S
170: $                         where first preconditioner is built from alpha I + S and second from
171: $                         alpha I + R

173:    Level: beginner

175: .seealso: PCCompositeSetType()
176: E*/
177: typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR,PC_COMPOSITE_GKB} PCCompositeType;

179: /*E
180:     PCFieldSplitSchurPreType - Determines how to precondition Schur complement

182:     Level: intermediate

184: .seealso: PCFieldSplitSetSchurPre()
185: E*/
186: 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;

188: /*E
189:     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use

191:     Level: intermediate

193: .seealso: PCFieldSplitSetSchurFactType()
194: E*/
195: typedef enum {
196:   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
197:   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
198:   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
199:   PC_FIELDSPLIT_SCHUR_FACT_FULL
200: } PCFieldSplitSchurFactType;

202: /*E
203:     PCPARMSGlobalType - Determines the global preconditioner method in PARMS

205:     Level: intermediate

207: .seealso: PCPARMSSetGlobal()
208: E*/
209: typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;

211: /*E
212:     PCPARMSLocalType - Determines the local preconditioner method in PARMS

214:     Level: intermediate

216: .seealso: PCPARMSSetLocal()
217: E*/
218: typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;

220: /*E
221:     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method

223:     Level: intermediate

225: .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseInterpolation()
226: E*/
227: typedef const char *PCGAMGType;
228: #define PCGAMGAGG         "agg"
229: #define PCGAMGGEO         "geo"
230: #define PCGAMGCLASSICAL   "classical"

232: typedef const char *PCGAMGClassicalType;
233: #define PCGAMGCLASSICALDIRECT   "direct"
234: #define PCGAMGCLASSICALSTANDARD "standard"

236: /*E
237:     PCMGType - Determines the type of multigrid method that is run.

239:    Level: beginner

241:    Values:
242: +  PC_MG_MULTIPLICATIVE (default) - traditional V or W cycle as determined by PCMGSetCycleType()
243: .  PC_MG_ADDITIVE - the additive multigrid preconditioner where all levels are
244:                 smoothed before updating the residual. This only uses the
245:                 down smoother, in the preconditioner the upper smoother is ignored
246: .  PC_MG_FULL - same as multiplicative except one also performs grid sequencing,
247:             that is starts on the coarsest grid, performs a cycle, interpolates
248:             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
249:             algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm
250:             calls the V-cycle only on the coarser level and has a post-smoother instead.
251: -  PC_MG_KASKADE - like full multigrid except one never goes back to a coarser level
252:                from a finer

254: .seealso: PCMGSetType(), PCMGSetCycleType(), PCMGSetCycleTypeOnLevel()

256: E*/
257: typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
258: #define PC_MG_CASCADE PC_MG_KASKADE;

260: /*E
261:     PCMGCycleType - Use V-cycle or W-cycle

263:    Level: beginner

265:    Values:
266: +  PC_MG_V_CYCLE
267: -  PC_MG_W_CYCLE

269: .seealso: PCMGSetCycleType()

271: E*/
272: typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;

274: /*E
275:     PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process

277:    Level: beginner

279:    Values:
280: +  PC_MG_GALERKIN_PMAT - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid
281: .  PC_MG_GALERKIN_MAT -  computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid
282: .  PC_MG_GALERKIN_BOTH - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once
283: -  PC_MG_GALERKIN_NONE - neither operator is computed via the Galerkin process, the user must provide the operator

285:    Users should never set PC_MG_GALERKIN_EXTERNAL, it is used by GAMG and ML

287: .seealso: PCMGSetCycleType()

289: E*/
290: typedef enum { PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, PC_MG_GALERKIN_EXTERNAL} PCMGGalerkinType;

292: /*E
293:     PCExoticType - Face based or wirebasket based coarse grid space

295:    Level: beginner

297: .seealso: PCExoticSetType(), PCEXOTIC
298: E*/
299: typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;

301: /*E
302:     PCPatchConstructType - The algorithm used to construct patches for the preconditioner

304:    Level: beginner

306: .seealso: PCPatchSetConstructType(), PCEXOTIC
307: E*/
308: typedef enum {PC_PATCH_STAR, PC_PATCH_VANKA, PC_PATCH_PARDECOMP, PC_PATCH_USER, PC_PATCH_PYTHON} PCPatchConstructType;

310: /*E
311:     PCFailedReason - indicates type of PC failure

313:     Level: beginner

315:     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
316: E*/
317: typedef enum {PC_NOERROR,PC_FACTOR_STRUCT_ZEROPIVOT,PC_FACTOR_NUMERIC_ZEROPIVOT,PC_FACTOR_OUTMEMORY,PC_FACTOR_OTHER,PC_SUBPC_ERROR} PCFailedReason;
318: #endif