pub trait PhaseSpaceRepr: Send + Sync {
Show 17 methods
// Required methods
fn compute_density(&self) -> DensityField;
fn advect_x(&mut self, displacement: &DisplacementField, dt: f64);
fn advect_v(&mut self, acceleration: &AccelerationField, dt: f64);
fn moment(&self, position: &[f64; 3], order: usize) -> Tensor;
fn total_mass(&self) -> f64;
fn casimir_c2(&self) -> f64;
fn entropy(&self) -> f64;
fn stream_count(&self) -> StreamCountField;
fn velocity_distribution(&self, position: &[f64; 3]) -> Vec<f64>;
fn as_any(&self) -> &dyn Any;
fn as_any_mut(&mut self) -> &mut dyn Any;
// Provided methods
fn total_kinetic_energy(&self) -> Option<f64> { ... }
fn to_snapshot(&self, time: f64) -> Option<PhaseSpaceSnapshot> { ... }
fn load_snapshot(
&mut self,
snap: PhaseSpaceSnapshot,
) -> Result<(), CausticError> { ... }
fn can_materialize(&self) -> bool { ... }
fn memory_bytes(&self) -> usize { ... }
fn set_progress(&mut self, _progress: Arc<StepProgress>) { ... }
}Expand description
Central trait for all phase-space storage and manipulation strategies.
Implementations differ in memory layout and algorithmic complexity:
UniformGrid6D: O(N⁶) brute-force gridTensorTrain: O(N³r³) low-rank decompositionSheetTracker: O(N³) Lagrangian cold sheet
§Examples
use caustic::{PhaseSpaceRepr, UniformGrid6D, PlummerIC, sample_on_grid};
use caustic::{Domain, DomainBuilder, SpatialBoundType, VelocityBoundType};
let domain = Domain::builder()
.spatial_extent(10.0)
.velocity_extent(5.0)
.spatial_resolution(8)
.velocity_resolution(8)
.spatial_bc(SpatialBoundType::Periodic)
.velocity_bc(VelocityBoundType::Open)
.build()
.unwrap();
let ic = PlummerIC::new(1.0, 1.0, 1.0);
let snap = sample_on_grid(&ic, &domain);
let repr = UniformGrid6D::from_snapshot(snap, domain);
// Query conserved quantities
let mass = repr.total_mass();
let c2 = repr.casimir_c2();
let ke = repr.total_kinetic_energy(); // Option<f64>
// Compute density for Poisson coupling
let density = repr.compute_density();Required Methods§
Sourcefn compute_density(&self) -> DensityField
fn compute_density(&self) -> DensityField
Integrate f over all velocities: ρ(x) = ∫f dv³. This is the coupling moment to the Poisson equation.
Sourcefn advect_x(&mut self, displacement: &DisplacementField, dt: f64)
fn advect_x(&mut self, displacement: &DisplacementField, dt: f64)
Drift sub-step: advect f in spatial coordinates by displacement Δx = v·dt. Pure translation in x at constant v.
Sourcefn advect_v(&mut self, acceleration: &AccelerationField, dt: f64)
fn advect_v(&mut self, acceleration: &AccelerationField, dt: f64)
Kick sub-step: advect f in velocity coordinates by Δv = g·dt. Pure translation in v at constant x.
Sourcefn moment(&self, position: &[f64; 3], order: usize) -> Tensor
fn moment(&self, position: &[f64; 3], order: usize) -> Tensor
Compute velocity moment of order n at given spatial position. Order 0 = density, 1 = mean velocity, 2 = dispersion tensor.
Sourcefn total_mass(&self) -> f64
fn total_mass(&self) -> f64
Total mass M = ∫f dx³dv³. Should be conserved to machine precision.
Sourcefn casimir_c2(&self) -> f64
fn casimir_c2(&self) -> f64
Casimir invariant C₂ = ∫f² dx³dv³. Increase over time indicates numerical diffusion.
Sourcefn entropy(&self) -> f64
fn entropy(&self) -> f64
Boltzmann entropy S = −∫f ln f dx³dv³. Should be exactly conserved; growth = numerical error.
Sourcefn stream_count(&self) -> StreamCountField
fn stream_count(&self) -> StreamCountField
Number of distinct velocity streams at each spatial point. Detects caustic surfaces (sheet folds).
Sourcefn velocity_distribution(&self, position: &[f64; 3]) -> Vec<f64>
fn velocity_distribution(&self, position: &[f64; 3]) -> Vec<f64>
Extract the local velocity distribution f(v|x) at a given spatial position. Used for dark matter detection predictions.
Sourcefn as_any(&self) -> &dyn Any
fn as_any(&self) -> &dyn Any
Downcast to concrete type for implementation-specific queries (e.g. HT rank data).
Sourcefn as_any_mut(&mut self) -> &mut dyn Any
fn as_any_mut(&mut self) -> &mut dyn Any
Mutable downcast for in-place modification (e.g. BUG integrator leaf updates).
Returns self as &mut dyn Any. Only HtTensor overrides this; the default
returns self but concrete downcasts will yield None.
Provided Methods§
Sourcefn total_kinetic_energy(&self) -> Option<f64>
fn total_kinetic_energy(&self) -> Option<f64>
Total kinetic energy T = ½∫fv² dx³dv³.
Returns None if this representation does not support direct kinetic
energy computation. All grid-based representations should implement this.
Sourcefn to_snapshot(&self, time: f64) -> Option<PhaseSpaceSnapshot>
fn to_snapshot(&self, time: f64) -> Option<PhaseSpaceSnapshot>
Extract a full 6D snapshot of the current state.
Returns None if this representation cannot produce a dense 6D snapshot
(e.g. because materialization would exceed memory). Check can_materialize
before calling if the result is optional in your context.
Sourcefn load_snapshot(
&mut self,
snap: PhaseSpaceSnapshot,
) -> Result<(), CausticError>
fn load_snapshot( &mut self, snap: PhaseSpaceSnapshot, ) -> Result<(), CausticError>
Replace the current state with data from a dense 6D snapshot.
Required for unsplit (method-of-lines) time integration, which manipulates
the distribution function directly rather than through drift/kick sub-steps.
Returns Err if this representation does not support snapshot loading.
Sourcefn can_materialize(&self) -> bool
fn can_materialize(&self) -> bool
Whether full 6D materialization fits in available memory. Default: true. Compressed representations (HT, TT) should override to check shape vs a memory threshold.
Sourcefn memory_bytes(&self) -> usize
fn memory_bytes(&self) -> usize
Approximate memory usage of this representation in bytes. Default returns 0; implementations should override for accurate tracking.
Sourcefn set_progress(&mut self, _progress: Arc<StepProgress>)
fn set_progress(&mut self, _progress: Arc<StepProgress>)
Attach shared progress state for intra-phase cell-level reporting.
Implementations can use this to report progress via StepProgress::set_intra_progress.