Skip to content

OEM — Orbit Ephemeris Message

An Orbit Ephemeris Message (OEM) carries time-ordered state vectors for spacecraft ephemeris exchange. The typical workflow is to parse an OEM file and convert it into an OrbitTrajectory for interpolation and analysis, or to generate an OEM from a propagator for distribution.

Parse and Access

Parse from file or string, then access header properties, segment metadata, and state vectors:

import brahe as bh
from brahe.ccsds import OEM

bh.initialize_eop()

# Parse from file (auto-detects KVN, XML, or JSON format)
oem = OEM.from_file("test_assets/ccsds/oem/OEMExample1.txt")

# Header properties
print(f"Format version: {oem.format_version}")
print(f"Originator:     {oem.originator}")
print(f"Classification: {oem.classification}")
print(f"Creation date:  {oem.creation_date}")

# Segments — OEM can contain multiple trajectory arcs
print(f"\nNumber of segments: {len(oem.segments)}")

# Access segment metadata
seg = oem.segments[0]
print("\nSegment 0:")
print(f"  Object name:   {seg.object_name}")
print(f"  Object ID:     {seg.object_id}")
print(f"  Center name:   {seg.center_name}")
print(f"  Ref frame:     {seg.ref_frame}")
print(f"  Time system:   {seg.time_system}")
print(f"  Start time:    {seg.start_time}")
print(f"  Stop time:     {seg.stop_time}")
print(f"  Interpolation: {seg.interpolation}")
print(f"  States:        {seg.num_states}")
print(f"  Covariances:   {seg.num_covariances}")

# Access individual state vectors
sv = seg.states[0]
print("\nFirst state vector:")
print(f"  Epoch:    {sv.epoch}")
print(
    f"  Position: [{sv.position[0]:.3f}, {sv.position[1]:.3f}, {sv.position[2]:.3f}] m"
)
print(
    f"  Velocity: [{sv.velocity[0]:.5f}, {sv.velocity[1]:.5f}, {sv.velocity[2]:.5f}] m/s"
)

# Iterate over all states in a segment
print("\nAll states in segment 0:")
for i, sv in enumerate(seg.states):
    print(
        f"  [{i}] {sv.epoch}  pos=({sv.position[0] / 1e3:.3f}, {sv.position[1] / 1e3:.3f}, {sv.position[2] / 1e3:.3f}) km"
    )

# Serialization
kvn = oem.to_string("KVN")
print(f"\nKVN output length: {len(kvn)} characters")
d = oem.to_dict()
print(f"Dict keys: {list(d.keys())}")
use brahe as bh;
use brahe::ccsds::OEM;

fn main() {
    bh::initialize_eop().unwrap();

    // Parse from file (auto-detects KVN, XML, or JSON format)
    let oem = OEM::from_file("test_assets/ccsds/oem/OEMExample1.txt").unwrap();

    // Header properties
    println!("Format version: {}", oem.header.format_version);
    println!("Originator:     {}", oem.header.originator);
    println!(
        "Classification: {}",
        oem.header.classification.as_deref().unwrap_or("None")
    );
    println!("Creation date:  {}", oem.header.creation_date);
    println!("\nNumber of segments: {}", oem.segments.len());

    // Access segment metadata
    let seg = &oem.segments[0];
    println!("\nSegment 0:");
    println!("  Object name:   {}", seg.metadata.object_name);
    println!("  Object ID:     {}", seg.metadata.object_id);
    println!("  Center name:   {}", seg.metadata.center_name);
    println!("  Ref frame:     {}", seg.metadata.ref_frame);
    println!("  Time system:   {}", seg.metadata.time_system);
    println!("  Start time:    {}", seg.metadata.start_time);
    println!("  Stop time:     {}", seg.metadata.stop_time);
    println!(
        "  Interpolation: {}",
        seg.metadata.interpolation.as_deref().unwrap_or("None")
    );
    println!("  States:        {}", seg.states.len());
    println!("  Covariances:   {}", seg.covariances.len());
    let sv = &seg.states[0];
    println!("\nFirst state vector:");
    println!("  Epoch:    {}", sv.epoch);
    println!(
        "  Position: [{:.3}, {:.3}, {:.3}] m",
        sv.position[0], sv.position[1], sv.position[2]
    );
    println!(
        "  Velocity: [{:.5}, {:.5}, {:.5}] m/s",
        sv.velocity[0], sv.velocity[1], sv.velocity[2]
    );
    println!("\nAll states in segment 0:");
    for (i, sv) in seg.states.iter().enumerate() {
        println!(
            "  [{}] {}  pos=({:.3}, {:.3}, {:.3}) km",
            i,
            sv.epoch,
            sv.position[0] / 1e3,
            sv.position[1] / 1e3,
            sv.position[2] / 1e3
        );
    }
    let kvn = oem.to_string(brahe::ccsds::CCSDSFormat::KVN).unwrap();
    println!("\nKVN output length: {} characters", kvn.len());
}
Output
Format version: 3.0
Originator:     NASA/JPL
Classification: public, test-data
Creation date:  1996-11-04 17:22:31.000 UTC

