tlm_adjoint.firedrake.block_system

Firedrake specific extensions to tlm_adjoint.block_system.

Module Contents

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

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:
alphascalar

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)

Nullspace and left nullspace spanned 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:
spacefiredrake.functionspaceimpl.WithGeometry

A scalar-valued function space containing the unity-valued function.

alphascalar

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)

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:
bcsfiredrake.bcs.DirichletBC or Sequence[firedrake.bcs.DirichletBC, …]

Homogeneous Dirichlet boundary conditions.

alphascalar

Defines the linear constraint matrix \(S = \alpha M\).

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.PETScMatrix(A)

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

Parameters:
Afiredrake.matrix.Matrix
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:
aufl.Form

Defines the sesquilinear form.

args, kwargs

Passed to the firedrake.assemble.assemble() function.

Returns:
PETScMatrix.

PETScMatrix defined by the assembled sesquilinear form.

class tlm_adjoint.firedrake.block_system.BlockMatrix(arg_space, action_space, blocks=None)

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

Parameters:
  • arg_space – Defines the space \(V\).

  • action_space – 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.block_system.Matrix, or None to indicate a zero block.

class tlm_adjoint.firedrake.block_system.LinearSolver(A, *args, **kwargs)

Solver for a linear system

\[A u = b,\]

using PETSc.

Parameters:
  • A – A tlm_adjoint.block_system.Matrix defining \(A\).

  • nullspace – A Nullspace or a Sequence of Nullspace objects defining the nullspace and left nullspace of \(A\). None indicates a NoneNullspace.

  • solver_parameters – A Mapping defining Krylov solver parameters.

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

  • comm – Communicator.

class tlm_adjoint.firedrake.block_system.WhiteNoiseSampler(space, rng, *, precondition=True, M=None, mfn_solver_parameters=None, ksp_solver_parameters=None)

White noise sampling.

Utility class for drawing independent spatial white noise samples. Generates a sample using

\[X = M^{-1} \Xi^{-T} \sqrt{ \Xi^T M \Xi } Z,\]

where

  • \(M\) is the mass matrix.

  • \(\Xi\) is a preconditioner.

  • \(Z\) is a vector whose elements are independent standard Gaussian samples.

The matrix square root is computed using SLEPc.

Parameters:
spacefiredrake.functionspaceimpl.WithGeometry

The function space.

rngnumpy.random.Generator

Pseudorandom number generator.

preconditionbool

If True then \(\Xi\) is set equal to the inverse of the (principal) square root of the diagonal of \(M\). Otherwise it is set equal to the identity.

Mfiredrake.matrix.Matrix

Mass matrix. Constructed by finite element assembly if not supplied.

mfn_solver_parametersMapping

slepc4py.SLEPc.MFN solver parameters, used for the matrix square root action.

ksp_solver_parametersMapping

Solver parameters, used for \(M^{-1}\).

Attributes:
spacefiredrake.functionspaceimpl.WithGeometry

The function space.

rngnumpy.random.Generator

Pseudorandom number generator.

dual_sample()

Generate a new sample in the dual space.

The result is given by

\[X = \Xi^{-T} \sqrt{ \Xi^T M \Xi } Z.\]
Returns:
firedrake.cofunction.Cofunction

The sample.

sample()

Generate a new sample.

The result is given by

\[X = M^{-1} \Xi^{-T} \sqrt{ \Xi^T M \Xi } Z.\]
Returns:
firedrake.function.Function

The sample.