Actual source code: petscpctypes.h

petsc-3.13.6 2020-09-29
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: .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 PCSHELL           "shell"
 31: #define PCBJACOBI         "bjacobi"
 32: #define PCMG              "mg"
 33: #define PCEISENSTAT       "eisenstat"
 34: #define PCILU             "ilu"
 35: #define PCICC             "icc"
 36: #define PCASM             "asm"
 37: #define PCGASM            "gasm"
 38: #define PCKSP             "ksp"
 39: #define PCCOMPOSITE       "composite"
 40: #define PCREDUNDANT       "redundant"
 41: #define PCSPAI            "spai"
 42: #define PCNN              "nn"
 43: #define PCCHOLESKY        "cholesky"
 44: #define PCPBJACOBI        "pbjacobi"
 45: #define PCVPBJACOBI       "vpbjacobi"
 46: #define PCMAT             "mat"
 47: #define PCHYPRE           "hypre"
 48: #define PCPARMS           "parms"
 49: #define PCFIELDSPLIT      "fieldsplit"
 50: #define PCTFS             "tfs"
 51: #define PCML              "ml"
 52: #define PCGALERKIN        "galerkin"
 53: #define PCEXOTIC          "exotic"
 54: #define PCCP              "cp"
 55: #define PCBFBT            "bfbt"
 56: #define PCLSC             "lsc"
 57: #define PCPYTHON          "python"
 58: #define PCPFMG            "pfmg"
 59: #define PCSYSPFMG         "syspfmg"
 60: #define PCREDISTRIBUTE    "redistribute"
 61: #define PCSVD             "svd"
 62: #define PCGAMG            "gamg"
 63: #define PCCHOWILUVIENNACL "chowiluviennacl"
 64: #define PCROWSCALINGVIENNACL "rowscalingviennacl"
 65: #define PCSAVIENNACL      "saviennacl"
 66: #define PCBDDC            "bddc"
 67: #define PCKACZMARZ        "kaczmarz"
 68: #define PCTELESCOPE       "telescope"
 69: #define PCPATCH           "patch"
 70: #define PCLMVM            "lmvm"
 71: #define PCHMG             "hmg"
 72: #define PCDEFLATION       "deflation"
 73: #define PCHPDDM           "hpddm"

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

 79:    Level: beginner

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

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

 89:    Level: advanced

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

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

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

105:    Level: intermediate

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

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

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.

126:    Level: beginner

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;

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

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

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

156:    Level: beginner

158: .seealso: PCGASMSetType()
159: E*/
160: typedef enum {PC_GASM_BASIC = 3,PC_GASM_RESTRICT = 1,PC_GASM_INTERPOLATE = 2,PC_GASM_NONE = 0} PCGASMType;

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

165: $  PC_COMPOSITE_ADDITIVE - results from Section 1.5 Writing Application Codes with PETSc of all preconditioners are added together
166: $  PC_COMPOSITE_MULTIPLICATIVE - preconditioners are applied sequentially to the residual freshly
167: $                                computed after the previous preconditioner Section 1.5 Writing Application Codes with PETSc
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,PC_COMPOSITE_GKB} PCCompositeType;

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

183:     Level: intermediate

185: .seealso: PCFieldSplitSetSchurPre()
186: E*/
187: 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: /*E
190:     PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use

192:     Level: intermediate

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

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

206:     Level: intermediate

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

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

215:     Level: intermediate

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

221: /*J
222:     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method

224:     Level: intermediate

226: $   PCGAMGAGG - (the default) smoothed aggregation algorithm, robust, very well tested
227: $   PCGAMGGEO - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not well tested
228: $   PCGAMGCLASSICAL - classical algebraic multigrid preconditioner, incomplete, poorly tested

230: .seealso: PCMG, PCSetType(), PCGAMGSetThreshold(), PCGAMGSetThreshold(), PCGAMGSetReuseInterpolation()
231: J*/
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: #define PC_MG_CASCADE PC_MG_KASKADE;

265: /*E
266:     PCMGCycleType - Use V-cycle or W-cycle

268:    Level: beginner

270:    Values:
271: +  PC_MG_V_CYCLE
272: -  PC_MG_W_CYCLE

274: .seealso: PCMGSetCycleType()

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

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

282:    Level: beginner

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

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

292: .seealso: PCMGSetCycleType()

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

297: /*E
298:     PCExoticType - Face based or wirebasket based coarse grid space

300:    Level: beginner

302: .seealso: PCExoticSetType(), PCEXOTIC
303: E*/
304: typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;