Number of segments: 3

Segment 0:
  Object name:   MARS GLOBAL SURVEYOR
  Object ID:     1996-062A
  Center name:   MARS BARYCENTER
  Ref frame:     J2000
  Time system:   UTC
  Start time:    1996-12-18 12:00:00.331 UTC
  Stop time:     1996-12-28 21:28:00.331 UTC
  Interpolation: HERMITE
  States:        4
  Covariances:   0

First state vector:
  Epoch:    1996-12-18 12:00:00.331 UTC
  Position: [2789619.000, -280045.000, -1746755.000] m
  Velocity: [4733.72000, -2495.86000, -1041.95000] m/s

All states in segment 0:
  [0] 1996-12-18 12:00:00.331 UTC  pos=(2789.619, -280.045, -1746.755) km
  [1] 1996-12-18 12:01:00.331 UTC  pos=(2783.419, -308.143, -1877.071) km
  [2] 1996-12-18 12:02:00.331 UTC  pos=(2776.033, -336.859, -2008.682) km
  [3] 1996-12-28 21:28:00.331 UTC  pos=(-3881.024, 563.959, -682.773) km

KVN output length: 4372 characters
Dict keys: ['header', 'segments']
Format version: 3
Originator:     NASA/JPL
Classification: public, test-data
Creation date:  1996-11-04 17:22:31.000 UTC

Number of segments: 3

Segment 0:
  Object name:   MARS GLOBAL SURVEYOR
  Object ID:     1996-062A
  Center name:   MARS BARYCENTER
  Ref frame:     J2000
  Time system:   UTC
  Start time:    1996-12-18 12:00:00.331 UTC
  Stop time:     1996-12-28 21:28:00.331 UTC
  Interpolation: HERMITE
  States:        4
  Covariances:   0

First state vector:
  Epoch:    1996-12-18 12:00:00.331 UTC
  Position: [2789619.000, -280045.000, -1746755.000] m
  Velocity: [4733.72000, -2495.86000, -1041.95000] m/s

All states in segment 0:
  [0] 1996-12-18 12:00:00.331 UTC  pos=(2789.619, -280.045, -1746.755) km
  [1] 1996-12-18 12:01:00.331 UTC  pos=(2783.419, -308.143, -1877.071) km
  [2] 1996-12-18 12:02:00.331 UTC  pos=(2776.033, -336.859, -2008.682) km
  [3] 1996-12-28 21:28:00.331 UTC  pos=(-3881.024, 563.959, -682.773) km

KVN output length: 4372 characters

Converting to OrbitTrajectory

The primary interoperability point for OEM data is conversion to brahe's OrbitTrajectory. Each OEM segment maps to a trajectory object, giving you Hermite interpolation at arbitrary epochs within the covered time span:

