dolfin-adjoint API reference

The entire dolfin-adjoint interface should be imported with a single call:

from dolfin import *
from dolfin_adjoint import *

It is essential that the importing of the dolfin_adjoint module happen after importing the dolfin module. dolfin-adjoint relies on overloading many of the key functions of dolfin to achieve its degree of automation.

Overloaded functions

dolfin_adjoint.assemble(*args, **kwargs)

When a form is assembled, the information about its nonlinear dependencies is lost, and it is no longer easy to manipulate. Therefore, dolfin_adjoint overloads the dolfin.assemble function to attach the form to the assembled object. This lets the automatic annotation work, even when the user calls the lower-level solve(A, x, b).

dolfin_adjoint.assemble_system(*args, **kwargs)

When a form is assembled, the information about its nonlinear dependencies is lost, and it is no longer easy to manipulate. Therefore, dolfin_adjoint overloads the dolfin.assemble_system function to attach the form to the assembled object. This lets the automatic annotation work, even when the user calls the lower-level solve(A, x, b).

dolfin_adjoint.solve(*args, **kwargs)

This solve routine wraps the real Dolfin solve call. Its purpose is to annotate the model, recording what solves occur and what forms are involved, so that the adjoint and tangent linear models may be constructed automatically by libadjoint.

To disable the annotation, just pass annotate=False to this routine, and it acts exactly like the Dolfin solve call. This is useful in cases where the solve is known to be irrelevant or diagnostic for the purposes of the adjoint computation (such as projecting fields to other function spaces for the purposes of visualisation).

dolfin_adjoint.project(v, V=None, bcs=None, mesh=None, solver_type='lu', preconditioner_type='default', form_compiler_parameters=None, annotate=None, name=None)

The project call performs an equation solve, and so it too must be annotated so that the adjoint and tangent linear models may be constructed automatically by libadjoint.

To disable the annotation of this function, just pass annotate=False. This is useful in cases where the solve is known to be irrelevant or diagnostic for the purposes of the adjoint computation (such as projecting fields to other function spaces for the purposes of visualisation).

dolfin_adjoint.interpolate(v, V, annotate=None, name=None)

The interpolate call changes Function data, and so it too must be annotated so that the adjoint and tangent linear models may be constructed automatically by libadjoint.

To disable the annotation of this function, just pass annotate=False. This is useful in cases where the interpolation is known to be irrelevant or diagnostic for the purposes of the adjoint computation (such as interpolating fields to other function spaces for the purposes of visualisation).

Overloaded objects

class dolfin_adjoint.LUSolver(*args)

This object is overloaded so that solves using this class are automatically annotated, so that libadjoint can automatically derive the adjoint and tangent linear models.

solve(*args, **kwargs)

To disable the annotation, just pass annotate=False to this routine, and it acts exactly like the Dolfin solve call. This is useful in cases where the solve is known to be irrelevant or diagnostic for the purposes of the adjoint computation (such as projecting fields to other function spaces for the purposes of visualisation).

class dolfin_adjoint.NewtonSolver(*args, **kwargs)

This object is overloaded so that solves using this class are automatically annotated, so that libadjoint can automatically derive the adjoint and tangent linear models.

solve(*args, **kwargs)

To disable the annotation, just pass annotate=False to this routine, and it acts exactly like the Dolfin solve call. This is useful in cases where the solve is known to be irrelevant or diagnostic for the purposes of the adjoint computation (such as projecting fields to other function spaces for the purposes of visualisation).

class dolfin_adjoint.KrylovSolver(*args)

This object is overloaded so that solves using this class are automatically annotated, so that libadjoint can automatically derive the adjoint and tangent linear models.

solve(*args, **kwargs)

To disable the annotation, just pass annotate=False to this routine, and it acts exactly like the Dolfin solve call. This is useful in cases where the solve is known to be irrelevant or diagnostic for the purposes of the adjoint computation (such as projecting fields to other function spaces for the purposes of visualisation).

class dolfin_adjoint.NonlinearVariationalSolver(problem, *args, **kwargs)

