tlm_adjoint.firedrake.solve
Finite element variational problem solution operations with Firedrake.
Module Contents
- class tlm_adjoint.firedrake.solve.EquationSolver(eq, x, bcs=None, *, J=None, form_compiler_parameters=None, solver_parameters=None, adjoint_solver_parameters=None, tlm_solver_parameters=None, cache_jacobian=None, cache_adjoint_jacobian=None, cache_tlm_jacobian=None, cache_rhs_assembly=None, match_quadrature=None)
Represents the solution of a finite element variational problem.
Caching is based on the approach described in
J. R. Maddison and P. E. Farrell, ‘Rapid development and adjoining of transient finite element models’, Computer Methods in Applied Mechanics and Engineering, 276, 95–121, 2014, doi: 10.1016/j.cma.2014.03.010
The arguments eq, x, bcs, J, form_compiler_parameters, and solver_parameters are based on the interface for the DOLFIN dolfin.solve function (see e.g. FEniCS 2017.1.0).
- Parameters:
eq – A
ufl.equation.Equationdefining the finite element variational problem.x – A
firedrake.function.Functiondefining the forward solution.bcs – Dirichlet boundary conditions.
J – A
ufl.Formdefining a Jacobian matrix approximation to use in a non-linear forward solve.form_compiler_parameters – Form compiler parameters.
solver_parameters – Linear or non-linear solver parameters.
adjoint_solver_parameters – Linear solver parameters to use in an adjoint solve.
tlm_solver_parameters – Linear solver parameters to use when solving tangent-linear problems.
cache_jacobian – Whether to cache the forward Jacobian matrix and linear solver data. Defaults to parameters[‘tlm_adjoint’][‘EquationSolver][‘cache_jacobian’]. If None then caching is autodetected.
cache_adjoint_jacobian – Whether to cache the adjoint Jacobian matrix and linear solver data. Defaults to cache_jacobian.
cache_tlm_jacobian – Whether to cache the Jacobian matrix and linear solver data when solving tangent-linear problems. Defaults to cache_jacobian.
cache_rhs_assembly – Whether to enable right-hand-side caching. If enabled then right-hand-side terms are divided into terms which are cached, terms which are converted into matrix multiplication by a cached matrix, and terms which are not cached. Defaults to parameters[‘tlm_adjoint’][‘EquationSolver’][‘cache_rhs_assembly’].
match_quadrature – Whether to set quadrature parameters consistently in the forward, adjoint, and tangent-linears. Defaults to parameters[‘tlm_adjoint’][‘EquationSolver’][‘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
EquationManageris paused.- Parameters:
X – A variable or a
Sequenceof variables storing the solution. May define an initial guess, and should be set by this method.deps – A
tupleof 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 or a
Sequenceof variables. Should not be modified.nl_deps – A
Sequenceof variables defining values for non-linear dependencies. Should not be modified.dep_Bs – A
Mappingwhose items are (dep_index, dep_B). Each dep_B is anAdjointRHSwhich 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 or
Sequenceof variables defining the initial guess for an iterative solve. May be modified or returned.nl_deps – A
Sequenceof variables defining values for non-linear dependencies. Should not be modified.B – The right-hand-side. A variable or
Sequenceof variables storing the value of the right-hand-side. May be modified or returned.
- Returns:
A variable or
Sequenceof variables storing the value of the adjoint solution. May return None to indicate a value of zero.
- tangent_linear(tlm_map)
Derive an
Equationcorresponding to a tangent-linear operation.- Parameters:
tlm_map – A
TangentLinearMapstoring values for tangent-linear variables.- Returns:
An
Equation, corresponding to the tangent-linear operation.
- class tlm_adjoint.firedrake.solve.LocalEquationSolver(eq, x, *, form_compiler_parameters=None, cache_jacobian=None, cache_adjoint_jacobian=None, cache_tlm_jacobian=None, cache_rhs_assembly=None, match_quadrature=None)
Represents the solution of a linear finite element variational problem, for the case where the matrix is element-wise local block diagonal.
Remaining arguments are passed to the
EquationSolverconstructor.- tangent_linear(tlm_map)
Derive an
Equationcorresponding to a tangent-linear operation.- Parameters:
tlm_map – A
TangentLinearMapstoring values for tangent-linear variables.- Returns:
An
Equation, corresponding to the tangent-linear operation.