Skip to content

documentation: A proposal for a tutorial targeted at first time users #362

@reumle

Description

@reumle

Hello, i am a new user of DiffOpt, with also minimal experience of JuMP, and I struggled a bit at start with DiffOpt. I feel that a tutorial might be handy for new users. The provided examples (i read only the first ones) seem more targeted to the functionnal optimization problems, than showing how to use DiffOpt per se.

If it helps here is a minimal session, targeted at showin a new user how to operate DiffOpt. One session for forward AD, and the next for reverse AD. Please look at it as a first attempt, that may need some rewriting before being really satisfactory.

I ran it in my session and works fine (Julia 11.4). I got to this thanks to Claude.

A minimal DiffOpt session using forward AD

# DiffOpt.jl — Differentiating an optimization problem w.r.t. its parameters
# Example: min β·x² subject to x ≥ α
# Forward mode: one differentiation call per parameter, retrieves all outputs at once.
#
# Analytical reference:
#   x*(α, β) = α        (constraint is active at optimum)
#   f*(α, β) = β·α²
#   dx*/dα = 1,   dx*/dβ = 0
#   df*/dα = 2·β·α = 12.0,   df*/dβ = α² = 4.0

using JuMP, DiffOpt, Ipopt
import MathOptInterface as MOI

# ── Step 1: standard JuMP model (no DiffOpt) ────────────────────────────────
# Build and solve the optimization problem as usual first.
# This step is independent of DiffOpt — verify your model works before
# making it differentiable.

α_val = 2.0
β_val = 3.0

model_base = Model(Ipopt.Optimizer)
set_silent(model_base)
@variable(model_base, x)
@constraint(model_base, x >= α_val)
@objective(model_base, Min, β_val * x^2)
optimize!(model_base)
println("Step 1 — x* = ", value(x))               # → 2.0
println("Step 1 — f* = ", objective_value(model_base))  # → 12.0

# ── Step 2: differentiable model ────────────────────────────────────────────
# Rebuild the model with two changes:
#   - wrap the solver with DiffOpt.nonlinear_diff_model()
#   - declare parameters with Parameter(...) instead of plain numbers
#
# Use nonlinear_diff_model when parameters appear in the objective function.
# Use quadratic_diff_model only when parameters appear in constraints only.

model = DiffOpt.nonlinear_diff_model(Ipopt.Optimizer)
set_silent(model)
@variable(model, x)
@variable(model, α in Parameter(α_val))   # α is a parameter, not a decision variable
@variable(model, β in Parameter(β_val))   # β is a parameter, not a decision variable
@constraint(model, x >= α)
@objective(model, Min, β * x^2)
optimize!(model)
println("Step 2 — x* = ", value(x))               # → 2.0  (same result as Step 1)
println("Step 2 — f* = ", objective_value(model))  # → 12.0

# ── Step 3: derivatives in FORWARD mode ─────────────────────────────────────
# Forward mode: set one parameter direction (perturbation), call
# forward_differentiate!, then read dx*/dθ and df*/dθ.
# Repeat for each parameter of interest.
# Best suited when the number of parameters is small.

# --- derivatives w.r.t. α ---
DiffOpt.empty_input_sensitivities!(model)       # always clear before each call
DiffOpt.set_forward_parameter(model, α, 1.0)   # perturbation direction: dα = 1
DiffOpt.forward_differentiate!(model)

dx_dα = DiffOpt.get_forward_variable(model, x)    # dx*/dα
df_dα = DiffOpt.get_forward_objective(model)       # df*/dα
println("dx*/dα = ", dx_dα)    # → 1.0
println("df*/dα = ", df_dα)    # → 12.0

# --- derivatives w.r.t. β ---
DiffOpt.empty_input_sensitivities!(model)       # clear before switching parameter
DiffOpt.set_forward_parameter(model, β, 1.0)   # perturbation direction: dβ = 1
DiffOpt.forward_differentiate!(model)

dx_dβ = DiffOpt.get_forward_variable(model, x)    # dx*/dβ
df_dβ = DiffOpt.get_forward_objective(model)       # df*/dβ
println("dx*/dβ = ", dx_dβ)    # → 0.0
println("df*/dβ = ", df_dβ)    # → 4.0

