Please check the build logs for more information.
See Builds for ideas on how to fix a failed build, or Metadata for how to configure docs.rs builds.
If you believe this is docs.rs' fault, open an issue.
libassist-sys: Raw FFI bindings and safe RAII wrappers for the ASSIST ephemeris-driven N-body integrator
A Rust crate by the Asteroid Institute, a program of the B612 Foundation
libassist-sys is a -sys crate exposing the C ABI of ASSIST, JPL's ephemeris-driven N-body force model layered on top of REBOUND. It vendors ASSIST as a git submodule, compiles it via the cc crate, and provides:
ffi— rawextern "C"bindings to ASSIST functions and types.Ephemeris— allocation-owning RAII wrapper aroundassist_ephem. Loadsde440.bsp+sb441-n16.bsp.Send + Sync— load once, share across threads.AssistSim— a REBOUNDSimulationwith ASSIST forces attached. Owns itsassist_extras, detaches and frees on drop. Supports non-grav coefficients (Marsden-Sekanina), variational equations for STM, andintegrate_or_interpolatefor cheap light-time iteration.Error— wraps REBOUND's typed integration-exit conditions plus anEphemerisErrorvariant for SPK / lookup failures.
Force model (ASSIST defaults): Sun + 8 planets + Moon + 16 massive asteroids + J2 Earth + general relativity. Optional non-gravitational acceleration.
No domain logic — no orbital elements, no observatories, no light-time iteration. Layer those on top via assist-rs.
Crate hierarchy
assist-rs ← domain types: Orbit, Origin, Observer, DataManager
└── libassist-sys ← raw ASSIST FFI + AssistSim/Ephemeris RAII (this crate)
└── librebound-sys ← raw REBOUND FFI + Simulation RAII
Installation
[]
= "1.2"
Requires a C compiler (clang recommended for ThinLTO, see Building). The crate depends on librebound-sys for the REBOUND ABI.
Usage
use ;
use Path;
let ephem = from_paths?;
let mut sim = new?;
sim.set_t;
default.apply;
let mut asim = new?;
asim.set_forces;
// Add a test particle in barycentric equatorial ICRF (AU, AU/day).
asim.sim_mut.add_test_particle;
let t_target = ephem.mjd_to_assist_time;
asim.integrate?;
let state = &asim.sim.particles;
println!;
# Ok::
For typed Orbit + heliocentric ecliptic API, light-time correction, observatories, and data-file management, use assist-rs instead — it re-exports everything from this crate.
Thread safety
Ephemeris is Send + Sync: the underlying assist_ephem is read-only after construction, the only mutator (set_jd_ref) takes &mut self, and as_ptr() returns *const to forbid back-door mutation. Multiple AssistSim instances on separate threads sharing one Arc<Ephemeris> is supported — REBOUND's hot paths use only static const tables, and the only process-shared mutable state is the SIGINT handler, which all simulations legitimately share. See wrappers.rs for the full audit.
AssistSim is Send but not Sync — one mutable simulation per thread.
Building
The crate vendors ASSIST under vendor/assist/ (git submodule). build.rs compiles vendor/assist/src/*.c against librebound-sys's REBOUND headers via the cc crate.
For maximum performance use clang (honours -flto=thin):
CC=clang
GCC builds are correct but ~5–7 % slower on hot ephemeris/force-evaluation paths.
Testing
Tests that need ephemeris data skip cleanly when the environment variables are unset. SPK files are available from the B612 Asteroid Institute data packages or directly from JPL.
Versioning
The crate's major.minor mirrors the vendored ASSIST release: libassist-sys 1.2.x wraps ASSIST 1.2.0. The patch component is reserved for Rust-side fixes that don't change the underlying C library. An ASSIST 1.3 release would become libassist-sys 1.3.0; an ASSIST 2.0 release would become libassist-sys 2.0.0.
The REBOUND version is whatever librebound-sys pulls in (currently 4.6.0). REBOUND-only upgrades that don't bump ASSIST will bump only librebound-sys; libassist-sys may follow with a patch release if any of its wrappers need adjustment.
Why this scheme
-sys crates are thin wrappers around a single upstream library. Mirroring upstream's major.minor makes the dependency obvious from the version alone — no need to consult a mapping table or [package.metadata.vendored]. The Rust side gets the patch slot for wrapper fixes (e.g. tightening Send/Sync bounds, fixing a leak in a wrapper, adding a typed Error variant), which is where ~all non-upstream changes land in practice.
This is not standard Rust semver — bumping the crate's major version isn't an intrinsic Rust-API break; it tracks the C ABI break in the upstream library. Callers depending on libassist-sys = "1" get any 1.x.y; callers needing strict pinning should pin to 1.2 or =1.2.0.
| libassist-sys | ASSIST | REBOUND |
|---|---|---|
| 1.2.x | 1.2.0 | 4.6.0 (via librebound-sys 4.6.x) |
License
GPL-3.0 — required by the vendored ASSIST and REBOUND sources. See LICENSE and vendor/assist/LICENSE.
References
- Holman et al. 2023, "ASSIST: An Ephemeris-Quality Test Particle Integrator", PSJ 4 69 (arXiv:2303.16246)
- Rein & Liu 2012, "REBOUND: An open-source multi-purpose N-body code for collisional dynamics", A&A 537 A128
- Rein & Spiegel 2015, "IAS15: a fast, adaptive, high-order integrator for gravitational dynamics", MNRAS 446 1424
Acknowledgments
This crate is a thin Rust wrapper. All credit for the underlying physics and integrator implementations belongs to the upstream projects:
- ASSIST — Matthew Holman and contributors. Source: https://github.com/matthewholman/assist. Vendored under
vendor/assist/with original LICENSE preserved atvendor/assist/LICENSE. - REBOUND — Hanno Rein and contributors. Source: https://github.com/hannorein/rebound. Used via the companion
librebound-syscrate, which vendors REBOUND with its LICENSE preserved.
If you use this crate in published work, please cite the ASSIST and REBOUND papers listed in References — not this crate.
The Rust wrapper layer is developed by the Asteroid Institute, a program of the B612 Foundation.