tlm_adjoint.firedrake.block_system

This module implements solvers for linear systems defined in mixed spaces.

The System class defines the block structure of the linear system, and solves the system using an outer Krylov solver. A custom preconditioner can be defined via the pc_fn callback to System.solve(), and this preconditioner can itself e.g. make use of further Krylov solvers. This provides a Python interface for custom block preconditioners.

Given a linear problem with a potentially singular matrix \(A\)

\[A u = b,\]

a System instead solves the linear problem

\[\left[ (I - M U (U^* M U)^{-1} U^*) A (I - V (V^* C V)^{-1} V^* C) + M U S V^* C \right] u = (I - M U (U^* M U)^{-1} U^*) b.\]

Here

  • \(U\) is a full rank matrix whose columns span the left nullspace for a modified system matrix \(\tilde{A}\).

  • \(V\) is a full rank matrix with the same number of columns as \(U\), whose columns span the nullspace for \(\tilde{A}\).

  • \(V^* C V\) and \(S\) are invertible matrices.

  • \(M\) is a Hermitian positive definite matrix.

Here the left nullspace for a matrix is defined to be the nullspace for its Hermitian transpose, and the modified system matrix \(\tilde{A}\) is defined

\[\tilde{A} = (I - M U (U^* M U)^{-1} U^*) A (I - V (V^* C V)^{-1} V^* C).\]

This has two primary use cases:

  1. Where a matrix \(A\) and right-hand-side \(b\) are constructed via finite element assembly on superspaces of the test space and trial space. The typical example is in the application of homogeneous essential Dirichlet boundary conditions.

  2. Where the matrix \(A\) is singular and \(b\) is orthogonal to the left nullspace of \(A\). Typically one would then choose \(U\) and \(V\) so that their columns respectively span the left nullspace and nullspace of \(A\), and the System then seeks a solution to the original problem subject to the linear constraints \(V^* C u = 0\).

Function spaces are defined via Firedrake function spaces, and Sequence objects containing Firedrake function spaces or similar Sequence objects. Similarly functions are defined via firedrake.function.Function or firedrake.cofunction.Cofunction objects, or Sequence objects containing firedrake.function.Function, firedrake.cofunction.Cofunction, or similar Sequence objects. This defines a basic tree structure which is useful e.g. when defining block matrices in terms of sub-block matrices.

Elements of the tree are accessed in a consistent order using a depth first search. Hence e.g.

((u_0, u_1), u_2)

and

(u_0, u_1, u_2)

where u_0, u_1, and u_2 are firedrake.function.Function or firedrake.cofunction.Cofunction objects, are both valid representations of a mixed space solution.

Module Contents

class tlm_adjoint.firedrake.block_system.MixedSpace(spaces)

Used to map between different versions of a mixed space.

This class defines two representations for the space:

  1. As a ‘split space’: A tree defining the mixed space. Stored using Firedrake function space and tuple objects, each corresponding to a node in the tree. Function spaces correspond to leaf nodes, and tuple objects to other nodes in the tree.

  2. As a ‘flattened space’: A Sequence containing leaf nodes of the split space with an ordering determined using a depth first search.

Provides methods to allow data to be copied to and from a compatible petsc4py.PETSc.Vec. This allows, for example, the construction:

u_0 = Function(space_0, name='u_0')
u_1 = Function(space_1, name='u_1')
u_2 = Function(space_2, name='u_2')

mixed_space = MixedSpace(((space_0, space_1), space_2))

and then data can be copied to a compatible petsc4py.PETSc.Vec via

mixed_space.to_petsc(u_petsc, ((u_0, u_1), u_2))

and from a compatible petsc4py.PETSc.Vec via

mixed_space.from_petsc(u_petsc, ((u_0, u_1), u_2))
Parameters:

spaces – The split space.

property comm

The communicator associated with the mixed space.

property split_space

The split space representation.

property flattened_space

The flattened space representation.

property local_size

The number of local degrees of freedom.

property global_size

The global number of degrees of freedom.

new_split()
Returns:

A new element in the split space.

from_petsc(u_petsc, u)

Copy data from a compatible petsc4py.PETSc.Vec.

Parameters:
to_petsc(u_petsc, u)

Copy data to a compatible petsc4py.PETSc.Vec. Does not update the ghost.

Parameters:
class tlm_adjoint.firedrake.block_system.Nullspace

Represents a matrix nullspace and left nullspace.

abstract apply_nullspace_transformation_lhs_right(x)

Apply the nullspace transformation associated with a matrix action on \(x\),

