Skip to content

Event Detection

The numerical propagator includes an event detection system that identifies specific orbital conditions during propagation. Events are defined by user-configurable detectors that monitor the spacecraft state and trigger when certain criteria are met. They can also be coupled with event callbacks to respond to detected events in real-time.

When an event is detected, the propagator uses a bisection algorithm to precisely locate the event time within a specified tolerance. The detected events are logged and can be accessed after propagation. Users can also configure how an event will affect the propagation, such as stopping propagation or continuing without interruption.

Events provide an extensible mechanism for implementing complex mission scenarios, such as maneuver execution, autonomous operations, and other condition-based actions.

The library also provides a set of premade event detectors for common scenarios, which can be used directly or serve as templates for custom detectors. You can find more details about premade events in the Premade Events documentation with a complete list of available types in the library API docuementation at Premade Event Detectors.

Event Types

Brahe provides three event fundamental event detector types:

Type Trigger Condition
TimeEvent Specific epoch reached
ValueEvent Computed quantity crosses a given value
BinaryEvent Boolean condition changes

Adding Event Detectors

Event detectors are added to the propagator before propagation:

1
2
3
4
prop = bh.NumericalOrbitPropagator(...)
prop.add_event_detector(event1)
prop.add_event_detector(event2)
prop.propagate_to(end_epoch)

Multiple detectors can be active simultaneously.

Accessing Event Results

After propagation, detected events are available via the event log:

1
2
3
4
events = prop.event_log()
for event in events:
    print(f"Event '{event.name}' at {event.window_open}")
    print(f"  State: {event.entry_state}")

Each event record contains:

Field Description
name Event detector name
window_open Epoch when event occurred
window_close Same as window_open for instantaneous events
entry_state State vector at event time

Time Events

Time events trigger at specific epochs. They're useful for scheduled operations like data collection windows, communication passes, or timed maneuvers.

import numpy as np
import brahe as bh

# Initialize EOP and space weather data (required for NRLMSISE-00 drag model)
bh.initialize_eop()
bh.initialize_sw()

