Library

Constraints

DirectTrajOpt.Constraints.EqualityConstraintType
struct EqualityConstraint

Represents a linear equality constraint.

Fields

  • ts::AbstractArray{Int}: the time steps at which the constraint is applied
  • js::AbstractArray{Int}: the components of the trajectory at which the constraint is applied
  • vals::Vector{R}: the values of the constraint
  • vardim::Int: the dimension of a single time step of the trajectory
  • label::String: a label for the constraint
source
DirectTrajOpt.Constraints.EqualityConstraintMethod
EqualityConstraint(
    name::Symbol,
    ts::Vector{Int},
    val::Vector{Float64},
    traj::NamedTrajectory;
    label="equality constraint on trajectory variable [name]"
)

Constructs equality constraint for trajectory variable in NamedTrajectory

source
DirectTrajOpt.Constraints.GlobalEqualityConstraintMethod
GlobalEqualityConstraint(
    name::Symbol,
    val::Vector{Float64},
    traj::NamedTrajectory;
    label="equality constraint on global variable [name]"
)::EqualityConstraint

Constructs equality constraint for global variable in NamedTrajectory

source

Integrators

DirectTrajOpt.Integrators.BilinearIntegratorType
BilinearIntegrator <: AbstractBilinearIntegrator

Integrator for control-linear dynamics of the form ẋ = G(u)x.

This integrator uses matrix exponential methods to compute accurate state transitions for bilinear systems where the system matrix depends linearly on the control input.

Fields

  • G::Function: Function mapping control u to system matrix G(u)
  • x_comps::Vector{Int}: Component indices for state variables
  • u_comps::Vector{Int}: Component indices for control variables
  • Δt_comp::Int: Component index for time step
  • z_dim::Int: Total dimension of knot point
  • x_dim::Int: State dimension
  • u_dim::Int: Control dimension

Constructors

BilinearIntegrator(G::Function, traj::NamedTrajectory, x::Symbol, u::Symbol)
BilinearIntegrator(G::Function, traj::NamedTrajectory, xs::Vector{Symbol}, u::Symbol)

Arguments

  • G: Function taking control u and returning state matrix (xdim × xdim)
  • traj: NamedTrajectory containing the optimization variables
  • x or xs: State variable name(s)
  • u: Control variable name

Dynamics

Computes the constraint: x{k+1} - exp(Δt * G(uk)) * x_k = 0

Example

# Linear dynamics: ẋ = (A + Σᵢ uᵢ Bᵢ) x
A = [-0.1 1.0; -1.0 -0.1]
B = [0.0 0.0; 0.0 1.0]
G = u -> A + u[1] * B

integrator = BilinearIntegrator(G, traj, :x, :u)
source
DirectTrajOpt.Integrators.denseMethod
dense(vals, structure, shape)

Convert sparse data to dense matrix.

Arguments

  • vals: vector of values
  • structure: vector of tuples of indices
  • shape: tuple of matrix dimensions
source

Objectives

DirectTrajOpt.Objectives.ObjectiveType
Objective

A structure representing an objective function for trajectory optimization.

An objective function consists of the cost function itself along with its first and second derivatives for gradient-based optimization. The structure supports automatic differentiation and sparse Hessian representations.

Fields

  • L::Function: Objective function that takes trajectory vector Z⃗ and returns a scalar cost
  • ∇L::Function: Gradient function returning ∂L/∂Z⃗
  • ∂²L::Union{Function, Nothing}: Hessian function returning non-zero Hessian entries
  • ∂²L_structure::Union{Function, Nothing}: Function returning sparsity structure of Hessian

Operators

Objectives support addition and scalar multiplication:

  • obj1 + obj2: Combine objectives by summing costs and derivatives
  • α * obj: Scale objective by constant α

Example

# Create a regularization objective
obj = QuadraticRegularizer(:u, traj, 1e-2)

# Combine multiple objectives  
total_obj = obj1 + obj2 + 0.1 * obj3
source
DirectTrajOpt.Objectives.GlobalObjectiveMethod
GlobalObjective(
    ℓ::Function,
    global_names::AbstractVector{Symbol},
    traj::NamedTrajectory;
    kwargs...
)
GlobalObjective(
    ℓ::Function,
    global_name::Symbol,
    traj::NamedTrajectory;
    kwargs...
)