This object is overloaded so that solves using this class are automatically annotated, so that libadjoint can automatically derive the adjoint and tangent linear models.

solve(annotate=None)

To disable the annotation, just pass annotate=False to this routine, and it acts exactly like the Dolfin solve call. This is useful in cases where the solve is known to be irrelevant or diagnostic for the purposes of the adjoint computation (such as projecting fields to other function spaces for the purposes of visualisation).

class dolfin_adjoint.NonlinearVariationalProblem(F, u, bcs=None, J=None, *args, **kwargs)

This object is overloaded so that solves using this class are automatically annotated, so that libadjoint can automatically derive the adjoint and tangent linear models.

class dolfin_adjoint.LinearVariationalSolver(problem, *args, **kwargs)

This object is overloaded so that solves using this class are automatically annotated, so that libadjoint can automatically derive the adjoint and tangent linear models.

solve(annotate=None)

To disable the annotation, just pass annotate=False to this routine, and it acts exactly like the Dolfin solve call. This is useful in cases where the solve is known to be irrelevant or diagnostic for the purposes of the adjoint computation (such as projecting fields to other function spaces for the purposes of visualisation).

class dolfin_adjoint.LinearVariationalProblem(a, L, u, bcs=None, *args, **kwargs)

This object is overloaded so that solves using this class are automatically annotated, so that libadjoint can automatically derive the adjoint and tangent linear models.

class dolfin_adjoint.Function(*args, **kwargs)

The Function class is overloaded so that you can give Functions names. For example,

u = Function(V, name="Velocity")

This allows you to refer to the Function by name throughout dolfin-adjoint, rather than needing to have the specific Function instance available.

For more details, see the dolfin-adjoint documentation.

assign(other, annotate=None, *args, **kwargs)

To disable the annotation, just pass annotate=False to this routine, and it acts exactly like the Dolfin assign call.

class dolfin_adjoint.Constant(value, cell=None, name=None)

The Constant class is overloaded so that you can give Constants names. For example,

nu = Constant(1.0e-4, name="Diffusivity")

This allows you to refer to the Constant by name throughout dolfin-adjoint, rather than needing to have the specific Constant instance available.

For more details, see the dolfin-adjoint documentation.

Driver functions

dolfin_adjoint.compute_gradient(J, param, forget=True, ignore=[], callback=<function <lambda>>, project=False)
dolfin_adjoint.compute_adjoint(functional, forget=True, ignore=[])
dolfin_adjoint.compute_tlm(parameter, forget=False)

Functional object

class dolfin_adjoint.Functional(timeform, verbose=False, name=None)

This class implements the libadjoint.Functional abstract base class for dolfin-adjoint. The core idea is that a functional is either

  1. an integral of a form over a certain time window, or
  2. a pointwise evaluation in time of a certain form, or
  3. a sum of terms like (a) and (b).

Some examples:

  • Integration over all time:

    J = Functional(inner(u, u)*dx*dt)
    
  • Integration over a certain time window:

    J = Functional(inner(u, u)*dx*dt[0:1])
    
  • Integration from a certain point until the end:

    J = Functional(inner(u, u)*dx*dt[0.5:])
    
  • Pointwise evaluation in time (does not need to line up with timesteps):

    J = Functional(inner(u, u)*dx*dt[0.5])
    
  • Pointwise evaluation at the start (e.g. for regularisation terms):

    J = Functional(inner(u, u)*dx*dt[START_TIME])
    
  • Pointwise evaluation at the end:

    J = Functional(inner(u, u)*dx*dt[FINISH_TIME])
    
  • And sums of these work too:

    J = Functional(inner(u, u)*dx*dt + inner(u, u)*dx*dt[FINISH_TIME])
    

If dt has been redefined, you can create your own time measure with TimeMeasure.

For anything but the evaluation at the final time to work, you need to annotate the timestepping of your model with adj_inc_timestep.

ReducedFunctional object

class dolfin_adjoint.ReducedFunctional(functional, controls, scale=1.0, eval_cb_pre=<function <lambda>>, eval_cb_post=<function <lambda>>, derivative_cb_pre=<function <lambda>>, derivative_cb_post=<function <lambda>>, replay_cb=<function <lambda>>, hessian_cb=<function <lambda>>, cache=None)