Same session using reverse differentiation

Note: step 1 and 2 identical to the previous tutorial !

# DiffOpt.jl — Differentiating an optimization problem w.r.t. its parameters
# Example: min β·x² subject to x ≥ α
# Reverse mode: one differentiation call retrieves derivatives w.r.t. ALL parameters at once.
#
# Analytical reference:
#   x*(α, β) = α        (constraint is active at optimum)
#   f*(α, β) = β·α²
#   dx*/dα = 1,   dx*/dβ = 0
#   df*/dα = 2·β·α = 12.0,   df*/dβ = α² = 4.0

using JuMP, DiffOpt, Ipopt
import MathOptInterface as MOI

# ── Step 1: standard JuMP model (no DiffOpt) ────────────────────────────────
# Build and solve the optimization problem as usual first.
# This step is independent of DiffOpt — verify your model works before
# making it differentiable.

α_val = 2.0
β_val = 3.0

model_base = Model(Ipopt.Optimizer)
set_silent(model_base)
@variable(model_base, x)
@constraint(model_base, x >= α_val)
@objective(model_base, Min, β_val * x^2)
optimize!(model_base)
println("Step 1 — x* = ", value(x))                    # → 2.0
println("Step 1 — f* = ", objective_value(model_base))  # → 12.0

# ── Step 2: differentiable model ────────────────────────────────────────────
# Rebuild the model with two changes:
#   - wrap the solver with DiffOpt.nonlinear_diff_model()
#   - declare parameters with Parameter(...) instead of plain numbers
#
# Use nonlinear_diff_model when parameters appear in the objective function.
# Use quadratic_diff_model only when parameters appear in constraints only.

model = DiffOpt.nonlinear_diff_model(Ipopt.Optimizer)
set_silent(model)
@variable(model, x)
@variable(model, α in Parameter(α_val))   # α is a parameter, not a decision variable
@variable(model, β in Parameter(β_val))   # β is a parameter, not a decision variable
@constraint(model, x >= α)
@objective(model, Min, β * x^2)
optimize!(model)
println("Step 2 — x* = ", value(x))                # → 2.0  (same result as Step 1)
println("Step 2 — f* = ", objective_value(model))   # → 12.0

# ── Step 3: derivatives in REVERSE mode ─────────────────────────────────────
# Reverse mode: set the sensitivity of the quantity of interest (x* or f*),
# call reverse_differentiate!, then read dq/dθ for ALL parameters at once.
# Best suited when the number of parameters is large.

# --- derivatives of x* w.r.t. α and β ---
# We declare that the quantity of interest is x* itself (sensitivity = 1.0).
DiffOpt.empty_input_sensitivities!(model)        # always clear before each call
DiffOpt.set_reverse_variable(model, x, 1.0)     # "differentiate x* w.r.t. everything"
DiffOpt.reverse_differentiate!(model)

dx_dα = DiffOpt.get_reverse_parameter(model, α)  # dx*/dα
dx_dβ = DiffOpt.get_reverse_parameter(model, β)  # dx*/dβ
println("dx*/dα = ", dx_dα)    # → 1.0
println("dx*/dβ = ", dx_dβ)    # → 0.0

# --- derivatives of f* = β·x*² w.r.t. α and β ---
# We declare that the quantity of interest is f* itself (sensitivity = 1.0).
DiffOpt.empty_input_sensitivities!(model)        # clear before switching quantity
DiffOpt.set_reverse_objective(model, 1.0)        # "differentiate f* w.r.t. everything"
DiffOpt.reverse_differentiate!(model)

df_dα = DiffOpt.get_reverse_parameter(model, α)  # df*/dα
df_dβ = DiffOpt.get_reverse_parameter(model, β)  # df*/dβ
println("df*/dα = ", df_dα)    # → 12.0
println("df*/dβ = ", df_dβ)    # → 4.0

Metadata

Metadata

Assignees

No one assigned

    Labels

    No labels
    No labels

    Type

    No type

    Projects

    No projects

    Milestone

    No milestone

    Relationships

    None yet

    Development

    No branches or pull requests

    Issue actions