The ForceModelConfig (Python) / ForceModelConfig (Rust) defines which physical forces affect the spacecraft during propagation. Brahe provides preset configurations for common scenarios and allows custom configurations for specific requirements.
importbraheasbh# Create a fully-configured force modelforce_config=bh.ForceModelConfig(# Gravity: Spherical harmonic model (EGM2008, 20x20 degree/order)gravity=bh.GravityConfiguration.spherical_harmonic(degree=20,order=20,model_type=bh.GravityModelType.EGM2008_360,),# Atmospheric drag: Harris-Priester model with parameter indicesdrag=bh.DragConfiguration(model=bh.AtmosphericModel.HARRIS_PRIESTER,area=bh.ParameterSource.parameter_index(1),# Index into parameter vectorcd=bh.ParameterSource.parameter_index(2),),# Solar radiation pressure: Conical eclipse modelsrp=bh.SolarRadiationPressureConfiguration(area=bh.ParameterSource.parameter_index(3),cr=bh.ParameterSource.parameter_index(4),eclipse_model=bh.EclipseModel.CONICAL,),# Third-body: Sun and Moon with DE440s ephemeristhird_body=bh.ThirdBodyConfiguration(ephemeris_source=bh.EphemerisSource.DE440s,bodies=[bh.ThirdBody.SUN,bh.ThirdBody.MOON],),# General relativistic correctionsrelativity=True,# Spacecraft mass (can also use parameter_index for estimation)mass=bh.ParameterSource.value(1000.0),# kg)print(f"Gravity: {force_config.gravity}")print(f"Drag: {force_config.drag}")print(f"SRP: {force_config.srp}")print(f"Third-body: {force_config.third_body}")print(f"Relativity: {force_config.relativity}")print(f"Mass: {force_config.mass}")# Gravity: SphericalHarmonic(source=EGM2008_360, degree=20, order=20)# Drag: DragConfiguration(model=HarrisPriester, area=ParameterIndex(1), cd=ParameterIndex(2))# SRP: SolarRadiationPressureConfiguration(area=ParameterIndex(3), cr=ParameterIndex(4), eclipse_model=Conical)# Third-body: ThirdBodyConfiguration(ephemeris_source=DE440s, bodies=[Sun, Moon])# Relativity: True# Mass: Value(1000.0)
Each sub-configuration is optional (None disables that force). The configuration is captured at propagator construction time and remains immutable during propagation.
Spacecraft parameters (mass \(m\), drag area \(A_d\), coefficient of drag \(C_d\), SRP area \(A_{SRP}\), coefficient of reflectivity \(C_r\)) can be specified in two ways via ParameterSource:
Value(f64) - Fixed constant embedded at construction. The value is baked into the dynamics function and cannot change during propagation.
ParameterIndex(usize) - Index into an parameter vector. This allows parameters to be varied or estimated as part of orbit determination or sensitivity analysis.
importbraheasbh# Point mass gravity is configured using GravityConfiguration.point_mass()# This uses only central body gravity (mu/r^2) - no spherical harmonics# Create point mass gravity configurationgravity=bh.GravityConfiguration.point_mass()# Use two_body() preset which includes point mass gravityforce_config=bh.ForceModelConfig.two_body()
usebraheasbh;fnmain(){// Point mass gravity configuration// Uses only central body gravity (mu/r^2)// No spherical harmonics, J2, or higher-order termslet_force_config=bh::ForceModelConfig{gravity:bh::GravityConfiguration::PointMass,drag:None,srp:None,third_body:None,relativity:false,mass:None,};}
Spherical Harmonics: High-fidelity gravity using EGM2008, GGM05S, or user-defined .gfc model. Degree and order control accuracy vs computation time.
importbraheasbh# ==============================================================================# Packaged Gravity Models# ==============================================================================# EGM2008 - High-fidelity NGA model (360x360 max)gravity_egm2008=bh.GravityConfiguration.spherical_harmonic(degree=20,order=20,model_type=bh.GravityModelType.EGM2008_360)# GGM05S - GRACE mission model (180x180 max)gravity_ggm05s=bh.GravityConfiguration.spherical_harmonic(degree=20,order=20,model_type=bh.GravityModelType.GGM05S)# JGM3 - Legacy model, fast computation (70x70 max)gravity_jgm3=bh.GravityConfiguration.spherical_harmonic(degree=20,order=20,model_type=bh.GravityModelType.JGM3)# ==============================================================================# Custom Model from File# ==============================================================================# Load custom gravity model from GFC format file# GravityModelType.from_file validates the path existscustom_model_type=bh.GravityModelType.from_file("data/gravity_models/EGM2008_360.gfc")gravity_custom=bh.GravityConfiguration.spherical_harmonic(degree=20,order=20,model_type=custom_model_type)
usebraheasbh;usebh::GravityModelType;fnmain(){// ==========================================================================// Packaged Gravity Models// ==========================================================================// EGM2008 - High-fidelity NGA model (360x360 max)let_gravity_egm2008=bh::GravityConfiguration::SphericalHarmonic{source:bh::GravityModelSource::ModelType(GravityModelType::EGM2008_360),degree:20,order:20,};// GGM05S - GRACE mission model (180x180 max)let_gravity_ggm05s=bh::GravityConfiguration::SphericalHarmonic{source:bh::GravityModelSource::ModelType(GravityModelType::GGM05S),degree:20,order:20,};// JGM3 - Legacy model, fast computation (70x70 max)let_gravity_jgm3=bh::GravityConfiguration::SphericalHarmonic{source:bh::GravityModelSource::ModelType(GravityModelType::JGM3),degree:20,order:20,};// ==========================================================================// Custom Model from File// ==========================================================================// Load custom gravity model from GFC format file// GravityModelType::from_file validates the path existsletcustom_model_type=GravityModelType::from_file("data/gravity_models/EGM2008_360.gfc").unwrap();let_gravity_custom=bh::GravityConfiguration::SphericalHarmonic{source:bh::GravityModelSource::ModelType(custom_model_type),degree:20,order:20,};}
where \(\rho\) is atmospheric density, \(v_{rel}\) is velocity relative to the atmosphere, \(C_D\) is drag coefficient, and \(A/m\) is area-to-mass ratio.
Three atmospheric models are available:
Harris-Priester: Fast model with diurnal density variations. Valid 100-1000 km altitude. No space weather data required.
importbraheasbh# Harris-Priester atmospheric drag configuration# - Valid for altitudes 100-1000 km# - Accounts for latitude-dependent diurnal bulge# - Does not require space weather data (F10.7, Ap)# Create drag configuration using parameter indices (default layout)drag_config=bh.DragConfiguration(model=bh.AtmosphericModel.HARRIS_PRIESTER,area=bh.ParameterSource.parameter_index(1),# drag_area from params[1]cd=bh.ParameterSource.parameter_index(2),# Cd from params[2])
usebraheasbh;fnmain(){// Harris-Priester atmospheric drag configuration// - Valid for altitudes 100-1000 km// - Accounts for latitude-dependent diurnal bulge// - Does not require space weather data (F10.7, Ap)// Using parameter indices (default layout)let_drag_config=bh::DragConfiguration{model:bh::AtmosphericModel::HarrisPriester,area:bh::ParameterSource::ParameterIndex(1),// drag_area from params[1]cd:bh::ParameterSource::ParameterIndex(2),// Cd from params[2]};}
NRLMSISE-00: High-fidelity empirical model using space weather data. Valid from ground to thermosphere (~1000 km). More computationally intensive.
importbraheasbh# Initialize space weather data providerbh.initialize_sw()# NRLMSISE-00 atmospheric drag configuration# - Naval Research Laboratory Mass Spectrometer and Incoherent Scatter Radar# - High-fidelity empirical model# - Valid from ground to thermospheric heights# - Uses space weather data (F10.7, Ap) when available# - More computationally expensive than Harris-Priester# Create drag configuration with NRLMSISE-00drag_config=bh.DragConfiguration(model=bh.AtmosphericModel.NRLMSISE00,area=bh.ParameterSource.parameter_index(1),# drag_area from params[1]cd=bh.ParameterSource.parameter_index(2),# Cd from params[2])
usebraheasbh;fnmain(){// Initialize space weather data providerbh::initialize_sw().unwrap();// NRLMSISE-00 atmospheric drag configuration// - Naval Research Laboratory Mass Spectrometer and Incoherent Scatter Radar// - High-fidelity empirical model// - Valid from ground to thermospheric heights// - Uses space weather data (F10.7, Ap) when available// - More computationally expensive than Harris-Priesterlet_drag_config=bh::DragConfiguration{model:bh::AtmosphericModel::NRLMSISE00,area:bh::ParameterSource::ParameterIndex(1),// drag_area from params[1]cd:bh::ParameterSource::ParameterIndex(2),// Cd from params[2]};}
Exponential: An expontential atmospheric density model defined by which provides a simple approximation that is fast for rough calculations:
\[ \rho(h) = \rho_0 e^{-\frac{h-h_0}{H}} \]
\(\rho_0\) is reference density at altitude \(h_0\) and \(H\) is scale height.
importbraheasbh# Create exponential atmospheric modelexp_model=bh.AtmosphericModel.exponential(scale_height=53000.0,# Scale height H in meters (53 km for ~300 km altitude)rho0=1.225e-11,# Reference density at h0 in kg/m^3h0=300000.0,# Reference altitude in meters (300 km))# Create drag configuration with exponential modeldrag_config=bh.DragConfiguration(model=exp_model,area=bh.ParameterSource.parameter_index(1),cd=bh.ParameterSource.parameter_index(2),)
usebraheasbh;usebh::GravityModelType;fnmain(){letdrag_config=bh::DragConfiguration{model:bh::AtmosphericModel::Exponential{scale_height:53000.0,// Scale height H in meters (53 km for ~300 km altitude)rho0:1.225e-11,// Reference density at h0 in kg/m^3h0:300000.0,// Reference altitude in meters (300 km)},area:bh::ParameterSource::ParameterIndex(1),cd:bh::ParameterSource::ParameterIndex(2),};// Create force model with exponential draglet_force_config=bh::ForceModelConfig{gravity:bh::GravityConfiguration::SphericalHarmonic{source:bh::GravityModelSource::ModelType(GravityModelType::EGM2008_360),degree:20,order:20,},drag:Some(drag_config),srp:None,third_body:None,relativity:false,mass:Some(bh::ParameterSource::ParameterIndex(0)),};}
where \(P_{\odot} \approx 4.56 \times 10^{-6}\) N/m² is solar pressure at 1 AU, \(C_R\) is reflectivity coefficient, \(\nu\) is shadow function (0-1), and \(\mathbf{r}_{\odot}\) is the Sun position vector.
Eclipse models determine shadow effects:
None: Always illuminated (fast, inaccurate in shadow)
Cylindrical: Sharp shadow boundary (simple, fast)
Conical: Penumbra and umbra regions (most accurate)
importbraheasbh# Solar Radiation Pressure configuration# Parameters:# - area: Cross-sectional area facing the Sun (m^2)# - cr: Coefficient of reflectivity (1.0=absorbing to 2.0=perfectly reflecting)# - eclipse_model: How to handle Earth's shadow# Option 1: No eclipse model (always illuminated)# Fast but inaccurate during eclipse periodssrp_cylindrical=bh.SolarRadiationPressureConfiguration(area=bh.ParameterSource.parameter_index(3),# srp_area from params[3]cr=bh.ParameterSource.parameter_index(4),# Cr from params[4]eclipse_model=bh.EclipseModel.NONE,)# Option 2: Cylindrical shadow model# Simple and fast, sharp shadow boundary (no penumbra)srp_cylindrical=bh.SolarRadiationPressureConfiguration(area=bh.ParameterSource.parameter_index(3),# srp_area from params[3]cr=bh.ParameterSource.parameter_index(4),# Cr from params[4]eclipse_model=bh.EclipseModel.CYLINDRICAL,)# Option 2: Conical shadow model (recommended)# Accounts for penumbra and umbra regionssrp_conical=bh.SolarRadiationPressureConfiguration(area=bh.ParameterSource.parameter_index(3),cr=bh.ParameterSource.parameter_index(4),eclipse_model=bh.EclipseModel.CONICAL,)
usebraheasbh;fnmain(){// Solar Radiation Pressure configuration// Parameters:// - area: Cross-sectional area facing the Sun (m^2)// - cr: Coefficient of reflectivity (1.0=absorbing to 2.0=perfectly reflecting)// - eclipse_model: How to handle Earth's shadow// Option 1: No eclipse model (always illuminated)// Fast but inaccurate during eclipse periodslet_srp_no_eclipse=bh::SolarRadiationPressureConfiguration{area:bh::ParameterSource::ParameterIndex(3),// srp_area from params[3]cr:bh::ParameterSource::ParameterIndex(4),// Cr from params[4]eclipse_model:bh::EclipseModel::None,};// Option 2: Cylindrical shadow model// Simple and fast, sharp shadow boundary (no penumbra)let_srp_cylindrical=bh::SolarRadiationPressureConfiguration{area:bh::ParameterSource::ParameterIndex(3),cr:bh::ParameterSource::ParameterIndex(4),eclipse_model:bh::EclipseModel::Cylindrical,};// Option 3: Conical shadow model (recommended)// Accounts for penumbra and umbra regionslet_srp_conical=bh::SolarRadiationPressureConfiguration{area:bh::ParameterSource::ParameterIndex(3),cr:bh::ParameterSource::ParameterIndex(4),eclipse_model:bh::EclipseModel::Conical,};}
importbraheasbh# Third-body perturbations configuration# Gravitational attraction from other celestial bodies# Option 1: Low-precision analytical ephemerides# Fast but less accurate (~km level errors for Sun/Moon)# Only Sun and Moon are availablethird_body_low=bh.ThirdBodyConfiguration(ephemeris_source=bh.EphemerisSource.LowPrecision,bodies=[bh.ThirdBody.SUN,bh.ThirdBody.MOON],)# Option 2: DE440s high-precision ephemerides (recommended)# Uses JPL Development Ephemeris 440 (small bodies version)# ~m level accuracy, valid 1550-2650 CE# All planets available, ~17 MB filethird_body_de440s=bh.ThirdBodyConfiguration(ephemeris_source=bh.EphemerisSource.DE440s,bodies=[bh.ThirdBody.SUN,bh.ThirdBody.MOON],)# Option 3: DE440 full-precision ephemerides# Highest accuracy (~mm level), valid 13200 BCE-17191 CE# All planets available, ~114 MB filethird_body_de440=bh.ThirdBodyConfiguration(ephemeris_source=bh.EphemerisSource.DE440,bodies=[bh.ThirdBody.SUN,bh.ThirdBody.MOON],)# Option 4: Include all major planets (high-fidelity)third_body_all_planets=bh.ThirdBodyConfiguration(ephemeris_source=bh.EphemerisSource.DE440s,bodies=[bh.ThirdBody.SUN,bh.ThirdBody.MOON,bh.ThirdBody.MERCURY,bh.ThirdBody.VENUS,bh.ThirdBody.MARS,bh.ThirdBody.JUPITER,bh.ThirdBody.SATURN,bh.ThirdBody.URANUS,bh.ThirdBody.NEPTUNE,],)
usebraheasbh;usebh::GravityModelType;fnmain(){// Third-body perturbations configuration// Gravitational attraction from other celestial bodies// Option 1: Low-precision analytical ephemerides// Fast but less accurate (~km level errors for Sun/Moon)// Only Sun and Moon are availablelet_third_body_low=bh::ThirdBodyConfiguration{ephemeris_source:bh::EphemerisSource::LowPrecision,bodies:vec![bh::ThirdBody::Sun,bh::ThirdBody::Moon],};// Option 2: DE440s high-precision ephemerides (recommended)// Uses JPL Development Ephemeris 440 (small bodies version)// ~m level accuracy, valid 1550-2650 CE// All planets available, ~17 MB fileletthird_body_de440s=bh::ThirdBodyConfiguration{ephemeris_source:bh::EphemerisSource::DE440s,bodies:vec![bh::ThirdBody::Sun,bh::ThirdBody::Moon],};// Option 3: DE440 full-precision ephemerides// Highest accuracy (~mm level), valid 13200 BCE-17191 CE// All planets available, ~114 MB filelet_third_body_de440=bh::ThirdBodyConfiguration{ephemeris_source:bh::EphemerisSource::DE440,bodies:vec![bh::ThirdBody::Sun,bh::ThirdBody::Moon],};// Option 4: Include all major planets (high-fidelity)let_third_body_all_planets=bh::ThirdBodyConfiguration{ephemeris_source:bh::EphemerisSource::DE440s,bodies:vec![bh::ThirdBody::Sun,bh::ThirdBody::Moon,bh::ThirdBody::Mercury,bh::ThirdBody::Venus,bh::ThirdBody::Mars,bh::ThirdBody::Jupiter,bh::ThirdBody::Saturn,bh::ThirdBody::Uranus,bh::ThirdBody::Neptune,],};// Create force model with Sun/Moon perturbations (common case)let_force_config=bh::ForceModelConfig{gravity:bh::GravityConfiguration::SphericalHarmonic{source:bh::GravityModelSource::ModelType(GravityModelType::EGM2008_360),degree:20,order:20,},drag:None,srp:None,third_body:Some(third_body_de440s),relativity:false,mass:None,};}
General relativistic corrections can be enabled via the relativity boolean flag. These effects are typically small but can be significant for precision orbit determination.
importbraheasbh# ParameterSource.value() creates a fixed constant parameter# Use when the parameter doesn't change and doesn't need to be estimated# Example: Fixed drag configuration# Mass, drag area, and Cd are all constantdrag_config=bh.DragConfiguration(model=bh.AtmosphericModel.HARRIS_PRIESTER,area=bh.ParameterSource.value(10.0),# Fixed 10 m^2 drag areacd=bh.ParameterSource.value(2.2),# Fixed Cd of 2.2)# Example: Fixed SRP configurationsrp_config=bh.SolarRadiationPressureConfiguration(area=bh.ParameterSource.value(15.0),# Fixed 15 m^2 SRP areacr=bh.ParameterSource.value(1.3),# Fixed Cr of 1.3eclipse_model=bh.EclipseModel.CONICAL,)# Create force model with all fixed parameters# Start from two_body preset and add componentsforce_config=bh.ForceModelConfig.two_body()force_config.gravity=bh.GravityConfiguration.spherical_harmonic(20,20)force_config.drag=drag_configforce_config.srp=srp_configforce_config.third_body=bh.ThirdBodyConfiguration(ephemeris_source=bh.EphemerisSource.LowPrecision,bodies=[bh.ThirdBody.SUN,bh.ThirdBody.MOON],)force_config.mass=bh.ParameterSource.value(500.0)# Fixed 500 kg mass
usebraheasbh;usebh::GravityModelType;fnmain(){// ParameterSource::Value creates a fixed constant parameter// Use when the parameter doesn't change and doesn't need to be estimated// Example: Fixed drag configuration// Mass, drag area, and Cd are all constantletdrag_config=bh::DragConfiguration{model:bh::AtmosphericModel::HarrisPriester,area:bh::ParameterSource::Value(10.0),// Fixed 10 m^2 drag areacd:bh::ParameterSource::Value(2.2),// Fixed Cd of 2.2};// Example: Fixed SRP configurationletsrp_config=bh::SolarRadiationPressureConfiguration{area:bh::ParameterSource::Value(15.0),// Fixed 15 m^2 SRP areacr:bh::ParameterSource::Value(1.3),// Fixed Cr of 1.3eclipse_model:bh::EclipseModel::Conical,};// Create force model with all fixed parameterslet_force_config=bh::ForceModelConfig{gravity:bh::GravityConfiguration::SphericalHarmonic{source:bh::GravityModelSource::ModelType(GravityModelType::EGM2008_360),degree:20,order:20,},drag:Some(drag_config),srp:Some(srp_config),third_body:Some(bh::ThirdBodyConfiguration{ephemeris_source:bh::EphemerisSource::LowPrecision,bodies:vec![bh::ThirdBody::Sun,bh::ThirdBody::Moon],}),relativity:false,mass:Some(bh::ParameterSource::Value(500.0)),// Fixed 500 kg mass};}
importbraheasbh# ParameterSource.parameter_index() references a value in the parameter vector# Use when parameters may change or need to be estimated# Default parameter layout:# Index 0: mass (kg)# Index 1: drag_area (m^2)# Index 2: Cd (dimensionless)# Index 3: srp_area (m^2)# Index 4: Cr (dimensionless)drag_config=bh.DragConfiguration(model=bh.AtmosphericModel.HARRIS_PRIESTER,area=bh.ParameterSource.parameter_index(1),# params[1] = drag_areacd=bh.ParameterSource.parameter_index(2),# params[2] = Cd)srp_config=bh.SolarRadiationPressureConfiguration(area=bh.ParameterSource.parameter_index(3),# params[3] = srp_areacr=bh.ParameterSource.parameter_index(4),# params[4] = Creclipse_model=bh.EclipseModel.CONICAL,)# Custom parameter layout examplecustom_drag=bh.DragConfiguration(model=bh.AtmosphericModel.HARRIS_PRIESTER,area=bh.ParameterSource.parameter_index(5),# Custom indexcd=bh.ParameterSource.parameter_index(10),# Custom index)
usebraheasbh;fnmain(){// ParameterSource::ParameterIndex references a value in the parameter vector// Use when parameters may change or need to be estimated// Default parameter layout:// Index 0: mass (kg)// Index 1: drag_area (m^2)// Index 2: Cd (dimensionless)// Index 3: srp_area (m^2)// Index 4: Cr (dimensionless)let_drag_config=bh::DragConfiguration{model:bh::AtmosphericModel::HarrisPriester,area:bh::ParameterSource::ParameterIndex(1),// params[1] = drag_areacd:bh::ParameterSource::ParameterIndex(2),// params[2] = Cd};let_srp_config=bh::SolarRadiationPressureConfiguration{area:bh::ParameterSource::ParameterIndex(3),// params[3] = srp_areacr:bh::ParameterSource::ParameterIndex(4),// params[4] = Creclipse_model:bh::EclipseModel::Conical,};// Custom parameter layout exampleprintln!("\nCustom parameter layout example:");println!(" You can map parameters to any indices:");let_custom_drag=bh::DragConfiguration{model:bh::AtmosphericModel::HarrisPriester,area:bh::ParameterSource::ParameterIndex(5),// Custom indexcd:bh::ParameterSource::ParameterIndex(10),// Custom index};}
importbraheasbh# Brahe provides several preset configurations for common scenarios# 1. two_body() - Point mass gravity only# Use for: Validation, comparison with Keplerian, quick estimatestwo_body=bh.ForceModelConfig.two_body()# 2. earth_gravity() - Spherical harmonic gravity only (20x20)# Use for: Studying gravity perturbations in isolationearth_gravity=bh.ForceModelConfig.earth_gravity()# 3. conservative_forces() - Gravity + third-body + relativity (no drag/SRP)# Use for: Long-term orbit evolution, conservative dynamics studiesconservative=bh.ForceModelConfig.conservative_forces()# 4. default() - Balanced configuration for LEO to GEO# Use for: General mission analysis, initial studiesdefault=bh.ForceModelConfig.default()# 5. leo_default() - Optimized for Low Earth Orbit# Use for: LEO missions where drag is dominantleo=bh.ForceModelConfig.leo_default()# 6. geo_default() - Optimized for Geostationary Orbit# Use for: GEO missions where SRP and third-body dominategeo=bh.ForceModelConfig.geo_default()# 7. high_fidelity() - Maximum precision# Use for: Precision orbit determination, research applicationshigh_fidelity=bh.ForceModelConfig.high_fidelity()
usebraheasbh;fnmain(){// Brahe provides several preset configurations for common scenarios// 1. two_body_gravity() - Point mass gravity only// Use for: Validation, comparison with Keplerian, quick estimateslet_two_body=bh::ForceModelConfig::two_body_gravity();// 2. earth_gravity() - Spherical harmonic gravity only (20x20)// Use for: Studying gravity perturbations in isolationlet_earth_gravity=bh::ForceModelConfig::earth_gravity();// 3. conservative_forces() - Gravity + third-body + relativity (no drag/SRP)// Use for: Long-term orbit evolution, conservative dynamics studieslet_conservative=bh::ForceModelConfig::conservative_forces();// 4. default() - Balanced configuration for LEO to GEO// Use for: General mission analysis, initial studieslet_default=bh::ForceModelConfig::default();// 5. leo_default() - Optimized for Low Earth Orbit// Use for: LEO missions where drag is dominantlet_leo=bh::ForceModelConfig::leo_default();// 6. geo_default() - Optimized for Geostationary Orbit// Use for: GEO missions where SRP and third-body dominatelet_geo=bh::ForceModelConfig::geo_default();// 7. high_fidelity() - Maximum precision// Use for: Precision orbit determination, research applicationslet_high_fidelity=bh::ForceModelConfig::high_fidelity();}