This class provides access to the reduced functional for given functional and controls. The reduced functional maps a point in control space to the associated functional value by implicitly solving the PDE that is annotated by dolfin-adjoint. The ReducedFunctional object can also compute functional derivatives with respect to the controls using the adjoint method.

__call__(**kwargs)
derivative(forget=True, project=False)
Evaluates the derivative of the reduced functional at the most
recently evaluated control value.
Args:
forget (Optional[bool]): Delete the forward state while solving the
adjoint equations. If you want to reevaluate derivative at the same point (or the Hessian) you will need to set this to False or None. Defaults to True.
project (Optional[bool]): If True, the returned value will be the L2
Riesz representer, if False it will be the l2 Riesz representative. The L2 projection requires one additional linear solve. Defaults to False.
Returns:
The functional derivative. The returned type is the same as the control type.
hessian(m_dot, project=False)

Evaluates the Hessian action at the most recently evaluated control value in direction m_dot.

Args:
m_dot: The direction in control space in which to compute the
Hessian. Must be of the same type as the Control (e.g. Function, Constant or lists of latter).
project (Optional[bool]): If True, the returned value will be the L2
Riesz representer, if False it will be the l2 Riesz representative. The L2 projection requires one additional linear solve. Defaults to False.
Returns:
The directional second derivative. The returned type is the same as the control type.

Note: Hessian evaluations never delete the forward state.

taylor_test(value=None, test_hessian=False, seed=None, perturbation_direction=None, size=None)

Run a Taylor test to check that the functional, gradient and (optionally) Hessian are consistent by running the Taylor test.

Args:
value (Optional): The point in control space where to perform the Taylor test.
Must be of the same type as the Control (e.g. Function, Constant or lists of latter). If value is None (default), the Taylor test will be performed at the control value of the latest evaluation.
test_hessian (Optional[boolean]): If True, the Taylor test also
includes the Hessian. Defaults to False.
seed (Optional[float]): The initial perturbation size for the Taylor
test.
perturbation_direction (Optional): The direction in which to perform
the Taylor test. Must be of the same type as the Control (e.g. Function, Constant or lists of latter). Defaults to a random direction.

size (Optional[int]): Number of perturbations for the test

Returns:
float: The minimum (higher-order) convergence rate of all performed tests.

The Taylor test also prints out detailed information about the convergence rate if the fenics.log_level is set INFO or higher.

class dolfin_adjoint.ReducedFunctionalNumPy(rf)

This class implements the reduced functional for given functional and controls based on numpy data structures.

This “NumPy version” of the dolfin_adjoint.ReducedFunctional is created from an existing ReducedFunctional object: rf_np = ReducedFunctionalNumPy(rf = rf)

__call__(m_array)

An implementation of the reduced functional evaluation that accepts the control values as an array of scalars

derivative(m_array=None, forget=True, project=False)

An implementation of the reduced functional derivative evaluation that accepts the controls as an array of scalars. If no control values are given, the result is derivative at the lastest forward run.

hessian(m_array, m_dot_array)

An implementation of the reduced functional hessian action evaluation that accepts the controls as an array of scalars. If m_array is None, the Hessian action at the latest forward run is returned.

pyopt_problem(constraints=None, bounds=None, name='Problem', ignore_model_errors=False)

Return a pyopt problem class that can be used with the PyOpt package, http://www.pyopt.org/

set_controls(array)
get_controls()

Control objects

dolfin_adjoint.Control(obj, *args, **kwargs)

Creates a dolfin-adjoint control.

class dolfin_adjoint.FunctionControl(coeff, value=None, perturbation=None)

This Parameter is used as input to the tangent linear model (TLM) when one wishes to compute dJ/d(initial condition) in a particular direction (perturbation).

class dolfin_adjoint.ConstantControl(a, coeff=1)

This Parameter is used as input to the tangent linear model (TLM) when one wishes to compute dJ/da, where a is a single scalar parameter.

