tlm_adjoint.firedrake.assembly

Finite element assembly operations with Firedrake.

Module Contents

class tlm_adjoint.firedrake.assembly.Assembly(x, rhs, *, form_compiler_parameters=None, match_quadrature=None)

Represents assignment to the result of finite element assembly:

x = assemble(rhs)

The forward residual \(\mathcal{F}\) is defined so that \(\partial \mathcal{F} / \partial x\) is the identity.

Parameters:
  • x – A variable defining the forward solution.

  • rhs – A ufl.form.BaseForm` to assemble. Should have arity 0 or 1, and should not depend on x.

  • form_compiler_parameters – Form compiler parameters.

  • match_quadrature – Whether to set quadrature parameters consistently in the forward, adjoint, and tangent-linears. Defaults to parameters[‘tlm_adjoint’][‘Assembly’][‘match_quadrature’].

drop_references()

Drop references to variables which store values.

forward_solve(x, deps=None)

Compute the forward solution.

Can assume that the currently active EquationManager is paused.

Parameters:
  • X – A variable if the forward solution has a single component, otherwise a Sequence of variables. May define an initial guess, and should be set by this method. Subclasses may replace this argument with x if the forward solution has a single component.

  • deps – A tuple of variables, defining values for dependencies. Only the elements corresponding to X may be modified. self.dependencies() should be used if not supplied.

subtract_adjoint_derivative_actions(adj_x, nl_deps, dep_Bs)

Subtract terms from other adjoint right-hand-sides.

Can be overridden for an optimized implementation, but otherwise uses Equation.adjoint_derivative_action().

Parameters:
  • adj_X – The adjoint solution. A variable if the adjoint solution has a single component, otherwise a Sequence of variables. Should not be modified. Subclasses may replace this argument with adj_x if the adjoint solution has a single component.

  • nl_deps – A Sequence of variables defining values for non-linear dependencies. Should not be modified.

  • dep_Bs – A Mapping whose items are (dep_index, dep_B). Each dep_B is an AdjointRHS which should be updated by subtracting adjoint derivative information computed by differentiating with respect to self.dependencies()[dep_index].

adjoint_jacobian_solve(adj_x, nl_deps, b)

Compute an adjoint solution.

Parameters:
  • adj_X – Either None, or a variable (if the adjoint solution has a single component) or Sequence of variables (otherwise) defining the initial guess for an iterative solve. May be modified or returned. Subclasses may replace this argument with adj_x if the adjoint solution has a single component.

  • nl_deps – A Sequence of variables defining values for non-linear dependencies. Should not be modified.

  • B – The right-hand-side. A variable (if the adjoint solution has a single component) or Sequence of variables (otherwise) storing the value of the right-hand-side. May be modified or returned. Subclasses may replace this argument with b if the adjoint solution has a single component.

Returns:

A variable or Sequence of variables storing the value of the adjoint solution. May return None to indicate a value of zero.

tangent_linear(tlm_map)

Derive an Equation corresponding to a tangent-linear operation.

Parameters:

tlm_map – A TangentLinearMap storing values for tangent-linear variables.

Returns:

An Equation, corresponding to the tangent-linear operation.