Actual source code: petscpctypes.h

petsc-3.8.4 2018-03-24
Report Typos and Errors
  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 PCSACUSP          "sacusp"        /* these four run on NVIDIA GPUs using CUSP */
 66: #define PCSACUSPPOLY      "sacusppoly"
 67: #define PCBICGSTABCUSP    "bicgstabcusp"
 68: #define PCAINVCUSP        "ainvcusp"
 69: #define PCCHOWILUVIENNACL "chowiluviennacl"
 70: #define PCROWSCALINGVIENNACL "rowscalingviennacl"
 71: #define PCSAVIENNACL      "saviennacl"
 72: #define PCBDDC            "bddc"
 73: #define PCKACZMARZ        "kaczmarz"
 74: #define PCTELESCOPE       "telescope"

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

 80:    Level: beginner

 82: .seealso:
 83: E*/
 84: typedef enum { PC_SIDE_DEFAULT=-1,PC_LEFT,PC_RIGHT,PC_SYMMETRIC} PCSide;
 85: #define PC_SIDE_MAX (PC_SYMMETRIC + 1)
 86: PETSC_EXTERN const char *const *const PCSides;

 88: /*E
 89:     PCRichardsonConvergedReason - reason a PCApplyRichardson method terminates

 91:    Level: advanced

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

 95: .seealso: PCApplyRichardson()
 96: E*/
 97: typedef enum {
 98:               PCRICHARDSON_CONVERGED_RTOL               =  2,
 99:               PCRICHARDSON_CONVERGED_ATOL               =  3,
100:               PCRICHARDSON_CONVERGED_ITS                =  4,
101:               PCRICHARDSON_DIVERGED_DTOL                = -4} PCRichardsonConvergedReason;

103: /*E
104:     PCJacobiType - What elements are used to form the Jacobi preconditioner

106:    Level: intermediate

108: .seealso:
109: E*/
110: typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType;
111: PETSC_EXTERN const char *const PCJacobiTypes[];

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;
133: PETSC_EXTERN const char *const PCASMTypes[];

135: /*E
136:     PCGASMType - Type of generalized additive Schwarz method to use (differs from ASM in allowing multiple processors per subdomain).

138:    Each subdomain has nested inner and outer parts.  The inner subdomains are assumed to form a non-overlapping covering of the computational
139:    domain, while the outer subdomains contain the inner subdomains and overlap with each other.  This preconditioner will compute
140:    a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in
141:    (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed.

143: $  PC_GASM_BASIC       - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied
144: $                        over the outer subdomains.  As a result, points in the overlap will receive the sum of the corrections
145: $                        from neighboring subdomains.
146: $                        Classical standard additive Schwarz.
147: $  PC_GASM_RESTRICT    - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only
148: $                        (i.e., zeroed out over the overlap portion of the outer subdomain before being applied).  As a result,
149: $                        each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering
150: $                        assumption).
151: $                        Default.
152: $  PC_GASM_INTERPOLATE - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is
153: $                        applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections
154: $                        from neighboring subdomains.
155: $
156: $  PC_GASM_NONE        - Residuals and corrections are zeroed out outside the local subdomains.
157: $                        Not very good.

159:    Level: beginner

161: .seealso: PCGASMSetType()
162: E*/
163: typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;
164: PETSC_EXTERN const char *const PCGASMTypes[];

166: /*E
167:     PCCompositeType - Determines how two or more preconditioner are composed

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

178:    Level: beginner

180: .seealso: PCCompositeSetType()
181: E*/
182: typedef enum {PC_COMPOSITE_ADDITIVE,PC_COMPOSITE_MULTIPLICATIVE,PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE,PC_COMPOSITE_SPECIAL,PC_COMPOSITE_SCHUR} PCCompositeType;
183: PETSC_EXTERN const char *const PCCompositeTypes[];

185: /*E
186:     PCFieldSplitSchurPreType - Determines how to precondition Schur complement

188:     Level: intermediate

190: .seealso: PCFieldSplitSetSchurPre()
191: E*/
192: 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;
193: PETSC_EXTERN const char *const PCFieldSplitSchurPreTypes[];

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