Constraint objects

class dolfin_adjoint.EqualityConstraint

This class represents equality constraints of the form

c_i(m) == 0

for 0 <= i < n, where m is the parameter.

function(m)

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.

jacobian(m)

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.

class dolfin_adjoint.InequalityConstraint

This class represents constraints of the form

c_i(m) >= 0

for 0 <= i < n, where m is the parameter.

function(m)

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.

jacobian(m)

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.

Annotation functions

dolfin_adjoint.adj_checkpointing(strategy, steps, snaps_on_disk, snaps_in_ram, verbose=False, replay=False, replay_comparison_tolerance=1e-10)
dolfin_adjoint.adj_start_timestep(time=0.0)

Dolfin does not supply us with information about timesteps, and so more information is required from the user for certain features. This function should be called at the start of the time loop with the initial time (defaults to 0).

See also: dolfin_adjoint.adj_inc_timestep

dolfin_adjoint.adj_inc_timestep(time=None, finished=False)

Dolfin does not supply us with information about timesteps, and so more information is required from the user for certain features. This function should be called at the end of the time loop with two arguments:

  • time – the time at the end of the timestep just computed
  • finished – whether this is the final timestep.

With this information, complex functional expressions using the Functional class can be used.

The finished argument is necessary because the final step of a functional integration must perform additional calculations.

See also: dolfin_adjoint.adj_start_timestep

Debugging functions

dolfin_adjoint.adj_html(*args, **kwargs)

This routine dumps the current state of the adjglobals.adjointer to a HTML visualisation. Use it like:

  • adj_html(“forward.html”, “forward”) # for the equations recorded on the forward run
  • adj_html(“adjoint.html”, “adjoint”) # for the equations to be assembled on the adjoint run
dolfin_adjoint.adj_check_checkpoints()
dolfin_adjoint.taylor_test(*args, **kwargs)
dolfin_adjoint.replay_dolfin(forget=False, tol=0.0, stop=False)

Generalised stability theory

dolfin_adjoint.compute_gst(ic, final, nsv, ic_norm='mass', final_norm='mass', which=1)

This function computes the generalised stability analysis of a simulation. Generalised stability theory computes the perturbations to a field (such as an initial condition, forcing term, etc.) that /grow the most/ over the finite time window of the simulation. For more details, see the mathematical documentation on the website.

  • ic – the input of the propagator
  • final – the output of the propagator
  • nsv – the number of optimal perturbations to compute
  • ic_norm – a symmetric positive-definite bilinear form that defines the norm on the input space
  • final_norm – a symmetric positive-definite bilinear form that defines the norm on the output space
  • which – which singular vectors to compute. Use e.g. slepc4py.SLEPc.EPS.Which.LARGEST_REAL

You can supply "mass" for ic_norm and final_norm to use the (default) mass matrices associated with these spaces.

For example:

gst = compute_gst("State", "State", nsv=10)
for i in range(gst.ncv): # number of converged vectors
  (sigma, u, v) = gst.get_gst(i, return_vectors=True)

Accessing tape

class dolfin_adjoint.DolfinAdjointVariable(coefficient, timestep=None, iteration=None)

A wrapper class for Dolfin objects to store additional information such as a time step, a iteration counter and the type of the variable (adjoint, forward or tangent linear).

__init__(coefficient, timestep=None, iteration=None)

Create a DolfinAdjointVariable associated with the provided coefficient.

If the coefficient is not known to dolfin_adjoint (i.e. if no equation for it was annotated), an Exception is thrown.

By default, the DolfinAdjointVariable references the latest timestep and iteration number, but may be overwritten with the timestep and the iteration parameters. Negative values may be used to reference the backwards.

tape_value(timestep=None, iteration=None)

Return the tape value associated with the variable (optionally for the given timestep and iteration).

iteration_count()

Return the annotated number of iterations at the variables timestep.

known_timesteps()

Return a list of timesteps for which this variable is annotated on the tape.

dolfin_adjoint.adj_reset()

Forget all annotation, and reset the entire dolfin-adjoint state.