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:

  • initialize()

  • update()

  • apply()

  • applyTranspose()

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:

  • initialize()

  • update()

  • apply()

  • applyTranspose()

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).

initialize(pc)[source]

Initialize any state in this preconditioner.

update(pc)[source]

Update any state in this preconditioner.

view(pc, viewer=None)[source]
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:

  • initialize()

  • update()

  • apply()

  • applyTranspose()

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:

  • initialize()

  • update()

  • apply()

  • applyTranspose()

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.

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.

form(pc, test, trial)[source]
initialize(pc)[source]

Initialize any state in this preconditioner.

update(pc)[source]

Update any state in this preconditioner.

view(pc, viewer=None)[source]
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:

  • initialize()

  • update()

  • apply()

  • applyTranspose()

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() and initialize().

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 get_appctx(pc)[source]
abstract initialize(pc)[source]

Initialize any state in this preconditioner.

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() and initialize().

abstract update(pc)[source]

Update any state in this preconditioner.

view(pc, viewer=None)[source]
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:

  • initialize()

  • update()

  • apply()

  • applyTranspose()

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^1 Riesz map is diagonalized 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 linear functions (tensor contractions). 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:

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_coef(J, quad_deg, discard_mixed=True, cell_average=True)[source]

Return the coefficients of the Jacobian form arguments and their gradient with respect to the reference coordinates.

Parameters:
  • J – the Jacobian bilinear form

  • quad_deg – the quadrature degree used for the coefficients

  • discard_mixed – discard entries in second order coefficient with mixed derivatives and mixed components

  • cell_average – to return the coefficients as DG_0 Functions

Returns:

a 2-tuple of coefficients: a dictionary mapping strings to firedrake.Functions with the coefficients of the form, assembly_callables: a list of assembly callables for each coefficient of the form

assemble_fdm_op(V, J, bcs, appctx)[source]

Assemble the sparse preconditioner with cell-wise constant coefficients.

Parameters:
  • V – the firedrake.FunctionSpace of the form arguments

  • J – the Jacobian bilinear form

  • bcs – an iterable of boundary conditions on V

  • appctx – the application context

Returns:

2-tuple with the preconditioner PETSc.Mat and its assembly callable

assemble_kron(A, V, bcs, eta, coefficients, Afdm, Dfdm, bcflags)[source]

Assemble the stiffness matrix in the FDM basis using Kronecker products of interval matrices

Parameters:
  • A – the PETSc.Mat to assemble

  • V – the firedrake.FunctionSpace of the form arguments

  • bcs – an iterable of firedrake.DirichletBCs

  • eta – a float penalty parameter for the symmetric interior penalty method

  • coefficients – a dict mapping strings to firedrake.Functions with the form coefficients

  • Afdm – the list with sparse interval matrices

  • Dfdm – the list with normal derivatives matrices

  • bcflags – the numpy.ndarray with BC facet flags returned by get_weak_bc_flags

initialize(pc)[source]

Initialize any state in this preconditioner.

update(pc)[source]

Update any state in this preconditioner.

view(pc, viewer=None)[source]

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.

initialize(pc)[source]

Initialize any state in this preconditioner.

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.

update(pc)[source]

Update any state in this preconditioner.

view(pc, viewer=None)[source]

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.

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.

initialize(obj)[source]

Initialize any state in this preconditioner.

update(pc)[source]

Update any state in this preconditioner.

view(pc, viewer=None)[source]

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.

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.

initialize(obj)[source]

Initialize any state in this preconditioner.

update(pc)[source]

Update any state in this preconditioner.

view(pc, viewer=None)[source]

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:

  • initialize()

  • update()

  • apply()

  • applyTranspose()

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:

  • initialize()

  • update()

  • apply()

  • applyTranspose()

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.

initialize(pc)[source]

Initialize any state in this preconditioner.

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".

update(pc)[source]

Update any state in this preconditioner.

view(pc, viewer=None)[source]

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.

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.

configure_patch(patch, pc)[source]
class firedrake.preconditioners.patch.PatchSNES[source]

Bases: SNESBase, PatchBase

Create a PC context suitable for PETSc.

Matrix free preconditioners should inherit from this class and implement:

  • initialize()

  • update()

  • apply()

  • applyTranspose()

configure_patch(patch, snes)[source]
step(snes, x, f, y)[source]
class firedrake.preconditioners.patch.PlaneSmoother[source]

Bases: object

static coords(dm, p, coordinates)[source]
sort_entities(dm, axis, dir, ndiv=None, divisions=None)[source]

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.

initialize(pc)[source]

Initialize any state in this preconditioner.

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\), \(M = \mathbb{I}\) and \(F_p = 1/\mathrm{Re} \\nabla^2 + u\cdot\\nabla\).

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 option pcd_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.

update(pc)[source]

Update any state in this preconditioner.

view(pc, viewer=None)[source]

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.

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.

coarsen_bc_value(bc, cV)[source]
configure_pmg(pc, pdm)[source]
class firedrake.preconditioners.pmg.PMGSNES[source]

Bases: SNESBase, PMGBase

Create a PC context suitable for PETSc.

Matrix free preconditioners should inherit from this class and implement:

  • initialize()

  • update()

  • apply()

  • applyTranspose()

coarsen_bc_value(bc, cV)[source]
configure_pmg(snes, pdm)[source]
step(snes, x, f, y)[source]

Module contents