firedrake.preconditioners package¶
Submodules¶
firedrake.preconditioners.asm module¶
- class firedrake.preconditioners.asm.ASMExtrudedStarPC[source]¶
Bases:
ASMStarPC
Patch-based PC using Star of mesh entities implmented as an
ASMPatchPC
.ASMExtrudedStarPC is an additive Schwarz preconditioner where each patch consists of all DoFs on the topological star of the mesh entity specified by pc_star_construct_dim.
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- get_patches(V)[source]¶
Get the patches used for PETSc PSASM
- Parameters:
V – the
FunctionSpace
.- Returns:
a list of index sets defining the ASM patches in local numbering (before lgmap.apply has been called).
- class firedrake.preconditioners.asm.ASMLinesmoothPC[source]¶
Bases:
ASMPatchPC
Linesmoother PC for extruded meshes implemented as an
ASMPatchPC
.ASMLinesmoothPC is an additive Schwarz preconditioner where each patch consists of all dofs associated with a vertical column (and hence extruded meshes are necessary). Three types of columns are possible: columns of horizontal faces (each column built over a face of the base mesh), columns of vertical faces (each column built over an edge of the base mesh), and columns of vertical edges (each column built over a vertex of the base mesh).
To select the column type or types for the patches, use ‘pc_linesmooth_codims’ to set integers giving the codimension of the base mesh entities for the columns. For example, ‘pc_linesmooth_codims 0,1’ creates patches for each cell and each facet of the base mesh.
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- get_patches(V)[source]¶
Get the patches used for PETSc PSASM
- Parameters:
V – the
FunctionSpace
.- Returns:
a list of index sets defining the ASM patches in local numbering (before lgmap.apply has been called).
- class firedrake.preconditioners.asm.ASMPatchPC[source]¶
Bases:
PCBase
PC for PETSc PCASM
should implement: -
get_patches()
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- apply(pc, x, y)[source]¶
Apply the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- applyTranspose(pc, x, y)[source]¶
Apply the transpose of the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- abstract get_patches(V)[source]¶
Get the patches used for PETSc PSASM
- Parameters:
V – the
FunctionSpace
.- Returns:
a list of index sets defining the ASM patches in local numbering (before lgmap.apply has been called).
- class firedrake.preconditioners.asm.ASMStarPC[source]¶
Bases:
ASMPatchPC
Patch-based PC using Star of mesh entities implmented as an
ASMPatchPC
.ASMStarPC is an additive Schwarz preconditioner where each patch consists of all DoFs on the topological star of the mesh entity specified by pc_star_construct_dim.
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- get_patches(V)[source]¶
Get the patches used for PETSc PSASM
- Parameters:
V – the
FunctionSpace
.- Returns:
a list of index sets defining the ASM patches in local numbering (before lgmap.apply has been called).
- class firedrake.preconditioners.asm.ASMVankaPC[source]¶
Bases:
ASMPatchPC
Patch-based PC using closure of star of mesh entities implmented as an
ASMPatchPC
.ASMVankaPC is an additive Schwarz preconditioner where each patch consists of all DoFs on the closure of the star of the mesh entity specified by pc_vanka_construct_dim (or codim).
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- get_patches(V)[source]¶
Get the patches used for PETSc PSASM
- Parameters:
V – the
FunctionSpace
.- Returns:
a list of index sets defining the ASM patches in local numbering (before lgmap.apply has been called).
firedrake.preconditioners.assembled module¶
- class firedrake.preconditioners.assembled.AssembledPC[source]¶
Bases:
PCBase
A matrix-free PC that assembles the operator.
Internally this makes a PETSc PC object that can be controlled by options using the extra options prefix
assembled_
.Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- apply(pc, x, y)[source]¶
Apply the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- class firedrake.preconditioners.assembled.AuxiliaryOperatorPC[source]¶
Bases:
AssembledPC
A preconditioner that builds a PC on a specified form. Mainly used for describing approximations to Schur complements.
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- abstract form(pc, test, trial)[source]¶
- Parameters:
pc – a PETSc.PC object. Use self.get_appctx(pc) to get the user-supplied application-context, if desired.
test – a TestFunction on this FunctionSpace.
trial – a TrialFunction on this FunctionSpace.
:returns (a, bcs), where a is a bilinear Form and bcs is a list of DirichletBC boundary conditions (possibly None).
firedrake.preconditioners.base module¶
- class firedrake.preconditioners.base.PCBase[source]¶
Bases:
PCSNESBase
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- abstract apply(pc, X, Y)[source]¶
Apply the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- abstract applyTranspose(pc, X, Y)[source]¶
Apply the transpose of the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- needs_python_amat = False¶
Set this to True if the A matrix needs to be Python (matfree).
- needs_python_pmat = False¶
Set this to False if the P matrix needs to be Python (matfree).
If the preconditioner also works with assembled matrices, then use False here.
- setUp(pc)[source]¶
Setup method called by PETSc.
Subclasses should probably not override this and instead implement
update()
andinitialize()
.
- class firedrake.preconditioners.base.PCSNESBase[source]¶
Bases:
object
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- static new_snes_ctx(pc, op, bcs, mat_type, fcp=None, options_prefix=None)[source]¶
Create a new SNES contex for nested preconditioning
- setUp(pc)[source]¶
Setup method called by PETSc.
Subclasses should probably not override this and instead implement
update()
andinitialize()
.
- class firedrake.preconditioners.base.SNESBase[source]¶
Bases:
PCSNESBase
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
firedrake.preconditioners.facet_split module¶
- class firedrake.preconditioners.facet_split.FacetSplitPC[source]¶
Bases:
PCBase
A preconditioner that splits a function into interior and facet DOFs.
Internally this creates a PETSc PC object that can be controlled by options using the extra options prefix
facet_
.This allows for statically-condensed preconditioners to be applied to linear systems involving the matrix applied to the full set of DOFs. Code generated for the matrix-free operator evaluation in the space with full DOFs will run faster than the one with interior-facet decoposition, since the full element has a simpler structure.
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- apply(pc, x, y)[source]¶
Apply the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- applyTranspose(pc, x, y)[source]¶
Apply the transpose of the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- needs_python_pmat = False¶
Set this to False if the P matrix needs to be Python (matfree).
If the preconditioner also works with assembled matrices, then use False here.
firedrake.preconditioners.fdm module¶
- class firedrake.preconditioners.fdm.FDMPC[source]¶
Bases:
PCBase
A preconditioner for tensor-product elements that changes the shape functions so that the H(d) (d in {grad, curl, div}) Riesz map is sparse on Cartesian cells, and assembles a global sparse matrix on which other preconditioners, such as ASMStarPC, can be applied.
Here we assume that the volume integrals in the Jacobian can be expressed as:
inner(d(v), alpha * d(u))*dx + inner(v, beta * u)*dx
where alpha and beta are possibly tensor-valued. The sparse matrix is obtained by approximating (v, alpha * u) and (v, beta * u) as diagonal mass matrices.
The PETSc options inspected by this class are: - ‘fdm_mat_type’: can be either ‘aij’ or ‘sbaij’ - ‘fdm_static_condensation’: are we assembling the Schur complement on facets?
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- allocate_matrix(V, J, bcs, fcp, pmat_type, use_static_condensation)[source]¶
Allocate the FDM sparse preconditioner.
- Parameters:
V – the
FunctionSpace
of the form argumentsJ – the Jacobian bilinear form
bcs – an iterable of boundary conditions on V
fcp – form compiler parameters to assemble coefficients
pmat_type – the PETSc.Mat.Type for the blocks in the diagonal
use_static_condensation – are we assembling the statically-condensed Schur complement on facets?
- Returns:
2-tuple with the preconditioner
PETSc.Mat
and a list of assembly callables
- apply(pc, x, y)[source]¶
Apply the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- applyTranspose(pc, x, y)[source]¶
Apply the transpose of the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- assemble_coefficients(J, fcp, block_diagonal=True)[source]¶
Obtain coefficients for the auxiliary operator as the diagonal of a weighted mass matrix in broken(V^k) * broken(V^{k+1}). See Section 3.2 of Brubeck2022b.
- Parameters:
J – the Jacobian bilinear
ufl.Form
,fcp – form compiler parameters to assemble the diagonal of the mass matrices.
block_diagonal – are we assembling the block diagonal of the mass matrices?
- Returns:
a 2-tuple of a dict with the zero-th order and second order coefficients keyed on
"beta"
and"alpha"
, and a list of assembly callables.
- assemble_reference_tensor(V, transpose=False, sort_interior=False)[source]¶
Return the reference tensor used in the diagonal factorisation of the sparse cell matrices. See Section 3.2 of Brubeck2022b.
- Parameters:
V – a
FunctionSpace
- Returns:
a
PETSc.Mat
interpolating V^k * d(V^k) onto broken(V^k) * broken(V^{k+1}) on the reference element.
- set_values(A, Vrow, Vcol, addv, triu=False)[source]¶
Assemble the stiffness matrix in the FDM basis using sparse reference tensors and diagonal mass matrices.
- Parameters:
A – the
PETSc.Mat
to assembleVrow – the
FunctionSpace
test spaceVcol – the
FunctionSpace
trial spaceaddv – a PETSc.Mat.InsertMode
triu – are we assembling only the upper triangular part?
- class firedrake.preconditioners.fdm.PoissonFDMPC[source]¶
Bases:
FDMPC
A preconditioner for tensor-product elements that changes the shape functions so that the H^1 Riesz map is sparse in the interior of a Cartesian cell, and assembles a global sparse matrix on which other preconditioners, such as ASMStarPC, can be applied.
Here we assume that the volume integrals in the Jacobian can be expressed as:
inner(grad(v), alpha(grad(u)))*dx + inner(v, beta(u))*dx
where alpha and beta are possibly tensor-valued. The sparse matrix is obtained by approximating alpha and beta by cell-wise constants and discarding the coefficients in alpha that couple together mixed derivatives and mixed components.
For spaces that are not H^1-conforming, this preconditioner will use the symmetric interior-penalty DG method. The penalty coefficient can be provided in the application context, keyed on
"eta"
.Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- assemble_coefficients(J, fcp)[source]¶
Obtain coefficients for the auxiliary operator as the diagonal of a weighted mass matrix in broken(V^k) * broken(V^{k+1}). See Section 3.2 of Brubeck2022b.
- Parameters:
J – the Jacobian bilinear
ufl.Form
,fcp – form compiler parameters to assemble the diagonal of the mass matrices.
block_diagonal – are we assembling the block diagonal of the mass matrices?
- Returns:
a 2-tuple of a dict with the zero-th order and second order coefficients keyed on
"beta"
and"alpha"
, and a list of assembly callables.
- assemble_reference_tensor(V)[source]¶
Return the reference tensor used in the diagonal factorisation of the sparse cell matrices. See Section 3.2 of Brubeck2022b.
- Parameters:
V – a
FunctionSpace
- Returns:
a
PETSc.Mat
interpolating V^k * d(V^k) onto broken(V^k) * broken(V^{k+1}) on the reference element.
- set_values(A, Vrow, Vcol, addv, triu=False)[source]¶
Assemble the stiffness matrix in the FDM basis using Kronecker products of interval matrices
- Parameters:
A – the
PETSc.Mat
to assembleVrow – the
FunctionSpace
test spaceVcol – the
FunctionSpace
trial spaceaddv – a PETSc.Mat.InsertMode
triu – are we assembling only the upper triangular part?
firedrake.preconditioners.gtmg module¶
- class firedrake.preconditioners.gtmg.GTMGPC[source]¶
Bases:
PCBase
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- apply(pc, X, Y)[source]¶
Apply the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- applyTranspose(pc, X, Y)[source]¶
Apply the transpose of the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- needs_python_pmat = False¶
Set this to False if the P matrix needs to be Python (matfree).
If the preconditioner also works with assembled matrices, then use False here.
firedrake.preconditioners.hiptmair module¶
- class firedrake.preconditioners.hiptmair.HiptmairPC[source]¶
Bases:
TwoLevelPC
A two-level method for H(curl) or H(div) problems with an auxiliary potential space in H^1 or H(curl), respectively.
Internally this creates a PETSc PCMG object that can be controlled by options using the extra options prefix
hiptmair_mg_
.This allows for effective multigrid relaxation methods with patch solves centered around vertices for H^1, edges for H(curl), or faces for H(div). For the lowest-order spaces this corresponds to point-Jacobi.
The H(div) auxiliary vector potential problem in H(curl) is singular for high-order. This can be overcome by pertubing the problem by a multiple of the mass matrix. The scaling factor can be provided (defaulting to 0) by providing a scalar in the application context, keyed on
"hiptmair_shift"
.Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- class firedrake.preconditioners.hiptmair.TwoLevelPC[source]¶
Bases:
PCBase
PC for two-level methods
should implement: -
coarsen()
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- apply(pc, X, Y)[source]¶
Apply the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- applyTranspose(pc, X, Y)[source]¶
Apply the transpose of the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- abstract coarsen(pc)[source]¶
Return a tuple with coarse bilinear form, coarse boundary conditions, and coarse-to-fine interpolation matrix
- needs_python_pmat = False¶
Set this to False if the P matrix needs to be Python (matfree).
If the preconditioner also works with assembled matrices, then use False here.
firedrake.preconditioners.hypre_ads module¶
- class firedrake.preconditioners.hypre_ads.HypreADS[source]¶
Bases:
PCBase
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- apply(pc, x, y)[source]¶
Apply the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
firedrake.preconditioners.hypre_ams module¶
- class firedrake.preconditioners.hypre_ams.HypreAMS[source]¶
Bases:
PCBase
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- apply(pc, x, y)[source]¶
Apply the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
firedrake.preconditioners.low_order module¶
- class firedrake.preconditioners.low_order.P1PC[source]¶
Bases:
PMGPC
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- coarsen_element(ele)[source]¶
Coarsen a given element to form the next problem down in the p-hierarchy.
If the supplied element should form the coarsest level of the p-hierarchy, raise ValueError. Otherwise, return a new
ufl.FiniteElement
.By default, this does power-of-2 coarsening in polynomial degree until we reach the coarse degree specified through PETSc options (1 by default).
- Parameters:
ele – a
ufl.FiniteElement
to coarsen.
- class firedrake.preconditioners.low_order.P1SNES[source]¶
Bases:
PMGSNES
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- coarsen_element(ele)[source]¶
Coarsen a given element to form the next problem down in the p-hierarchy.
If the supplied element should form the coarsest level of the p-hierarchy, raise ValueError. Otherwise, return a new
ufl.FiniteElement
.By default, this does power-of-2 coarsening in polynomial degree until we reach the coarse degree specified through PETSc options (1 by default).
- Parameters:
ele – a
ufl.FiniteElement
to coarsen.
firedrake.preconditioners.massinv module¶
- class firedrake.preconditioners.massinv.MassInvPC[source]¶
Bases:
PCBase
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- apply(pc, X, Y)[source]¶
Apply the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- applyTranspose(pc, X, Y)¶
Apply the transpose of the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- needs_python_pmat = True¶
A matrix free operator that inverts the mass matrix in the provided space.
Internally this creates a PETSc KSP object that can be controlled by options using the extra options prefix
Mp_
.For Stokes problems, to be spectrally equivalent to the Schur complement, the mass matrix should be weighted by the viscosity. This can be provided (defaulting to constant viscosity) by providing a field defining the viscosity in the application context, keyed on
"mu"
.
firedrake.preconditioners.patch module¶
- class firedrake.preconditioners.patch.PatchPC[source]¶
Bases:
PCBase
,PatchBase
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- apply(pc, x, y)[source]¶
Apply the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
firedrake.preconditioners.pcd module¶
- class firedrake.preconditioners.pcd.PCDPC[source]¶
Bases:
PCBase
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- apply(pc, x, y)[source]¶
Apply the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- applyTranspose(pc, x, y)[source]¶
Apply the transpose of the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.
- needs_python_pmat = True¶
A Pressure-Convection-Diffusion preconditioner for Navier-Stokes.
This preconditioner approximates the inverse of the pressure schur complement for the Navier-Stokes equations by.
\[S^{-1} \sim K^{-1} F_p M^{-1}\]Where \(K = \nabla^2\), \(F_p = (1/\mathrm{Re}) \nabla^2 + u\cdot\nabla\) and \(M = \mathbb{I}\).
The inverse of \(K\) is approximated by a KSP which can be controlled using the options prefix
pcd_Kp_
.The inverse of \(M\) is similarly approximated by a KSP which can be controlled using the options prefix
pcd_Mp_
.\(F_p\) requires both the Reynolds number and the current velocity. You must provide these with options using the glbaol option
Re
for the Reynolds number and the prefixed optionpcd_velocity_space
which should be the index into the full space that gives velocity field.Note
Currently, the boundary conditions applied to the PCD operator are correct for characteristic velocity boundary conditions, but sub-optimal for in and outflow boundaries.
firedrake.preconditioners.pmg module¶
- class firedrake.preconditioners.pmg.PMGPC[source]¶
Bases:
PCBase
,PMGBase
Create a PC context suitable for PETSc.
Matrix free preconditioners should inherit from this class and implement:
- apply(pc, x, y)[source]¶
Apply the preconditioner to X, putting the result in Y.
Both X and Y are PETSc Vecs, Y is not guaranteed to be zero on entry.