# Create initial epoch and state
epoch = bh.Epoch.from_datetime(2024, 1, 1, 12, 0, 0.0, 0.0, bh.TimeSystem.UTC)
oe = np.array([bh.R_EARTH + 500e3, 0.001, 97.8, 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
prop = bh.NumericalOrbitPropagator(
    epoch,
    state,
    bh.NumericalPropagationConfig.default(),
    bh.ForceModelConfig.default(),
    params,
)

# Add time events at specific epochs
event_30min = bh.TimeEvent(epoch + 1800.0, "30-minute mark")
event_1hr = bh.TimeEvent(epoch + 3600.0, "1-hour mark")

# Add a terminal event that stops propagation
event_terminal = bh.TimeEvent(epoch + 5400.0, "90-minute stop").set_terminal()

prop.add_event_detector(event_30min)
prop.add_event_detector(event_1hr)
prop.add_event_detector(event_terminal)

# Propagate for 2 hours (will stop at 90 minutes due to terminal event)
prop.propagate_to(epoch + 7200.0)

# Check detected events
events = prop.event_log()
print(f"Detected {len(events)} events:")
for event in events:
    dt = event.window_open - epoch
    print(f"  '{event.name}' at t+{dt:.1f}s")

# Verify propagation stopped at terminal event
final_time = prop.current_epoch - epoch
print(f"\nPropagation stopped at: t+{final_time:.1f}s (requested: t+7200s)")

# Validate
assert len(events) == 3  # All three events detected
assert abs(final_time - 5400.0) < 1.0  # Stopped at 90 min

print("\nExample validated successfully!")
use brahe as bh;
use bh::events::DTimeEvent;
use bh::traits::DStatePropagator;
use nalgebra as na;

fn main() {
    // Initialize EOP and space weather data (required for NRLMSISE-00 drag model)
    bh::initialize_eop().unwrap();
    bh::initialize_sw().unwrap();

    // Create initial epoch and state
    let epoch = bh::Epoch::from_datetime(2024, 1, 1, 12, 0, 0.0, 0.0, bh::TimeSystem::UTC);
    let oe = na::SVector::<f64, 6>::new(
        bh::R_EARTH + 500e3,
        0.001,
        97.8,
        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
    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();

    // Add time events at specific epochs
    let event_30min = DTimeEvent::new(epoch + 1800.0, "30-minute mark".to_string());
    let event_1hr = DTimeEvent::new(epoch + 3600.0, "1-hour mark".to_string());

    // Add a terminal event that stops propagation
    let event_terminal =
        DTimeEvent::new(epoch + 5400.0, "90-minute stop".to_string()).set_terminal();

    prop.add_event_detector(Box::new(event_30min));
    prop.add_event_detector(Box::new(event_1hr));
    prop.add_event_detector(Box::new(event_terminal));

    // Propagate for 2 hours (will stop at 90 minutes due to terminal event)
    prop.propagate_to(epoch + 7200.0);

    // Check detected events
    let events = prop.event_log();
    println!("Detected {} events:", events.len());
    for event in events {
        let dt = event.window_open - epoch;
        println!("  '{}' at t+{:.1}s", event.name, dt);
    }

    // Verify propagation stopped at terminal event
    let final_time = prop.current_epoch() - epoch;
    println!(
        "\nPropagation stopped at: t+{:.1}s (requested: t+7200s)",
        final_time
    );

    // Validate
    assert_eq!(events.len(), 3); // All three events detected
    assert!((final_time - 5400.0).abs() < 1.0); // Stopped at 90 min

    println!("\nExample validated successfully!");
}

Value Events

Value events trigger when a user-defined function crosses a value value. This is the most flexible event type, enabling detection of arbitrary orbital conditions.

Value events are defined with a value function which accepts the current epoch and state vector, returning a scalar value.

Event Direction

Value events can detect:

  • INCREASING - Value crosses value from below
  • DECREASING - Value crosses value from above
  • ANY - Any value crossing

Custom Value Functions

The value function receives the current epoch and state vector, returning a scalar value:

import numpy as np
import brahe as bh

# Initialize EOP and space weather data (required for NRLMSISE-00 drag model)
bh.initialize_eop()
bh.initialize_sw()

# Create initial epoch and state
epoch = bh.Epoch.from_datetime(2024, 1, 1, 12, 0, 0.0, 0.0, bh.TimeSystem.UTC)
oe = np.array([bh.R_EARTH + 500e3, 0.01, 45.0, 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
prop = bh.NumericalOrbitPropagator(
    epoch,
    state,
    bh.NumericalPropagationConfig.default(),
    bh.ForceModelConfig.default(),
    params,
)


# Define custom value function: detect when z-component crosses zero
# This detects equator crossings (ascending and descending node)
def z_position(epoch, state):
    """Return z-component of position (meters)."""
    return state[2]


# Create ValueEvent: detect when z crosses 0 (equator crossing)
# Ascending node: z goes from negative to positive (INCREASING)
ascending_node = bh.ValueEvent(
    "Ascending Node",
    z_position,
    0.0,  # target value
    bh.EventDirection.INCREASING,
)

# Descending node: z goes from positive to negative (DECREASING)
descending_node = bh.ValueEvent(
    "Descending Node",
    z_position,
    0.0,
    bh.EventDirection.DECREASING,
)

prop.add_event_detector(ascending_node)
prop.add_event_detector(descending_node)

# Propagate for 3 orbits
orbital_period = bh.orbital_period(oe[0])
prop.propagate_to(epoch + 3 * orbital_period)

# Check detected events
events = prop.event_log()
ascending = [e for e in events if "Ascending" in e.name]
descending = [e for e in events if "Descending" in e.name]

print("Equator crossings over 3 orbits:")
print(f"  Ascending nodes: {len(ascending)}")
print(f"  Descending nodes: {len(descending)}")

for event in events[:6]:  # Show first 6
    dt = event.window_open - epoch
    z = event.entry_state[2]
    print(f"  '{event.name}' at t+{dt:.1f}s (z={z:.1f} m)")

# Validate
assert len(ascending) == 3  # One per orbit
assert len(descending) == 3  # One per orbit

print("\nExample validated successfully!")
use brahe as bh;
use bh::events::{DValueEvent, EventDirection};
use bh::traits::DStatePropagator;
use nalgebra as na;
use std::f64::consts::PI;

fn main() {
    // Initialize EOP and space weather data (required for NRLMSISE-00 drag model)
    bh::initialize_eop().unwrap();
    bh::initialize_sw().unwrap();

    // Create initial epoch and state
    let epoch = bh::Epoch::from_datetime(2024, 1, 1, 12, 0, 0.0, 0.0, bh::TimeSystem::UTC);
    let oe = na::SVector::<f64, 6>::new(
        bh::R_EARTH + 500e3,
        0.01,
        45.0,
        0.0,
        0.0,
        0.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
    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();

    // Create ValueEvent: detect when z crosses 0 (equator crossing)
    // Ascending node: z goes from negative to positive (INCREASING)
    let ascending_node_fn =
        |_t: bh::Epoch, state: &na::DVector<f64>, _params: Option<&na::DVector<f64>>| state[2];
    let ascending_node = DValueEvent::new(
        "Ascending Node",
        ascending_node_fn, // z-component
        0.0,               // target value
        EventDirection::Increasing,
    );

    // Descending node: z goes from positive to negative (DECREASING)
    let descending_node_fn =
        |_t: bh::Epoch, state: &na::DVector<f64>, _params: Option<&na::DVector<f64>>| state[2];
    let descending_node = DValueEvent::new(
        "Descending Node",
        descending_node_fn,
        0.0,
        EventDirection::Decreasing,
    );

    prop.add_event_detector(Box::new(ascending_node));
    prop.add_event_detector(Box::new(descending_node));

    // Propagate for 3 orbits
    let orbital_period = 2.0 * PI * (oe[0].powi(3) / bh::GM_EARTH).sqrt();
    prop.propagate_to(epoch + 3.0 * orbital_period);

    // Check detected events
    let events = prop.event_log();
    let ascending: Vec<_> = events
        .iter()
        .filter(|e| e.name.contains("Ascending"))
        .collect();
    let descending: Vec<_> = events
        .iter()
        .filter(|e| e.name.contains("Descending"))
        .collect();

    println!("Equator crossings over 3 orbits:");
    println!("  Ascending nodes: {}", ascending.len());
    println!("  Descending nodes: {}", descending.len());

    for event in events.iter().take(6) {
        let dt = event.window_open - epoch;
        let z = event.entry_state[2];
        println!("  '{}' at t+{:.1}s (z={:.1} m)", event.name, dt, z);
    }

    // Validate
    assert_eq!(ascending.len(), 3); // One per orbit
    assert_eq!(descending.len(), 3); // One per orbit

    println!("\nExample validated successfully!");
}

Binary Events

Binary events detect when a boolean condition transitions between true and false. They use EdgeType to specify which transition to detect:

  • RISING_EDGE - Condition becomes true (false → true)
  • FALLING_EDGE - Condition becomes false (true → false)
  • ANY_EDGE - Either transition

The binary condition function receives the current epoch and state vector, returning a boolean value.

import numpy as np
import brahe as bh

# Initialize EOP data
bh.initialize_eop()

# Create initial epoch and state
epoch = bh.Epoch.from_datetime(2024, 1, 1, 12, 0, 0.0, 0.0, bh.TimeSystem.UTC)
oe = np.array([bh.R_EARTH + 500e3, 0.01, 45.0, 15.0, 30.0, 45.0])
state = bh.state_koe_to_eci(oe, bh.AngleFormat.DEGREES)


# Define condition function: check if spacecraft is in northern hemisphere
def in_northern_hemisphere(epoch, state):
    """Returns True if z-position is positive (northern hemisphere)."""
    return state[2] > 0


# Create binary events for hemisphere crossings
# Rising edge: false → true (entering northern hemisphere)
enter_north = bh.BinaryEvent(
    "Enter Northern",
    in_northern_hemisphere,
    bh.EdgeType.RISING_EDGE,
)

# Falling edge: true → false (exiting northern hemisphere)
exit_north = bh.BinaryEvent(
    "Exit Northern",
    in_northern_hemisphere,
    bh.EdgeType.FALLING_EDGE,
)

# Create propagator
prop = bh.NumericalOrbitPropagator(
    epoch,
    state,
    bh.NumericalPropagationConfig.default(),
    bh.ForceModelConfig.two_body(),
    None,
)

prop.add_event_detector(enter_north)
prop.add_event_detector(exit_north)

# Propagate for 2 orbits
orbital_period = bh.orbital_period(oe[0])
prop.propagate_to(epoch + 2 * orbital_period)

# Check detected events
events = prop.event_log()
enters = [e for e in events if "Enter" in e.name]
exits = [e for e in events if "Exit" in e.name]

print("Hemisphere crossings over 2 orbits:")
print(f"  Entered northern: {len(enters)} times")
print(f"  Exited northern:  {len(exits)} times")

print("\nEvent timeline:")
for event in events[:8]:  # First 8 events
    dt = event.window_open - epoch
    z_km = event.entry_state[2] / 1e3
    print(f"  t+{dt:7.1f}s: {event.name:16} (z = {z_km:+.1f} km)")

# Validate
assert len(enters) == 2  # Once per orbit
assert len(exits) == 2  # Once per orbit

print("\nExample validated successfully!")
use brahe as bh;
use bh::events::{DBinaryEvent, EdgeType};
use bh::traits::DStatePropagator;
use nalgebra as na;
use std::f64::consts::PI;

fn main() {
    // Initialize EOP data
    bh::initialize_eop().unwrap();

    // Create initial epoch and state
    let epoch = bh::Epoch::from_datetime(2024, 1, 1, 12, 0, 0.0, 0.0, bh::TimeSystem::UTC);
    let oe = na::SVector::<f64, 6>::new(bh::R_EARTH + 500e3, 0.01, 45.0, 0.0, 0.0, 0.0);
    let state = bh::state_koe_to_eci(oe, bh::AngleFormat::Degrees);

    // Create binary events for hemisphere crossings
    // Rising edge: false → true (entering northern hemisphere)
    let enter_north = DBinaryEvent::new(
        "Enter Northern",
        |_t, state: &na::DVector<f64>, _params| {
            // Returns True if z-position is positive (northern hemisphere)
            state[2] > 0.0
        },
        EdgeType::RisingEdge,
    );

    // Falling edge: true → false (exiting northern hemisphere)
    let exit_north = DBinaryEvent::new(
        "Exit Northern",
        |_t, state: &na::DVector<f64>, _params| state[2] > 0.0,
        EdgeType::FallingEdge,
    );

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

    prop.add_event_detector(Box::new(enter_north));
    prop.add_event_detector(Box::new(exit_north));

    // Propagate for 2 orbits
    let orbital_period = 2.0 * PI * (oe[0].powi(3) / bh::GM_EARTH).sqrt();
    prop.propagate_to(epoch + 2.0 * orbital_period);

    // Check detected events
    let events = prop.event_log();
    let enters: Vec<_> = events.iter().filter(|e| e.name.contains("Enter")).collect();
    let exits: Vec<_> = events.iter().filter(|e| e.name.contains("Exit")).collect();

    println!("Hemisphere crossings over 2 orbits:");
    println!("  Entered northern: {} times", enters.len());
    println!("  Exited northern:  {} times", exits.len());

    println!("\nEvent timeline:");
    for event in events.iter().take(8) {
        let dt = event.window_open - epoch;
        let z_km = event.entry_state[2] / 1e3;
        println!(
            "  t+{:7.1}s: {:16} (z = {:+.1} km)",
            dt, event.name, z_km
        );
    }

    // Validate
    assert_eq!(enters.len(), 2); // Once per orbit
    assert_eq!(exits.len(), 2); // Once per orbit

    println!("\nExample validated successfully!");
}

See Also