306: /*E
307:    PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains.

309:    Level: intermediate

311:    Values:
312: +  PC_BDDC_INTERFACE_EXT_DIRICHLET - solves Dirichlet interior problem; this is the standard BDDC algorithm
313: -  PC_BDDC_INTERFACE_EXT_LUMP - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP"

315: E*/
316: typedef enum {
317:   PC_BDDC_INTERFACE_EXT_DIRICHLET,
318:   PC_BDDC_INTERFACE_EXT_LUMP
319: } PCBDDCInterfaceExtType;

321: /*E
322:     PCPatchConstructType - The algorithm used to construct patches for the preconditioner

324:    Level: beginner

326: .seealso: PCPatchSetConstructType(), PCEXOTIC
327: E*/
328: typedef enum {PC_PATCH_STAR, PC_PATCH_VANKA, PC_PATCH_PARDECOMP, PC_PATCH_USER, PC_PATCH_PYTHON} PCPatchConstructType;

330: /*E
331:     PCDeflationSpaceType - Type of deflation

333:     Values:
334: +   PC_DEFLATION_SPACE_HAAR        - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off
335: .   PC_DEFLATION_SPACE_DB2         - MATCOMPOSITE of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet)
336: .   PC_DEFLATION_SPACE_DB4         - same as above, but with db4 (4 coefficient Daubechies)
337: .   PC_DEFLATION_SPACE_DB8         - same as above, but with db8 (8 coefficient Daubechies)
338: .   PC_DEFLATION_SPACE_DB16        - same as above, but with db16 (16 coefficient Daubechies)
339: .   PC_DEFLATION_SPACE_BIORTH22    - same as above, but with biorthogonal 2.2 (6 coefficients)
340: .   PC_DEFLATION_SPACE_MEYER       - same as above, but with Meyer/FIR (62 coefficients)
341: .   PC_DEFLATION_SPACE_AGGREGATION - aggregates local indices (given by operator matix distribution) into a subdomain
342: -   PC_DEFLATION_SPACE_USER        - indicates space set by user

344:     Notes:
345:       Wavelet-based space (except Haar) can be used in multilevel deflation.

347:     Level: intermediate

349: .seealso: PCDeflationSetSpaceToCompute(), PCDEFLATION
350: E*/
351: typedef enum {
352:   PC_DEFLATION_SPACE_HAAR,
353:   PC_DEFLATION_SPACE_DB2,
354:   PC_DEFLATION_SPACE_DB4,
355:   PC_DEFLATION_SPACE_DB8,
356:   PC_DEFLATION_SPACE_DB16,
357:   PC_DEFLATION_SPACE_BIORTH22,
358:   PC_DEFLATION_SPACE_MEYER,
359:   PC_DEFLATION_SPACE_AGGREGATION,
360:   PC_DEFLATION_SPACE_USER
361: } PCDeflationSpaceType;

363: /*E
364:     PCHPDDMCoarseCorrectionType - Type of coarse correction used by PCHPDDM

366:     Level: intermediate

368:     Values:
369: +   PC_HPDDM_COARSE_CORRECTION_DEFLATED (default) - eq. (1) in PCHPDDMShellApply()
370: .   PC_HPDDM_COARSE_CORRECTION_ADDITIVE - eq. (2)
371: -   PC_HPDDM_COARSE_CORRECTION_BALANCED - eq. (3)

373: .seealso: PCHPDDM, PCSetType(), PCHPDDMShellApply()
374: E*/
375: typedef enum { PC_HPDDM_COARSE_CORRECTION_DEFLATED, PC_HPDDM_COARSE_CORRECTION_ADDITIVE, PC_HPDDM_COARSE_CORRECTION_BALANCED } PCHPDDMCoarseCorrectionType;

377: /*E
378:     PCFailedReason - indicates type of PC failure

380:     Level: beginner

382:     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
383: E*/
384: typedef enum {PC_NOERROR,PC_FACTOR_STRUCT_ZEROPIVOT,PC_FACTOR_NUMERIC_ZEROPIVOT,PC_FACTOR_OUTMEMORY,PC_FACTOR_OTHER,PC_SUBPC_ERROR} PCFailedReason;

386: /*E
387:     PCGAMGLayoutType - Layout for reduced grids

389:     Level: intermediate

391: .seealso: PCGAMGSetCoarseGridLayoutType()
392:     Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h
393: E*/
394: typedef enum {PCGAMG_LAYOUT_COMPACT,PCGAMG_LAYOUT_SPREAD} PCGAMGLayoutType;

396: #endif