import brahe as bh
from brahe.ccsds import OEM

bh.initialize_eop()

# Parse an OEM file
oem = OEM.from_file("test_assets/ccsds/oem/OEMExample5.txt")
seg = oem.segments[0]
print(f"Segment: {seg.object_name}, {seg.num_states} states, frame={seg.ref_frame}")

# Convert segment 0 to an OrbitTrajectory
traj = oem.segment_to_trajectory(0)
print(f"\nTrajectory: {len(traj)} states")
print(f"  Frame: {traj.frame}")
print(f"  Start: {traj.start_epoch()}")
print(f"  End:   {traj.end_epoch()}")
print(f"  Span:  {traj.timespan():.0f} seconds")

# Access states by index
epc, state = traj.get(0)
print("\nFirst state:")
print(f"  Epoch: {epc}")
print(
    f"  Position: [{state[0] / 1e3:.3f}, {state[1] / 1e3:.3f}, {state[2] / 1e3:.3f}] km"
)
print(f"  Velocity: [{state[3]:.3f}, {state[4]:.3f}, {state[5]:.3f}] m/s")

# Interpolate at an arbitrary epoch between states
epc0, _ = traj.get(0)
epc1, _ = traj.get(1)
mid_epoch = epc0 + (epc1 - epc0) / 2.0
interp_state = traj.interpolate(mid_epoch)
print(f"\nInterpolated state at {mid_epoch}:")
print(
    f"  Position: [{interp_state[0] / 1e3:.3f}, {interp_state[1] / 1e3:.3f}, {interp_state[2] / 1e3:.3f}] km"
)

# Convert all segments at once
oem_multi = OEM.from_file("test_assets/ccsds/oem/OEMExample1.txt")
trajs = oem_multi.to_trajectories()
print(f"\nMulti-segment OEM: {len(trajs)} trajectories")
for i, t in enumerate(trajs):
    print(f"  [{i}] {len(t)} states, span={t.timespan():.0f}s")
use brahe as bh;
use brahe::ccsds::OEM;
use brahe::traits::Trajectory;

fn main() {
    bh::initialize_eop().unwrap();

    // Parse an OEM file
    let oem = OEM::from_file("test_assets/ccsds/oem/OEMExample5.txt").unwrap();
    let seg = &oem.segments[0];
    println!(
        "Segment: {}, {} states, frame={}",
        seg.metadata.object_name,
        seg.states.len(),
        seg.metadata.ref_frame
    );
    let traj = oem.segment_to_trajectory(0).unwrap();
    println!("\nTrajectory: {} states", traj.len());
    println!("  Frame: {:?}", traj.frame);
    println!("  Start: {}", traj.start_epoch().unwrap());
    println!("  End:   {}", traj.end_epoch().unwrap());
    println!("  Span:  {:.0} seconds", traj.timespan().unwrap());

    // Access states by index
    let (epoch, state) = traj.first().unwrap();
    println!("\nFirst state:");
    println!("  Epoch: {}", epoch);
    println!(
        "  Position: [{:.3}, {:.3}, {:.3}] km",
        state[0] / 1e3,
        state[1] / 1e3,
        state[2] / 1e3
    );
    println!(
        "  Velocity: [{:.3}, {:.3}, {:.3}] m/s",
        state[3], state[4], state[5]
    );

    // Convert all segments from a multi-segment OEM
    let oem_multi = OEM::from_file("test_assets/ccsds/oem/OEMExample1.txt").unwrap();
    let trajs = oem_multi.to_trajectories().unwrap();
    println!("\nMulti-segment OEM: {} trajectories", trajs.len());
    for (i, t) in trajs.iter().enumerate() {
        println!("  [{}] {} states, span={:.0}s", i, t.len(), t.timespan().unwrap());
    }
}
Output
Segment: ISS, 49 states, frame=GCRF

