tlm_adjoint.fenics.interpolation

Interpolation operations with FEniCS.

Module Contents

class tlm_adjoint.fenics.interpolation.ExprInterpolation(x, rhs)

Represents interpolation of rhs onto the space for x.

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

Parameters:
  • x – The forward solution.

  • rhs – A ufl.core.expr.Expr defining the expression to interpolate onto the space for x. Should not depend on x.

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.

adjoint_derivative_action(nl_deps, dep_index, adj_x)

Return the action of the adjoint of a derivative of the forward residual on the adjoint solution. This is the negative of an adjoint right-hand-side term.

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

  • dep_index – An int. The derivative is defined by differentiation of the forward residual with respect to self.dependencies()[dep_index].

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

Returns:

The action of the adjoint of a derivative on the adjoint solution. Will be passed to subtract_adjoint_derivative_action(), and valid types depend upon the adjoint variable type. Typically this will be a variable, or a two element tuple (alpha, F), where alpha is a numbers.Complex and F a variable, with the value defined by the product of alpha and F.

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.

class tlm_adjoint.fenics.interpolation.Interpolation(x, y, *, x_coords=None, P=None, tolerance=0.0)

Represents interpolation of the scalar-valued function y onto the space for x.

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

Internally this builds (or uses a supplied) interpolation matrix for the local process only. This behaves correctly if the there are no edges between owned and non-owned nodes in the degree of freedom graph associated with the discrete function space for y.

Parameters:
  • x – A scalar-valued DOLFIN Function defining the forward solution.

  • y – A scalar-valued DOLFIN Function to interpolate onto the space for x.

  • X_coords – A numpy.ndarray defining the coordinates at which to interpolate y. Shape is (n, d) where n is the number of process local degrees of freedom for x and d is the geometric dimension. Defaults to the process local degree of freedom locations for x. Ignored if P is supplied.

  • P – The interpolation matrix. A scipy.sparse.spmatrix.

  • tolerance – Maximum permitted distance (as returned by the DOLFIN BoundingBoxTree.compute_closest_entity method) of an interpolation point from a cell in the mesh for y. Ignored if P is supplied.

class tlm_adjoint.fenics.interpolation.PointInterpolation(X, y, X_coords=None, *, P=None, tolerance=0.0)

Represents interpolation of a scalar-valued function at given points.

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

Internally this builds (or uses a supplied) interpolation matrix for the local process only. This behaves correctly if the there are no edges between owned and non-owned nodes in the degree of freedom graph associated with the discrete function space for y.

Parameters:
  • X – A scalar variable, or a Sequence of scalar variables, defining the forward solution.

  • y – A scalar-valued DOLFIN Function to interpolate.

  • X_coords – A numpy.ndarray defining the coordinates at which to interpolate y. Shape is (n, d) where n is the number of interpolation points and d is the geometric dimension. Ignored if P is supplied.

  • P – The interpolation matrix. A scipy.sparse.spmatrix.

  • tolerance – Maximum permitted distance (as returned by the DOLFIN BoundingBoxTree.compute_closest_entity method) of an interpolation point from a cell in the mesh for y. Ignored if P is supplied.

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.

adjoint_derivative_action(nl_deps, dep_index, adj_X)

Return the action of the adjoint of a derivative of the forward residual on the adjoint solution. This is the negative of an adjoint right-hand-side term.

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

  • dep_index – An int. The derivative is defined by differentiation of the forward residual with respect to self.dependencies()[dep_index].

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

Returns:

The action of the adjoint of a derivative on the adjoint solution. Will be passed to subtract_adjoint_derivative_action(), and valid types depend upon the adjoint variable type. Typically this will be a variable, or a two element tuple (alpha, F), where alpha is a numbers.Complex and F a variable, with the value defined by the product of alpha and F.

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.