\[x \rightarrow (I - V (V^* C V)^{-1} V^* C) x.\]
Parameters:

x – Defines \(x\).

abstract apply_nullspace_transformation_lhs_left(y)

Apply the left nullspace transformation associated with a matrix action,

\[y \rightarrow (I - M U (U^* M U)^{-1} U^*) y.\]
Parameters:

y – Defines \(y\).

abstract constraint_correct_lhs(x, y)

Add the linear constraint term to \(y\),

\[y \rightarrow y + M U S V^* C x.\]
Parameters:
  • x – Defines \(x\).

  • y – Defines \(y\).

abstract pc_constraint_correct_soln(u, b)

Add the preconditioner linear constraint term to \(u\),

\[u \rightarrow u + V \tilde{S}^{-1} U^* b,\]

with

\[\tilde{S}^{-1} = \left( V^* C V \right)^{-1} S^{-1} \left( U^* M U \right)^{-1}.\]
Parameters:
  • u – Defines \(u\).

  • b – Defines \(b\).

correct_soln(x)

Correct the linear system solution so that it is orthogonal to space spanned by the columns of \(V\).

Parameters:

x – The linear system solution, to be corrected.

pre_mult_correct_lhs(x)

Apply the pre-left-multiplication nullspace transformation.

Parameters:

x – Defines the vector on which the matrix action is computed.

post_mult_correct_lhs(x, y)

Apply the post-left-multiplication nullspace transformation, and add the linear constraint term.

Parameters:
  • x – Defines the vector on which the matrix action is computed, and used to add the linear constraint term. If None is supplied then the linear constraint term is not added.

  • y – Defines the result of the matrix action on x.

correct_rhs(b)

Correct the linear system right-hand-side so that it is orthogonal to the space spanned by the columns of \(U\).

Parameters:

b – The linear system right-hand-side, to be corrected.

pc_pre_mult_correct(b)

Apply the pre-preconditioner-application nullspace transformation.

Parameters:

b – Defines the vector on which the preconditioner action is computed.

pc_post_mult_correct(u, b)

Apply the post-preconditioner-application left nullspace transformation, and add the linear constraint term.

Parameters:
  • u – Defines the result of the preconditioner action on b.

  • b – Defines the vector on which the preconditioner action is computed, and used to add the linear constraint term. If None is supplied then the linear constraint term is not added.

class tlm_adjoint.firedrake.block_system.NoneNullspace

An empty nullspace and left nullspace.

apply_nullspace_transformation_lhs_right(x)

Apply the nullspace transformation associated with a matrix action on \(x\),

\[x \rightarrow (I - V (V^* C V)^{-1} V^* C) x.\]
Parameters:

x – Defines \(x\).

apply_nullspace_transformation_lhs_left(y)

Apply the left nullspace transformation associated with a matrix action,

\[y \rightarrow (I - M U (U^* M U)^{-1} U^*) y.\]
Parameters:

y – Defines \(y\).

constraint_correct_lhs(x, y)

Add the linear constraint term to \(y\),

\[y \rightarrow y + M U S V^* C x.\]
Parameters:
  • x – Defines \(x\).

  • y – Defines \(y\).

pc_constraint_correct_soln(u, b)

Add the preconditioner linear constraint term to \(u\),

\[u \rightarrow u + V \tilde{S}^{-1} U^* b,\]

with

\[\tilde{S}^{-1} = \left( V^* C V \right)^{-1} S^{-1} \left( U^* M U \right)^{-1}.\]
Parameters:
  • u – Defines \(u\).

  • b – Defines \(b\).

class tlm_adjoint.firedrake.block_system.ConstantNullspace(*, alpha=1.0)

A nullspace and left nullspace spanned by the vector of ones.

Here \(V = U\), \(U\) is a single column matrix whose elements are ones, \(C = M\), and \(M\) is an identity matrix.

Parameters:

alpha – Defines the linear constraint matrix \(S = \left( \alpha / N \right)\) where \(N\) is the length of the vector of ones.

apply_nullspace_transformation_lhs_right(x)

Apply the nullspace transformation associated with a matrix action on \(x\),

\[x \rightarrow (I - V (V^* C V)^{-1} V^* C) x.\]
Parameters:

x – Defines \(x\).

apply_nullspace_transformation_lhs_left(y)

Apply the left nullspace transformation associated with a matrix action,

\[y \rightarrow (I - M U (U^* M U)^{-1} U^*) y.\]
Parameters:

y – Defines \(y\).

constraint_correct_lhs(x, y)

Add the linear constraint term to \(y\),