Trajectory: 49 states
  Frame: GCRF
  Start: 2017-04-11 22:31:43.122 UTC
  End:   2017-04-12 22:31:43.122 UTC
  Span:  86400 seconds

First state:
  Epoch: 2017-04-11 22:31:43.122 UTC
  Position: [2906.275, 4076.358, 4561.364] km
  Velocity: [-6879.497, 1449.531, 3081.318] m/s

Interpolated state at 2017-04-11 22:46:43.122 UTC:
  Position: [-1917.595, 1696.341, 2465.366] km

Multi-segment OEM: 3 trajectories
  [0] 4 states, span=898080s
  [1] 4 states, span=100735s
  [2] 5 states, span=100735s
Segment: ISS, 49 states, frame=GCRF

Trajectory: 49 states
  Frame: OrbitFrame(Geocentric Celestial Reference Frame)
  Start: 2017-04-11 22:31:43.122 UTC
  End:   2017-04-12 22:31:43.122 UTC
  Span:  86400 seconds

First state:
  Epoch: 2017-04-11 22:31:43.122 UTC
  Position: [2906.275, 4076.358, 4561.364] km
  Velocity: [-6879.497, 1449.531, 3081.318] m/s

Multi-segment OEM: 3 trajectories
  [0] 4 states, span=898080s
  [1] 4 states, span=100735s
  [2] 5 states, span=100735s

How OEM Messages Are Organized

An OEM message begins with a header that records the format version, creation date, and originator. The bulk of the data lives in one or more segments, each of which has its own metadata block and a sequence of state vectors.

Multiple segments exist because a single file may need to cover different trajectory arcs. A maneuver boundary, a change in reference frame, or a gap in tracking data each warrant a new segment. Within a segment, the metadata block records the object identity, center body, reference frame, time system, time span, and interpolation settings. The state vectors follow — each line provides an epoch plus position and velocity (and optionally acceleration). If covariance data is available, it appears as one or more 6\(\times\)6 symmetric matrices attached to the segment, each with its own epoch and optional reference frame override.

Creating and Writing OEMs

Build an OEM programmatically by defining a header, adding segments with metadata, and populating state vectors. The resulting message can be serialized to KVN, XML, or JSON:

import brahe as bh
import numpy as np
from brahe.ccsds import OEM

bh.initialize_eop()

# Create a new OEM with header info
oem = OEM(originator="BRAHE_EXAMPLE")
oem.classification = "unclassified"
oem.message_id = "OEM-2024-001"

# Define a LEO orbit and propagate with KeplerianPropagator (two-body)
epoch = bh.Epoch.from_datetime(2024, 6, 15, 0, 0, 0.0, 0.0, bh.TimeSystem.UTC)
oe = np.array([bh.R_EARTH + 500e3, 0.001, 51.6, 15.0, 30.0, 0.0])
prop = bh.KeplerianPropagator.from_keplerian(epoch, oe, bh.AngleFormat.DEGREES, 60.0)

# Add a segment with metadata
step = 60.0  # 60-second spacing
n_states = 5
stop_epoch = epoch + step * (n_states - 1)

seg_idx = oem.add_segment(
    object_name="LEO SAT",
    object_id="2024-100A",
    center_name="EARTH",
    ref_frame="EME2000",
    time_system="UTC",
    start_time=epoch,
    stop_time=stop_epoch,
    interpolation="LAGRANGE",
    interpolation_degree=7,
)

# Propagate to build trajectory, then bulk-add states to segment
prop.propagate_to(stop_epoch)
seg = oem.segments[seg_idx]
seg.add_trajectory(prop.trajectory)

print(f"Created OEM with {len(oem.segments)} segment, {seg.num_states} states")

# Write to KVN string
kvn = oem.to_string("KVN")
print(f"\nKVN output ({len(kvn)} chars):")
print(kvn[:500])

# Write to file
oem.to_file("/tmp/brahe_example_oem.txt", "KVN")
print("\nWritten to /tmp/brahe_example_oem.txt")