198:     Level: intermediate

200: .seealso: PCFieldSplitSetSchurFactType()
201: E*/
202: typedef enum {
203:   PC_FIELDSPLIT_SCHUR_FACT_DIAG,
204:   PC_FIELDSPLIT_SCHUR_FACT_LOWER,
205:   PC_FIELDSPLIT_SCHUR_FACT_UPPER,
206:   PC_FIELDSPLIT_SCHUR_FACT_FULL
207: } PCFieldSplitSchurFactType;
208: PETSC_EXTERN const char *const PCFieldSplitSchurFactTypes[];

210: /*E
211:     PCPARMSGlobalType - Determines the global preconditioner method in PARMS

213:     Level: intermediate

215: .seealso: PCPARMSSetGlobal()
216: E*/
217: typedef enum {PC_PARMS_GLOBAL_RAS,PC_PARMS_GLOBAL_SCHUR,PC_PARMS_GLOBAL_BJ} PCPARMSGlobalType;
218: PETSC_EXTERN const char *const PCPARMSGlobalTypes[];
219: /*E
220:     PCPARMSLocalType - Determines the local preconditioner method in PARMS

222:     Level: intermediate

224: .seealso: PCPARMSSetLocal()
225: E*/
226: typedef enum {PC_PARMS_LOCAL_ILU0,PC_PARMS_LOCAL_ILUK,PC_PARMS_LOCAL_ILUT,PC_PARMS_LOCAL_ARMS} PCPARMSLocalType;
227: PETSC_EXTERN const char *const PCPARMSLocalTypes[];

229: /*E
230:     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method

232:     Level: intermediate

234: .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseInterpolation()
235: E*/
236: typedef const char *PCGAMGType;
237: #define PCGAMGAGG         "agg"
238: #define PCGAMGGEO         "geo"
239: #define PCGAMGCLASSICAL   "classical"

241: typedef const char *PCGAMGClassicalType;
242: #define PCGAMGCLASSICALDIRECT   "direct"
243: #define PCGAMGCLASSICALSTANDARD "standard"

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

248:    Level: beginner

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

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

265: E*/
266: typedef enum { PC_MG_MULTIPLICATIVE,PC_MG_ADDITIVE,PC_MG_FULL,PC_MG_KASKADE } PCMGType;
267: PETSC_EXTERN const char *const PCMGTypes[];
268: #define PC_MG_CASCADE PC_MG_KASKADE;

270: /*E
271:     PCMGCycleType - Use V-cycle or W-cycle

273:    Level: beginner

275:    Values:
276: +  PC_MG_V_CYCLE
277: -  PC_MG_W_CYCLE

279: .seealso: PCMGSetCycleType()

281: E*/
282: typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
283: PETSC_EXTERN const char *const PCMGCycleTypes[];

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

288:    Level: beginner

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

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

298: .seealso: PCMGSetCycleType()

300: E*/
301: typedef enum { PC_MG_GALERKIN_BOTH,PC_MG_GALERKIN_PMAT,PC_MG_GALERKIN_MAT, PC_MG_GALERKIN_NONE, PC_MG_GALERKIN_EXTERNAL} PCMGGalerkinType;
302: PETSC_EXTERN const char *const PCMGGalerkinTypes[];

304: /*E
305:     PCExoticType - Face based or wirebasket based coarse grid space

307:    Level: beginner

309: .seealso: PCExoticSetType(), PCEXOTIC
310: E*/
311: typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
312: PETSC_EXTERN const char *const PCExoticTypes[];
313: PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);

315: /*E
316:     PCFailedReason - indicates type of PC failure

318:     Level: beginner

320:     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
321: E*/
322: typedef enum {PC_NOERROR,PC_FACTOR_STRUCT_ZEROPIVOT,PC_FACTOR_NUMERIC_ZEROPIVOT,PC_FACTOR_OUTMEMORY,PC_FACTOR_OTHER,PC_SUBPC_ERROR} PCFailedReason;
323: PETSC_EXTERN const char *const PCFailedReasons[];
324: #endif