1 #if !defined(PETSCPCTYPES_H) 2 #define PETSCPCTYPES_H 3 4 /* SUBMANSEC = PC */ 5 6 /*S 7 PC - Abstract PETSc object that manages all preconditioners including direct solvers such as `PCLU` 8 9 Level: beginner 10 11 .seealso: `PCCreate()`, `PCSetType()`, `PCType` 12 S*/ 13 typedef struct _p_PC *PC; 14 15 /*J 16 PCType - String with the name of a PETSc preconditioner method. 17 18 Level: beginner 19 20 Note: 21 `PCRegister()` is used to register preconditioners that are then accessible via `PCSetType()` 22 23 .seealso: `PCSetType()`, `PC`, `PCCreate()`, `PCRegister()`, `PCSetFromOptions()`, `PCLU`, `PCJACOBI`, `PCBJACOBI` 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 PCQR "qr" 31 #define PCSHELL "shell" 32 #define PCAMGX "amgx" 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 PCBJKOKKOS "bjkokkos" 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 PCVPBJACOBI "vpbjacobi" 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 PCSMG "smg" 63 #define PCSYSPFMG "syspfmg" 64 #define PCREDISTRIBUTE "redistribute" 65 #define PCSVD "svd" 66 #define PCGAMG "gamg" 67 #define PCCHOWILUVIENNACL "chowiluviennacl" 68 #define PCROWSCALINGVIENNACL "rowscalingviennacl" 69 #define PCSAVIENNACL "saviennacl" 70 #define PCBDDC "bddc" 71 #define PCKACZMARZ "kaczmarz" 72 #define PCTELESCOPE "telescope" 73 #define PCPATCH "patch" 74 #define PCLMVM "lmvm" 75 #define PCHMG "hmg" 76 #define PCDEFLATION "deflation" 77 #define PCHPDDM "hpddm" 78 #define PCH2OPUS "h2opus" 79 #define PCMPI "mpi" 80 81 /*E 82 PCSide - If the preconditioner is to be applied to the left, right 83 or symmetrically around the operator. 84 85 Level: beginner 86 87 .seealso: 88 E*/ 89 typedef enum { 90 PC_SIDE_DEFAULT = -1, 91 PC_LEFT, 92 PC_RIGHT, 93 PC_SYMMETRIC 94 } PCSide; 95 #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 96 97 /*E 98 PCRichardsonConvergedReason - reason a `PCApplyRichardson() method terminated 99 100 Level: advanced 101 102 Developer Note: 103 this must match petsc/finclude/petscpc.h and the `KSPConvergedReason` values in petscksp.h 104 105 .seealso: `PCApplyRichardson()` 106 E*/ 107 typedef enum { 108 PCRICHARDSON_CONVERGED_RTOL = 2, 109 PCRICHARDSON_CONVERGED_ATOL = 3, 110 PCRICHARDSON_CONVERGED_ITS = 4, 111 PCRICHARDSON_DIVERGED_DTOL = -4 112 } PCRichardsonConvergedReason; 113 114 /*E 115 PCJacobiType - What elements are used to form the Jacobi preconditioner 116 117 Level: intermediate 118 119 .seealso: 120 E*/ 121 typedef enum { 122 PC_JACOBI_DIAGONAL, 123 PC_JACOBI_ROWMAX, 124 PC_JACOBI_ROWSUM 125 } PCJacobiType; 126 127 /*E 128 PCASMType - Type of additive Schwarz method to use 129 130 $ `PC_ASM_BASIC` - Symmetric version where residuals from the ghost points are used 131 $ and computed values in ghost regions are added together. 132 $ Classical standard additive Schwarz. 133 $ `PC_ASM_RESTRICT` - Residuals from ghost points are used but computed values in ghost 134 $ region are discarded. 135 $ Default. 136 $ `PC_ASM_INTERPOLATE` - Residuals from ghost points are not used, computed values in ghost 137 $ region are added back in. 138 $ `PC_ASM_NONE` - Residuals from ghost points are not used, computed ghost values are 139 $ discarded. 140 $ Not very good. 141 142 Level: beginner 143 144 .seealso: `PCASMSetType()` 145 E*/ 146 typedef enum { 147 PC_ASM_BASIC = 3, 148 PC_ASM_RESTRICT = 1, 149 PC_ASM_INTERPOLATE = 2, 150 PC_ASM_NONE = 0 151 } PCASMType; 152 153 /*E 154 PCGASMType - Type of generalized additive Schwarz method to use (differs from `PCASM` in allowing multiple processors per subdomain). 155 156 Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 157 domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 158 a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 159 (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 160 161 $ `PC_GASM_BASIC` - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 162 $ over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 163 $ from neighboring subdomains. 164 $ Classical standard additive Schwarz. 165 $ `PC_GASM_RESTRICT` - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 166 $ (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 167 $ each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 168 $ assumption). 169 $ Default. 170 $ `PC_GASM_INTERPOLATE` - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 171 $ applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 172 $ from neighboring subdomains. 173 $ 174 $ `PC_GASM_NONE` - Residuals and corrections are zeroed out outside the local subdomains. 175 $ Not very good. 176 177 Level: beginner 178 179 .seealso: `PCGASMSetType()` 180 E*/ 181 typedef enum { 182 PC_GASM_BASIC = 3, 183 PC_GASM_RESTRICT = 1, 184 PC_GASM_INTERPOLATE = 2, 185 PC_GASM_NONE = 0 186 } PCGASMType; 187 188 /*E 189 PCCompositeType - Determines how two or more preconditioner are composed 190 191 $ `PC_COMPOSITE_ADDITIVE` - results from application of all preconditioners are added together 192 $ `PC_COMPOSITE_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly 193 $ computed after the previous preconditioner application 194 $ `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly 195 $ computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners) 196 $ `PC_COMPOSITE_SPECIAL` - This is very special for a matrix of the form alpha I + R + S 197 $ where first preconditioner is built from alpha I + S and second from 198 $ alpha I + R 199 200 Level: beginner 201 202 .seealso: `PCCompositeSetType()` 203 E*/ 204 typedef enum { 205 PC_COMPOSITE_ADDITIVE, 206 PC_COMPOSITE_MULTIPLICATIVE, 207 PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, 208 PC_COMPOSITE_SPECIAL, 209 PC_COMPOSITE_SCHUR, 210 PC_COMPOSITE_GKB 211 } PCCompositeType; 212 213 /*E 214 PCFieldSplitSchurPreType - Determines how to precondition a Schur complement 215 216 Level: intermediate 217 218 .seealso: `PCFieldSplitSetSchurPre()` 219 E*/ 220 typedef enum { 221 PC_FIELDSPLIT_SCHUR_PRE_SELF, 222 PC_FIELDSPLIT_SCHUR_PRE_SELFP, 223 PC_FIELDSPLIT_SCHUR_PRE_A11, 224 PC_FIELDSPLIT_SCHUR_PRE_USER, 225 PC_FIELDSPLIT_SCHUR_PRE_FULL 226 } PCFieldSplitSchurPreType; 227 228 /*E 229 PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 230 231 Level: intermediate 232 233 .seealso: `PCFieldSplitSetSchurFactType()` 234 E*/ 235 typedef enum { 236 PC_FIELDSPLIT_SCHUR_FACT_DIAG, 237 PC_FIELDSPLIT_SCHUR_FACT_LOWER, 238 PC_FIELDSPLIT_SCHUR_FACT_UPPER, 239 PC_FIELDSPLIT_SCHUR_FACT_FULL 240 } PCFieldSplitSchurFactType; 241 242 /*E 243 PCPARMSGlobalType - Determines the global preconditioner method in `PCPARMS` 244 245 Level: intermediate 246 247 .seealso: `PCPARMSSetGlobal()` 248 E*/ 249 typedef enum { 250 PC_PARMS_GLOBAL_RAS, 251 PC_PARMS_GLOBAL_SCHUR, 252 PC_PARMS_GLOBAL_BJ 253 } PCPARMSGlobalType; 254 255 /*E 256 PCPARMSLocalType - Determines the local preconditioner method in `PCPARMS` 257 258 Level: intermediate 259 260 .seealso: `PCPARMSSetLocal()` 261 E*/ 262 typedef enum { 263 PC_PARMS_LOCAL_ILU0, 264 PC_PARMS_LOCAL_ILUK, 265 PC_PARMS_LOCAL_ILUT, 266 PC_PARMS_LOCAL_ARMS 267 } PCPARMSLocalType; 268 269 /*J 270 PCGAMGType - type of generalized algebraic multigrid `PCGAMG` method 271 272 Level: intermediate 273 274 $ `PCGAMGAGG` - (the default) smoothed aggregation algorithm, robust, very well tested 275 $ `PCGAMGGEO` - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not well tested 276 $ `PCGAMGCLASSICAL` - classical algebraic multigrid preconditioner, incomplete, poorly tested 277 278 .seealso: `PCMG`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()` 279 J*/ 280 typedef const char *PCGAMGType; 281 #define PCGAMGAGG "agg" 282 #define PCGAMGGEO "geo" 283 #define PCGAMGCLASSICAL "classical" 284 285 typedef const char *PCGAMGClassicalType; 286 #define PCGAMGCLASSICALDIRECT "direct" 287 #define PCGAMGCLASSICALSTANDARD "standard" 288 289 /*E 290 PCMGType - Determines the type of multigrid method that is run. 291 292 Level: beginner 293 294 Values: 295 + `PC_MG_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `PCMGSetCycleType()` 296 . `PC_MG_ADDITIVE` - the additive multigrid preconditioner where all levels are 297 smoothed before updating the residual. This only uses the 298 down smoother, in the preconditioner the upper smoother is ignored 299 . `PC_MG_FULL` - same as multiplicative except one also performs grid sequencing, 300 that is starts on the coarsest grid, performs a cycle, interpolates 301 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 302 algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 303 calls the V-cycle only on the coarser level and has a post-smoother instead. 304 - `PC_MG_KASKADE` - like full multigrid except one never goes back to a coarser level 305 from a finer 306 307 .seealso: `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()` 308 309 E*/ 310 typedef enum { 311 PC_MG_MULTIPLICATIVE, 312 PC_MG_ADDITIVE, 313 PC_MG_FULL, 314 PC_MG_KASKADE 315 } PCMGType; 316 #define PC_MG_CASCADE PC_MG_KASKADE; 317 318 /*E 319 PCMGCycleType - Use V-cycle or W-cycle 320 321 Level: beginner 322 323 Values: 324 + `PC_MG_V_CYCLE` - use the v cycle 325 - `PC_MG_W_CYCLE` - use the w cycle 326 327 .seealso: `PCMGSetCycleType()` 328 329 E*/ 330 typedef enum { 331 PC_MG_CYCLE_V = 1, 332 PC_MG_CYCLE_W = 2 333 } PCMGCycleType; 334 335 /*E 336 PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process 337 338 Level: beginner 339 340 Values: 341 + `PC_MG_GALERKIN_PMAT` - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid 342 . `PC_MG_GALERKIN_MAT` - computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid 343 . `PC_MG_GALERKIN_BOTH` - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once 344 - `PC_MG_GALERKIN_NONE` - neither operator is computed via the Galerkin process, the user must provide the operator 345 346 Users should never set `PC_MG_GALERKIN_EXTERNAL`, it is used by `PCHYPRE` and `PCML` 347 348 .seealso: `PCMGSetCycleType()` 349 350 E*/ 351 typedef enum { 352 PC_MG_GALERKIN_BOTH, 353 PC_MG_GALERKIN_PMAT, 354 PC_MG_GALERKIN_MAT, 355 PC_MG_GALERKIN_NONE, 356 PC_MG_GALERKIN_EXTERNAL 357 } PCMGGalerkinType; 358 359 /*E 360 PCExoticType - Face based or wirebasket based coarse grid space 361 362 Level: beginner 363 364 .seealso: `PCExoticSetType()`, `PCEXOTIC` 365 E*/ 366 typedef enum { 367 PC_EXOTIC_FACE, 368 PC_EXOTIC_WIREBASKET 369 } PCExoticType; 370 371 /*E 372 PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains. 373 374 Level: intermediate 375 376 Values: 377 + `PC_BDDC_INTERFACE_EXT_DIRICHLET` - solves Dirichlet interior problem; this is the standard BDDC algorithm 378 - `PC_BDDC_INTERFACE_EXT_LUMP` - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP" 379 380 E*/ 381 typedef enum { 382 PC_BDDC_INTERFACE_EXT_DIRICHLET, 383 PC_BDDC_INTERFACE_EXT_LUMP 384 } PCBDDCInterfaceExtType; 385 386 /*E 387 PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation 388 389 Level: beginner 390 391 .seealso: `PCMGSetAdaptCoarseSpaceType()`, `PCMG` 392 E*/ 393 typedef enum { 394 PCMG_ADAPT_NONE, 395 PCMG_ADAPT_POLYNOMIAL, 396 PCMG_ADAPT_HARMONIC, 397 PCMG_ADAPT_EIGENVECTOR, 398 PCMG_ADAPT_GENERALIZED_EIGENVECTOR, 399 PCMG_ADAPT_GDSW 400 } PCMGCoarseSpaceType; 401 402 /*E 403 PCPatchConstructType - The algorithm used to construct patches for the preconditioner 404 405 Level: beginner 406 407 .seealso: `PCPatchSetConstructType()`, `PCEXOTIC` 408 E*/ 409 typedef enum { 410 PC_PATCH_STAR, 411 PC_PATCH_VANKA, 412 PC_PATCH_PARDECOMP, 413 PC_PATCH_USER, 414 PC_PATCH_PYTHON 415 } PCPatchConstructType; 416 417 /*E 418 PCDeflationSpaceType - Type of deflation 419 420 Values: 421 + `PC_DEFLATION_SPACE_HAAR` - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off 422 . `PC_DEFLATION_SPACE_DB2` - MATCOMPOSITE of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet) 423 . `PC_DEFLATION_SPACE_DB4` - same as above, but with db4 (4 coefficient Daubechies) 424 . `PC_DEFLATION_SPACE_DB8` - same as above, but with db8 (8 coefficient Daubechies) 425 . `PC_DEFLATION_SPACE_DB16` - same as above, but with db16 (16 coefficient Daubechies) 426 . `PC_DEFLATION_SPACE_BIORTH22` - same as above, but with biorthogonal 2.2 (6 coefficients) 427 . `PC_DEFLATION_SPACE_MEYER` - same as above, but with Meyer/FIR (62 coefficients) 428 . `PC_DEFLATION_SPACE_AGGREGATION` - aggregates local indices (given by operator matix distribution) into a subdomain 429 - `PC_DEFLATION_SPACE_USER` - indicates space set by user 430 431 Notes: 432 Wavelet-based space (except Haar) can be used in multilevel deflation. 433 434 Level: intermediate 435 436 .seealso: `PCDeflationSetSpaceToCompute()`, `PCDEFLATION` 437 E*/ 438 typedef enum { 439 PC_DEFLATION_SPACE_HAAR, 440 PC_DEFLATION_SPACE_DB2, 441 PC_DEFLATION_SPACE_DB4, 442 PC_DEFLATION_SPACE_DB8, 443 PC_DEFLATION_SPACE_DB16, 444 PC_DEFLATION_SPACE_BIORTH22, 445 PC_DEFLATION_SPACE_MEYER, 446 PC_DEFLATION_SPACE_AGGREGATION, 447 PC_DEFLATION_SPACE_USER 448 } PCDeflationSpaceType; 449 450 /*E 451 PCHPDDMCoarseCorrectionType - Type of coarse correction used by `PCHPDDM` 452 453 Level: intermediate 454 455 Values: 456 + `PC_HPDDM_COARSE_CORRECTION_DEFLATED` (default) - eq. (1) in PCHPDDMShellApply() 457 . `PC_HPDDM_COARSE_CORRECTION_ADDITIVE` - eq. (2) 458 - `PC_HPDDM_COARSE_CORRECTION_BALANCED` - eq. (3) 459 460 .seealso: `PCHPDDM`, `PCSetType()`, `PCHPDDMShellApply()` 461 E*/ 462 typedef enum { 463 PC_HPDDM_COARSE_CORRECTION_DEFLATED, 464 PC_HPDDM_COARSE_CORRECTION_ADDITIVE, 465 PC_HPDDM_COARSE_CORRECTION_BALANCED 466 } PCHPDDMCoarseCorrectionType; 467 468 /*E 469 PCFailedReason - indicates type of `PC` failure 470 471 Level: beginner 472 473 Developer Note: 474 Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h 475 E*/ 476 typedef enum { 477 PC_SETUP_ERROR = -1, 478 PC_NOERROR, 479 PC_FACTOR_STRUCT_ZEROPIVOT, 480 PC_FACTOR_NUMERIC_ZEROPIVOT, 481 PC_FACTOR_OUTMEMORY, 482 PC_FACTOR_OTHER, 483 PC_SUBPC_ERROR 484 } PCFailedReason; 485 486 /*E 487 PCGAMGLayoutType - Layout for reduced grids 488 489 Level: intermediate 490 491 .seealso: `PCGAMGSetCoarseGridLayoutType()` 492 493 Any additions/changes here MUST also be made in include/petsc/finclude/petscpc.h 494 E*/ 495 typedef enum { 496 PCGAMG_LAYOUT_COMPACT, 497 PCGAMG_LAYOUT_SPREAD 498 } PCGAMGLayoutType; 499 500 #endif 501