Create an objective that only involves the global components.

source
DirectTrajOpt.Objectives.KnotPointObjectiveMethod
KnotPointObjective(
    ℓ::Function,
    names::AbstractVector{Symbol},
    traj::NamedTrajectory,
    params::AbstractVector;
    kwargs...
)
KnotPointObjective(
    ℓ::Function,
    names::AbstractVector{Symbol},
    traj::NamedTrajectory;
    kwargs...
)
KnotPointObjective(
    ℓ::Function,
    name::Symbol,
    args...;
    kwargs...
)

Create a knot point summed objective function for trajectory optimization, where ℓ(x, p) on trajectory knot point variables x with parameters p. If the parameters argument is omitted, ℓ(x) is assumed to be a function of x only.

For multiple variables, the function can accept either:

  • Separate arguments: ℓ(x, u) for variables [:x, :u]
  • Concatenated vector: ℓ(xu) where xu = [x; u]

The constructor automatically detects which form expects.

Arguments

  • ℓ::Function: Function that defines the objective, ℓ(x, p) or ℓ(x).
    • For single variable: ℓ(x) where x is the variable values at a knot point
    • For multiple variables: ℓ(x, u) or ℓ(xu) depending on preference
  • names::AbstractVector{Symbol}: Names of the trajectory variables to be optimized.
  • traj::NamedTrajectory: The trajectory on which the objective is defined.
  • params::AbstractVector: Parameters p for the objective function ℓ, for each time.

Keyword Arguments

  • times::AbstractVector{Int}=1:traj.N: Time indices at which the objective is evaluated.
  • Qs::AbstractVector{Float64}=ones(traj.N): Weights for the objective function at each time.

Examples

# Single variable objective
obj = KnotPointObjective(
    x -> norm(x)^2,
    [:x], traj
)

# Multiple variables with separate arguments (recommended)
obj = KnotPointObjective(
    (x, u) -> x[1]^2 + u[1]^2,
    [:x, :u], traj
)

# Multiple variables with concatenated vector
obj = KnotPointObjective(
    xu -> xu[1]^2 + xu[3]^2,  # xu = [x[1], x[2], u[1]]
    [:x, :u], traj
)
source
DirectTrajOpt.Objectives.QuadraticRegularizerMethod
QuadraticRegularizer(
    name::Symbol,
    traj::NamedTrajectory,
    R::Union{Real, AbstractVector{<:Real}};
    baseline::AbstractMatrix{<:Real}=zeros(traj.dims[name], traj.N),
    times::AbstractVector{Int}=1:traj.N
)

Create a quadratic regularization objective for a trajectory component.

Minimizes the weighted squared deviation from a baseline trajectory, integrated over time:

\[J = \sum_{k \in \text{times}} \frac{1}{2} (v_k - v_\text{baseline})^T R (v_k - v_\text{baseline}) \Delta t\]

where v_k is the trajectory component at knot point k.

Arguments

  • name::Symbol: Name of the trajectory component to regularize
  • traj::NamedTrajectory: The trajectory containing the component
  • R: Regularization weight(s). Can be:
    • Scalar: same weight for all components
    • Vector: individual weights for each component dimension
  • baseline::AbstractMatrix: Target values (default: zeros). Size: (component_dim × N)
  • times::AbstractVector{Int}: Time indices to include in regularization (default: all)

Returns

  • Objective: Regularization objective with gradient and Hessian

Examples

# Regularize control with uniform weight
obj = QuadraticRegularizer(:u, traj, 1e-2)

# Regularize with different weights per component
obj = QuadraticRegularizer(:u, traj, [1e-2, 1e-3])

# Regularize around a reference trajectory
obj = QuadraticRegularizer(:x, traj, 1.0, baseline=x_ref)

# Only regularize middle time steps
obj = QuadraticRegularizer(:u, traj, 1e-2, times=2:traj.N-1)
source

Dynamics

DirectTrajOpt.Dynamics.TrajectoryDynamicsType
TrajectoryDynamics

Represents the dynamics of a trajectory optimization problem through integrators.

This structure encapsulates the dynamics constraints that enforce consistency between consecutive time steps in the trajectory. It uses integrators to define how the state evolves from one time step to the next and provides automatic differentiation support through Jacobian and Hessian computations.

