Actual source code: petscpctypes.h

petsc-3.7.5 2017-01-01
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 below 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 PCBDDC            "bddc"
 70: #define PCKACZMARZ        "kaczmarz"
 71: #define PCTELESCOPE       "telescope"

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

 77:    Level: beginner

 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;

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

 88:    Level: advanced

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

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

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

103:    Level: intermediate

105: .seealso:
106: E*/
107: typedef enum { PC_JACOBI_DIAGONAL,PC_JACOBI_ROWMAX,PC_JACOBI_ROWSUM} PCJacobiType;
108: PETSC_EXTERN const char *const PCJacobiTypes[];

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

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

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

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

175:    Level: beginner

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

182: /*E
183:     PCFieldSplitSchurPreType - Determines how to precondition Schur complement

185:     Level: intermediate

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

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

195:     Level: intermediate

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

207: /*E
208:     PCPARMSGlobalType - Determines the global preconditioner method in PARMS

210:     Level: intermediate

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

219:     Level: intermediate

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

226: /*E
227:     PCGAMGType - type of generalized algebraic multigrid (PCGAMG) method

229:     Level: intermediate

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

238: typedef const char *PCGAMGClassicalType;
239: #define PCGAMGCLASSICALDIRECT   "direct"
240: #define PCGAMGCLASSICALSTANDARD "standard"

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

245:    Level: beginner

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

260: .seealso: PCMGSetType()

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

267: /*E
268:     PCMGCycleType - Use V-cycle or W-cycle

270:    Level: beginner

272:    Values:
273: +  PC_MG_V_CYCLE
274: -  PC_MG_W_CYCLE

276: .seealso: PCMGSetCycleType()

278: E*/
279: typedef enum { PC_MG_CYCLE_V = 1,PC_MG_CYCLE_W = 2 } PCMGCycleType;
280: PETSC_EXTERN const char *const PCMGCycleTypes[];

282: /*E
283:     PCExoticType - Face based or wirebasket based coarse grid space

285:    Level: beginner

287: .seealso: PCExoticSetType(), PCEXOTIC
288: E*/
289: typedef enum { PC_EXOTIC_FACE,PC_EXOTIC_WIREBASKET } PCExoticType;
290: PETSC_EXTERN const char *const PCExoticTypes[];
291: PETSC_EXTERN PetscErrorCode PCExoticSetType(PC,PCExoticType);

293: /*E
294:     PCFailedReason - indicates type of PC failure

296:     Level: beginner

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