\[y \rightarrow y + M U S V^* C x.\]
Parameters:
  • x – Defines \(x\).

  • y – Defines \(y\).

pc_constraint_correct_soln(u, b)

Add the preconditioner linear constraint term to \(u\),

\[u \rightarrow u + V \tilde{S}^{-1} U^* b,\]

with

\[\tilde{S}^{-1} = \left( V^* C V \right)^{-1} S^{-1} \left( U^* M U \right)^{-1}.\]
Parameters:
  • u – Defines \(u\).

  • b – Defines \(b\).

class tlm_adjoint.firedrake.block_system.UnityNullspace(space, *, alpha=1.0)

A nullspace and left nullspace defined by the unity-valued function.

Here \(V = U\), \(U\) is a single column matrix containing the degree-of-freedom vector for the unity-valued function, \(C = M\), and \(M\) is the mass matrix.

Parameters:
  • space – A scalar-valued function space containing the unity-valued function.

  • alpha – Defines the linear constraint matrix \(S = \alpha \left( U^* M U \right)^{-1}\).

apply_nullspace_transformation_lhs_right(x)

Apply the nullspace transformation associated with a matrix action on \(x\),

\[x \rightarrow (I - V (V^* C V)^{-1} V^* C) x.\]
Parameters:

x – Defines \(x\).

apply_nullspace_transformation_lhs_left(y)

Apply the left nullspace transformation associated with a matrix action,

\[y \rightarrow (I - M U (U^* M U)^{-1} U^*) y.\]
Parameters:

y – Defines \(y\).

constraint_correct_lhs(x, y)

Add the linear constraint term to \(y\),

\[y \rightarrow y + M U S V^* C x.\]
Parameters:
  • x – Defines \(x\).

  • y – Defines \(y\).

pc_constraint_correct_soln(u, b)

Add the preconditioner linear constraint term to \(u\),

\[u \rightarrow u + V \tilde{S}^{-1} U^* b,\]

with

\[\tilde{S}^{-1} = \left( V^* C V \right)^{-1} S^{-1} \left( U^* M U \right)^{-1}.\]
Parameters:
  • u – Defines \(u\).

  • b – Defines \(b\).

class tlm_adjoint.firedrake.block_system.DirichletBCNullspace(bcs, *, alpha=1.0)

A nullspace and left nullspace associated with homogeneous Dirichlet boundary conditions.

Here \(V = U\), \(U\) is a zero-one matrix with exactly one non-zero per column corresponding to one boundary condition degree-of-freedom, \(C = M\), and \(M\) is an identity matrix.

Parameters:
apply_nullspace_transformation_lhs_right(x)

Apply the nullspace transformation associated with a matrix action on \(x\),

\[x \rightarrow (I - V (V^* C V)^{-1} V^* C) x.\]
Parameters:

x – Defines \(x\).

apply_nullspace_transformation_lhs_left(y)

Apply the left nullspace transformation associated with a matrix action,

\[y \rightarrow (I - M U (U^* M U)^{-1} U^*) y.\]
Parameters:

y – Defines \(y\).

constraint_correct_lhs(x, y)

Add the linear constraint term to \(y\),

\[y \rightarrow y + M U S V^* C x.\]
Parameters:
  • x – Defines \(x\).

  • y – Defines \(y\).

pc_constraint_correct_soln(u, b)

Add the preconditioner linear constraint term to \(u\),

\[u \rightarrow u + V \tilde{S}^{-1} U^* b,\]

with

\[\tilde{S}^{-1} = \left( V^* C V \right)^{-1} S^{-1} \left( U^* M U \right)^{-1}.\]
Parameters:
  • u – Defines \(u\).

  • b – Defines \(b\).

class tlm_adjoint.firedrake.block_system.BlockNullspace(nullspaces)

Nullspaces for a mixed space.

Parameters:

nullspaces – A Nullspace or a Sequence of Nullspace objects defining the nullspace. None indicates a NoneNullspace.

apply_nullspace_transformation_lhs_right(x)

Apply the nullspace transformation associated with a matrix action on \(x\),

\[x \rightarrow (I - V (V^* C V)^{-1} V^* C) x.\]
Parameters:

x – Defines \(x\).

apply_nullspace_transformation_lhs_left(y)

Apply the left nullspace transformation associated with a matrix action,

\[y \rightarrow (I - M U (U^* M U)^{-1} U^*) y.\]
Parameters:

y – Defines \(y\).

constraint_correct_lhs(x, y)

Add the linear constraint term to \(y\),