Fields

  • F!::Function: In-place function computing dynamics violations δ = f(zₖ, zₖ₊₁)
  • ∂F!::Function: In-place function computing Jacobian of dynamics
  • ∂fs::Vector{SparseMatrixCSC}: Cached Jacobian matrices for each time step
  • μ∂²F!::Function: In-place function computing Hessian of Lagrangian
  • μ∂²fs::Vector{SparseMatrixCSC}: Cached Hessian matrices for each time step
  • μ∂²F_structure::SparseMatrixCSC: Sparsity structure of full trajectory Hessian
  • dim::Int: Total dimension of dynamics (sum of all integrator state dimensions)

Constructor

TrajectoryDynamics(
    integrators::Vector{<:AbstractIntegrator},
    traj::NamedTrajectory;
    verbose=false
)

Create trajectory dynamics from integrators and a trajectory structure.

Example

G = rand(2, 2)
integrator = BilinearIntegrator(G, traj, :x, :u)
dynamics = TrajectoryDynamics(integrator, traj)
source

Problems

DirectTrajOpt.Problems.DirectTrajOptProblemType
mutable struct DirectTrajOptProblem

A direct trajectory optimization problem containing all information needed for setup and solution.

Fields

  • trajectory::NamedTrajectory: The trajectory containing optimization variables and data
  • objective::Objective: The objective function to minimize
  • dynamics::TrajectoryDynamics: The system dynamics (integrators)
  • constraints::Vector{<:AbstractConstraint}: Constraints on the trajectory

Constructors

DirectTrajOptProblem(
    traj::NamedTrajectory,
    obj::Objective,
    integrators::Vector{<:AbstractIntegrator};
    constraints::Vector{<:AbstractConstraint}=AbstractConstraint[]
)

Create a problem from a trajectory, objective, and integrators. Trajectory constraints (initial, final, bounds) are automatically extracted and added to the constraint list.

Example

traj = NamedTrajectory((x = rand(2, 10), u = rand(1, 10)), timestep=:Δt)
obj = QuadraticRegularizer(:u, traj, 1.0)
integrator = BilinearIntegrator(G, traj, :x, :u)
prob = DirectTrajOptProblem(traj, obj, integrator)
source
DirectTrajOpt.Problems.get_trajectory_constraintsMethod
get_trajectory_constraints(traj::NamedTrajectory)

Extract and create constraints from a NamedTrajectory's initial, final, and bounds specifications.

Arguments

  • traj::NamedTrajectory: Trajectory with specified initial conditions, final conditions, and/or bounds

Returns

  • Vector{AbstractConstraint}: Vector of constraints including:
    • Initial value equality constraints (from traj.initial)
    • Final value equality constraints (from traj.final)
    • Bounds constraints (from traj.bounds)

Details

The function automatically handles time indices based on which constraints are specified:

  • If both initial and final constraints exist for a component, bounds apply to interior points (2:N-1)
  • If only initial exists, bounds apply from second point onward (2:N)
  • If only final exists, bounds apply up to second-to-last point (1:N-1)
  • If neither exist, bounds apply to all time points (1:N)
source

Problem Solvers

Problem Solvers

DirectTrajOpt.Solvers.solve!Method
solve!(
    prob::DirectTrajOptProblem;
    options::IpoptOptions=IpoptOptions(),
    max_iter::Int=options.max_iter,
    verbose::Bool=true,
    linear_solver::String=options.linear_solver,
    print_level::Int=options.print_level,
    callback=nothing
)

Solve a direct trajectory optimization problem using Ipopt.

Arguments

  • prob::DirectTrajOptProblem: The trajectory optimization problem to solve.
  • options::IpoptOptions: Ipopt solver options. Default is IpoptOptions().
  • max_iter::Int: Maximum number of iterations for the optimization solver.
  • verbose::Bool: If true, print solver progress information.
  • linear_solver::String: Linear solver to use (e.g., "mumps", "pardiso", "ma27", "ma57", "ma77", "ma86", "ma97").
  • print_level::Int: Ipopt print level (0-12). Higher values provide more detailed output.
  • callback::Function: Optional callback function to execute during optimization.

Returns

  • nothing: The problem's trajectory is updated in place with the optimized solution.

Example

prob = DirectTrajOptProblem(trajectory, objective, dynamics)
solve!(prob; max_iter=100, verbose=true)
source