AstroSolve.jl
AstroSolve provides a framework for solving astrodynamics design problems through constrained optimization and event-driven architecture. AstroSolve builds an event sequence as a Directed Acyclic Graph (DAG), which will soon be fully differentiable using AD. The DAG model for optimization is inspired by the architecture used in NASA's Copernicus software.
Reference: Jacob Williams, Robert Falck, and Izaak Beekman. "Application of Modern Fortran to Spacecraft Trajectory Design and Optimization." AIAA/AAS Space Flight Mechanics Meeting, Kissimmee, FL, 2018. Available online
AstroSolve's architecture is based on a few key user-facing components:
| Component | Description |
|---|---|
| Models | Spacecraft, maneuvers, propagators, etc. that model the physics of the system |
| SolverVariables | Optimization variables that the solver can adjust during trajectory design |
| Constraints | Equality and inequality constraints that must be satisfied in the solution |
| Events | Discrete mission phases that execute actions (propagate, maneuver, etc.) and apply constraints |
| Sequence | DAG structure that defines event dependencies and execution order |
Quick Start
The example below uses AstroSolve's Event Sequence architecture to solve a simple orbit raising problem using IPOPT.
using Epicycle
# ============================================================================
# Shared Resources
# ============================================================================
# Create spacecraft with default orbital state
sat = Spacecraft()
# Create force models, integrator, and propagator
gravity = PointMassGravity(earth, (moon, sun))
forces = ForceModel(gravity)
integ = IntegratorConfig(Tsit5(); dt=10.0, reltol=1e-9, abstol=1e-9)
prop = OrbitPropagator(forces, integ)
# ============================================================================
# Event 1: TOI - Transfer Orbit Insertion
# ============================================================================
# Define TOI maneuver
toi = ImpulsiveManeuver(
axes = VNB(),
element1 = 0.1,
)
# Define solver variable for TOI delta-V
toi_var = SolverVariable(
calc = ManeuverCalc(toi, sat, DeltaVVector()),
name = "toi",
lower_bound = [-10.0, 0.0, 0.0],
upper_bound = [10.0, 0.0, 0.0],
)
# Define TOI event struct with event function, solver variables, and constraints
toi_fun() = maneuver!(sat, toi)
toi_event = Event(
name = "TOI Maneuver",
event = toi_fun,
vars = [toi_var],
funcs = []
)
# ============================================================================
# Event 2: Propagate to Apoapsis
# ============================================================================
# Define constraint on position magnitude at apoapsis
pos_target = 55000.0
pos_con = Constraint(
calc = OrbitCalc(sat, PosMag()),
lower_bounds = [pos_target],
upper_bounds = [pos_target],
scale = [1.0],
)
# Define propagation event to apoapsis with position constraint
prop_fun() = propagate!(prop, sat, StopAt(sat, PosDotVel(), 0.0; direction=-1))
prop_event = Event(
name = "Propagate to Apoapsis",
event = prop_fun,
funcs = [pos_con]
)
# ============================================================================
# Trajectory Optimization
# ============================================================================
# Create sequence and add events
seq = Sequence()
add_sequence!(seq, toi_event, prop_event)
# Solve trajectory optimization using default settings (finite differences, IPOPT)
result = solve_trajectory!(seq)
# Write a report documenting sequence and solution
report_sequence(seq)
report_solution(seq, result)Core Functions
This section contains reference material for the key user facing elements of AstroSolve. Note: Constraint is currently in the AstroCallbacks module but will be moved to AstroSolve in a future release.
AstroSolve.SolverVariable — Type
SolverVariable{C<:AbstractCalc}Represents a solver-controlled variable defined by a calculation container for trajectory optimization.
Fields
calc::C: Calculation container (e.g.,OrbitCalc,ManeuverCalc,BodyCalc)numvars::Int: Number of scalar variables for this calculationlower_bound::Vector{Float64}: Lower bounds for optimization (lengthnumvars, defaults to-Inf)upper_bound::Vector{Float64}: Upper bounds for optimization (lengthnumvars, defaults to+Inf)shift::Vector{Float64}: Variable shifting for numerical conditioning (lengthnumvars, defaults to0.0)scale::Vector{Float64}: Variable scaling for numerical conditioning (lengthnumvars, defaults to1.0)name::String: Optional human-readable identifier
Constructor
SolverVariable(; calc, lower_bound=nothing, upper_bound=nothing, shift=nothing,
scale=nothing, name="")Arguments
calc: Calculation container implementingAbstractCalcinterfacelower_bound: Scalar or vector of lower bounds (broadcast tonumvarsif scalar)upper_bound: Scalar or vector of upper bounds (broadcast tonumvarsif scalar)shift: Scalar or vector of shift values for numerical conditioningscale: Scalar or vector of scale factors for numerical conditioningname: Optional descriptive name for the variable
Notes
- See
AbstractCalcfor supported calculation types numvarsis automatically determined from the calculation variable type- Bounds, shift, and scale vectors are stored as
Float64for optimizer compatibility - The calculation must support setting values (
calc_is_settable(calc.var) == true)
Examples
using Epicycle
# Spacecraft and maneuver objects
sc = Spacecraft(); toi = ImpulsiveManeuver()
# 3D delta-V vector variable with individual component bounds
var_toi = SolverVariable(
calc = ManeuverCalc(toi, sc, DeltaVVector()),
lower_bound = [-2.0, -2.0, -2.0],
upper_bound = [2.0, 2.0, 2.0],
scale = [1.0, 1.0, 1.0],
name = "Departure DV Vector"
)
# Orbital element variable
var_sma = SolverVariable(
calc = OrbitCalc(sc, SMA()),
lower_bound = 6700.0,
upper_bound = 42000.0,
name = "Target Semi-Major Axis"
)AstroSolve.Event — Type
EventDefines an event node for trajectory sequence execution in mission design workflows.
Fields
name::String: Human-readable identifier for the eventevent::Function: Function closure to execute when this node runsvars::Vector{Any}: Vector of solver variables (SolverVariable) associated with this eventfuncs::Vector{Any}: Vector of constraint/objective functions for this event
Notes
Events are the fundamental building blocks of trajectory sequences, encapsulating both the action to be performed and the optimization variables/constraints associated with that action.
Examples
# Simple propagation event
prop_event = Event(
name = "Propagate to Apoapsis",
event = () -> propagate!(dynsys, integ, StopAtApoapsis(sat)),
)
# Event with solver variables and constraints
sc = Spacecraft(); toi = ImpulsiveManeuver();
var_dv = SolverVariable(
calc = ManeuverCalc(toi, sc, DeltaVVector()),
)
fun_dv = Constraint(
calc = ManeuverCalc(toi, sc, DeltaVMag()),
lower_bounds = [0.1],
upper_bounds = [0.1],
scale = [1.0],
)
maneuver_event = Event(
event = () -> maneuver!(sc, toi),
vars = [var_dv],
funcs = [fun_dv],
name = "Departure Maneuver",
)AstroSolve.Sequence — Type
SequenceA directed acyclic graph (DAG) representation of trajectory events for optimization.
In Epicycle, trajectories are modeled as sequences of events (maneuvers, propagations, etc.) that must occur in a specific order. The Sequence struct represents this as a DAG where events are nodes and dependencies are edges.
Fields
events::Vector{Event}: Collection of all events in the sequenceadj_map::Dict{Event, Vector{Event}}: Adjacency map defining event dependencies
Examples
seq = Sequence()
add_events!(seq, event1, Event[]) # Add root event
add_events!(seq, event2, [event1]) # event2 depends on event1
add_events!(seq, event3, [event1, event2]) # event3 depends on bothAstroSolve.add_events! — Function
add_events!(seq::Sequence, event::Event, dependencies::Vector{Event})Add an event to the sequence with specified dependencies.
This function adds event to the sequence DAG, ensuring that all events in dependencies must occur before event during execution. The internal adjacency map is updated to reflect these dependencies.
Arguments
seq::Sequence: Sequence to modifyevent::Event: Event to add to the sequencedependencies::Vector{Event}: Events directly linked to and precedingevent.
Examples
seq = Sequence()
add_events!(seq, maneuver_event, [prop_event]) # maneuver after propagation
add_events!(seq, final_event, [maneuver_event, other_event]) # final after bothAstroSolve.solve_trajectory! — Function
solve_trajectory!(seq::Sequence; record_iterations::Bool=false)Solve trajectory sequence using default SNOW/IPOPT configuration.
Keyword Arguments
record_iterations::Bool=false: Record solver iteration history for diagnostics. Whentrue, iteration segments are stored inspacecraft.history.iterations. Whenfalse(default), only the final solution is recorded inspacecraft.history.segments.
Notes
- Spacecraft history flags are automatically managed during optimization
- Original flag values are restored after optimization completes
- Final solution respects original recording settings (no forced recording)
solve_trajectory!(seq::Sequence, options::SNOW.Options; record_iterations::Bool=false)Solve trajectory sequence using specified SNOW optimization options.
Arguments
seq::Sequence: Event sequence defining the trajectory optimization problemoptions::SNOW.Options: SNOW/IPOPT solver options
Keyword Arguments
record_iterations::Bool=false: Record solver iteration history for diagnostics. Whentrue, iteration segments are stored inspacecraft.history.iterations. Whenfalse(default), only the final solution is recorded inspacecraft.history.segments.
Returns
Named tuple with:
variables: Optimal variable valuesobjective: Optimal objective function valueconstraints: Constraint values at optimal solutioninfo: Solver convergence information
Notes
- Automatically manages spacecraft history recording flags during optimization
- During iterations: records to
iterationsvector ifrecord_iterations=true - After convergence: restores original flags and records final solution to
segments - Original flag values are preserved and restored after optimization
AstroSolve.report_sequence — Function
report_sequence(seq::Sequence)Generate a comprehensive report summarizing the trajectory sequence configuration.
This function creates a human-readable summary of the sequence structure showing:
- Sequence overview (total events, variables, constraints)
- Event details with dependencies and execution order
- Variable information (names, sizes, bounds)
- Constraint summary (types, sizes, bounds)
Arguments
seq::Sequence: Trajectory sequence to analyze
Returns
- Nothing (prints directly to stdout for better formatting)
AstroSolve.report_solution — Function
report_solution(seq::Sequence, result)Generate a comprehensive report showing the optimized trajectory solution.
This function displays the optimization results including variable values, constraint satisfaction, and solution summary. It works with the result from solve_trajectory!().
Arguments
seq::Sequence: Original trajectory sequenceresult: Result tuple from solve_trajectory! with (variables, objective, info)
Returns
- Nothing (prints directly to stdout for better formatting)
API Reference
AstroSolve.EventAstroSolve.EventAstroSolve.SequenceAstroSolve.SequenceAstroSolve.SequenceManagerAstroSolve.SequenceManagerAstroSolve.SolverVariableAstroSolve.SolverVariableAstroSolve._push_stateful!AstroSolve.add_events!AstroSolve.add_events!AstroSolve.add_sequence!AstroSolve.apply_eventAstroSolve.copy_struct_fields!AstroSolve.default_snow_optionsAstroSolve.find_all_stateful_structsAstroSolve.format_boundsAstroSolve.get_bounds_descriptionAstroSolve.get_calc_descriptionAstroSolve.get_enhanced_constraint_descriptionAstroSolve.get_enhanced_variable_descriptionAstroSolve.get_event_constraint_mappingAstroSolve.get_fun_lower_boundsAstroSolve.get_fun_upper_boundsAstroSolve.get_fun_valuesAstroSolve.get_simplified_type_nameAstroSolve.get_sol_varAstroSolve.get_subject_descriptionAstroSolve.get_var_lower_boundsAstroSolve.get_var_lower_boundsAstroSolve.get_var_scalesAstroSolve.get_var_shiftsAstroSolve.get_var_upper_boundsAstroSolve.get_var_upper_boundsAstroSolve.get_var_valuesAstroSolve.get_var_valuesAstroSolve.order_unique_varsAstroSolve.report_sequenceAstroSolve.report_sequenceAstroSolve.report_solutionAstroSolve.report_solutionAstroSolve.reset_stateful_structs!AstroSolve.set_sol_varAstroSolve.set_var_valuesAstroSolve.set_var_valuesAstroSolve.show_constraint_componentsAstroSolve.show_variable_componentsAstroSolve.solve_trajectory!AstroSolve.solve_trajectory!AstroSolve.solve_trajectory!AstroSolve.solver_fun!AstroSolve.topo_sortBase.showBase.show
AstroSolve.Event — Method
Event(; name::String = "", event::Function = () -> nothing, vars = [], funcs = [])Outer constructor for Event with kwargs.
AstroSolve.Sequence — Method
Sequence()Create an empty trajectory sequence with no events or dependencies.
AstroSolve.SequenceManager — Type
SequenceManagerOrchestrator for trajectory sequence execution and optimization.
The SequenceManager prepares and manages a Sequence for iterative optimization. It topologically sorts events, orders variables and constraints, extracts bounds and scaling parameters, and manages stateful object restoration.
Fields
sequence::Sequence: Original sequence definitionsorted_events::Vector{Event}: Topologically sorted events for executionordered_vars::Vector{SolverVariable}: Unique variables in dependency orderordered_funcs::Vector{Constraint}: All constraints in execution orderfun_sizes::Vector{Int}: Number of scalar constraints per constraint functionvar_shift::Vector: Concatenated variable shift parameters for scalingvar_scale::Vector: Concatenated variable scale parameters for scalingvar_lower_bounds::Vector: Concatenated lower bounds for all variablesvar_upper_bounds::Vector: Concatenated upper bounds for all variablesstateful_structs::Vector{Any}: All stateful objects referenced in sequenceinitial_stateful_structs::Vector{Any}: Deep copies for state restoration
Examples
seq = Sequence()
add_events!(seq, event1, Event[])
add_events!(seq, event2, [event1])
sm = SequenceManager(seq)AstroSolve.SequenceManager — Method
SequenceManager(seq::Sequence)Construct a sequence manager from a trajectory sequence.
This constructor performs several optimization steps:
- Topologically sorts events to ensure dependency ordering
- Extracts unique variables while preserving dependency order
- Collects all constraints with their sizes
- Identifies and snapshots all stateful objects for reset capability
- Pre-computes variable bounds and scaling parameters
Notes
ErrorException: If sequence contains cycles (not a valid DAG)
AstroSolve.SolverVariable — Method
SolverVariable(; calc, lower_bound=nothing, upper_bound=nothing, shift=nothing, scale=nothing, name="")Convenience constructor that automatically infers the calculation type parameter C from the calc argument.
This eliminates the need to explicitly specify the type parameter when creating SolverVariable instances.
Usage Comparison
# Without convenience constructor (explicit type parameter):
var1 = SolverVariable{OrbitCalc}(calc=OrbitCalc(sat, SemiMajorAxis()))
# With convenience constructor (type inferred):
var2 = SolverVariable(calc=OrbitCalc(sat, SemiMajorAxis())) # PreferredThis constructor delegates to SolverVariable{C}(...) with the inferred type.
AstroSolve._push_stateful! — Method
_push_stateful!(out::Vector{Any}, seen::Set{UInt}, objs)Push unique stateful objects to output vector, preserving discovery order.
Only objects that are marked as stateful (via is_astrosolve_stateful) and not already seen are added to the output collection.
Arguments
out::Vector{Any}: Output collection to append toseen::Set{UInt}: Set tracking object IDs of already discovered objectsobjs: Collection of objects to check and potentially add
Returns
Vector{Any}: The modified output vector (for chaining)
AstroSolve.add_events! — Method
add_events!(seq::Sequence, event::Event, dependencies::Vector{Event})Add an event to the sequence with specified dependencies.
This function adds event to the sequence DAG, ensuring that all events in dependencies must occur before event during execution. The internal adjacency map is updated to reflect these dependencies.
Arguments
seq::Sequence: Sequence to modifyevent::Event: Event to add to the sequencedependencies::Vector{Event}: Events directly linked to and precedingevent.
Examples
seq = Sequence()
add_events!(seq, maneuver_event, [prop_event]) # maneuver after propagation
add_events!(seq, final_event, [maneuver_event, other_event]) # final after bothAstroSolve.add_sequence! — Method
add_sequence!(seq::Sequence, events::Event...)Add a linear sequence of events where each event depends on the previous one.
This is a convenience function for the common case of a linear event chain. Each event is added only with dependencies from the event immediately before it.
Arguments
seq::Sequence: Sequence to add events toevents::Event...: Events in execution order (first executes first)
Examples
seq = Sequence()
# These two are equivalent:
add_sequence!(seq, toi_event, prop_event, moi_event)
# Equivalent to:
add_events!(seq, toi_event, Event[])
add_events!(seq, prop_event, [toi_event])
add_events!(seq, moi_event, [prop_event])AstroSolve.apply_event — Method
apply_event(event::Event)Execute the event's function closure.
AstroSolve.copy_struct_fields! — Method
copy_struct_fields!(dest, src)Recursively copy all fields from src to dest struct, preserving mutable references.
This function performs deep copying for stateful struct restoration. For arrays, it overwrites contents rather than replacing references. For nested mutable structs, it recursively copies fields.
Arguments
dest: Destination struct to copy intosrc: Source struct to copy from
Returns
- Modified
deststruct
AstroSolve.default_snow_options — Method
default_snow_options()Create default SNOW optimization options.
AstroSolve.find_all_stateful_structs — Method
find_all_stateful_structs(ordered_vars, sorted_events)Find all unique stateful structs referenced by SolverVariables and events.
Stateful structs are objects that maintain state during trajectory propagation and need to be reset between optimization iterations. This function discovers all such objects referenced throughout the sequence.
Arguments
ordered_vars::Vector{SolverVariable}: Ordered optimization variablessorted_events::Vector{Event}: Topologically sorted events
Returns
Vector{Any}: Collection of unique stateful objects in discovery order
AstroSolve.format_bounds — Function
format_bounds(lower, upper, size=1)Format numeric bounds for display with consistent precision.
AstroSolve.get_bounds_description — Method
get_bounds_description(lower, upper)Describe bounds as equality, inequality, or range constraints.
AstroSolve.get_calc_description — Method
get_calc_description(calc)Extract calculation description from calc objects.
AstroSolve.get_enhanced_constraint_description — Method
get_enhanced_constraint_description(func)Get enhanced constraint description with calc info for display.
AstroSolve.get_enhanced_variable_description — Method
get_enhanced_variable_description(var::SolverVariable)Get enhanced variable description with calc info for display.
AstroSolve.get_event_constraint_mapping — Method
get_event_constraint_mapping(sm::SequenceManager)Create mapping from events to their constraints for reporting.
AstroSolve.get_fun_lower_bounds — Method
get_fun_lower_bounds(sm::SequenceManager)Extract lower bounds from all constraint functions in the sequence manager.
Arguments
sm::SequenceManager: Sequence manager containing constraint functions
Returns
Vector: Concatenated lower bounds for all constraints
AstroSolve.get_fun_upper_bounds — Method
get_fun_upper_bounds(sm::SequenceManager)Extract upper bounds from all constraint functions in the sequence manager.
Arguments
sm::SequenceManager: Sequence manager containing constraint functions
Returns
Vector: Concatenated upper bounds for all constraints
AstroSolve.get_fun_values — Method
get_fun_values(sm::SequenceManager)Evaluate all constraint functions in the sequence manager and return concatenated function values.
This function evaluates all constraints in the current state and returns their values as a flat vector. The ordering matches the constraint order established during SequenceManager construction.
Arguments
sm::SequenceManager: Sequence manager containing constraints to evaluate
Returns
Vector{Float64}: Concatenated constraint function values
Notes
This function assumes all stateful objects are in the correct state for constraint evaluation (typically after event sequence execution).
AstroSolve.get_simplified_type_name — Method
get_simplified_type_name(obj)Simplify type names by removing module prefixes and parameters.
AstroSolve.get_sol_var — Method
get_sol_var(var::SolverVariable)Get the solver variable values from the struct.
AstroSolve.get_subject_description — Method
get_subject_description(subject)Get subject description with type and name information.
AstroSolve.get_var_lower_bounds — Method
get_var_lower_bounds(ordered_vars)Extract lower bounds from all solver variables into a flat vector.
Arguments
ordered_vars::Vector{SolverVariable}: Variables to extract bounds from
Returns
Vector: Concatenated lower bounds for optimization constraints
AstroSolve.get_var_lower_bounds — Method
get_var_lower_bounds(sm::SequenceManager)Extract lower bounds from all variables in the sequence manager.
Arguments
sm::SequenceManager: Sequence manager containing variables
Returns
Vector: Concatenated lower bounds for optimization constraints
AstroSolve.get_var_scales — Method
get_var_scales(ordered_vars)Extract scale parameters from all solver variables into a flat vector.
Arguments
ordered_vars::Vector{SolverVariable}: Variables to extract scales from
Returns
Vector: Concatenated scale values for optimization scaling
AstroSolve.get_var_shifts — Method
get_var_shifts(ordered_vars)Extract shift parameters from all solver variables into a flat vector.
Arguments
ordered_vars::Vector{SolverVariable}: Variables to extract shifts from
Returns
Vector: Concatenated shift values for optimization scaling
AstroSolve.get_var_upper_bounds — Method
get_var_upper_bounds(ordered_vars)Extract upper bounds from all solver variables into a flat vector.
Arguments
ordered_vars::Vector{SolverVariable}: Variables to extract bounds from
Returns
Vector: Concatenated upper bounds for optimization constraints
AstroSolve.get_var_upper_bounds — Method
get_var_upper_bounds(sm::SequenceManager)Extract upper bounds from all variables in the sequence manager.
Arguments
sm::SequenceManager: Sequence manager containing variables
Returns
Vector: Concatenated upper bounds for optimization constraints
AstroSolve.get_var_values — Method
get_var_values(ordered_vars)Extract current values from all SolverVariables into a flat vector.
Arguments
ordered_vars::Vector{SolverVariable}: Variables to extract values from
Returns
Vector: Concatenated values from all variables
AstroSolve.get_var_values — Method
get_var_values(sm::SequenceManager)Extract current values from all variables in the sequence manager.
Arguments
sm::SequenceManager: Sequence manager containing variables
Returns
Vector: Concatenated values from all variables
AstroSolve.order_unique_vars — Method
order_unique_vars(sorted_events::Vector{Event})Extract unique solver variables from events in dependency order.
This function processes topologically sorted events and extracts all unique SolverVariables while preserving the order of first appearance. This ensures variables are ordered according to their dependency relationships.
Arguments
sorted_events::Vector{Event}: Events in topological order
Returns
Vector{SolverVariable}: Unique variables in dependency order
AstroSolve.report_sequence — Method
report_sequence(seq::Sequence)Generate a comprehensive report summarizing the trajectory sequence configuration.
This function creates a human-readable summary of the sequence structure showing:
- Sequence overview (total events, variables, constraints)
- Event details with dependencies and execution order
- Variable information (names, sizes, bounds)
- Constraint summary (types, sizes, bounds)
Arguments
seq::Sequence: Trajectory sequence to analyze
Returns
- Nothing (prints directly to stdout for better formatting)
AstroSolve.report_solution — Method
report_solution(seq::Sequence, result)Generate a comprehensive report showing the optimized trajectory solution.
This function displays the optimization results including variable values, constraint satisfaction, and solution summary. It works with the result from solve_trajectory!().
Arguments
seq::Sequence: Original trajectory sequenceresult: Result tuple from solve_trajectory! with (variables, objective, info)
Returns
- Nothing (prints directly to stdout for better formatting)
AstroSolve.reset_stateful_structs! — Method
reset_stateful_structs!(sm::SequenceManager)Reset all stateful structs in the sequence manager to their initial state.
This function is called before each optimization iteration to ensure a clean starting state for trajectory propagation. It restores all spacecraft, maneuvers, and other stateful objects to their initial conditions.
Arguments
sm::SequenceManager: Sequence manager containing stateful objects to reset
AstroSolve.set_sol_var — Method
set_sol_var(var::SolverVariable, val::Vector)Set the value(s) of the solver variable struct
AstroSolve.set_var_values — Method
set_var_values(sm::SequenceManager, x::AbstractVector)Set all solver variables in the sequence manager to values from vector x.
Arguments
sm::SequenceManager: Sequence manager containing ordered variablesx::AbstractVector: Flat vector of variable values
AstroSolve.set_var_values — Method
set_var_values(vec::Vector, ordered_vars::Vector{SolverVariable})Set values for an ordered collection of solver variables from a flat vector.
This function distributes values from a flat optimization vector to the appropriate SolverVariables, handling the different sizes of vector vs scalar variables correctly.
Arguments
vec::Vector: Flat vector of values to distributeordered_vars::Vector{SolverVariable}: Variables to set in order
Examples
# Set 3 variables: 2 scalars + 1 3-vector
x = [1.0, 2.0, 0.1, 0.2, 0.3] # 5 total values
set_var_values(x, [scalar1, scalar2, vector3])AstroSolve.show_constraint_components — Method
show_constraint_components(func, is_last_constraint::Bool)Show constraint component breakdown with proper indentation.
AstroSolve.show_variable_components — Method
show_variable_components(var::SolverVariable, is_last_var::Bool)Show variable component breakdown with proper indentation.
AstroSolve.solve_trajectory! — Method
solve_trajectory!(seq::Sequence, options::SNOW.Options; record_iterations::Bool=false)Solve trajectory sequence using specified SNOW optimization options.
Arguments
seq::Sequence: Event sequence defining the trajectory optimization problemoptions::SNOW.Options: SNOW/IPOPT solver options
Keyword Arguments
record_iterations::Bool=false: Record solver iteration history for diagnostics. Whentrue, iteration segments are stored inspacecraft.history.iterations. Whenfalse(default), only the final solution is recorded inspacecraft.history.segments.
Returns
Named tuple with:
variables: Optimal variable valuesobjective: Optimal objective function valueconstraints: Constraint values at optimal solutioninfo: Solver convergence information
Notes
- Automatically manages spacecraft history recording flags during optimization
- During iterations: records to
iterationsvector ifrecord_iterations=true - After convergence: restores original flags and records final solution to
segments - Original flag values are preserved and restored after optimization
AstroSolve.solve_trajectory! — Method
solve_trajectory!(seq::Sequence; record_iterations::Bool=false)Solve trajectory sequence using default SNOW/IPOPT configuration.
Keyword Arguments
record_iterations::Bool=false: Record solver iteration history for diagnostics. Whentrue, iteration segments are stored inspacecraft.history.iterations. Whenfalse(default), only the final solution is recorded inspacecraft.history.segments.
Notes
- Spacecraft history flags are automatically managed during optimization
- Original flag values are restored after optimization completes
- Final solution respects original recording settings (no forced recording)
AstroSolve.solver_fun! — Method
solver_fun!(F::AbstractVector, x::AbstractVector, sm::SequenceManager)Core optimization function: set variables, execute sequence, and evaluate constraints.
This is the primary interface between optimization solvers and Epicycle trajectory sequences. It performs a complete trajectory execution cycle:
- Reset: Restore all stateful objects to initial conditions
- Set Variables: Apply optimization variables to events
- Execute: Run events in topological order, applying maneuvers and propagations, etc.
- Evaluate: Collect function values at appropriate times during execution
Arguments
F::AbstractVector: Output vector to fill with constraint valuesx::AbstractVector: Input optimization variables (flat vector)sm::SequenceManager: Configured sequence manager
Returns
Int: Status code (0 for success)
Notes
- Function values are evaluated immediately after their associated events to capture the correct spacecraft state at that point in the trajectory
AstroSolve.topo_sort — Method
topo_sort(seq::Sequence)Perform topological sorting of events using Kahn's algorithm.
This function sorts the events in the sequence. If event A depends on event B, then B will appear before A in the sorted output.
Arguments
seq::Sequence: Sequence containing events and their dependencies
Returns
Vector{Event}: Events sorted in dependency order (dependencies first)
Throws
ErrorException: If the sequence contains cycles (not a valid DAG)
Algorithm
Uses Kahn's algorithm:
- Calculate in-degree for all events
- Start with events having no dependencies (in-degree 0)
- Remove events and update in-degrees until all processed
- Detect cycles if any events remain unprocessed
A More Complex Example
This example demonstrates GEO orbit insertion including plane change, using 8 sequential events (3 maneuvers, 5 propagations) with constraints applied at multiple trajectory points. The solution employes a bi-elliptic transfer.
using Epicycle
# ============================================================================
# Shared Resources
# ============================================================================
# Create spacecraft
sat = Spacecraft(
state = CartesianState([3737.792, -4607.692, -2845.644, 5.411, 5.367, -1.566]),
time = Time("2000-01-01T11:59:28.000", UTC(), ISOT()),
name = "GeoSat-1"
)
# Create force models, integrator, and propagator
gravity = PointMassGravity(earth, ()) # Only Earth gravity
forces = ForceModel(gravity)
integ = IntegratorConfig(DP8(); abstol=1e-12, reltol=1e-12, dt=60.0)
prop = OrbitPropagator(forces, integ)
# ============================================================================
# Event 1: Propagate to Equatorial Plane Crossing
# ============================================================================
# Define propagation event to equatorial plane crossing
prop_to_z_crossing_1_fun() = propagate!(prop, sat, StopAt(sat, PosZ(), 0.0))
prop_to_z_crossing_1_event = Event(
name = "Prop to Z 1",
event = prop_to_z_crossing_1_fun,
)
# ============================================================================
# Event 2: TOI - Transfer Orbit Insertion
# ============================================================================
# Define TOI maneuver
toi = ImpulsiveManeuver(
axes = VNB(),
element1 = 2.518,
element2 = 0.0,
element3 = 0.0,
)
# Define solver variable for TOI delta-V
toi_var = SolverVariable(
calc = ManeuverCalc(toi, sat, DeltaVVector()),
name = "toi_v",
lower_bound = [0.0, 0.0, 0.0],
upper_bound = [8.0, 0.0, 0.0],
)
# Define TOI event struct with event function, solver variables, and constraints
toi_fun() = maneuver!(sat, toi)
toi_event = Event(
name = "TOI",
event = toi_fun,
vars = [toi_var],
)
# ============================================================================
# Event 3: Propagate to Apoapsis
# ============================================================================
# Define constraint on radius at apoapsis
apogee_radius_con = Constraint(
calc = OrbitCalc(sat, PosMag()),
lower_bounds = [85000.0],
upper_bounds = [85000.0],
scale = [1.0],
)
# Define propagation event to apoapsis with radius constraint
prop_to_apogee_fun() = propagate!(prop, sat, StopAt(sat, PosDotVel(), 0.0; direction=-1))
prop_to_apogee_event = Event(
name = "Prop to Apoapsis",
event = prop_to_apogee_fun,
funcs = [apogee_radius_con],
)
# ============================================================================
# Event 4: Propagate to Perigee
# ============================================================================
prop_to_perigee_1_fun() = propagate!(prop, sat, StopAt(sat, PosDotVel(), 0.0; direction=1))
prop_to_perigee_1_event = Event(
name = "Prop to Perigee 1",
event = prop_to_perigee_1_fun,
)
# ============================================================================
# Event 5: Propagate to Equatorial Plane Crossing Again
# ============================================================================
prop_to_z_crossing_2_fun() = propagate!(prop, sat, StopAt(sat, PosZ(), 0.0))
prop_to_z_crossing_2_event = Event(
name = "Prop to Z 2",
event = prop_to_z_crossing_2_fun,
)
# ============================================================================
# Event 6: MCC - Mid-Course Correction
# ============================================================================
# Define MCC maneuver
mcc = ImpulsiveManeuver(
axes = VNB(),
element1 = 0.559,
element2 = 0.588,
element3 = 0.0,
)
# Define solver variable for MCC delta-V
mcc_var = SolverVariable(
calc = ManeuverCalc(mcc, sat, DeltaVVector()),
name = "mcc_vn",
lower_bound = [-1.0, -1.0, -0.001],
upper_bound = [4.0, 1.0, 0.001],
)
# Define MCC event struct with event function and solver variables
mcc_fun() = maneuver!(sat, mcc)
mcc_event = Event(
name = "MCC",
event = mcc_fun,
vars = [mcc_var],
)
# ============================================================================
# Event 7: Propagate to Perigee and Check Constraints
# ============================================================================
# Define constraints on inclination and perigee radius
inclination_con = Constraint(
calc = OrbitCalc(sat, Inc()),
lower_bounds = [deg2rad(2.0)],
upper_bounds = [deg2rad(2.0)],
scale = [1.0],
)
perigee_radius_con = Constraint(
calc = OrbitCalc(sat, PosMag()),
lower_bounds = [42195.0],
upper_bounds = [42195.0],
scale = [1.0],
)
# Define propagation event to perigee with inclination and perigee radius constraints
prop_to_perigee_2_fun() = propagate!(prop, sat, StopAt(sat, PosDotVel(), 0.0; direction=1))
prop_to_perigee_2_event = Event(
name = "Prop to Perigee 2",
event = prop_to_perigee_2_fun,
funcs = [inclination_con, perigee_radius_con],
)
# ============================================================================
# Event 8: MOI - GEO Orbit Insertion
# ============================================================================
# Define MOI maneuver
moi = ImpulsiveManeuver(
axes = VNB(),
element1 = 0.282,
element2 = 0.0,
element3 = 0.0,
)
# Define solver variable for MOI delta-V
moi_var = SolverVariable(
calc = ManeuverCalc(moi, sat, DeltaVVector()),
name = "moi_v",
lower_bound = [-1.0, -0.001, -0.001],
upper_bound = [4.0, 0.001, 0.001],
)
# Define constraint on semi-major axis at GEO
final_sma_con = Constraint(
calc = OrbitCalc(sat, SMA()),
lower_bounds = [42166.90],
upper_bounds = [42166.90],
scale = [1.0],
)
# Define MOI event struct with event function, solver variables, and constraints
moi_fun() = maneuver!(sat, moi)
moi_event = Event(
name = "MOI",
event = moi_fun,
vars = [moi_var],
funcs = [final_sma_con],
)
# ============================================================================
# Trajectory Optimization
# ============================================================================
seq = Sequence()
add_sequence!(seq, prop_to_z_crossing_1_event, toi_event, prop_to_apogee_event,
prop_to_perigee_1_event, prop_to_z_crossing_2_event, mcc_event,
prop_to_perigee_2_event, moi_event)
# Set up IPOPT options
ipopt_options = Dict(
"max_iter" => 1000,
"tol" => 1e-6,
"output_file" => "ipopt_geo_transfer$(rand(UInt)).out",
"file_print_level" => 5,
"print_level" => 5,
)
snow_options = Options(derivatives=ForwardFD(), solver=IPOPT(ipopt_options))
# Solve trajectory optimization with iteration recording enabled
result = solve_trajectory!(seq, snow_options; record_iterations=true)
report_sequence(seq)
report_solution(seq, result)
# Propagate about one day to see final orbit
propagate!(prop, sat, StopAt(sat, PropDurationDays(), 1.1))
# Visualize with iterations
view = View3D()
add_spacecraft!(view, sat; show_iterations=true)
display_view(view)