# Verify round-trip
oem2 = OEM.from_file("/tmp/brahe_example_oem.txt")
print(f"Round-trip: {len(oem2.segments)} segment, {oem2.segments[0].num_states} states")
use brahe as bh;
use brahe::ccsds::{
    CCSDSFormat, CCSDSRefFrame, CCSDSTimeSystem, OEM, OEMMetadata, OEMSegment, OEMStateVector,
};
use brahe::traits::DStateProvider;
use nalgebra as na;

fn main() {
    bh::initialize_eop().unwrap();

    // Create a new OEM with header info
    let mut oem = OEM::new("BRAHE_EXAMPLE".to_string());

    // Define a LEO orbit and propagate with KeplerianPropagator (two-body)
    let epoch =
        bh::Epoch::from_datetime(2024, 6, 15, 0, 0, 0.0, 0.0, bh::TimeSystem::UTC);
    let oe = na::SVector::<f64, 6>::new(bh::R_EARTH + 500e3, 0.001, 51.6, 15.0, 30.0, 0.0);
    let prop = bh::KeplerianPropagator::from_keplerian(epoch, oe, bh::AngleFormat::Degrees, 60.0);

    // Create segment metadata
    let step = 60.0_f64; // 60-second spacing
    let n_states = 5usize;
    let stop_epoch = epoch + step * (n_states - 1) as f64;

    let metadata = OEMMetadata::new(
        "LEO SAT".to_string(),
        "2024-100A".to_string(),
        "EARTH".to_string(),
        CCSDSRefFrame::EME2000,
        CCSDSTimeSystem::UTC,
        epoch,
        stop_epoch,
    )
    .with_interpolation("LAGRANGE".to_string(), Some(7));

    let mut seg = OEMSegment::new(metadata);

    // Populate states from the Keplerian propagator
    for i in 0..n_states {
        let t = epoch + i as f64 * step;
        let s = prop.state(t).unwrap();
        seg.push_state(OEMStateVector::new(
            t,
            [s[0], s[1], s[2]],
            [s[3], s[4], s[5]],
        ));
    }

    oem.push_segment(seg);

    println!(
        "Created OEM with {} segment, {} states",
        oem.segments.len(),
        oem.segments[0].states.len()
    );
    let kvn = oem.to_string(CCSDSFormat::KVN).unwrap();
    println!("\nKVN output ({} chars):", kvn.len());
    let preview: String = kvn.chars().take(500).collect();
    println!("{}", preview);

    // Write to file
    oem.to_file("/tmp/brahe_example_oem.txt", CCSDSFormat::KVN)
        .unwrap();
    println!("\nWritten to /tmp/brahe_example_oem.txt");

    // Verify round-trip
    let oem2 = OEM::from_file("/tmp/brahe_example_oem.txt").unwrap();
    println!(
        "Round-trip: {} segment, {} states",
        oem2.segments.len(),
        oem2.segments[0].states.len()
    );
}
Output
Created OEM with 1 segment, 5 states

KVN output (998 chars):
CCSDS_OEM_VERS = 3.0
CLASSIFICATION = unclassified
CREATION_DATE = 2026-03-23T02:20:37.1934732633
ORIGINATOR = BRAHE_EXAMPLE
MESSAGE_ID = OEM-2024-001

META_START
OBJECT_NAME = LEO SAT
OBJECT_ID = 2024-100A
CENTER_NAME = EARTH
REF_FRAME = EME2000
TIME_SYSTEM = UTC
START_TIME = 2024-06-15T00:00:00.000
STOP_TIME = 2024-06-15T00:04:00.000
INTERPOLATION = LAGRANGE
INTERPOLATION_DEGREE = 7
META_STOP
2024-06-15T00:00:00.000     5195.590962     3601.468794     2692.479496    -4.741222842     2.97334092

