Preconditioning infrastructure¶
Firedrake has tight coupling with the PETSc library which provides support for a wide range of preconditioning strategies, see the relevant PETSc documentation for an overview.
In addition to these algebraic approaches, Firedrake offers a flexible
framework for defining preconditioners that need to construct or apply
auxiliary operators. The basic approach is described in
[KM18]. Here we provide a brief overview of the
preconditioners available in Firedrake that use this approach. To use
these preconditioners, one sets "pc_type": "python" and
"pc_python_type": "fully_qualified.NameOfPC" in the
solver_parameters.
Additive Schwarz methods¶
Small-block overlapping additive Schwarz preconditioners built on top of PCASM that can be used as components of robust multigrid schemes when using geometric multigrid.
- ASMPatchPC
- Abstract base class for which one must implement - ASMPatchPC.get_patches()to extract sets of degrees of freedom. Needs to be used with assembled sparse matrices (- "mat_type": "aij").
- ASMStarPC
- Constructs patches by gathering degrees of freedom in the star of specified mesh entities. 
- ASMVankaPC
- Constructs patches using the Vanka scheme by gathering degrees of freedom in the closure of the star of specified mesh entities. 
- ASMLinesmoothPC
- Constructs patches gathering degrees of freedom in vertical columns on - extruded meshes.
- ASMExtrudedStarPC
- Like - ASMStarPCbut on extruded meshes.
In addition to these algebraic approaches to constructing patches, Firedrake also interfaces with PCPATCH for both linear and nonlinear overlapping Schwarz methods. The approach is described in detail in [FKMW21]. These preconditioners can be used with both sparse matrices and Firedrake’s matrix-free operators, and can be applied either additively or multiplicatively within an MPI rank and additively between ranks.
- PatchPC
- Small-block overlapping Schwarz smoother with topological definition of patches. Does not support extruded meshes. 
- PatchSNES
- Nonlinear overlapping Schwarz smoother with topological definition of patches. Does not support extruded meshes. 
- PlaneSmoother
- A Python construction class for - PatchPCand- PatchSNESthat approximately groups mesh entities into lines or planes (useful for advection-dominated problems).
Note
The additive Schwarz preconditioners listed here construct patches around
mesh entities.  Crucially, the mesh must have an overlapping parallel domain
decomposition that supports the patches. This is set via the
distribution_parameters kwarg of the Mesh() constructor.  For
instance, vertex-star patches require
distribution_parameters["overlap_type"] = (DistributedMeshOverlapType.VERTEX, 1)
while Vanka patches require
distribution_parameters["overlap_type"] = (DistributedMeshOverlapType.VERTEX, 2)
Multigrid methods¶
Firedrake has support for rediscretised geometric multigrid on both
normal and extruded meshes, with regular refinement. This is obtained
by constructing a mesh hierarchy
and then using "pc_type": "mg". In addition to this basic support,
it also has out of the box support for a number of problem-specific
preconditioners.
- HypreADS
- An interface to Hypre’s auxiliary space divergence solver. Currently only implemented for lowest-order Raviart-Thomas elements. 
- HypreAMS
- An interface to Hypre’s auxiliary space Maxwell solver. Currently only implemented for lowest order Nedelec elements of the first kind. 
- PMGPC
- Generic p-coarsening rediscretised linear multigrid. If the problem is built on a mesh hierarchy then the coarse grid can do further h-multigrid with geometric coarsening. 
- P1PC
- Coarsening directly to linear elements. 
- PMGSNES
- Generic p-coarsening nonlinear multigrid. If the problem is built on a mesh hierarchy then the coarse grid can do further h-multigrid with geometric coarsening. 
- P1SNES
- Coarsening directly to linear elements. 
- GTMGPC
- A two-level non-nested multigrid scheme for the hybridised mixed method that operates on the trace variable, using the approach of [GT09]. The Firedrake implementation is described in [BGGMuller21]. 
Assembled and auxiliary operators¶
Many preconditioning schemes call for auxiliary operators, these are
facilitated by variations on Firedrake’s
AssembledPC which can be used to deliver an
assembled operator inside a nested solver where the outer matrix is a
matrix-free operator. Matrix-free operators can be used “natively”
with PETSc’s "jacobi" preconditioner, since they can provide their
diagonal cheaply. For more complicated things, one must assemble an
operator instead.
- AssembledPC
- Assemble an operator as a sparse matrix and then apply an inner preconditioner. For example, this might be used to assemble a coarse grid in an (otherwise matrix-free) multigrid solver. 
- AuxiliaryOperatorPC
- Abstract base class for preconditioners built from assembled auxiliary operators. One should subclass this preconditioner and override the - PCSNESBase.form()method. This can be used to provide bilinear forms to the solver that were not there in the original problem, for example, the pressure mass matrix for block preconditioners of the Stokes equations.
- FDMPC
- An auxiliary operator that uses piecewise-constant coefficients that is assembled in the basis of shape functions that diagonalize separable problems in the interior of each cell. Currently implemented for quadrilateral and hexahedral cells. The assembled matrix becomes as sparse as a low-order refined preconditioner, to which one may apply other preconditioners such as - ASMStarPCor- ASMExtrudedStarPC. See details in [BF22] and [BF24].
- MassInvPC
- Preconditioner for applying an inverse mass matrix. 
- PCDPC
- A preconditioner providing the Pressure-Convection-Diffusion approximation to the Schur complement for the Navier-Stokes equations. Note that this implementation only treats problems with characteristic velocity boundary conditions correctly. 
Hybridisation and element-wise static condensation¶
Firedrake has a number of preconditioners that use the Slate facility for element-wise linear algebra on
assembled tensors. These are described in detail in [GMHC20].
- HybridizationPC
- A preconditioner for hybridisable H(div) mixed methods that breaks the vector-valued space, and enforces continuity through introduction of a trace variable. The (now-broken) problem is eliminated element-wise onto the trace space to leave a single-variable global problem, whose solver can be configured. 
- SCPC
- A preconditioner that performs element-wise static condensation onto a single field. 
Jack Betteridge, Thomas H. Gibson, Ivan G. Graham, and Eike H. Müller. Multigrid preconditioners for the hybridised discontinuous Galerkin discretisation of the shallow water equations. Journal of Computational Physics, 426:109948, 2021. URL: https://www.sciencedirect.com/science/article/pii/S0021999120307221, doi:10.1016/j.jcp.2020.109948.
Patrick E. Farrell, Matthew G. Knepley, Lawrence Mitchell, and Florian Wechsung. PCPATCH: software for the topological construction of multigrid relaxation methods. ACM Transactions on Mathematical Software, 47(25):1–22, 2021. URL: https://arxiv.org/abs/1912.08516, arXiv:1912.08516, doi:10.1145/3445791.
Jayadeep Gopalakrishnan and Shuguang Tan. A convergent multigrid cycle for the hybridized mixed method. Numerical Linear Algebra with Applications, 16(9):689–714, sep 2009. URL: https://doi.org/10.1002/nla.636, doi:10.1002/nla.636.
