time-integrated functional depending on constant-in-time control

Issue #75 new
Stephan Kramer created an issue

In the attached example (based on test_firedrake/identity/test_identity) I have a functional with a time integral that explicitly depends on the control, a function s. Because s doesn't actually get touched during the forward model this seems to confuse dolfin-adjoint with the following error:

LibadjointErrorInvalidInputs: Error in adj_iteration_count: No iteration found for supplied variable w_2:1:0:Forward

This can be worked around by adding a spurious projection to itself of s every time step (commented out in the example test). The same issue seems to be present in firedrake-adjoint and dolfin-adjoint.

Comments (1)

  1. Miguel Salazar

    The following code works for me:

    from dolfin import *
    from dolfin_adjoint import *
    
    mesh = UnitIntervalMesh(mpi_comm_world(), 2)
    
    W = FunctionSpace(mesh, "CG", 1)
    rho_const = interpolate(Constant(5.0), W, annotate=False)
    rho = interpolate(rho_const, W, name="Control")
    
    ##################### BRIEF EXPLANATION #############################
    # It seems that we need the state_1 and state_2 to make it work
    #
    
    g_const = interpolate(Constant(5.0), W)
    g = interpolate(g_const, W, name="Control2")
    
    
    u = Function(W, name="State")
    u1 = Function(W, name="State_1")
    u0 = Function(W, name="State_0")
    
    u_ = TrialFunction(W)
    v = TestFunction(W)
    
    F = rho*u_*v * dx - Constant(1)*g*v*dx
    adj_start_timestep(0.0)
    
    solve(lhs(F) == rhs(F), u)
    u0.assign(u1, annotate=True)
    u1.assign(u, annotate=True)
    rho.assign(rho_const, annotate=True)
    
    adj_inc_timestep(time=100.0, finished=False)
    
    solve(lhs(F) == rhs(F), u)
    u0.assign(u1, annotate=True)
    u1.assign(u, annotate=True)
    rho.assign(rho_const, annotate=True)
    
    adj_inc_timestep(time=200.0, finished=True)
    
    adj_html("forward.html", "forward")
    adj_html("adjoint.html", "adjoint")
    
    success = replay_dolfin(tol=0.0, stop=True)
    print success
    J = Functional((rho*0.5 * inner(u1, u1) * dx)*dt)
    
    # Reduced functional with single control
    m = FunctionControl(rho)
    
    Jhat = ReducedFunctional(J, m)
    Jhat.derivative()
    Jhat(rho)
    Jhat.hessian(rho)
    
    direction = interpolate(Constant(1), W)
    assert Jhat.taylor_test(rho, test_hessian=False, perturbation_direction=direction) > 1.9
    

    I had to include the functions u0 and u1 to pass the taylor test. I am not entirely sure why. I also wonder what's the cost of doing the assign for the control variable. Now you need to save also the design variable even though it is a constant.

  2. Log in to comment