Written to /tmp/brahe_example_oem.txt
Round-trip: 1 segment, 5 states
Created OEM with 1 segment, 5 states

KVN output (942 chars):
CCSDS_OEM_VERS = 3.0
CREATION_DATE = 2026-03-23T02:15:49.1149560477
ORIGINATOR = BRAHE_EXAMPLE

META_START
OBJECT_NAME = LEO SAT
OBJECT_ID = 2024-100A
CENTER_NAME = EARTH
REF_FRAME = EME2000
TIME_SYSTEM = UTC
START_TIME = 2024-06-15T00:00:00.000
STOP_TIME = 2024-06-15T00:04:00.000
INTERPOLATION = LAGRANGE
INTERPOLATION_DEGREE = 7
META_STOP
2024-06-15T00:00:00.000     6878.136300        0.000001        0.051600     0.015000000     0.030000000     0.000000000
2024-06-15T00:01:00.000     6878.13630

Written to /tmp/brahe_example_oem.txt
Round-trip: 1 segment, 5 states

Round-Trip Fidelity

Writing and re-parsing an OEM preserves all metadata, state vectors, and covariance data. Numeric precision may vary slightly due to floating-point formatting, but values are preserved within the precision of the output format.

Generating from a Propagator

Propagate an orbit numerically, extract the trajectory, and build an OEM for distribution:

import brahe as bh
import numpy as np
from brahe.ccsds import OEM

bh.initialize_eop()
bh.initialize_sw()

# Define initial state
epoch = bh.Epoch.from_datetime(2024, 6, 15, 0, 0, 0.0, 0.0, bh.TimeSystem.UTC)
oe = np.array([bh.R_EARTH + 500e3, 0.001, 51.6, 15.0, 30.0, 45.0])
state = bh.state_koe_to_eci(oe, bh.AngleFormat.DEGREES)
params = np.array([500.0, 2.0, 2.2, 2.0, 1.3])

# Create propagator with default force model
prop = bh.NumericalOrbitPropagator(
    epoch,
    state,
    bh.NumericalPropagationConfig.default(),
    bh.ForceModelConfig.default(),
    params,
)

# Propagate for 90 minutes
target_epoch = epoch + 5400.0
prop.propagate_to(target_epoch)
print(f"Propagated from {epoch} to {prop.current_epoch()}")

# Get the accumulated trajectory
traj = prop.trajectory
print(f"Trajectory: {len(traj)} states, span={traj.timespan():.0f}s")

# Build an OEM from the trajectory states using the trajectory kwarg
oem = OEM(originator="BRAHE_PROP")
stop_epoch = prop.current_epoch()
seg_idx = oem.add_segment(
    object_name="LEO SAT",
    object_id="2024-100A",
    center_name="EARTH",
    ref_frame="EME2000",
    time_system="UTC",
    start_time=epoch,
    stop_time=stop_epoch,
    interpolation="LAGRANGE",
    interpolation_degree=7,
    trajectory=traj,
)
seg = oem.segments[seg_idx]

print(f"\nOEM: {len(oem.segments)} segment, {seg.num_states} states")

# Write to KVN
kvn = oem.to_string("KVN")
print(f"KVN output: {len(kvn)} characters")

# Verify by re-parsing
oem2 = OEM.from_str(kvn)
print(f"Round-trip: {oem2.segments[0].num_states} states")
use brahe as bh;
use bh::ccsds::{
    CCSDSFormat, CCSDSRefFrame, CCSDSTimeSystem, OEM, OEMMetadata, OEMSegment, OEMStateVector,
};
use bh::traits::{DStatePropagator, Trajectory};
use nalgebra as na;

