Skip to main content

Crate ionotrace

Crate ionotrace 

Source
Expand description

ionotrace — High-performance ionospheric ray tracing engine

Simulates HF radio wave propagation through a magnetized, collisional ionosphere using Hamilton’s equations of motion, based on the OT 75-76 algorithm.

Compiles to native (with parallel ray tracing via Rayon) or WebAssembly.

§Quick Start

Trace a single 10 MHz ray at 20° elevation:

use ionotrace::TraceConfig;

let result = TraceConfig::new(10.0, 20.0).trace().unwrap();

println!("Max height: {:.1} km", result.max_height);
println!("Ground range: {:.1} km", result.ground_range_km);
println!("Returned: {}", result.returned_to_ground);

§Custom Physics

Override ionospheric models via ModelParams:

use ionotrace::{TraceConfig, ModelParams};
use ionotrace::params::{ElectronDensityModel, MagneticFieldModel};

let mut config = TraceConfig::new(15.0, 30.0);
config.params = ModelParams::builder()
    .ed_model(ElectronDensityModel::DualChapman)
    .mag_model(MagneticFieldModel::Igrf14)
    .fc(8.0)    // critical frequency (MHz)
    .hm(300.0)  // peak height (km)
    .build()
    .unwrap();

§Fan Traces

Sweep a range of elevation angles (parallelized on native via Rayon):

use ionotrace::{fan_trace, FanTraceConfig};

let config = FanTraceConfig {
    freq_mhz: 10.0,
    elev_min: 5.0,
    elev_max: 80.0,
    elev_step: 5.0,
    ..FanTraceConfig::default()
};

let result = fan_trace(&config).unwrap();
println!("Traced {} rays", result.n_rays);

for ray in &result.rays {
    if ray.ground {
        println!("  {:.0}° → {:.0} km", ray.elev, ray.range_km);
    }
}

§Target Solver

Find the elevation/azimuth to hit a specific location:

use ionotrace::{solve_target, TargetConfig, SearchSpec};

let config = TargetConfig {
    target_lat_deg: 50.0,   // target latitude
    target_lon_deg: 5.0,    // target longitude
    tx_lat_deg: 40.0,       // transmitter latitude
    freq_mhz: SearchSpec::Fixed(10.0),
    error_limit_km: 20.0,
    ..TargetConfig::default()
};

let result = solve_target(&config).unwrap();
if let Some(best) = &result.best {
    println!("Aim: {:.1}° elev, {:.1}° az → {:.1} km error",
        best.elevation_deg, best.azimuth_deg, best.error_km);
}

§Exporting

Save results to CSV or JSON:

use ionotrace::{TraceConfig, export_trace_csv, export_json};

let result = TraceConfig::new(10.0, 20.0).trace().unwrap();
let csv = export_trace_csv(&result).unwrap();
let json = export_json(&result).unwrap();

§Physics Models

The engine is modular — swap models independently:

ComponentOptions
Electron densityChapman, ELECT1, Linear, Quasi-Parabolic, Variable Chapman, Dual Chapman
Magnetic fieldDipole, Constant, Cubic, IGRF-14
Refractive indexFull Appleton-Hartree, with/without field and collisions
CollisionsDouble-exponential, Single-exponential, Constant
PerturbationsTorus, Trough, Shock, Bulge, Exponential

See params for all configuration options.

Re-exports§

pub use error::TraceError;
pub use export::export_fan_trace_csv;
pub use export::export_json;
pub use export::export_trace_csv;
pub use fan::fan_trace;
pub use fan::FanRay;
pub use fan::FanRayPoint;
pub use fan::FanTraceConfig;
pub use fan::FanTraceResult;
pub use params::ModelParams;
pub use target::solve_target;
pub use target::SearchSpec;
pub use target::TargetConfig;
pub use target::TargetResult;
pub use target::TargetSolution;
pub use tracer::TraceConfig;
pub use tracer::TracePoint;
pub use tracer::TraceResult;

Modules§

error
Error types for ionospheric ray tracing.
export
Export utilities for trace results — CSV and JSON.
fan
Fan tracing — sweep a range of elevation angles in one call.
models
Physics models for ionospheric ray tracing.
params
Physical constants and model parameters for ionospheric ray tracing.
target
Target location solver — find ray launch parameters to hit a geographic target.
tracer
Ray tracing execution and result types.