API Reference
Complete API documentation for AstroModels.
AstroModels.AstroModels — Module
Module containing physical models such as spacecraft.
AstroModels.CADModel — Type
CADModel(; file_path="", scale=1.0, visible=false)Represents a 3D CAD model for visualization of spacecraft or other objects.
Keyword Arguments
file_path::String = ""— Path to CAD model file (empty string means no model). Ifvisible=true, must not be empty.scale::Float64 = 1.0— Scaling factor for the model. Must be positive (> 0.0).visible::Bool = false— Whether the model should be displayed (default hidden until model is set)
Examples
model = CADModel(file_path="data/SpacecraftCADModel.obj", scale=10.0, visible=true)AstroModels.CADModel — Method
CADModel(; file_path="", scale=1.0, visible=false)Kwarg constructor with defaults for all fields.
AstroModels.HistorySegment — Type
HistorySegmentRepresents a continuous segment of spacecraft trajectory history.
A history segment stores time-ordered state data for a portion of a spacecraft's trajectory. Segments are used to delineate discontinuities (e.g., maneuvers) and to organize mission events (e.g., propagation phases).
Fields
times::Vector{Time}: Ordered vector of time pointsstates::Vector{CartesianState{Float64}}: Corresponding state vectors (always Float64)coordinate_system::CoordinateSystem: Reference frame for this segmentname::String: Optional segment identifier (empty string if unnamed)metadata::Dict{String,Any}: Extensible metadata storage for custom attributes
Notes
- States are always stored as Float64 regardless of input type (e.g., Dual numbers)
- Immutable struct, but internal vectors can be mutated for efficiency
timesandstatesare parallel vectors of equal length
Examples
using AstroModels, AstroUniverse, AstroEpochs, AstroFrames, AstroStates
# Create empty segment
coord_sys = CoordinateSystem(earth, ICRFAxes())
segment = HistorySegment(coord_sys, name="initial_orbit")
# Create segment from existing data
times = [Time("2024-01-01T00:00:00", TAI(), ISOT()),
Time("2024-01-01T01:00:00", TAI(), ISOT())]
states = [CartesianState([7000.0, 0.0, 0.0, 0.0, 7.5, 0.0]),
CartesianState([7100.0, 100.0, 0.0, 0.0, 7.4, 0.1])]
segment = HistorySegment(times, states, coord_sys, name="propagation_1")AstroModels.HistorySegment — Method
HistorySegment(coord_system::CoordinateSystem;
name::String="",
metadata::Dict{String,Any}=Dict{String,Any}())Create an empty history segment with specified coordinate system.
AstroModels.HistorySegment — Method
HistorySegment(times::Vector{<:Time},
states::Vector{CartesianState{Float64}},
coord_system::CoordinateSystem;
name::String="",
metadata::Dict{String,Any}=Dict{String,Any}())Create a history segment from existing time and state data.
AstroModels.Spacecraft — Type
mutable struct Spacecraft{S<:OrbitState, TT<:Time, CS<:AbstractCoordinateSystem, T<:Real}Spacecraft struct with state, time, mass, name, history, and coordinate system and other data
Fields
- state::S — orbital state as an OrbitState struct
- time::TT — epoch as a Time struct
- mass::T — total mass
- name::String — user label.
- history::SpacecraftHistory — trajectory history organized into segments
- coord_sys::CS — coordinate system (origin and axes) associated with the spacecraft.
- cad_model::CADModel — 3D model for visualization
Notes:
- Use the keyword constructor to create spacecraft and only define the fields you want to change from the defaults.
- State can be provided two ways as shown in the example below
- history stores trajectory data in segments; use
history.segmentsto access individual HistorySegments - Numeric parameter T is chosen by promotion: T = promote_type(eltype(state), typeof(time.jd1), typeof(mass)).
Examples
using AstroModels, AstroStates, AstroEpochs, AstroFrames, AstroUniverse
sc = Spacecraft(
state = CartesianState([7000.0, 300.0, 0.0, 0.0, 7.5, 0.03]),
time = Time("2015-09-21T00:00:00", TAI(), ISOT()),
mass = 1000.0,
name = "Demo"
)
# alternative epoch definition
sc = Spacecraft(state = OrbitState([7000.0, 300.0, 0.0, 0.0, 7.5, 0.03],Cartesian()))AstroModels.Spacecraft — Method
Spacecraft(; state = CartesianState([7000.0, 0.0, 0.0, 0.0, 7.5, 0.0]),
time = Time("2015-09-21T12:23:12", UTC(), ISOT()),
mass = 1000.0,
name = "unnamed",
history = nothing,
coord_sys = CoordinateSystem(earth, ICRFAxes()))Kwarg outer constructor for Spacecraft with defaults for all fields.
AstroModels.Spacecraft — Method
Spacecraft(state::Union{AbstractState,OrbitState}, time::TT; mass=1000.0, name="unnamed",
history=nothing, coord_sys=CoordinateSystem(earth, ICRFAxes())) where {TT<:Time}Outer positional constructor for Spacecraft that promotes numeric types as needed.
AstroModels.SpacecraftHistory — Type
SpacecraftHistoryContainer for spacecraft trajectory history, organized into segments.
A SpacecraftHistory stores the complete trajectory of a spacecraft as a collection of segments. Each segment represents a continuous portion of the trajectory, with segments typically separated by discontinuities (e.g., maneuvers) or mission phase boundaries.
Fields
segments::Vector{HistorySegment}: Solution trajectory segmentsiterations::Vector{HistorySegment}: Solver iteration segments (diagnostic data)record_segments::Bool: Enable solution segment recording (default: true)record_iterations::Bool: Enable iteration segment recording (default: false)
Notes
- Mutable struct to allow solver to control recording flags
- Empty history is valid (zero segments)
- Provides iteration, indexing, and standard collection interfaces
- Iteration and indexing operate on
segmentsonly (solution trajectory) - Access
iterationsdirectly for diagnostic data:history.iterations
Examples
using AstroModels, AstroEpochs, AstroUniverse, AstroFrames
# Basic usage - record solution trajectory (default)
history = SpacecraftHistory()
coord_sys = CoordinateSystem(earth, ICRFAxes())
segment1 = HistorySegment(coord_sys, name="orbit_1")
push_segment!(history, segment1)
# Iterate over solution segments
for (i, segment) in enumerate(history)
println("Segment ", i, ": ", segment.name, " (", length(segment.times), " points)")
end
# Access by index (operates on solution segments only)
first_segment = history[1]
println("History contains ", length(history), " solution segments")
# Advanced: Solver iteration recording for diagnostics
# Disable solution recording, enable iteration recording
history.record_segments = false
history.record_iterations = true
# During optimization, segments go to iterations vector
for iter in 1:10
segment = HistorySegment(coord_sys, name="iteration_$iter")
push_segment!(history, segment) # Routes to history.iterations
end
# Switch back for final solution
history.record_segments = true
history.record_iterations = false
final_segment = HistorySegment(coord_sys, name="solution")
push_segment!(history, final_segment) # Routes to history.segments
# Result: 1 solution segment, 10 diagnostic iterations
println(length(history)) # 1 (solution only)
println(length(history.iterations)) # 10 (diagnostic data)AstroModels.SpacecraftHistory — Method
SpacecraftHistory(segments::Vector{HistorySegment})Create a SpacecraftHistory from an existing vector of segments with default recording settings.
Defaults
iterations: Empty vector (no iteration data)record_segments = true: Solution trajectory recording enabledrecord_iterations = false: Iteration recording disabled
AstroModels.SpacecraftHistory — Method
SpacecraftHistory()Create an empty SpacecraftHistory with default settings.
Defaults
record_segments = true: Solution trajectory recording enabledrecord_iterations = false: Iteration recording disabled (diagnostic mode off)
AstroModels._indent_and_print — Method
_indent_and_print(io::IO, obj, prefix::AbstractString)Indent and print composed objects using their own show methods
AstroModels._value_to_float64 — Method
_value_to_float64(x::Real)Helper to safely convert any Real to Float64, handling Dual numbers from AD.
AstroModels.get_state — Method
get_state(sc::Spacecraft, target::AbstractOrbitStateType) -> AbstractOrbitStateReturn the spacecraft's orbital state as the concrete type specified by target, converting if needed. Does not mutate sc.
Notes
- If conversion requires the gravitational parameter μ, it is taken from
sc.coord_sys.origin.muwhen available; otherwise an error is thrown. - If the current state already matches
target, the existing state is returned without conversion to the new type. - State is returned in the coordinate system of the spacecraft; no coordinate transformations are performed.
Examples
using AstroModels
sc = Spacecraft(
state=CartesianState([7000.0, 300.0, 0.0, 0.0, 7.5, 0.03]),
time=Time("2015-09-21T12:23:12", TAI(), ISOT())
);AstroModels.push_segment! — Method
push_segment!(history::SpacecraftHistory, segment::HistorySegment)Add a segment to the history, routing based on recording flags.
Routing Logic
- If
record_iterations=true: adds tohistory.iterations(diagnostic data) - Else if
record_segments=true: adds tohistory.segments(solution trajectory) - Else: no-op (segment not recorded)
Notes
This routing is controlled by the solver during trajectory optimization:
- During optimization iterations:
record_iterations=true, record_segments=false - Final solution run:
record_iterations=false, record_segments=true
AstroModels.set_posvel! — Method
set_posvel!(sc::Spacecraft, x::AbstractVector{<:Real})Set the Cartesian position-velocity vector [x, y, z, vx, vy, vz] for sc in-place.
Supported
OrbitStatewithstatetype == Cartesian().
Limitations
- Other state representations currently throw an ArgumentError. Use
get_state(sc, Cartesian())first.
Examples
using AstroModels, AstroStates, AstroEpochs
sc = Spacecraft(
state=CartesianState([7000.0, 300.0, 0.0, 0.0, 7.5, 0.03]),
time=Time("2015-09-21T12:23:12", TAI(), ISOT())
);
set_posvel!(sc, [7050.0, 0.0, 0.0, 0.0, 7.6, 0.0]);
to_posvel(sc)
# output
6-element Vector{Float64}:
7050.0
0.0
0.0
0.0
7.6
0.0AstroModels.state_eltype — Method
state_eltype(os::OrbitState) = eltype(os.state)Returns the type of the elements in the state vector of an OrbitState.
AstroModels.to_float64 — Method
to_float64(state::CartesianState{Float64})
to_float64(state::CartesianState)Convert a CartesianState of any numeric type to CartesianState{Float64}.
This function handles automatic conversion from automatic differentiation types (Dual numbers), BigFloat, and other numeric types to Float64 for history storage.
Arguments
state::CartesianState: State to convert (any numeric type parameter)
Returns
CartesianState{Float64}: State with Float64 components
Notes
- No-op for CartesianState{Float64} (returns input unchanged)
- Uses
to_vector()and broadcasts Float64 conversion
Examples
# From Dual numbers (automatic differentiation)
state_dual = CartesianState{Dual}([...])
state_f64 = to_float64(state_dual)
# From BigFloat
state_big = CartesianState{BigFloat}([...])
state_f64 = to_float64(state_big)
# From Float64 (no-op)
state = CartesianState([7000.0, 0.0, 0.0, 0.0, 7.5, 0.0])
same_state = to_float64(state) AstroModels.to_float64 — Method
to_float64(time::Time{Float64})
to_float64(time::Time)Convert a Time of any numeric type to Time{Float64}.
This function handles automatic conversion from automatic differentiation types (Dual numbers), BigFloat, and other numeric types to Float64 for history storage.
Arguments
time::Time: Time to convert (any numeric type parameter)
Returns
Time{Float64}: Time with Float64 components
Notes
- No-op for Time{Float64} (returns input unchanged)
- Preserves time scale and format during conversion
Examples
# From Dual numbers (automatic differentiation)
time_dual = Time{Dual}(...)
time_f64 = to_float64(time_dual)
# From Float64 (no-op)
time = Time("2024-01-01T00:00:00", TAI(), ISOT())
same_time = to_float64(time) # Returns input unchangedAstroModels.to_posvel — Method
to_posvel(sc::Spacecraft) -> Vector{<:Real}Return the Cartesian position-velocity vector [x, y, z, vx, vy, vz] for sc in its current coordinate system, converting the stored orbital state if needed.
Currently supported
- OrbitState with
statetype == Cartesian()— returns the internal 6-vector (fast path). - CartesianState — returns
to_vector(state).
Limitations
- Other state types currently throw an ArgumentError. Use
get_state(sc, Cartesian())first or extend conversions.
Examples
using AstroModels, AstroStates, AstroEpochs, AstroFrames, AstroUniverse
sc = Spacecraft(
state=CartesianState([7000.0, 300.0, 0.0, 0.0, 7.5, 0.03]),
time=Time("2015-09-21T12:23:12", TAI(), ISOT()));
to_posvel(sc)
# output
6-element Vector{Float64}:
7000.0
300.0
0.0
0.0
7.5
0.03Base.deepcopy_internal — Method
Base.deepcopy_internal(sc::Spacecraft, dict::IdDict)Deep copy a spacecraft to ensure no aliasing of inner mutable fields
Base.deepcopy_internal — Method
Base.deepcopy_internal(history::SpacecraftHistory, dict::IdDict)Deep copy implementation for SpacecraftHistory to ensure proper copying of mutable fields.
Base.promote — Method
Base.promote(sc::Spacecraft{S,TT,CS,T}, ::Type{Tnew}) where {S,TT,CS,T,Tnew<:Real}Promotes a Spacecraft to a new numeric type Tnew for automatic differentiation support. The state, time, and mass are promoted to Tnew, while history remains as Float64 for efficiency.
This enables AD workflows where computation types (e.g., ForwardDiff.Dual) are promoted while preserving Float64 ephemeris storage.
Arguments
sc::Spacecraft: The spacecraft to promote::Type{Tnew}: Target numeric type (e.g., ForwardDiff.Dual{Nothing,Float64,3})
Returns
Spacecraft{S_new, TT_new, CS, Tnew}: Promoted spacecraft
Example
using ForwardDiff
sc = Spacecraft(state=CartesianState([7000.0, 0.0, 0.0, 0.0, 7.5, 0.0]))
sc_dual = promote(sc, ForwardDiff.Dual{Nothing,Float64,3})Index
AstroModels.AstroModelsAstroModels.CADModelAstroModels.CADModelAstroModels.HistorySegmentAstroModels.HistorySegmentAstroModels.HistorySegmentAstroModels.SpacecraftAstroModels.SpacecraftAstroModels.SpacecraftAstroModels.SpacecraftHistoryAstroModels.SpacecraftHistoryAstroModels.SpacecraftHistoryAstroModels._indent_and_printAstroModels._value_to_float64AstroModels.get_stateAstroModels.push_segment!AstroModels.set_posvel!AstroModels.state_eltypeAstroModels.to_float64AstroModels.to_float64AstroModels.to_posvelBase.copyBase.deepcopy_internalBase.deepcopy_internalBase.promoteBase.showBase.showBase.showBase.show