firedrake.preconditioners package

Submodules

firedrake.preconditioners.asm module

class firedrake.preconditioners.asm.ASMLinesmoothPC[source]

Bases: firedrake.preconditioners.asm.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: firedrake.preconditioners.base.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: firedrake.preconditioners.asm.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: firedrake.preconditioners.asm.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: firedrake.preconditioners.base.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: firedrake.preconditioners.assembled.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: firedrake.preconditioners.base.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: firedrake.preconditioners.base.PCSNESBase

Create a PC context suitable for PETSc.

Matrix free preconditioners should inherit from this class and implement:

  • initialize()

  • update()

  • apply()

  • applyTranspose()

firedrake.preconditioners.gtmg module

class firedrake.preconditioners.gtmg.GTMGPC[source]

Bases: firedrake.preconditioners.base.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: firedrake.preconditioners.base.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: firedrake.preconditioners.base.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: firedrake.preconditioners.pmg.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).

It raises a NotImplementedError for ufl.MixedElement`s, as I don't know if there's a sensible default strategy to implement here. It is intended that the user subclass `PMGPC to override this method for their problem.

Parameters

ele – a ufl.FiniteElement to coarsen.

firedrake.preconditioners.massinv module

class firedrake.preconditioners.massinv.MassInvPC[source]

Bases: firedrake.preconditioners.base.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 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: firedrake.preconditioners.base.PCBase, firedrake.preconditioners.patch.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: firedrake.preconditioners.base.SNESBase, firedrake.preconditioners.patch.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)[source]

firedrake.preconditioners.pcd module

class firedrake.preconditioners.pcd.PCDPC[source]

Bases: firedrake.preconditioners.base.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.MixedInterpolationMatrix(Vf, Vc, Vf_bcs, Vc_bcs)[source]

Bases: object

Interpolation matrix for a mixed finite element space.

mult(mat, resf, resc)[source]
multTranspose(mat, xc, xf, inc=False)[source]
multTransposeAdd(mat, x, y, w)[source]
class firedrake.preconditioners.pmg.PMGBase[source]

Bases: firedrake.preconditioners.base.PCSNESBase

A class for implementing p-multigrid Internally, this sets up a DM with a custom coarsen routine that p-coarsens the problem. This DM is passed to an internal PETSc PC of type MG and with options prefix pmg_. The relaxation to apply on every p-level is described by pmg_mg_levels_, and the coarse solve by pmg_mg_coarse_. Geometric multigrid or any other solver in firedrake may be applied to the coarse problem.

Other PETSc options inspected by this class in particular are: - ‘pmg_mg_coarse_degree’: to specify the degree of the coarse level - ‘pmg_mg_levels_transfer_mat_type’: can be either ‘aij’ or ‘matfree’

The p-coarsening is implemented in the coarsen_element routine. This takes in a ufl.FiniteElement and either returns a new, coarser element, or raises a ValueError (if the supplied element should be the coarsest one of the hierarchy).

The default coarsen_element is to perform power-of-2 reduction of the polynomial degree. For mixed systems a NotImplementedError is raised, as I don’t know how to make a sensible default for this. It is expected that many (most?) applications of this preconditioner will subclass PMGBase to override coarsen_element.

Create a PC context suitable for PETSc.

Matrix free preconditioners should inherit from this class and implement:

coarsen(fdm, comm)[source]
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).

It raises a NotImplementedError for ufl.MixedElement`s, as I don't know if there's a sensible default strategy to implement here. It is intended that the user subclass `PMGPC to override this method for their problem.

Parameters

ele – a ufl.FiniteElement to coarsen.

create_injection(dmc, dmf)[source]
create_interpolation(dmc, dmf)[source]
static create_transfer(cctx, fctx, mat_type, cbcs, fbcs)[source]
initialize(pc)[source]

Initialize any state in this preconditioner.

static reconstruct_degree(ele, N)[source]

Reconstruct a given element, modifying its polynomial degree.

Can only set a single degree along all axes of a TensorProductElement.

Parameters
update(pc)[source]

Update any state in this preconditioner.

view(pc, viewer=None)[source]
class firedrake.preconditioners.pmg.PMGPC[source]

Bases: firedrake.preconditioners.base.PCBase, firedrake.preconditioners.pmg.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: firedrake.preconditioners.base.SNESBase, firedrake.preconditioners.pmg.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]
class firedrake.preconditioners.pmg.StandaloneInterpolationMatrix(Vf, Vc, Vf_bcs, Vc_bcs)[source]

Bases: object

Interpolation matrix for a single standalone space.

static make_blas_kernels(Vf, Vc)[source]

Interpolation and restriction kernels between CG / DG tensor product spaces on quads and hexes.

Works by tabulating the coarse 1D Lagrange basis functions as the (Nf+1)-by-(Nc+1) matrix Jhat, and using the fact that the 2D / 3D tabulation is the tensor product J = kron(Jhat, kron(Jhat, Jhat))

make_kernels(Vf, Vc)[source]

Interpolation and restriction kernels between arbitrary elements.

This is temporary while we wait for dual evaluation in FInAT.

mult(mat, resf, resc)[source]

Implement restriction: restrict residual on fine grid resf to coarse grid resc.

multTranspose(mat, xc, xf, inc=False)[source]

Implement prolongation: prolong correction on coarse grid xc to fine grid xf.

multTransposeAdd(mat, x, y, w)[source]
static multiplicity(V)[source]
static prolongation_transfer_kernel_action(Vf, expr)[source]
firedrake.preconditioners.pmg.get_line_element(V)[source]
firedrake.preconditioners.pmg.get_line_nodes(V)[source]
firedrake.preconditioners.pmg.prolongation_matrix_aij(Pk, P1, Pk_bcs, P1_bcs)[source]
firedrake.preconditioners.pmg.prolongation_matrix_matfree(Vf, Vc, Vf_bcs, Vc_bcs)[source]
firedrake.preconditioners.pmg.prolongation_transfer_kernel_aij(Pk, P1)[source]
firedrake.preconditioners.pmg.tensor_product_space_query(V)[source]

Checks whether the custom transfer kernels support the FunctionSpace V.

V must be either CG(N) or DG(N) on quads or hexes (same N along every direction).

Parameters

V – FunctionSpace

Returns

4-tuple of (use_tensorproduct, degree, family, variant)

Module contents