fn main() {
    bh::initialize_eop().unwrap();
    bh::initialize_sw().unwrap();

    // Define initial state
    let epoch = bh::Epoch::from_datetime(2024, 6, 15, 0, 0, 0.0, 0.0, bh::TimeSystem::UTC);
    let oe = na::SVector::<f64, 6>::new(bh::R_EARTH + 500e3, 0.001, 51.6, 15.0, 30.0, 45.0);
    let state = bh::state_koe_to_eci(oe, bh::AngleFormat::Degrees);
    let params = na::DVector::from_vec(vec![500.0, 2.0, 2.2, 2.0, 1.3]);

    // Create propagator with default force model
    let mut prop = bh::DNumericalOrbitPropagator::new(
        epoch,
        na::DVector::from_column_slice(state.as_slice()),
        bh::NumericalPropagationConfig::default(),
        bh::ForceModelConfig::default(),
        Some(params),
        None,
        None,
        None,
    )
    .unwrap();

    // Propagate for 90 minutes
    let target_epoch = epoch + 5400.0;
    prop.propagate_to(target_epoch);
    println!("Propagated from {} to {}", epoch, prop.current_epoch());

    // Get the accumulated trajectory
    let traj = prop.trajectory();
    println!(
        "Trajectory: {} states, span={:.0}s",
        traj.len(),
        traj.timespan().unwrap()
    );

    // Build an OEM from the trajectory states
    let mut oem = OEM::new("BRAHE_PROP".to_string());
    let stop_epoch = prop.current_epoch();

    let metadata = OEMMetadata::new(
        "LEO SAT".to_string(),
        "2024-100A".to_string(),
        "EARTH".to_string(),
        CCSDSRefFrame::EME2000,
        CCSDSTimeSystem::UTC,
        epoch,
        stop_epoch,
    )
    .with_interpolation("LAGRANGE".to_string(), Some(7));

    let mut seg = OEMSegment::new(metadata);

    // Extract states from trajectory and add to OEM
    for i in 0..traj.len() {
        let (epc, s) = traj.get(i).unwrap();
        seg.push_state(OEMStateVector::new(
            epc,
            [s[0], s[1], s[2]],
            [s[3], s[4], s[5]],
        ));
    }

    let num_states = seg.states.len();
    oem.push_segment(seg);
    println!("\nOEM: {} segment, {} states", oem.segments.len(), num_states);

    // Write to KVN
    let kvn = oem.to_string(CCSDSFormat::KVN).unwrap();
    println!("KVN output: {} characters", kvn.len());

    // Verify by re-parsing
    let oem2 = OEM::from_str(&kvn).unwrap();
    println!("Round-trip: {} states", oem2.segments[0].states.len());
}
Output
1
2
3
4
5
6
Propagated from 2024-06-15 00:00:00.000 UTC to 2024-06-15 01:30:00.000 UTC
Trajectory: 91 states, span=5400s

OEM: 1 segment, 91 states
KVN output: 11259 characters
Round-trip: 91 states
1
2
3
4
5
6
Propagated from 2024-06-15 00:00:00.000 UTC to 2024-06-15 01:30:00.000 UTC
Trajectory: 91 states, span=5400s

OEM: 1 segment, 91 states
KVN output: 11259 characters
Round-trip: 91 states

KVN Format Example

A minimal OEM KVN file looks like:

CCSDS_OEM_VERS = 3.0
CREATION_DATE = 2024-01-15T00:00:00
ORIGINATOR = BRAHE

META_START
OBJECT_NAME = MY SATELLITE
OBJECT_ID = 2024-001A
CENTER_NAME = EARTH
REF_FRAME = EME2000
TIME_SYSTEM = UTC
START_TIME = 2024-01-15T00:00:00
STOP_TIME = 2024-01-15T01:00:00
META_STOP

2024-01-15T00:00:00  6878.137  0.000  0.000  0.000  7.612  0.000
2024-01-15T00:30:00  -3439.068  5957.355  0.000  -6.593  -3.806  0.000
2024-01-15T01:00:00  -3439.068  -5957.355  0.000  6.593  -3.806  0.000

The data lines contain epoch followed by position (km) and velocity (km/s), space-separated.


See Also