Library
Common Interface
DirectTrajOpt.CommonInterface — Module
Common interface for components (integrators and constraints).
This module defines the generic interface functions that both integrators and constraints implement through multiple dispatch. This avoids naming conflicts when both modules are imported.
DirectTrajOpt.CommonInterface.eval_hessian_of_lagrangian — Function
eval_hessian_of_lagrangian(component, traj::NamedTrajectory, μ::AbstractVector)High-level method to evaluate and return the full Hessian of the Lagrangian for the component.
For integrators: Computes the Hessian using ForwardDiff across all timesteps. For constraints: Calls hessianoflagrangian to fill compact storage, then assembles the full sparse Hessian.
Arguments
component: An integrator or constrainttraj: The trajectory providing values and structureμ: Lagrange multipliers
Returns
- A sparse matrix representing the full Hessian μ'∇²f
DirectTrajOpt.CommonInterface.eval_jacobian — Function
eval_jacobian(component, traj::NamedTrajectory)High-level method to evaluate and return the full Jacobian for the component.
For integrators: Computes the Jacobian using ForwardDiff across all timesteps. For constraints: Calls jacobian! to fill compact storage, then assembles the full sparse Jacobian.
Arguments
component: An integrator or constrainttraj: The trajectory providing values and structure
Returns
- A sparse matrix representing the full Jacobian
DirectTrajOpt.CommonInterface.evaluate! — Function
evaluate!(values, component, traj::NamedTrajectory)Evaluate the component (constraint or dynamics) and store the result in-place in values.
For integrators: Computes dynamics violations δ = f(x{k+1}, xk, u_k, ...) for all timesteps. For constraints: Computes constraint violations g(x) for all applicable timesteps/variables.
Arguments
values: Pre-allocated vector to store the evaluation resultscomponent: An integrator or constrainttraj: The trajectory providing values and structure
Returns
- Nothing (modifies
valuesin-place)
DirectTrajOpt.CommonInterface.hessian_of_lagrangian — Function
hessian_of_lagrangian(component, μ, args...)Compute the Hessian of the Lagrangian (weighted by multipliers μ) for the component.
Arguments
component: An integrator, constraint, or other componentμ: Lagrange multipliersargs...: Component-specific arguments (e.g., knot points, trajectory values)
Returns
- A sparse matrix representing μ'∇²f
DirectTrajOpt.CommonInterface.hessian_of_lagrangian! — Function
hessian_of_lagrangian!(μ∂²f, component, μ, args...)Compute the Hessian of the Lagrangian in-place, storing the result in μ∂²f.
Arguments
μ∂²f: Pre-allocated sparse matrix for the Hessiancomponent: An integrator, constraint, or other componentμ: Lagrange multipliersargs...: Component-specific arguments (e.g., knot points, trajectory values)
DirectTrajOpt.CommonInterface.hessian_structure — Function
hessian_structure(component, traj::NamedTrajectory)Return the sparsity structure of the Hessian of the Lagrangian for the given component.
Arguments
component: An integrator, constraint, or other componenttraj: The trajectory providing dimension information
Returns
- A sparse matrix representing the Hessian structure
DirectTrajOpt.CommonInterface.jacobian! — Function
jacobian!(∂f, component, args...)Compute the Jacobian of the component in-place, storing the result in ∂f.
Arguments
∂f: Pre-allocated sparse matrix for the Jacobiancomponent: An integrator, constraint, or other componentargs...: Component-specific arguments (e.g., knot points, trajectory values)
DirectTrajOpt.CommonInterface.jacobian_structure — Function
jacobian_structure(component, traj::NamedTrajectory)Return the sparsity structure of the Jacobian for the given component. This should return a sparse matrix with the same structure as the Jacobian, where non-zero entries indicate where partial derivatives exist.
Arguments
component: An integrator, constraint, or other componenttraj: The trajectory providing dimension information
Returns
- A sparse matrix representing the Jacobian structure
Constraints
DirectTrajOpt.Constraints.AllEqualConstraint — Type
struct AllEqualConstraint <: AbstractLinearConstraintConstraint that all components of a variable should be equal to each other. Commonly used for fixed timesteps.
Fields
var_name::Symbol: Variable name to constraincomponent_index::Int: Which component of the variable (1 for scalar variables)label::String: Constraint label
DirectTrajOpt.Constraints.BoundsConstraint — Type
struct BoundsConstraint <: AbstractLinearConstraintRepresents a box constraint defined by variable names. Indices and concrete bounds are computed in constrain!.
Fields
var_names::Union{Symbol, Vector{Symbol}}: Variable name(s) to constraintimes::Union{Nothing, Vector{Int}}: Time indices (nothing for global variables)bounds_values::Union{Float64, Vector{Float64}, Tuple{Vector{Float64}, Vector{Float64}}}: Bound specificationis_global::Bool: Whether this is a global variable constraintsubcomponents::Union{Nothing, UnitRange{Int}}: Optional subcomponent selectionlabel::String: Constraint label
DirectTrajOpt.Constraints.BoundsConstraint — Method
BoundsConstraint(
name::Symbol,
ts::Vector{Int},
bounds::Union{Float64, Vector{Float64}, Tuple{Vector{Float64}, Vector{Float64}}};
subcomponents=nothing,
label="bounds constraint on trajectory variable $name"
)Constructs box constraint for trajectory variable. Indices are computed when applied to a trajectory.
Arguments
name: Variable namets: Time indicesbounds: Can be:- Scalar: symmetric bounds [-bounds, bounds]
- Vector: symmetric bounds [-bounds, bounds] element-wise
- Tuple: (lowerbounds, upperbounds)
DirectTrajOpt.Constraints.EqualityConstraint — Type
struct EqualityConstraint <: AbstractLinearConstraintRepresents a linear equality constraint defined by variable names. Indices are computed when constraint is applied in constrain!.
Fields
var_names::Union{Symbol, Vector{Symbol}}: Variable name(s) to constraintimes::Union{Nothing, Vector{Int}}: Time indices (nothing for global variables)values::Vector{Float64}: Constraint valuesis_global::Bool: Whether this is a global variable constraintlabel::String: Constraint label
DirectTrajOpt.Constraints.EqualityConstraint — Method
EqualityConstraint(
name::Symbol,
ts::Vector{Int},
val::Union{Float64, Vector{Float64}};
label="equality constraint on trajectory variable $name"
)Constructs equality constraint for trajectory variable. Indices are computed when applied to a trajectory.
DirectTrajOpt.Constraints.L1SlackConstraint — Type
struct L1SlackConstraint <: AbstractLinearConstraintLinear constraint tying a slack variable to the absolute value of another variable.
For each timestep k and component i, enforces:
\[v_{k,i} \leq s_{k,i}, \quad -v_{k,i} \leq s_{k,i} \quad \Longleftrightarrow \quad |v_{k,i}| \leq s_{k,i}\]
The bound s ≥ 0 is expected to come from the trajectory's bounds on the slack component. When combined with a LinearRegularizer on the slack variable, this yields an exact L1 penalty on v.
Fields
var_name::Symbol: Variable to penalize (e.g.:du)slack_name::Symbol: Slack variable name (e.g.:s_du)times::Vector{Int}: Time indices where constraint is appliedlabel::String: Constraint label
DirectTrajOpt.Constraints.L1SlackConstraint — Method
L1SlackConstraint(
var_name::Symbol,
slack_name::Symbol,
traj::NamedTrajectory;
times::AbstractVector{Int}=1:traj.N,
label="L1 slack constraint: |var_name| ≤ slack_name"
)Construct an L1 slack constraint tying |var_name| to slack_name.
DirectTrajOpt.Constraints.NonlinearGlobalConstraint — Type
NonlinearGlobalConstraint{F} <: AbstractNonlinearConstraintConstraint applied to global (trajectory-wide) parameters only.
Computes Jacobians and Hessians on-the-fly using automatic differentiation. For pre-allocated optimization, see Piccolissimo.OptimizedNonlinearGlobalConstraint.
Fields
g::F: Constraint function mapping global variables -> constraint valuesglobal_names::Vector{Symbol}: Names of global variables the constraint depends onequality::Bool: If true, g(globals) = 0; if false, g(globals) ≤ 0dim::Int: Dimension of constraint outputglobal_dim::Int: Combined dimension of all constrained global variables
DirectTrajOpt.Constraints.NonlinearGlobalKnotPointConstraint — Type
NonlinearGlobalKnotPointConstraint{F} <: AbstractNonlinearConstraintConstraint applied at individual knot points that also depends on global parameters.
Computes Jacobians and Hessians on-the-fly using automatic differentiation. For pre-allocated optimization, see Piccolissimo.OptimizedNonlinearGlobalKnotPointConstraint.
Fields
g::F: Constraint function mapping (knotpointvars..., global_vars..., params) -> constraint valuesvar_names::Vector{Symbol}: Names of knot point variables the constraint depends onglobal_names::Vector{Symbol}: Names of global variables the constraint depends ontimes::Vector{Int}: Time indices where constraint is appliedequality::Bool: If true, g(x, globals) = 0; if false, g(x, globals) ≤ 0params::Vector: Parameters for each time indexg_dim::Int: Dimension of constraint output at each time stepvar_dim::Int: Combined dimension of knot point variablesglobal_dim::Int: Combined dimension of global variablescombined_dim::Int: vardim + globaldimdim::Int: Total constraint dimension (g_dim * length(times))
DirectTrajOpt.Constraints.NonlinearKnotPointConstraint — Type
NonlinearKnotPointConstraint{F} <: AbstractNonlinearConstraintConstraint applied at individual knot points over a trajectory.
Computes Jacobians and Hessians on-the-fly using automatic differentiation. For pre-allocated optimization, see Piccolissimo.OptimizedNonlinearKnotPointConstraint.
Fields
g::F: Constraint function mapping (variables..., params) -> constraint valuesvar_names::Vector{Symbol}: Names of trajectory variables the constraint depends onequality::Bool: If true, g(x) = 0; if false, g(x) ≤ 0times::Vector{Int}: Time indices where constraint is appliedparams::Vector: Parameters for each time index (e.g., time-varying targets)g_dim::Int: Dimension of constraint output at each time stepvar_dim::Int: Combined dimension of all constrained variablesdim::Int: Total constraint dimension (g_dim * length(times))
DirectTrajOpt.Constraints.NonlinearKnotPointConstraint — Method
(constraint::NonlinearKnotPointConstraint)(δ, zₖ::KnotPoint, k::Int)Evaluate the constraint at a single knot point.
DirectTrajOpt.Constraints.SymmetryConstraint — Type
struct SymmetryConstraint <: AbstractLinearConstraintConstraint enforcing symmetry in trajectory variables across time. Even symmetry: x[t] = x[N-t+1] Odd symmetry: x[t] = -x[N-t+1]
Fields
var_name::Symbol: Variable name to constraincomponent_indices::Vector{Int}: Which components of the variableeven::Bool: True for even symmetry (x[t] = x[N-t+1]), false for odd (-x[t] = x[N-t+1])include_timestep::Bool: Whether to also enforce even symmetry on timestepslabel::String: Constraint label
DirectTrajOpt.Constraints.TimeConsistencyConstraint — Type
struct TimeConsistencyConstraint <: AbstractLinearConstraintConstraint that enforces consistency between time values and timesteps: t{k+1} = tk + Δt_k for k = 1, ..., T-1
This is used when both absolute times (:t) and timesteps (:Δt) are stored in the trajectory and need to remain consistent during optimization.
Fields
time_name::Symbol: Name of the time variable (default:t)timestep_name::Symbol: Name of the timestep variable (default:Δt)label::String: Constraint label
DirectTrajOpt.Constraints.TimeConsistencyConstraint — Method
TimeConsistencyConstraint(;
time_name::Symbol=:t,
timestep_name::Symbol=:Δt,
label="time consistency constraint (t_{k+1} = t_k + Δt_k)"
)Construct a constraint enforcing t{k+1} = tk + Δt_k for all k.
Arguments
time_name: Name of the time variable in the trajectory (default:t)timestep_name: Name of the timestep variable in the trajectory (default:Δt)label: Constraint label for logging/debugging
DirectTrajOpt.Constraints.TotalConstraint — Type
struct TotalConstraint <: AbstractLinearConstraintConstraint that the sum of a variable's components equals a target value. Commonly used for trajectory duration constraints.
Fields
var_name::Symbol: Variable name to sumcomponent_index::Int: Which component of the variable (1 for scalar variables)value::Float64: Target sum valuelabel::String: Constraint label
Note
When applied to the trajectory's timestep variable, only the first N-1 timesteps are summed (the last knot point has no duration after it). For other variables, all N values are summed.
DirectTrajOpt.CommonInterface.eval_hessian_of_lagrangian — Method
eval_hessian_of_lagrangian(constraint::NonlinearGlobalConstraint, traj::NamedTrajectory, μ::AbstractVector)Compute and return the full Hessian of the Lagrangian using automatic differentiation.
DirectTrajOpt.CommonInterface.eval_hessian_of_lagrangian — Method
eval_hessian_of_lagrangian(constraint::NonlinearGlobalKnotPointConstraint, traj::NamedTrajectory, μ::AbstractVector)Compute and return the full Hessian of the Lagrangian using automatic differentiation.
DirectTrajOpt.CommonInterface.eval_hessian_of_lagrangian — Method
eval_hessian_of_lagrangian(constraint::NonlinearKnotPointConstraint, traj::NamedTrajectory, μ::AbstractVector)Compute and return the full Hessian of the Lagrangian using automatic differentiation.
DirectTrajOpt.CommonInterface.eval_jacobian — Method
eval_jacobian(constraint::NonlinearGlobalConstraint, traj::NamedTrajectory)Compute and return the full Jacobian using automatic differentiation.
DirectTrajOpt.CommonInterface.eval_jacobian — Method
eval_jacobian(constraint::NonlinearGlobalKnotPointConstraint, traj::NamedTrajectory)Compute and return the full Jacobian using automatic differentiation.
DirectTrajOpt.CommonInterface.eval_jacobian — Method
eval_jacobian(constraint::NonlinearKnotPointConstraint, traj::NamedTrajectory)Compute and return the full Jacobian using automatic differentiation.
DirectTrajOpt.CommonInterface.evaluate! — Method
evaluate!(values::AbstractVector, constraint::NonlinearGlobalConstraint, traj::NamedTrajectory)Evaluate the global constraint, storing results in-place in values. This is part of the common interface with integrators.
DirectTrajOpt.CommonInterface.evaluate! — Method
evaluate!(values::AbstractVector, constraint::NonlinearGlobalKnotPointConstraint, traj::NamedTrajectory)Evaluate the global knot point constraint, storing results in-place in values. This is part of the common interface with integrators.
DirectTrajOpt.CommonInterface.evaluate! — Method
evaluate!(values::AbstractVector, constraint::NonlinearKnotPointConstraint, traj::NamedTrajectory)Evaluate the constraint at all specified time indices, storing results in-place in values. This is part of the common interface with integrators.
DirectTrajOpt.Constraints.DurationConstraint — Method
DurationConstraint(value::Float64; label="duration constraint of $value")Constraint that the total trajectory duration equals a target value. The trajectory's timestep variable is inferred when applied.
Note
Duration is computed as the sum of the first N-1 timesteps, since the final knot point represents the end state and has no duration after it.
DirectTrajOpt.Constraints.GlobalBoundsConstraint — Method
GlobalBoundsConstraint(
name::Symbol,
bounds::Union{Float64, Vector{Float64}, Tuple{Vector{Float64}, Vector{Float64}}};
label="bounds constraint on global variable $name"
)Constructs box constraint for global variable. Indices are computed when applied to a trajectory.
DirectTrajOpt.Constraints.GlobalEqualityConstraint — Method
GlobalEqualityConstraint(
name::Symbol,
val::Union{Float64, Vector{Float64}};
label="equality constraint on global variable $name"
)::EqualityConstraintConstructs equality constraint for global variable. Indices are computed when applied to a trajectory.
DirectTrajOpt.Constraints.SymmetricControlConstraint — Method
SymmetricControlConstraint(
name::Symbol,
idx::Vector{Int};
even=true,
include_timestep=true,
label="symmetry constraint on $name"
)Constraint enforcing symmetry on control variables. Indices are computed when applied to a trajectory.
DirectTrajOpt.Constraints.TimeStepsAllEqualConstraint — Method
TimeStepsAllEqualConstraint(;label="timesteps all equal constraint")Constraint that all timesteps are equal (for fixed-timestep trajectories). The trajectory's timestep variable is inferred when applied.
DirectTrajOpt.Constraints.get_full_hessian — Function
get_full_hessian(constraint, traj::NamedTrajectory)Assemble the full sparse Hessian matrix from compact per-timestep blocks.
DirectTrajOpt.Constraints.get_full_jacobian — Function
get_full_jacobian(constraint, traj::NamedTrajectory)Assemble the full sparse Jacobian matrix from compact per-timestep blocks.
DirectTrajOpt.Constraints.hessian_of_lagrangian! — Function
hessian_of_lagrangian!(constraint, traj::NamedTrajectory, μ::AbstractVector)Compute the Hessian of the Lagrangian (μ'g) for the constraint in-place.
DirectTrajOpt.Constraints.test_constraint — Method
test_constraint(
constraint::AbstractNonlinearConstraint,
traj::NamedTrajectory;
show_jacobian_diff=false,
show_hessian_diff=false,
test_equality=true,
atol=1e-5,
rtol=1e-5
)Test that constraint Jacobian and Hessian match finite difference approximations.
Arguments
constraint: Constraint to testtraj: Trajectory to evaluate constraint on
Keyword Arguments
show_jacobian_diff=false: Print detailed Jacobian differencesshow_hessian_diff=false: Print detailed Hessian differencestest_equality=true: Test element-wise equality (vs norm-based test)atol=1e-5: Absolute tolerancertol=1e-5: Relative tolerance
Returns
Tuple of (∂g, ∂gfinitediff, μ∂²g, μ∂²gfinitediff) for inspection
Example
g(x) = [norm(x) - 1.0]
constraint = NonlinearKnotPointConstraint(g, :x, traj)
test_constraint(constraint, traj)Integrators
DirectTrajOpt.Integrators.BilinearIntegrator — Type
BilinearIntegrator <: AbstractBilinearIntegratorIntegrator 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_name::Symbol: State variable nameu_name::Symbol: Control variable namex_dim::Int: Dimension of state variablevar_dim::Int: Combined dimension of all variables this integrator depends on (2*xdim + udim + 1)dim::Int: Total constraint dimension (x_dim * (N-1))∂fs::Vector{SparseMatrixCSC{Float64, Int}}: Pre-allocated compact Jacobian storage (xdim × vardim per timestep)μ∂²fs::Vector{SparseMatrixCSC{Float64, Int}}: Pre-allocated compact Hessian storage (vardim × vardim per timestep)
Constructors
BilinearIntegrator(G::Function, x::Symbol, u::Symbol, traj::NamedTrajectory)Arguments
G: Function taking control u and returning state matrix (xdim × xdim)x: State variable nameu: Control variable nametraj: Trajectory structure used to determine dimensions and pre-allocate storage
Dynamics
Computes the constraint: x{k+1} - exp(Δt * G(uk)) * x_k = 0 Dependencies: xₖ, uₖ, Δtₖ, xₖ₊₁
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, :x, :u, traj)DirectTrajOpt.Integrators.DerivativeIntegrator — Type
DerivativeIntegrator <: AbstractIntegratorIntegrator for derivative constraints of the form xₖ₊₁ - xₖ - Δt * ẋₖ = 0.
This enforces smoothness by relating a variable to its derivative.
Fields
f::Function: Constraint function f(xₖ₊₁, xₖ, ẋₖ, Δtₖ) = xₖ₊₁ - xₖ - Δtₖ * ẋₖx_name::Symbol: Variable nameẋ_name::Symbol: Derivative variable namex_dim::Int: Dimension of variablevar_dim::Int: Combined dimension (2*x_dim + 1 for xₖ, ẋₖ, Δtₖ, xₖ₊₁)dim::Int: Total constraint dimension (x_dim * (N-1))∂fs::Vector{SparseMatrixCSC{Float64, Int}}: Compact Jacobian storageμ∂²fs::Vector{SparseMatrixCSC{Float64, Int}}: Compact Hessian storage
Example
# Enforce velocity smoothness: vₖ₊₁ - vₖ - Δt * aₖ = 0
integrator = DerivativeIntegrator(:v, :a, traj)DirectTrajOpt.Integrators.dense — Method
dense(vals, structure, shape)Convert sparse data to dense matrix.
Arguments
vals: vector of valuesstructure: vector of tuples of indicesshape: tuple of matrix dimensions
DirectTrajOpt.Integrators.show_diffs — Method
show_diffs(A::Matrix, B::Matrix)Show differences between matrices.
Objectives
DirectTrajOpt.Objectives.AbstractObjective — Type
AbstractObjectiveAbstract type for all objective functions in trajectory optimization.
Concrete objective types must implement:
objective_value(obj, traj): Evaluate the objective at trajectorygradient!(∇, obj, traj): Compute gradient in-place (gradient is always dense)hessian_structure(obj, traj): Return sparsity structure as sparse matrixget_full_hessian(obj, traj): Return the full Hessian matrix
Objectives support addition and scalar multiplication through CompositeObjective.
Note: Unlike constraints and integrators, objective gradients are always dense, so no gradient_structure method is needed. The gradient! method fills the entire ∇ vector.
DirectTrajOpt.Objectives.CompositeObjective — Type
CompositeObjective <: AbstractObjectiveRepresents a weighted sum or composition of multiple objectives.
Fields
objectives::Vector{AbstractObjective}: Individual objectives to combineweights::Vector{Float64}: Weight for each objective
DirectTrajOpt.Objectives.GlobalKnotPointObjective — Type
GlobalKnotPointObjective <: AbstractObjectiveKnot point objective that includes both time-varying and global trajectory components.
Objective function ℓ operates on extracted variable values:
\[J = \sum_{k \in \text{times}} Q_k \ell([x_k; g], p_k)\]
where ℓ receives both knot point variables and global variables concatenated.
Fields
ℓ::Function: Objective function mapping (knotvars + globalvars, params) → scalar costvar_names::Vector{Symbol}: Names of trajectory variables at knot pointsglobal_names::Vector{Symbol}: Names of global trajectory variablestimes::Vector{Int}: Time indices where objective is evaluatedparams::Vector: Parameters for each time indexQs::Vector{Float64}: Weights for each time index
DirectTrajOpt.Objectives.GlobalObjective — Type
GlobalObjective <: AbstractObjectiveObjective that only involves global (non-time-varying) trajectory components.
Objective function ℓ operates on extracted global variable values:
\[J = Q \cdot \ell(\text{global\_vars})\]
Fields
ℓ::Function: Objective function mapping global variables → scalar costglobal_names::Vector{Symbol}: Names of global trajectory variablesQ::Float64: Weight for the objective
Constructor
GlobalObjective(
ℓ::Function,
global_names::Union{Symbol, AbstractVector{Symbol}},
traj::NamedTrajectory;
Q::Float64=1.0
)DirectTrajOpt.Objectives.KnotPointObjective — Type
KnotPointObjective <: AbstractObjectiveKnot point summed objective function for trajectory optimization.
Stores the objective function ℓ that operates on extracted variable values:
\[J = \sum_{k \in \text{times}} Q_k \ell(x_k, p_k)\]
where ℓ is evaluated on trajectory variables at each knot point.
Fields
ℓ::Function: Objective function mapping (variables..., params) -> scalar costvar_names::Vector{Symbol}: Names of trajectory variables the objective depends ontimes::Vector{Int}: Time indices where objective is evaluatedparams::Vector: Parameters for each time indexQs::Vector{Float64}: Weights for each time index∂²Ls::Vector{SparseMatrixCSC{Float64, Int}}: Preallocated sparse Hessian storage (one per timestep)
Constructor
KnotPointObjective(
ℓ::Function,
names::Union{Symbol, AbstractVector{Symbol}},
traj::NamedTrajectory,
params::AbstractVector;
times::AbstractVector{Int}=1:traj.N,
Qs::AbstractVector{Float64}=ones(length(times))
)For single variable: ℓ(x, p) where x is the variable values at a knot point For multiple variables: ℓ(x, u, p) where each argument corresponds to a variable in names
Examples
# Single variable
obj = KnotPointObjective((x, _) -> norm(x)^2, :x, traj, fill(nothing, traj.N))
# Multiple variables - concatenated
obj = KnotPointObjective((xu, _) -> xu[1]^2 + xu[2]^2, [:x, :u], traj, fill(nothing, traj.N))
# With parameters and weights
obj = KnotPointObjective(
(x, p) -> norm(x - p)^2, :x, traj, [x_targets...];
times=1:10, Qs=[1.0, 2.0, ...]
)DirectTrajOpt.Objectives.LinearRegularizer — Type
LinearRegularizer <: AbstractObjectiveLinear regularization objective for a trajectory component.
Computes:
\[J = \sum_{k \in \text{times}} \sum_i R_i \cdot v_{k,i} \cdot \Delta t_k\]
Used for L1 penalty via slack variables: when applied to a non-negative slack variable s ≥ 0 satisfying |du| ≤ s, minimizing Σ R_i s_i Δt yields the exact L1 norm of du.
Gradients and Hessians are computed analytically. The Hessian has only cross-terms ∂²J/∂v∂Δt = R_i (no diagonal).
Fields
name::Symbol: Name of the variable to regularizeR::Vector{Float64}: Per-component weightstimes::Vector{Int}: Time indices where regularization is applied
Constructor
LinearRegularizer(
name::Symbol,
traj::NamedTrajectory,
R::Union{Real, AbstractVector{<:Real}};
times::AbstractVector{Int}=1:traj.N
)DirectTrajOpt.Objectives.MinimumTimeObjective — Type
MinimumTimeObjective <: AbstractObjectiveObjective that minimizes total trajectory duration.
Computes:
\[J = D \sum_{k=1}^{N-1} \Delta t_k\]
Fields
D::Float64: Scaling factor for minimum time objective
Constructor
MinimumTimeObjective(traj::NamedTrajectory; D::Float64=1.0)
MinimumTimeObjective(traj::NamedTrajectory, D::Real)DirectTrajOpt.Objectives.NullObjective — Type
NullObjective <: AbstractObjectiveA zero objective that contributes nothing to the cost. Useful as a placeholder or when only constraints matter.
DirectTrajOpt.Objectives.QuadraticRegularizer — Type
QuadraticRegularizer <: AbstractObjectiveQuadratic regularization objective for a trajectory component.
Computes:
\[J = \sum_{k \in \text{times}} \frac{1}{2} (v_k - v_\text{baseline})^T R (v_k - v_\text{baseline}) \Delta t\]
Gradients and Hessians are computed analytically.
Fields
name::Symbol: Name of the variable to regularizeR::Vector{Float64}: Diagonal weight matrixbaseline::Matrix{Float64}: Baseline values (column per timestep)times::Vector{Int}: Time indices where regularization is applied
Constructor
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
)DirectTrajOpt.Objectives.TerminalObjective — Method
TerminalObjective(ℓ, names::Vector{Symbol}, traj; Q=1.0)Create a terminal objective that operates on multiple variables concatenated together.
This is useful for objectives that need to access multiple state variables at the final timestep, such as coherent fidelity objectives.
Arguments
ℓ::Function: Loss function taking concatenated values from all named variablesnames::Vector{Symbol}: Names of variables to concatenatetraj::NamedTrajectory: The trajectory
Keyword Arguments
Q::Float64=1.0: Weight on the objective
DirectTrajOpt.Objectives.TerminalObjective — Method
TerminalObjective(
ℓ::Function,
name::Symbol,
global_names::Union{Symbol, AbstractVector{Symbol}},
traj::NamedTrajectory;
Q::Float64=1.0
)Create a terminal (final time) objective that includes both knot point and global variables. This is a convenience wrapper around GlobalKnotPointObjective with times=[traj.N] and Qs=[Q].
Arguments
ℓ::Function: Objective function mapping concatenated [knotvars; globalvars] → scalarname::Symbol: Name of the knot point variableglobal_names: Name(s) of global variable(s)traj::NamedTrajectory: The trajectory
Example
# Terminal objective with knot point state and global parameter
TerminalObjective(
xg -> norm(xg[1:2] - xg[3:4])^2, # Distance from state to goal
:x, :x_goal, traj; Q=100.0
)DirectTrajOpt.Objectives.get_full_hessian — Function
get_full_hessian(obj::AbstractObjective, traj::NamedTrajectory)Compute and return the full Hessian matrix of the objective.
DirectTrajOpt.Objectives.gradient! — Function
gradient!(∇::AbstractVector, obj::AbstractObjective, traj::NamedTrajectory)Compute the gradient of the objective in-place. The gradient is always dense.
DirectTrajOpt.Objectives.hessian_structure — Function
hessian_structure(obj::AbstractObjective, traj::NamedTrajectory)Return the sparsity structure of the Hessian as a sparse matrix with non-zero entries where the Hessian has non-zero values.
DirectTrajOpt.Objectives.objective_value — Function
objective_value(obj::AbstractObjective, traj::NamedTrajectory)Evaluate the objective function at the given trajectory.
DirectTrajOpt.Objectives.test_objective — Method
test_objective(
obj::AbstractObjective,
traj::NamedTrajectory;
show_gradient_diff=false,
show_hessian_diff=false,
test_equality=true,
atol=1e-5,
rtol=1e-5
)Test an objective's gradient and Hessian implementations against finite differences.
Similar to test_integrator, this validates that computed derivatives match finite differences.
Arguments
obj::AbstractObjective: The objective to testtraj::NamedTrajectory: Trajectory defining the problem structure
Keyword Arguments
show_gradient_diff: Print element-wise differences in gradientshow_hessian_diff: Print element-wise differences in Hessiantest_equality: Test element-wise equality (vs. norm-based)atol: Absolute tolerance for comparisonsrtol: Relative tolerance for comparisons
Problems
DirectTrajOpt.Problems.DirectTrajOptProblem — Type
mutable struct DirectTrajOptProblemA direct trajectory optimization problem containing all information needed for setup and solution.
Fields
trajectory::NamedTrajectory: The trajectory containing optimization variables and dataobjective::AbstractObjective: The objective function to minimizeintegrators::Vector{<:AbstractIntegrator}: The integrators defining system dynamicsconstraints::Vector{<:AbstractConstraint}: Constraints on the trajectory
Constructors
DirectTrajOptProblem(
traj::NamedTrajectory,
obj::AbstractObjective,
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. The dynamics object is created by the evaluator at solve time.
Example
traj = NamedTrajectory((x = rand(2, 10), u = rand(1, 10)), timestep=:Δt)
obj = QuadraticRegularizer(:u, traj, 1.0)
integrator = BilinearIntegrator(G, :x, :u)
prob = DirectTrajOptProblem(traj, obj, integrator)DirectTrajOpt.Problems.get_trajectory_constraints — Method
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)
- Initial value equality constraints (from
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)
Problem Solvers
Problem Solvers
DirectTrajOpt.IpoptSolverExt.IpoptOptions — Type
Solver options for Ipopt
https://coin-or.github.io/Ipopt/OPTIONS.html#OPT_print_options_documentationDirectTrajOpt.IpoptSolverExt._fill_hessian_values! — Method
_fill_hessian_values!(H, evaluator, Z, σ, μ)Fill Hessian of Lagrangian values directly into output vector without building intermediate sparse matrices. Uses linear index map and direct SparseArrays access to eliminate allocations.
DirectTrajOpt.IpoptSolverExt._fill_jacobian_values! — Method
_fill_jacobian_values!(∂, evaluator, Z)Fill Jacobian values directly into the output vector without building intermediate sparse matrices. Uses pre-computed linear index map and direct SparseArrays access to eliminate allocations.
DirectTrajOpt.IpoptSolverExt._update_trajectory_cache! — Method
_update_trajectory_cache!(evaluator, Z⃗)Update the cached trajectory in-place with new data from Z⃗. Avoids repeated allocation of NamedTrajectory wrappers.
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 isIpoptOptions().max_iter::Int: Maximum number of iterations for the optimization solver.verbose::Bool: Iftrue, 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)