tlm_adjoint.fenics.assembly
Finite element assembly operations with FEniCS.
Module Contents
- class tlm_adjoint.fenics.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
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 or a
Sequence
of variables storing the solution. May define an initial guess, and should be set by this method.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 or a
Sequence
of variables. Should not be modified.
- 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 elementtuple
(alpha, F), where alpha is anumbers.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 or
Sequence
of variables defining the initial guess for an iterative solve. May be modified or returned.nl_deps – A
Sequence
of variables defining values for non-linear dependencies. Should not be modified.B – The right-hand-side. A variable or
Sequence
of variables storing the value of the right-hand-side. May be modified or returned.
- 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.