1a4963045SJacob Faibussowitsch #pragma once 2b0753f9dSMatthew G. Knepley 3ac09b921SBarry Smith /* SUBMANSEC = PC */ 4ac09b921SBarry Smith 5b0753f9dSMatthew G. Knepley /*S 687497f52SBarry Smith PC - Abstract PETSc object that manages all preconditioners including direct solvers such as `PCLU` 7b0753f9dSMatthew G. Knepley 8b0753f9dSMatthew G. Knepley Level: beginner 9b0753f9dSMatthew G. Knepley 1016a05f60SBarry Smith .seealso: [](doc_linsolve), [](sec_pc), `PCCreate()`, `PCSetType()`, `PCType` 11b0753f9dSMatthew G. Knepley S*/ 12b0753f9dSMatthew G. Knepley typedef struct _p_PC *PC; 13b0753f9dSMatthew G. Knepley 14b0753f9dSMatthew G. Knepley /*J 1516a05f60SBarry Smith PCType - String with the name of a PETSc preconditioner 16b0753f9dSMatthew G. Knepley 17b0753f9dSMatthew G. Knepley Level: beginner 18b0753f9dSMatthew G. Knepley 1987497f52SBarry Smith Note: 2087497f52SBarry Smith `PCRegister()` is used to register preconditioners that are then accessible via `PCSetType()` 21b0753f9dSMatthew G. Knepley 2216a05f60SBarry Smith .seealso: [](doc_linsolve), [](sec_pc), `PCSetType()`, `PC`, `PCCreate()`, `PCRegister()`, `PCSetFromOptions()`, `PCLU`, `PCJACOBI`, `PCBJACOBI` 23b0753f9dSMatthew G. Knepley J*/ 24b0753f9dSMatthew G. Knepley typedef const char *PCType; 25b0753f9dSMatthew G. Knepley #define PCNONE "none" 26b0753f9dSMatthew G. Knepley #define PCJACOBI "jacobi" 27b0753f9dSMatthew G. Knepley #define PCSOR "sor" 28b0753f9dSMatthew G. Knepley #define PCLU "lu" 29a2fc1e05SToby Isaac #define PCQR "qr" 30b0753f9dSMatthew G. Knepley #define PCSHELL "shell" 31e6f8f311SMark Adams #define PCAMGX "amgx" 32b0753f9dSMatthew G. Knepley #define PCBJACOBI "bjacobi" 33b0753f9dSMatthew G. Knepley #define PCMG "mg" 34b0753f9dSMatthew G. Knepley #define PCEISENSTAT "eisenstat" 35b0753f9dSMatthew G. Knepley #define PCILU "ilu" 36b0753f9dSMatthew G. Knepley #define PCICC "icc" 37b0753f9dSMatthew G. Knepley #define PCASM "asm" 38b0753f9dSMatthew G. Knepley #define PCGASM "gasm" 39b0753f9dSMatthew G. Knepley #define PCKSP "ksp" 40e607c864SMark Adams #define PCBJKOKKOS "bjkokkos" 41b0753f9dSMatthew G. Knepley #define PCCOMPOSITE "composite" 42b0753f9dSMatthew G. Knepley #define PCREDUNDANT "redundant" 43b0753f9dSMatthew G. Knepley #define PCSPAI "spai" 44b0753f9dSMatthew G. Knepley #define PCNN "nn" 45b0753f9dSMatthew G. Knepley #define PCCHOLESKY "cholesky" 46b0753f9dSMatthew G. Knepley #define PCPBJACOBI "pbjacobi" 470da83c2eSBarry Smith #define PCVPBJACOBI "vpbjacobi" 48b0753f9dSMatthew G. Knepley #define PCMAT "mat" 49b0753f9dSMatthew G. Knepley #define PCHYPRE "hypre" 50b0753f9dSMatthew G. Knepley #define PCPARMS "parms" 51b0753f9dSMatthew G. Knepley #define PCFIELDSPLIT "fieldsplit" 52b0753f9dSMatthew G. Knepley #define PCTFS "tfs" 53b0753f9dSMatthew G. Knepley #define PCML "ml" 54b0753f9dSMatthew G. Knepley #define PCGALERKIN "galerkin" 55b0753f9dSMatthew G. Knepley #define PCEXOTIC "exotic" 56b0753f9dSMatthew G. Knepley #define PCCP "cp" 57b0753f9dSMatthew G. Knepley #define PCBFBT "bfbt" 58b0753f9dSMatthew G. Knepley #define PCLSC "lsc" 59b0753f9dSMatthew G. Knepley #define PCPYTHON "python" 60b0753f9dSMatthew G. Knepley #define PCPFMG "pfmg" 611c188c59Sftrigaux #define PCSMG "smg" 62b0753f9dSMatthew G. Knepley #define PCSYSPFMG "syspfmg" 63b0753f9dSMatthew G. Knepley #define PCREDISTRIBUTE "redistribute" 64b0753f9dSMatthew G. Knepley #define PCSVD "svd" 65b0753f9dSMatthew G. Knepley #define PCGAMG "gamg" 664b3f184cSKarl Rupp #define PCCHOWILUVIENNACL "chowiluviennacl" 6770baa948SKarl Rupp #define PCROWSCALINGVIENNACL "rowscalingviennacl" 6807598726SKarl Rupp #define PCSAVIENNACL "saviennacl" 69b0753f9dSMatthew G. Knepley #define PCBDDC "bddc" 70b0753f9dSMatthew G. Knepley #define PCKACZMARZ "kaczmarz" 7168ddcbeaSBarry Smith #define PCTELESCOPE "telescope" 724bbf5ea8SMatthew G. Knepley #define PCPATCH "patch" 73b9ac7092SAlp Dener #define PCLMVM "lmvm" 74360ee056SFande Kong #define PCHMG "hmg" 7537eeb815SJakub Kruzik #define PCDEFLATION "deflation" 7638cfc46eSPierre Jolivet #define PCHPDDM "hpddm" 7753022affSStefano Zampini #define PCH2OPUS "h2opus" 78f1f2ae84SBarry Smith #define PCMPI "mpi" 79b0753f9dSMatthew G. Knepley 80b0753f9dSMatthew G. Knepley /*E 81b0753f9dSMatthew G. Knepley PCSide - If the preconditioner is to be applied to the left, right 82b0753f9dSMatthew G. Knepley or symmetrically around the operator. 83b0753f9dSMatthew G. Knepley 8416a05f60SBarry Smith Values: 8516a05f60SBarry Smith + `PC_LEFT` - applied after the operator is applied 8616a05f60SBarry Smith . `PC_RIGHT` - applied before the operator is applied 8716a05f60SBarry Smith - `PC_SYMMETRIC` - a portion of the preconditioner is applied before the operator and the transpose of this portion is applied after the operator is applied. 8816a05f60SBarry Smith 89b0753f9dSMatthew G. Knepley Level: beginner 90b0753f9dSMatthew G. Knepley 9116a05f60SBarry Smith Note: 9216a05f60SBarry Smith Certain `KSPType` support only a subset of `PCSide` values 9316a05f60SBarry Smith 9416a05f60SBarry Smith .seealso: [](sec_pc), `PC`, `KSPSetPCSide()` 95b0753f9dSMatthew G. Knepley E*/ 969371c9d4SSatish Balay typedef enum { 979371c9d4SSatish Balay PC_SIDE_DEFAULT = -1, 989371c9d4SSatish Balay PC_LEFT, 999371c9d4SSatish Balay PC_RIGHT, 1009371c9d4SSatish Balay PC_SYMMETRIC 1019371c9d4SSatish Balay } PCSide; 102b0753f9dSMatthew G. Knepley #define PC_SIDE_MAX (PC_SYMMETRIC + 1) 103b0753f9dSMatthew G. Knepley 104b0753f9dSMatthew G. Knepley /*E 10509e53591SBarry Smith PCRichardsonConvergedReason - reason a `PCRICHARDSON` `PCApplyRichardson()` method terminated 106b0753f9dSMatthew G. Knepley 107b0753f9dSMatthew G. Knepley Level: advanced 108b0753f9dSMatthew G. Knepley 10987497f52SBarry Smith Developer Note: 11016a05f60SBarry Smith This must match `include/petsc/finclude/petscpc.h` and the `KSPConvergedReason` values in `include/petscksp.h 111b0753f9dSMatthew G. Knepley 11209e53591SBarry Smith .seealso: [](sec_pc), `PCRICHARDSON`, `PC`, `PCApplyRichardson()` 113b0753f9dSMatthew G. Knepley E*/ 114b0753f9dSMatthew G. Knepley typedef enum { 115b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_RTOL = 2, 116b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ATOL = 3, 117b0753f9dSMatthew G. Knepley PCRICHARDSON_CONVERGED_ITS = 4, 1189371c9d4SSatish Balay PCRICHARDSON_DIVERGED_DTOL = -4 1199371c9d4SSatish Balay } PCRichardsonConvergedReason; 120b0753f9dSMatthew G. Knepley 121b0753f9dSMatthew G. Knepley /*E 12216a05f60SBarry Smith PCJacobiType - What elements of the matrix are used to form the Jacobi preconditioner 123b0753f9dSMatthew G. Knepley 12409e53591SBarry Smith Values: 12509e53591SBarry Smith + `PC_JACOBI_DIAGONAL` - use the diagonal entry, if it is zero use one 12609e53591SBarry Smith . `PC_JACOBI_ROWMAX` - use the maximum absolute value in the row 12709e53591SBarry Smith - `PC_JACOBI_ROWSUM` - use the sum of the values in the row (not the absolute values) 12809e53591SBarry Smith 129b0753f9dSMatthew G. Knepley Level: intermediate 130b0753f9dSMatthew G. Knepley 13109e53591SBarry Smith .seealso: [](sec_pc), `PCJACOBI`, `PC` 132b0753f9dSMatthew G. Knepley E*/ 1339371c9d4SSatish Balay typedef enum { 1349371c9d4SSatish Balay PC_JACOBI_DIAGONAL, 1359371c9d4SSatish Balay PC_JACOBI_ROWMAX, 1369371c9d4SSatish Balay PC_JACOBI_ROWSUM 1379371c9d4SSatish Balay } PCJacobiType; 138b0753f9dSMatthew G. Knepley 139b0753f9dSMatthew G. Knepley /*E 140b0753f9dSMatthew G. Knepley PCASMType - Type of additive Schwarz method to use 141b0753f9dSMatthew G. Knepley 14209e53591SBarry Smith Values: 14309e53591SBarry Smith + `PC_ASM_BASIC` - Symmetric version where residuals from the ghost points are used 14409e53591SBarry Smith and computed values in ghost regions are added together. 14509e53591SBarry Smith Classical standard additive Schwarz. 14609e53591SBarry Smith . `PC_ASM_RESTRICT` - Residuals from ghost points are used but computed values in ghost 14709e53591SBarry Smith region are discarded. 14809e53591SBarry Smith Default. 14909e53591SBarry Smith . `PC_ASM_INTERPOLATE` - Residuals from ghost points are not used, computed values in ghost 15009e53591SBarry Smith region are added back in. 15109e53591SBarry Smith - `PC_ASM_NONE` - Residuals from ghost points are not used, computed ghost values are 15209e53591SBarry Smith discarded. 15309e53591SBarry Smith Not very good. 154b0753f9dSMatthew G. Knepley 155b0753f9dSMatthew G. Knepley Level: beginner 156b0753f9dSMatthew G. Knepley 15716a05f60SBarry Smith .seealso: [](sec_pc), `PC`, `PCASM`, `PCASMSetType()`, `PCGASMType` 158b0753f9dSMatthew G. Knepley E*/ 1599371c9d4SSatish Balay typedef enum { 1609371c9d4SSatish Balay PC_ASM_BASIC = 3, 1619371c9d4SSatish Balay PC_ASM_RESTRICT = 1, 1629371c9d4SSatish Balay PC_ASM_INTERPOLATE = 2, 1639371c9d4SSatish Balay PC_ASM_NONE = 0 1649371c9d4SSatish Balay } PCASMType; 165b0753f9dSMatthew G. Knepley 166b0753f9dSMatthew G. Knepley /*E 16787497f52SBarry Smith PCGASMType - Type of generalized additive Schwarz method to use (differs from `PCASM` in allowing multiple processors per subdomain). 168b0753f9dSMatthew G. Knepley 16909e53591SBarry Smith Values: 17009e53591SBarry Smith + `PC_GASM_BASIC` - Symmetric version where the full from the outer subdomain is used, and the resulting correction is applied 17109e53591SBarry Smith over the outer subdomains. As a result, points in the overlap will receive the sum of the corrections 17209e53591SBarry Smith from neighboring subdomains. 17309e53591SBarry Smith Classical standard additive Schwarz. 17409e53591SBarry Smith . `PC_GASM_RESTRICT` - Residual from the outer subdomain is used but the correction is restricted to the inner subdomain only 17509e53591SBarry Smith (i.e., zeroed out over the overlap portion of the outer subdomain before being applied). As a result, 17609e53591SBarry Smith each point will receive a correction only from the unique inner subdomain containing it (nonoverlapping covering 17709e53591SBarry Smith assumption). 17809e53591SBarry Smith Default. 17909e53591SBarry Smith . `PC_GASM_INTERPOLATE` - Residual is zeroed out over the overlap portion of the outer subdomain, but the resulting correction is 18009e53591SBarry Smith applied over the outer subdomain. As a result, points in the overlap will receive the sum of the corrections 18109e53591SBarry Smith from neighboring subdomains. 18209e53591SBarry Smith - `PC_GASM_NONE` - Residuals and corrections are zeroed out outside the local subdomains. 18309e53591SBarry Smith Not very good. 184b0753f9dSMatthew G. Knepley 185b0753f9dSMatthew G. Knepley Level: beginner 186b0753f9dSMatthew G. Knepley 18716a05f60SBarry Smith Note: 18816a05f60SBarry Smith Each subdomain has nested inner and outer parts. The inner subdomains are assumed to form a non-overlapping covering of the computational 18916a05f60SBarry Smith domain, while the outer subdomains contain the inner subdomains and overlap with each other. This preconditioner will compute 19016a05f60SBarry Smith a subdomain correction over each *outer* subdomain from a residual computed there, but its different variants will differ in 19116a05f60SBarry Smith (a) how the outer subdomain residual is computed, and (b) how the outer subdomain correction is computed. 19216a05f60SBarry Smith 19316a05f60SBarry Smith .seealso: [](sec_pc), `PCGASM`, `PCASM`, `PC`, `PCGASMSetType()`, `PCASMType` 194b0753f9dSMatthew G. Knepley E*/ 1959371c9d4SSatish Balay typedef enum { 1969371c9d4SSatish Balay PC_GASM_BASIC = 3, 1979371c9d4SSatish Balay PC_GASM_RESTRICT = 1, 1989371c9d4SSatish Balay PC_GASM_INTERPOLATE = 2, 1999371c9d4SSatish Balay PC_GASM_NONE = 0 2009371c9d4SSatish Balay } PCGASMType; 201b0753f9dSMatthew G. Knepley 202b0753f9dSMatthew G. Knepley /*E 20316a05f60SBarry Smith PCCompositeType - Determines how two or more preconditioner are composed with the `PCType` of `PCCOMPOSITE` 204b0753f9dSMatthew G. Knepley 20509e53591SBarry Smith Values: 20609e53591SBarry Smith + `PC_COMPOSITE_ADDITIVE` - results from application of all preconditioners are added together 20709e53591SBarry Smith . `PC_COMPOSITE_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly 20809e53591SBarry Smith computed after the previous preconditioner application 20909e53591SBarry Smith . `PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE` - preconditioners are applied sequentially to the residual freshly 21009e53591SBarry Smith computed from first preconditioner to last and then back (Use only for symmetric matrices and preconditioners) 21109e53591SBarry Smith . `PC_COMPOSITE_SPECIAL` - This is very special for a matrix of the form alpha I + R + S 21209e53591SBarry Smith where first preconditioner is built from alpha I + S and second from 21309e53591SBarry Smith alpha I + R 21409e53591SBarry Smith . `PC_COMPOSITE_SCHUR` - composes the Schur complement of the matrix from two blocks, see `PCFIELDSPLIT` 21509e53591SBarry Smith - `PC_COMPOSITE_GKB` - the generalized Golub-Kahan bidiagonalization preconditioner, see `PCFIELDSPLIT` 216b0753f9dSMatthew G. Knepley 217b0753f9dSMatthew G. Knepley Level: beginner 218b0753f9dSMatthew G. Knepley 21909e53591SBarry Smith .seealso: [](sec_pc), `PCCOMPOSITE`, `PCFIELDSPLIT`, `PC`, `PCCompositeSetType()` 220b0753f9dSMatthew G. Knepley E*/ 2219371c9d4SSatish Balay typedef enum { 2229371c9d4SSatish Balay PC_COMPOSITE_ADDITIVE, 2239371c9d4SSatish Balay PC_COMPOSITE_MULTIPLICATIVE, 2249371c9d4SSatish Balay PC_COMPOSITE_SYMMETRIC_MULTIPLICATIVE, 2259371c9d4SSatish Balay PC_COMPOSITE_SPECIAL, 2269371c9d4SSatish Balay PC_COMPOSITE_SCHUR, 2279371c9d4SSatish Balay PC_COMPOSITE_GKB 2289371c9d4SSatish Balay } PCCompositeType; 229b0753f9dSMatthew G. Knepley 230b0753f9dSMatthew G. Knepley /*E 23187497f52SBarry Smith PCFieldSplitSchurPreType - Determines how to precondition a Schur complement 232b0753f9dSMatthew G. Knepley 23316a05f60SBarry Smith Values: 234*13044ca3SPierre Jolivet + `PC_FIELDSPLIT_SCHUR_PRE_SELF` - the preconditioner for the Schur complement is generated from the symbolic representation of the Schur complement matrix. 235*13044ca3SPierre Jolivet The only preconditioners that currently work with this symbolic representation matrix object are `PCLSC` and `PCHPDDM` 236*13044ca3SPierre Jolivet . `PC_FIELDSPLIT_SCHUR_PRE_SELFP` - the preconditioning for the Schur complement is generated from an explicitly-assembled approximation Sp = A11 - A10 inv(diag(A00)) A01. 23716a05f60SBarry Smith This is only a good preconditioner when diag(A00) is a good preconditioner for A00. Optionally, A00 can be 23816a05f60SBarry Smith lumped before extracting the diagonal using the additional option `-fieldsplit_1_mat_schur_complement_ainv_type lump` 23916a05f60SBarry Smith . `PC_FIELDSPLIT_SCHUR_PRE_A11` - the preconditioner for the Schur complement is generated from the block diagonal part of the matrix used to define the preconditioner, 24016a05f60SBarry Smith associated with the Schur complement (i.e. A11), not the Schur complement matrix 24116a05f60SBarry Smith . `PC_FIELDSPLIT_SCHUR_PRE_USER` - the preconditioner for the Schur complement is generated from the user provided matrix (pre argument 24216a05f60SBarry Smith to this function). 24316a05f60SBarry Smith - `PC_FIELDSPLIT_SCHUR_PRE_FULL` - the preconditioner for the Schur complement is generated from the exact Schur complement matrix representation 24416a05f60SBarry Smith computed internally by `PCFIELDSPLIT` (this is expensive) useful mostly as a test that the Schur complement approach can work for your problem 24516a05f60SBarry Smith 246b0753f9dSMatthew G. Knepley Level: intermediate 247b0753f9dSMatthew G. Knepley 24809e53591SBarry Smith .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurPre()`, `PC` 249b0753f9dSMatthew G. Knepley E*/ 2509371c9d4SSatish Balay typedef enum { 2519371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_SELF, 2529371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_SELFP, 2539371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_A11, 2549371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_USER, 2559371c9d4SSatish Balay PC_FIELDSPLIT_SCHUR_PRE_FULL 2569371c9d4SSatish Balay } PCFieldSplitSchurPreType; 257b0753f9dSMatthew G. Knepley 258b0753f9dSMatthew G. Knepley /*E 259b0753f9dSMatthew G. Knepley PCFieldSplitSchurFactType - determines which off-diagonal parts of the approximate block factorization to use 260b0753f9dSMatthew G. Knepley 26116a05f60SBarry Smith Values: 26216a05f60SBarry Smith + `PC_FIELDSPLIT_SCHUR_FACT_DIAG` - the preconditioner is solving `D` 26316a05f60SBarry Smith . `PC_FIELDSPLIT_SCHUR_FACT_LOWER` - the preconditioner is solving `L D` 26416a05f60SBarry Smith . `PC_FIELDSPLIT_SCHUR_FACT_UPPER` - the preconditioner is solving `D U` 26516a05f60SBarry Smith - `PC_FIELDSPLIT_SCHUR_FACT_FULL` - the preconditioner is solving `L(D U)` 26616a05f60SBarry Smith 26716a05f60SBarry Smith where the matrix is factorized as 26816a05f60SBarry Smith .vb 26916a05f60SBarry Smith (A B) = (1 0) (A 0) (1 Ainv*B) = L D U 27016a05f60SBarry Smith (C E) (C*Ainv 1) (0 S) (0 1) 27116a05f60SBarry Smith .ve 27216a05f60SBarry Smith 273b0753f9dSMatthew G. Knepley Level: intermediate 274b0753f9dSMatthew G. Knepley 27509e53591SBarry Smith .seealso: [](sec_pc), `PCFIELDSPLIT`, `PCFieldSplitSetSchurFactType()`, `PC` 276b0753f9dSMatthew G. Knepley E*/ 277b0753f9dSMatthew G. Knepley typedef enum { 278b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_DIAG, 279b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_LOWER, 280b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_UPPER, 281b0753f9dSMatthew G. Knepley PC_FIELDSPLIT_SCHUR_FACT_FULL 282b0753f9dSMatthew G. Knepley } PCFieldSplitSchurFactType; 283b0753f9dSMatthew G. Knepley 284b0753f9dSMatthew G. Knepley /*E 28587497f52SBarry Smith PCPARMSGlobalType - Determines the global preconditioner method in `PCPARMS` 286b0753f9dSMatthew G. Knepley 287b0753f9dSMatthew G. Knepley Level: intermediate 288b0753f9dSMatthew G. Knepley 28909e53591SBarry Smith .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetGlobal()`, `PC` 290b0753f9dSMatthew G. Knepley E*/ 2919371c9d4SSatish Balay typedef enum { 2929371c9d4SSatish Balay PC_PARMS_GLOBAL_RAS, 2939371c9d4SSatish Balay PC_PARMS_GLOBAL_SCHUR, 2949371c9d4SSatish Balay PC_PARMS_GLOBAL_BJ 2959371c9d4SSatish Balay } PCPARMSGlobalType; 2969d8081ecSMatthew G. Knepley 297b0753f9dSMatthew G. Knepley /*E 29887497f52SBarry Smith PCPARMSLocalType - Determines the local preconditioner method in `PCPARMS` 299b0753f9dSMatthew G. Knepley 300b0753f9dSMatthew G. Knepley Level: intermediate 301b0753f9dSMatthew G. Knepley 30209e53591SBarry Smith .seealso: [](sec_pc), `PCPARMS`, `PCPARMSSetLocal()`, `PC` 303b0753f9dSMatthew G. Knepley E*/ 3049371c9d4SSatish Balay typedef enum { 3059371c9d4SSatish Balay PC_PARMS_LOCAL_ILU0, 3069371c9d4SSatish Balay PC_PARMS_LOCAL_ILUK, 3079371c9d4SSatish Balay PC_PARMS_LOCAL_ILUT, 3089371c9d4SSatish Balay PC_PARMS_LOCAL_ARMS 3099371c9d4SSatish Balay } PCPARMSLocalType; 310b0753f9dSMatthew G. Knepley 311f0fc11ceSJed Brown /*J 31287497f52SBarry Smith PCGAMGType - type of generalized algebraic multigrid `PCGAMG` method 313b0753f9dSMatthew G. Knepley 31409e53591SBarry Smith Values: 31509e53591SBarry Smith + `PCGAMGAGG` - (the default) smoothed aggregation algorithm, robust, very well tested 3165e4ac8c8Smarkadams4 . `PCGAMGGEO` - geometric coarsening, uses mesh generator to produce coarser meshes, limited to triangles, not supported, referance implementation (2D) 3175e4ac8c8Smarkadams4 - `PCGAMGCLASSICAL` - classical algebraic multigrid preconditioner, incomplete, not supported, referance implementation 31809e53591SBarry Smith 319b0753f9dSMatthew G. Knepley Level: intermediate 320b0753f9dSMatthew G. Knepley 32109e53591SBarry Smith .seealso: [](sec_pc), `PCGAMG`, `PCMG`, `PC`, `PCSetType()`, `PCGAMGSetThreshold()`, `PCGAMGSetThreshold()`, `PCGAMGSetReuseInterpolation()` 322f0fc11ceSJed Brown J*/ 323b0753f9dSMatthew G. Knepley typedef const char *PCGAMGType; 324b0753f9dSMatthew G. Knepley #define PCGAMGAGG "agg" 325b0753f9dSMatthew G. Knepley #define PCGAMGGEO "geo" 326b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICAL "classical" 327b0753f9dSMatthew G. Knepley 328b0753f9dSMatthew G. Knepley typedef const char *PCGAMGClassicalType; 329b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALDIRECT "direct" 330b0753f9dSMatthew G. Knepley #define PCGAMGCLASSICALSTANDARD "standard" 331b0753f9dSMatthew G. Knepley 332b0753f9dSMatthew G. Knepley /*E 333b0753f9dSMatthew G. Knepley PCMGType - Determines the type of multigrid method that is run. 334b0753f9dSMatthew G. Knepley 335b0753f9dSMatthew G. Knepley Values: 33687497f52SBarry Smith + `PC_MG_MULTIPLICATIVE` (default) - traditional V or W cycle as determined by `PCMGSetCycleType()` 33787497f52SBarry Smith . `PC_MG_ADDITIVE` - the additive multigrid preconditioner where all levels are 338b0753f9dSMatthew G. Knepley smoothed before updating the residual. This only uses the 339b0753f9dSMatthew G. Knepley down smoother, in the preconditioner the upper smoother is ignored 34087497f52SBarry Smith . `PC_MG_FULL` - same as multiplicative except one also performs grid sequencing, 341b0753f9dSMatthew G. Knepley that is starts on the coarsest grid, performs a cycle, interpolates 342b0753f9dSMatthew G. Knepley 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 343b0753f9dSMatthew G. Knepley algorithm supports smoothing on before the restriction on each level in the initial restriction to the coarsest stage. In addition that algorithm 344b0753f9dSMatthew G. Knepley calls the V-cycle only on the coarser level and has a post-smoother instead. 34516a05f60SBarry Smith - `PC_MG_KASKADE` - like full multigrid except one never goes back to a coarser level from a finer 34616a05f60SBarry Smith 34716a05f60SBarry Smith Level: beginner 348b0753f9dSMatthew G. Knepley 34909e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetType()`, `PCMGSetCycleType()`, `PCMGSetCycleTypeOnLevel()` 350b0753f9dSMatthew G. Knepley E*/ 3519371c9d4SSatish Balay typedef enum { 3529371c9d4SSatish Balay PC_MG_MULTIPLICATIVE, 3539371c9d4SSatish Balay PC_MG_ADDITIVE, 3549371c9d4SSatish Balay PC_MG_FULL, 3559371c9d4SSatish Balay PC_MG_KASKADE 3569371c9d4SSatish Balay } PCMGType; 357b0753f9dSMatthew G. Knepley #define PC_MG_CASCADE PC_MG_KASKADE; 358b0753f9dSMatthew G. Knepley 359b0753f9dSMatthew G. Knepley /*E 360b0753f9dSMatthew G. Knepley PCMGCycleType - Use V-cycle or W-cycle 361b0753f9dSMatthew G. Knepley 362b0753f9dSMatthew G. Knepley Values: 36387497f52SBarry Smith + `PC_MG_V_CYCLE` - use the v cycle 36487497f52SBarry Smith - `PC_MG_W_CYCLE` - use the w cycle 365b0753f9dSMatthew G. Knepley 36616a05f60SBarry Smith Level: beginner 36716a05f60SBarry Smith 36809e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()` 369b0753f9dSMatthew G. Knepley E*/ 3709371c9d4SSatish Balay typedef enum { 3719371c9d4SSatish Balay PC_MG_CYCLE_V = 1, 3729371c9d4SSatish Balay PC_MG_CYCLE_W = 2 3739371c9d4SSatish Balay } PCMGCycleType; 374b0753f9dSMatthew G. Knepley 375b0753f9dSMatthew G. Knepley /*E 3762134b1e4SBarry Smith PCMGalerkinType - Determines if the coarse grid operators are computed via the Galerkin process 3772134b1e4SBarry Smith 3782134b1e4SBarry Smith Values: 37987497f52SBarry Smith + `PC_MG_GALERKIN_PMAT` - computes the pmat (matrix from which the preconditioner is built) via the Galerkin process from the finest grid 38087497f52SBarry Smith . `PC_MG_GALERKIN_MAT` - computes the mat (matrix used to apply the operator) via the Galerkin process from the finest grid 38187497f52SBarry Smith . `PC_MG_GALERKIN_BOTH` - computes both the mat and pmat via the Galerkin process (if pmat == mat the construction is only done once 38287497f52SBarry Smith - `PC_MG_GALERKIN_NONE` - neither operator is computed via the Galerkin process, the user must provide the operator 3832134b1e4SBarry Smith 38409e53591SBarry Smith Level: beginner 38509e53591SBarry Smith 38609e53591SBarry Smith Note: 38787497f52SBarry Smith Users should never set `PC_MG_GALERKIN_EXTERNAL`, it is used by `PCHYPRE` and `PCML` 3882134b1e4SBarry Smith 38909e53591SBarry Smith .seealso: [](sec_pc), `PCMG`, `PC`, `PCMGSetCycleType()` 3902134b1e4SBarry Smith E*/ 3919371c9d4SSatish Balay typedef enum { 3929371c9d4SSatish Balay PC_MG_GALERKIN_BOTH, 3939371c9d4SSatish Balay PC_MG_GALERKIN_PMAT, 3949371c9d4SSatish Balay PC_MG_GALERKIN_MAT, 3959371c9d4SSatish Balay PC_MG_GALERKIN_NONE, 3969371c9d4SSatish Balay PC_MG_GALERKIN_EXTERNAL 3979371c9d4SSatish Balay } PCMGGalerkinType; 3982134b1e4SBarry Smith 3992134b1e4SBarry Smith /*E 400b0753f9dSMatthew G. Knepley PCExoticType - Face based or wirebasket based coarse grid space 401b0753f9dSMatthew G. Knepley 402b0753f9dSMatthew G. Knepley Level: beginner 403b0753f9dSMatthew G. Knepley 40409e53591SBarry Smith .seealso: [](sec_pc), `PCExoticSetType()`, `PCEXOTIC` 405b0753f9dSMatthew G. Knepley E*/ 4069371c9d4SSatish Balay typedef enum { 4079371c9d4SSatish Balay PC_EXOTIC_FACE, 4089371c9d4SSatish Balay PC_EXOTIC_WIREBASKET 4099371c9d4SSatish Balay } PCExoticType; 410b0753f9dSMatthew G. Knepley 4118c1cd74cSHong Zhang /*E 412bc960bbfSJed Brown PCBDDCInterfaceExtType - Defines how interface balancing is extended into the interior of subdomains. 413bc960bbfSJed Brown 414bc960bbfSJed Brown Values: 41587497f52SBarry Smith + `PC_BDDC_INTERFACE_EXT_DIRICHLET` - solves Dirichlet interior problem; this is the standard BDDC algorithm 41687497f52SBarry Smith - `PC_BDDC_INTERFACE_EXT_LUMP` - skips interior solve; sometimes called M_1 and associated with "lumped FETI-DP" 417bc960bbfSJed Brown 41816a05f60SBarry Smith Level: intermediate 41916a05f60SBarry Smith 42009e53591SBarry Smith .seealso: [](sec_pc), `PCBDDC`, `PC` 421bc960bbfSJed Brown E*/ 422bc960bbfSJed Brown typedef enum { 423bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_DIRICHLET, 424bc960bbfSJed Brown PC_BDDC_INTERFACE_EXT_LUMP 425bc960bbfSJed Brown } PCBDDCInterfaceExtType; 426bc960bbfSJed Brown 427bc960bbfSJed Brown /*E 428f3b08a26SMatthew G. Knepley PCMGCoarseSpaceType - Function space for coarse space for adaptive interpolation 429f3b08a26SMatthew G. Knepley 430f3b08a26SMatthew G. Knepley Level: beginner 431f3b08a26SMatthew G. Knepley 43209e53591SBarry Smith .seealso: [](sec_pc), `PCMGSetAdaptCoarseSpaceType()`, `PCMG`, `PC` 433f3b08a26SMatthew G. Knepley E*/ 4349371c9d4SSatish Balay typedef enum { 4359371c9d4SSatish Balay PCMG_ADAPT_NONE, 4369371c9d4SSatish Balay PCMG_ADAPT_POLYNOMIAL, 4379371c9d4SSatish Balay PCMG_ADAPT_HARMONIC, 4389371c9d4SSatish Balay PCMG_ADAPT_EIGENVECTOR, 4399371c9d4SSatish Balay PCMG_ADAPT_GENERALIZED_EIGENVECTOR, 4409371c9d4SSatish Balay PCMG_ADAPT_GDSW 4419371c9d4SSatish Balay } PCMGCoarseSpaceType; 442f3b08a26SMatthew G. Knepley 443f3b08a26SMatthew G. Knepley /*E 44416a05f60SBarry Smith PCPatchConstructType - The algorithm used to construct patches for the `PCPATCH` preconditioner 4454bbf5ea8SMatthew G. Knepley 4464bbf5ea8SMatthew G. Knepley Level: beginner 4474bbf5ea8SMatthew G. Knepley 44809e53591SBarry Smith .seealso: [](sec_pc), `PCPatchSetConstructType()`, `PCPATCH`, `PC` 4494bbf5ea8SMatthew G. Knepley E*/ 4509371c9d4SSatish Balay typedef enum { 4519371c9d4SSatish Balay PC_PATCH_STAR, 4529371c9d4SSatish Balay PC_PATCH_VANKA, 4539371c9d4SSatish Balay PC_PATCH_PARDECOMP, 4549371c9d4SSatish Balay PC_PATCH_USER, 4559371c9d4SSatish Balay PC_PATCH_PYTHON 4569371c9d4SSatish Balay } PCPatchConstructType; 4574bbf5ea8SMatthew G. Knepley 4584bbf5ea8SMatthew G. Knepley /*E 459e26ad46dSJakub Kruzik PCDeflationSpaceType - Type of deflation 460e662bc50SJakub Kruzik 461e662bc50SJakub Kruzik Values: 46287497f52SBarry Smith + `PC_DEFLATION_SPACE_HAAR` - directly assembled based on Haar (db2) wavelet with overflowed filter cuted-off 46309e53591SBarry Smith . `PC_DEFLATION_SPACE_DB2` - `MATCOMPOSITE` of 1-lvl matices based on db2 (2 coefficient Daubechies / Haar wavelet) 46487497f52SBarry Smith . `PC_DEFLATION_SPACE_DB4` - same as above, but with db4 (4 coefficient Daubechies) 46587497f52SBarry Smith . `PC_DEFLATION_SPACE_DB8` - same as above, but with db8 (8 coefficient Daubechies) 46687497f52SBarry Smith . `PC_DEFLATION_SPACE_DB16` - same as above, but with db16 (16 coefficient Daubechies) 46787497f52SBarry Smith . `PC_DEFLATION_SPACE_BIORTH22` - same as above, but with biorthogonal 2.2 (6 coefficients) 46887497f52SBarry Smith . `PC_DEFLATION_SPACE_MEYER` - same as above, but with Meyer/FIR (62 coefficients) 469d5b43468SJose E. Roman . `PC_DEFLATION_SPACE_AGGREGATION` - aggregates local indices (given by operator matrix distribution) into a subdomain 47087497f52SBarry Smith - `PC_DEFLATION_SPACE_USER` - indicates space set by user 471e662bc50SJakub Kruzik 4721fdb00f9SJakub Kruzik Level: intermediate 4731fdb00f9SJakub Kruzik 47409e53591SBarry Smith Note: 47509e53591SBarry Smith Wavelet-based space (except Haar) can be used in multilevel deflation. 47609e53591SBarry Smith 47709e53591SBarry Smith .seealso: [](sec_pc), `PCDeflationSetSpaceToCompute()`, `PCDEFLATION`, `PC` 478e662bc50SJakub Kruzik E*/ 479e662bc50SJakub Kruzik typedef enum { 480e662bc50SJakub Kruzik PC_DEFLATION_SPACE_HAAR, 481e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB2, 482e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB4, 483e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB8, 484e662bc50SJakub Kruzik PC_DEFLATION_SPACE_DB16, 485e662bc50SJakub Kruzik PC_DEFLATION_SPACE_BIORTH22, 486e662bc50SJakub Kruzik PC_DEFLATION_SPACE_MEYER, 487e662bc50SJakub Kruzik PC_DEFLATION_SPACE_AGGREGATION, 488e662bc50SJakub Kruzik PC_DEFLATION_SPACE_USER 489e662bc50SJakub Kruzik } PCDeflationSpaceType; 490e662bc50SJakub Kruzik 491e662bc50SJakub Kruzik /*E 49287497f52SBarry Smith PCHPDDMCoarseCorrectionType - Type of coarse correction used by `PCHPDDM` 49338cfc46eSPierre Jolivet 49438cfc46eSPierre Jolivet Values: 49509e53591SBarry Smith + `PC_HPDDM_COARSE_CORRECTION_DEFLATED` (default) - eq. (1) in `PCHPDDMShellApply()` 49687497f52SBarry Smith . `PC_HPDDM_COARSE_CORRECTION_ADDITIVE` - eq. (2) 49787497f52SBarry Smith - `PC_HPDDM_COARSE_CORRECTION_BALANCED` - eq. (3) 49838cfc46eSPierre Jolivet 49916a05f60SBarry Smith Level: intermediate 50016a05f60SBarry Smith 50109e53591SBarry Smith .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCSetType()`, `PCHPDDMShellApply()` 50238cfc46eSPierre Jolivet E*/ 5039371c9d4SSatish Balay typedef enum { 5049371c9d4SSatish Balay PC_HPDDM_COARSE_CORRECTION_DEFLATED, 5059371c9d4SSatish Balay PC_HPDDM_COARSE_CORRECTION_ADDITIVE, 5069371c9d4SSatish Balay PC_HPDDM_COARSE_CORRECTION_BALANCED 5079371c9d4SSatish Balay } PCHPDDMCoarseCorrectionType; 50838cfc46eSPierre Jolivet 50938cfc46eSPierre Jolivet /*E 510*13044ca3SPierre Jolivet PCHPDDMSchurPreType - Type of `PCHPDDM` preconditioner for a `MATSCHURCOMPLEMENT` generated by `PCFIELDSPLIT` with `PCFieldSplitSchurPreType` set to `PC_FIELDSPLIT_SCHUR_PRE_SELF` 511*13044ca3SPierre Jolivet 512*13044ca3SPierre Jolivet Values: 513*13044ca3SPierre Jolivet + `PC_HPDDM_SCHUR_PRE_LEAST_SQUARES` (default) - only with a near-zero A11 block and A10 = A01^T; a preconditioner for solving A01^T A00^-1 A01 x = b is built by approximating the Schur complement with (inv(sqrt(diag(A00))) A01)^T (inv(sqrt(diag(A00))) A01) and by considering the associated linear least squares problem 514*13044ca3SPierre Jolivet - `PC_HPDDM_SCHUR_PRE_GENEO` - only with A10 = A01^T, `PCHPDDMSetAuxiliaryMat()` called on the `PC` of the A00 block, and if A11 is nonzero, then `PCHPDDMSetAuxiliaryMat()` must be called on the associated `PC` as well (it is built automatically for the user otherwise); the Schur complement `PC` is set internally to `PCKSP`, with the prefix `-fieldsplit_1_pc_hpddm_`; the operator associated to the `PC` is spectrally equivalent to the original Schur complement 515*13044ca3SPierre Jolivet 516*13044ca3SPierre Jolivet Level: advanced 517*13044ca3SPierre Jolivet 518*13044ca3SPierre Jolivet .seealso: [](sec_pc), `PCHPDDM`, `PC`, `PCFIELDSPLIT`, `PC_FIELDSPLIT_SCHUR_PRE_SELF`, `PCFieldSplitSetSchurPre()`, `PCHPDDMSetAuxiliaryMat()` 519*13044ca3SPierre Jolivet E*/ 520*13044ca3SPierre Jolivet typedef enum { 521*13044ca3SPierre Jolivet PC_HPDDM_SCHUR_PRE_LEAST_SQUARES, 522*13044ca3SPierre Jolivet PC_HPDDM_SCHUR_PRE_GENEO, 523*13044ca3SPierre Jolivet } PCHPDDMSchurPreType; 524*13044ca3SPierre Jolivet 525*13044ca3SPierre Jolivet /*E 52687497f52SBarry Smith PCFailedReason - indicates type of `PC` failure 5278c1cd74cSHong Zhang 5288c1cd74cSHong Zhang Level: beginner 5298c1cd74cSHong Zhang 53087497f52SBarry Smith Developer Note: 53116a05f60SBarry Smith Any additions/changes here MUST also be made in `include/petsc/finclude/petscpc.h` 53209e53591SBarry Smith 53309e53591SBarry Smith .seealso: [](sec_pc), `PC` 5348c1cd74cSHong Zhang E*/ 5359371c9d4SSatish Balay typedef enum { 5369371c9d4SSatish Balay PC_SETUP_ERROR = -1, 5379371c9d4SSatish Balay PC_NOERROR, 5389371c9d4SSatish Balay PC_FACTOR_STRUCT_ZEROPIVOT, 5399371c9d4SSatish Balay PC_FACTOR_NUMERIC_ZEROPIVOT, 5409371c9d4SSatish Balay PC_FACTOR_OUTMEMORY, 5419371c9d4SSatish Balay PC_FACTOR_OTHER, 54287b47708SBarry Smith PC_INCONSISTENT_RHS, 5439371c9d4SSatish Balay PC_SUBPC_ERROR 5449371c9d4SSatish Balay } PCFailedReason; 545ce7c7f2fSMark Adams 546ce7c7f2fSMark Adams /*E 547ce7c7f2fSMark Adams PCGAMGLayoutType - Layout for reduced grids 548ce7c7f2fSMark Adams 549ce7c7f2fSMark Adams Level: intermediate 550ce7c7f2fSMark Adams 55109e53591SBarry Smith Developer Note: 55216a05f60SBarry Smith Any additions/changes here MUST also be made in `include/petsc/finclude/petscpc.h` 55309e53591SBarry Smith 55409e53591SBarry Smith .seealso: [](sec_pc), `PCGAMG`, `PC`, `PCGAMGSetCoarseGridLayoutType()` 555ce7c7f2fSMark Adams E*/ 5569371c9d4SSatish Balay typedef enum { 5579371c9d4SSatish Balay PCGAMG_LAYOUT_COMPACT, 5589371c9d4SSatish Balay PCGAMG_LAYOUT_SPREAD 5599371c9d4SSatish Balay } PCGAMGLayoutType; 560