\[y \rightarrow y + M U S V^* C x.\]
Parameters:
  • x – Defines \(x\).

  • y – Defines \(y\).

pc_constraint_correct_soln(u, b)

Add the preconditioner linear constraint term to \(u\),

\[u \rightarrow u + V \tilde{S}^{-1} U^* b,\]

with

\[\tilde{S}^{-1} = \left( V^* C V \right)^{-1} S^{-1} \left( U^* M U \right)^{-1}.\]
Parameters:
  • u – Defines \(u\).

  • b – Defines \(b\).

class tlm_adjoint.firedrake.block_system.Matrix(arg_space, action_space)

Represents a matrix \(A\) mapping \(V \rightarrow W\).

Parameters:
  • arg_space – Defines the space V.

  • action_space – Defines the space W.

property arg_space

The space defining \(V\).

property action_space

The space defining \(W\).

abstract mult_add(x, y)

Add \(A x\) to \(y\).

Parameters:
  • x – Defines \(x\). Should not be modified.

  • y – Defines \(y\).

class tlm_adjoint.firedrake.block_system.PETScMatrix(arg_space, action_space, a)

A tlm_adjoint.firedrake.block_system.Matrix associated with a firedrake.matrix.Matrix \(A\) mapping \(V \rightarrow W\).

Parameters:
mult_add(x, y)

Add \(A x\) to \(y\).

Parameters:
  • x – Defines \(x\). Should not be modified.

  • y – Defines \(y\).

tlm_adjoint.firedrake.block_system.form_matrix(a, *args, **kwargs)

Construct a PETScMatrix associated with a given sesquilinear form.

Parameters:

a – A ufl.Form defining the sesquilinear form.

Returns:

The PETScMatrix.

Remaining arguments are passed to the firedrake.assemble.assemble() function.

class tlm_adjoint.firedrake.block_system.BlockMatrix(arg_spaces, action_spaces, blocks=None)

A matrix \(A\) mapping \(V \rightarrow W\), where \(V\) and \(W\) are defined by mixed spaces.

Parameters:
  • arg_spaces – Defines the space V.

  • action_spaces – Defines the space W.

  • block – A Mapping defining the blocks of the matrix. Items are ((i, j), block) where the block in the i th and j th column is defined by block. Each block is a tlm_adjoint.firedrake.block_system.Matrix or ufl.Form, or None to indicate a zero block.

mult_add(x, y)

Add \(A x\) to \(y\).

Parameters:
  • x – Defines \(x\). Should not be modified.

  • y – Defines \(y\).

class tlm_adjoint.firedrake.block_system.System(arg_spaces, action_spaces, blocks, *, nullspaces=None, comm=None)

A linear system

\[A u = b.\]
Parameters:
solve(u, b, *, solver_parameters=None, pc_fn=None, pre_callback=None, post_callback=None, correct_initial_guess=True, correct_solution=True)

Solve the linear system.

Parameters:
  • u – Defines the solution \(u\).

  • b – Defines the right-hand-side \(b\).

  • solver_parameters

    A Mapping defining outer Krylov solver parameters. Parameters (a number of which are based on FEniCS solver parameters) are:

    • ’linear_solver’: The Krylov solver type, default ‘fgmres’.

    • ’pc_side’: Overrides the PETSc default preconditioning side.

    • ’relative_tolerance’: Relative tolerance. Required.

    • ’absolute_tolerance’: Absolute tolerance. Required.

    • ’divergence_limit’: Overrides the default divergence limit.

    • ’maximum_iterations’: Maximum number of iterations. Default 1000.

    • ’norm_type’: Overrides the default convergence norm definition.

    • ’nonzero_initial_guess’: Whether to use a non-zero initial guess, defined by the input u. Default True.

    • ’gmres_restart’: Overrides the default GMRES restart parameter.

  • pc_fn

    Defines the application of a preconditioner. A callable

    def pc_fn(u, b):
    

    The preconditioner is applied to b, and the result stored in u. Defaults to an identity.

  • pre_callback – A callable accepting a single petsc4py.PETSc.KSP argument. Used for detailed manual configuration. Called after all other configuration options are set, but before the petsc4py.PETSc.KSP.setUp() method is called.

  • post_callback – A callable accepting a single petsc4py.PETSc.KSP argument. Called after the petsc4py.PETSc.KSP.solve() method has been called.

  • correct_initial_guess – Whether to apply a nullspace correction to the initial guess.

  • correct_solution – Whether to apply a nullspace correction to the solution.

Returns:

The number of Krylov iterations.