Warning
You are reading a version of the website built against the unstable main
branch. This content is liable to change without notice and may be inappropriate for your use case.
You can find the documentation for the current stable release here.
firedrake.adjoint package¶
Submodules¶
firedrake.adjoint.ensemble_reduced_functional module¶
- class firedrake.adjoint.ensemble_reduced_functional.EnsembleReducedFunctional(functional, control, ensemble, scatter_control=True, gather_functional=None, derivative_components=None, scale=1.0, tape=None, eval_cb_pre=<function EnsembleReducedFunctional.<lambda>>, eval_cb_post=<function EnsembleReducedFunctional.<lambda>>, derivative_cb_pre=<function EnsembleReducedFunctional.<lambda>>, derivative_cb_post=<function EnsembleReducedFunctional.<lambda>>, hessian_cb_pre=<function EnsembleReducedFunctional.<lambda>>, hessian_cb_post=<function EnsembleReducedFunctional.<lambda>>)[source]¶
Bases:
AbstractReducedFunctional
Enable solving simultaneously reduced functionals in parallel.
Consider a functional \(J\) and its gradient \(\dfrac{dJ}{dm}\), where \(m\) is the control parameter. Let us assume that \(J\) is the sum of \(N\) functionals \(J_i(m)\), i.e.,
\[J = \sum_{i=1}^{N} J_i(m).\]The gradient over a summation is a linear operation. Therefore, we can write the gradient \(\dfrac{dJ}{dm}\) as
\[\frac{dJ}{dm} = \sum_{i=1}^{N} \frac{dJ_i}{dm},\]The
EnsembleReducedFunctional
allows simultaneous evaluation of \(J_i\) and \(\dfrac{dJ_i}{dm}\). After that, the allreduceEnsemble
operation is employed to sum the functionals and their gradients over an ensemble communicator.If gather_functional is present, then all the values of J are communicated to all ensemble ranks, and passed in a list to gather_functional, which is a reduced functional that expects a list of that size of the relevant types.
- Parameters:
functional (pyadjoint.OverloadedType) – An instance of an OverloadedType, usually
pyadjoint.AdjFloat
. This should be the functional that we want to reduce.control (pyadjoint.Control or list of pyadjoint.Control) – A single or a list of Control instances, which you want to map to the functional.
ensemble (Ensemble) – An instance of the
Ensemble
. It is used to communicate the functionals and their derivatives between the ensemble members.scatter_control (bool) – Whether scattering a control (or a list of controls) over the ensemble communicator
Ensemble.ensemble comm
.gather_functional (An instance of the
pyadjoint.ReducedFunctional
.) – that takes in all of the Js.derivative_components (list of int) – The indices of the controls that the derivative should be computed with respect to. If present, it overwrites
derivative_cb_pre
andderivative_cb_post
.scale (float) – A scaling factor applied to the functional and its gradient(with respect to the control).
tape (pyadjoint.Tape) – A tape object that the reduced functional will use to evaluate the functional and its gradients (or derivatives).
eval_cb_pre – Callback function before evaluating the functional. Input is a list of Controls.
eval_cb_pos – Callback function after evaluating the functional. Inputs are the functional value and a list of Controls.
derivative_cb_pre – Callback function before evaluating gradients (or derivatives). Input is a list of gradients (or derivatives). Should return a list of Controls (usually the same list as the input) to be passed to
pyadjoint.compute_gradient()
.derivative_cb_post – Callback function after evaluating derivatives. Inputs are the functional, a list of gradients (or derivatives), and controls. All of them are the checkpointed versions. Should return a list of gradients (or derivatives) (usually the same list as the input) to be returned from
self.derivative
.hessian_cb_pre – Callback function before evaluating the Hessian. Input is a list of Controls.
hessian_cb_post – Callback function after evaluating the Hessian. Inputs are the functional, a list of Hessian, and controls.
See also
Notes
The functionals \(J_i\) and the control must be defined over a common ensemble.comm communicator. To understand more about how ensemble parallelism works, please refer to the Firedrake manual.
- __call__(values)[source]¶
Computes the reduced functional with supplied control value.
- Parameters:
values (pyadjoint.OverloadedType) – If you have multiple controls this should be a list of new values for each control in the order you listed the controls to the constructor. If you have a single control it can either be a list or a single object. Each new value should have the same type as the corresponding control.
- Returns:
The computed value. Typically of instance of
pyadjoint.AdjFloat
.- Return type:
- property controls¶
Return the list of controls on which this functional depends.
- derivative(adj_input=1.0, apply_riesz=False)[source]¶
Compute derivatives of a functional with respect to the control parameters.
- Parameters:
- Returns:
dJdm_total (pyadjoint.OverloadedType)
The result of Allreduce operations of
dJdm_local
intodJdm_total
over the`Ensemble.ensemble_comm`.
- hessian(m_dot, hessian_input=None, evaluate_tlm=True, apply_riesz=False)[source]¶
The Hessian is not yet implemented for ensemble reduced functional.
- Raises:
NotImplementedError – This method is not yet implemented for ensemble reduced functional.
- tlm(m_dot)[source]¶
Return the action of the tangent linear model of the functional.
The tangent linear model is evaluated w.r.t. the control on a vector m_dot, around the last supplied value of the control.
- Parameters:
m_dot (pyadjoint.OverloadedType) – The direction in which to compute the action of the tangent linear model.
- Returns:
pyadjoint.OverloadedType (The action of the tangent linear model in the)
direction m_dot. Should be an instance of the same type as the functional.
firedrake.adjoint.ufl_constraints module¶
- class firedrake.adjoint.ufl_constraints.UFLConstraint(form, control)[source]¶
Bases:
Constraint
Easily implement scalar constraints using UFL.
The form must be a 0-form that depends on a Function control.
- function(m)[source]¶
Evaluate c(m), where c(m) == 0 for equality constraints and c(m) >= 0 for inequality constraints.
c(m) must return a numpy array or a dolfin Function or Constant.
- hessian_action(m, dm, dp, result)[source]¶
Computes the Hessian action of c(m) in direction dm and dp.
Stores the result in result.
- jacobian(m)[source]¶
Returns the full Jacobian matrix as a list of vector-like objects representing the gradient of the constraint function with respect to the parameter m.
The objects returned must be of the same type as m’s data.
- jacobian_action(m, dm, result)[source]¶
Computes the Jacobian action of c(m) in direction dm.
Stores the result in result.
- class firedrake.adjoint.ufl_constraints.UFLEqualityConstraint(form, control)[source]¶
Bases:
UFLConstraint
,EqualityConstraint
- class firedrake.adjoint.ufl_constraints.UFLInequalityConstraint(form, control)[source]¶
Bases:
UFLConstraint
,InequalityConstraint
Module contents¶
The public interface to Firedrake’s adjoint.
To start taping, run:
from firedrake.adjoint import *
continue_annotation()