#![allow(clippy::needless_range_loop)]
#![allow(missing_docs)]
use serde::{Deserialize, Serialize};
#[derive(Debug, Clone, Serialize, Deserialize)]
pub struct PySphConfig {
pub h: f64,
pub rest_density: f64,
pub stiffness: f64,
pub viscosity: f64,
pub gravity: [f64; 3],
pub particle_mass: f64,
pub floor_enabled: bool,
pub floor_y: f64,
pub floor_restitution: f64,
}
impl PySphConfig {
pub fn water() -> Self {
Self {
h: 0.1,
rest_density: 1000.0,
stiffness: 200.0,
viscosity: 0.01,
gravity: [0.0, -9.81, 0.0],
particle_mass: 0.02,
floor_enabled: false,
floor_y: 0.0,
floor_restitution: 0.3,
}
}
pub fn gas() -> Self {
Self {
h: 0.2,
rest_density: 1.2,
stiffness: 50.0,
viscosity: 1.81e-5,
gravity: [0.0, 0.0, 0.0],
particle_mass: 0.001,
floor_enabled: false,
floor_y: 0.0,
floor_restitution: 1.0,
}
}
}
impl Default for PySphConfig {
fn default() -> Self {
Self::water()
}
}
#[derive(Debug, Clone)]
pub struct PySphSimulation {
positions: Vec<[f64; 3]>,
velocities: Vec<[f64; 3]>,
densities: Vec<f64>,
pressures: Vec<f64>,
config: PySphConfig,
time: f64,
step_count: u64,
}
impl PySphSimulation {
pub fn new(config: PySphConfig) -> Self {
Self {
positions: Vec::new(),
velocities: Vec::new(),
densities: Vec::new(),
pressures: Vec::new(),
config,
time: 0.0,
step_count: 0,
}
}
pub fn add_particle(&mut self, pos: [f64; 3]) -> usize {
let idx = self.positions.len();
self.positions.push(pos);
self.velocities.push([0.0; 3]);
self.densities.push(self.config.rest_density);
self.pressures.push(0.0);
idx
}
pub fn add_particle_with_velocity(&mut self, pos: [f64; 3], vel: [f64; 3]) -> usize {
let idx = self.add_particle(pos);
self.velocities[idx] = vel;
idx
}
pub fn add_particle_block(
&mut self,
origin: [f64; 3],
nx: usize,
ny: usize,
nz: usize,
dx: f64,
) {
for iz in 0..nz {
for iy in 0..ny {
for ix in 0..nx {
let pos = [
origin[0] + ix as f64 * dx,
origin[1] + iy as f64 * dx,
origin[2] + iz as f64 * dx,
];
self.add_particle(pos);
}
}
}
}
pub fn particle_count(&self) -> usize {
self.positions.len()
}
pub fn position(&self, i: usize) -> Option<[f64; 3]> {
self.positions.get(i).copied()
}
pub fn velocity(&self, i: usize) -> Option<[f64; 3]> {
self.velocities.get(i).copied()
}
pub fn density(&self, i: usize) -> Option<f64> {
self.densities.get(i).copied()
}
pub fn time(&self) -> f64 {
self.time
}
pub fn step_count(&self) -> u64 {
self.step_count
}
pub fn all_positions(&self) -> Vec<f64> {
self.positions
.iter()
.flat_map(|p| p.iter().copied())
.collect()
}
pub fn all_velocities(&self) -> Vec<f64> {
self.velocities
.iter()
.flat_map(|v| v.iter().copied())
.collect()
}
pub fn get_density_field(&self) -> Vec<f64> {
self.densities.clone()
}
pub fn get_pressure_field(&self) -> Vec<f64> {
self.pressures.clone()
}
pub fn mean_density(&self) -> f64 {
if self.densities.is_empty() {
return 0.0;
}
self.densities.iter().sum::<f64>() / self.densities.len() as f64
}
pub fn step(&mut self, dt: f64) {
let n = self.positions.len();
if n == 0 {
self.time += dt;
self.step_count += 1;
return;
}
let h = self.config.h;
let h2 = h * h;
let m = self.config.particle_mass;
let rho0 = self.config.rest_density;
let k = self.config.stiffness;
let mu = self.config.viscosity;
let g = self.config.gravity;
let poly6_coeff = 315.0 / (64.0 * std::f64::consts::PI * h.powi(9));
let spiky_coeff = -45.0 / (std::f64::consts::PI * h.powi(6));
let visc_coeff = 45.0 / (std::f64::consts::PI * h.powi(6));
for i in 0..n {
let mut rho = 0.0f64;
for j in 0..n {
let dx = self.positions[i][0] - self.positions[j][0];
let dy = self.positions[i][1] - self.positions[j][1];
let dz = self.positions[i][2] - self.positions[j][2];
let r2 = dx * dx + dy * dy + dz * dz;
if r2 < h2 {
let diff = h2 - r2;
rho += m * poly6_coeff * diff * diff * diff;
}
}
self.densities[i] = rho.max(1e-3);
self.pressures[i] = k * (self.densities[i] / rho0 - 1.0);
}
let mut forces = vec![[0.0f64; 3]; n];
for i in 0..n {
let mut fp = [0.0f64; 3];
let mut fv = [0.0f64; 3];
for j in 0..n {
if i == j {
continue;
}
let dx = self.positions[i][0] - self.positions[j][0];
let dy = self.positions[i][1] - self.positions[j][1];
let dz = self.positions[i][2] - self.positions[j][2];
let r2 = dx * dx + dy * dy + dz * dz;
if r2 >= h2 || r2 < 1e-20 {
continue;
}
let r = r2.sqrt();
let hr = h - r;
let pf = -m * (self.pressures[i] + self.pressures[j])
/ (2.0 * self.densities[j].max(1e-6))
* spiky_coeff
* hr
* hr
/ r;
fp[0] += pf * dx;
fp[1] += pf * dy;
fp[2] += pf * dz;
let vf = mu * m / self.densities[j].max(1e-6) * visc_coeff * hr;
fv[0] += vf * (self.velocities[j][0] - self.velocities[i][0]);
fv[1] += vf * (self.velocities[j][1] - self.velocities[i][1]);
fv[2] += vf * (self.velocities[j][2] - self.velocities[i][2]);
}
let rho_i = self.densities[i];
forces[i][0] = (fp[0] + fv[0]) / rho_i + g[0];
forces[i][1] = (fp[1] + fv[1]) / rho_i + g[1];
forces[i][2] = (fp[2] + fv[2]) / rho_i + g[2];
}
for i in 0..n {
self.velocities[i][0] += forces[i][0] * dt;
self.velocities[i][1] += forces[i][1] * dt;
self.velocities[i][2] += forces[i][2] * dt;
self.positions[i][0] += self.velocities[i][0] * dt;
self.positions[i][1] += self.velocities[i][1] * dt;
self.positions[i][2] += self.velocities[i][2] * dt;
}
if self.config.floor_enabled {
let floor_y = self.config.floor_y;
let rest_coeff = self.config.floor_restitution;
for i in 0..n {
if self.positions[i][1] < floor_y {
self.positions[i][1] = floor_y;
if self.velocities[i][1] < 0.0 {
self.velocities[i][1] = -self.velocities[i][1] * rest_coeff;
}
}
}
}
self.time += dt;
self.step_count += 1;
}
}
#[cfg(test)]
mod tests {
use crate::PySphConfig;
use crate::PySphSimulation;
fn water_sim() -> PySphSimulation {
PySphSimulation::new(PySphConfig::water())
}
#[test]
fn test_sph_creation_empty() {
let sim = water_sim();
assert_eq!(sim.particle_count(), 0);
assert!((sim.time()).abs() < 1e-15);
}
#[test]
fn test_sph_add_particle() {
let mut sim = water_sim();
let idx = sim.add_particle([1.0, 2.0, 3.0]);
assert_eq!(idx, 0);
assert_eq!(sim.particle_count(), 1);
let p = sim.position(0).unwrap();
assert!((p[0] - 1.0).abs() < 1e-12);
}
#[test]
fn test_sph_add_particle_with_velocity() {
let mut sim = water_sim();
sim.add_particle_with_velocity([0.0; 3], [1.0, 2.0, 3.0]);
let v = sim.velocity(0).unwrap();
assert!((v[0] - 1.0).abs() < 1e-12);
assert!((v[1] - 2.0).abs() < 1e-12);
}
#[test]
fn test_sph_particle_block() {
let mut sim = water_sim();
sim.add_particle_block([0.0; 3], 2, 3, 4, 0.1);
assert_eq!(sim.particle_count(), 24);
}
#[test]
fn test_sph_step_moves_particle_under_gravity() {
let mut sim = water_sim();
sim.add_particle([0.0, 1.0, 0.0]);
let y0 = sim.position(0).unwrap()[1];
sim.step(0.001);
let y1 = sim.position(0).unwrap()[1];
assert!(y1 < y0, "particle should fall: y0={} y1={}", y0, y1);
}
#[test]
fn test_sph_step_empty_no_panic() {
let mut sim = water_sim();
sim.step(0.01);
assert!((sim.time() - 0.01).abs() < 1e-15);
}
#[test]
fn test_sph_time_advances() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
sim.step(0.01);
sim.step(0.01);
assert!((sim.time() - 0.02).abs() < 1e-12);
assert_eq!(sim.step_count(), 2);
}
#[test]
fn test_sph_density_field_length() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
sim.add_particle([1.0, 0.0, 0.0]);
let d = sim.get_density_field();
assert_eq!(d.len(), 2);
}
#[test]
fn test_sph_pressure_field_length() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
let p = sim.get_pressure_field();
assert_eq!(p.len(), 1);
}
#[test]
fn test_sph_all_positions_length() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
sim.add_particle([1.0, 0.0, 0.0]);
assert_eq!(sim.all_positions().len(), 6);
}
#[test]
fn test_sph_all_velocities_length() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
assert_eq!(sim.all_velocities().len(), 3);
}
#[test]
fn test_sph_floor_boundary() {
let cfg = PySphConfig {
floor_enabled: true,
floor_y: 0.0,
floor_restitution: 0.5,
..PySphConfig::water()
};
let mut sim = PySphSimulation::new(cfg);
sim.add_particle([0.0, 0.001, 0.0]);
sim.velocities[0] = [0.0, -5.0, 0.0];
sim.step(0.01);
let y = sim.position(0).unwrap()[1];
assert!(y >= 0.0, "particle below floor: y={}", y);
}
#[test]
fn test_sph_mean_density_nonzero_after_step() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
sim.step(0.001);
assert!(sim.mean_density() > 0.0);
}
#[test]
fn test_sph_gas_config() {
let cfg = PySphConfig::gas();
assert!(cfg.rest_density < 10.0);
assert!(cfg.particle_mass < 0.01);
}
#[test]
fn test_sph_position_out_of_bounds() {
let sim = water_sim();
assert!(sim.position(99).is_none());
}
}
impl PySphSimulation {
pub fn neighbors(&self, i: usize, r: f64) -> Option<Vec<usize>> {
let pos_i = self.positions.get(i)?;
let r2 = r * r;
let neighbors: Vec<usize> = self
.positions
.iter()
.enumerate()
.filter(|(j, pos_j)| {
if *j == i {
return false;
}
let dx = pos_i[0] - pos_j[0];
let dy = pos_i[1] - pos_j[1];
let dz = pos_i[2] - pos_j[2];
dx * dx + dy * dy + dz * dz <= r2
})
.map(|(j, _)| j)
.collect();
Some(neighbors)
}
pub fn neighbor_count(&self, i: usize, r: f64) -> Option<usize> {
self.neighbors(i, r).map(|v| v.len())
}
pub fn closest_neighbor(&self, i: usize) -> Option<(usize, f64)> {
let pos_i = self.positions.get(i)?;
let mut best = None::<(usize, f64)>;
for (j, pos_j) in self.positions.iter().enumerate() {
if j == i {
continue;
}
let dx = pos_i[0] - pos_j[0];
let dy = pos_i[1] - pos_j[1];
let dz = pos_i[2] - pos_j[2];
let dist = (dx * dx + dy * dy + dz * dz).sqrt();
match best {
None => best = Some((j, dist)),
Some((_, d)) if dist < d => best = Some((j, dist)),
_ => {}
}
}
best
}
pub fn center_of_mass(&self) -> [f64; 3] {
let n = self.positions.len();
if n == 0 {
return [0.0; 3];
}
let mut cx = 0.0f64;
let mut cy = 0.0f64;
let mut cz = 0.0f64;
for pos in &self.positions {
cx += pos[0];
cy += pos[1];
cz += pos[2];
}
let inv_n = 1.0 / n as f64;
[cx * inv_n, cy * inv_n, cz * inv_n]
}
pub fn max_speed(&self) -> f64 {
self.velocities
.iter()
.map(|v| (v[0] * v[0] + v[1] * v[1] + v[2] * v[2]).sqrt())
.fold(0.0f64, f64::max)
}
pub fn kinetic_energy(&self) -> f64 {
let m = self.config.particle_mass;
self.velocities
.iter()
.map(|v| 0.5 * m * (v[0] * v[0] + v[1] * v[1] + v[2] * v[2]))
.sum()
}
pub fn max_density(&self) -> f64 {
self.densities
.iter()
.cloned()
.fold(f64::NEG_INFINITY, f64::max)
}
pub fn min_density(&self) -> f64 {
self.densities.iter().cloned().fold(f64::INFINITY, f64::min)
}
pub fn particle_aabb(&self) -> Option<[f64; 6]> {
if self.positions.is_empty() {
return None;
}
let mut mn = [f64::INFINITY; 3];
let mut mx = [f64::NEG_INFINITY; 3];
for pos in &self.positions {
for k in 0..3 {
if pos[k] < mn[k] {
mn[k] = pos[k];
}
if pos[k] > mx[k] {
mx[k] = pos[k];
}
}
}
Some([mn[0], mn[1], mn[2], mx[0], mx[1], mx[2]])
}
pub fn remove_below_y(&mut self, y_min: f64) -> usize {
let before = self.positions.len();
let keep: Vec<bool> = self.positions.iter().map(|p| p[1] >= y_min).collect();
let mut new_pos = Vec::new();
let mut new_vel = Vec::new();
let mut new_den = Vec::new();
let mut new_prs = Vec::new();
for (i, &keep_i) in keep.iter().enumerate() {
if keep_i {
new_pos.push(self.positions[i]);
new_vel.push(self.velocities[i]);
new_den.push(self.densities[i]);
new_prs.push(self.pressures[i]);
}
}
self.positions = new_pos;
self.velocities = new_vel;
self.densities = new_den;
self.pressures = new_prs;
before - self.positions.len()
}
}
#[derive(Debug, Clone)]
#[allow(dead_code)]
pub struct SphStats {
pub particle_count: usize,
pub mean_density: f64,
pub max_density: f64,
pub min_density: f64,
pub max_speed: f64,
pub kinetic_energy: f64,
pub time: f64,
pub step_count: u64,
}
impl PySphSimulation {
pub fn collect_stats(&self) -> SphStats {
SphStats {
particle_count: self.particle_count(),
mean_density: self.mean_density(),
max_density: self.max_density(),
min_density: if self.densities.is_empty() {
0.0
} else {
self.min_density()
},
max_speed: self.max_speed(),
kinetic_energy: self.kinetic_energy(),
time: self.time,
step_count: self.step_count,
}
}
}
#[cfg(test)]
mod sph_ext_tests {
use crate::PySphConfig;
use crate::PySphSimulation;
fn water_sim() -> PySphSimulation {
PySphSimulation::new(PySphConfig::water())
}
#[test]
fn test_sph_neighbors_empty() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
let n = sim.neighbors(0, 1.0).unwrap();
assert!(n.is_empty());
}
#[test]
fn test_sph_neighbors_finds_close() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
sim.add_particle([0.05, 0.0, 0.0]); sim.add_particle([5.0, 0.0, 0.0]); let n = sim.neighbors(0, 0.1).unwrap();
assert_eq!(n.len(), 1);
assert_eq!(n[0], 1);
}
#[test]
fn test_sph_neighbor_count() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
sim.add_particle([0.05, 0.0, 0.0]);
sim.add_particle([0.05, 0.05, 0.0]);
let count = sim.neighbor_count(0, 0.15).unwrap();
assert_eq!(count, 2);
}
#[test]
fn test_sph_closest_neighbor() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
sim.add_particle([1.0, 0.0, 0.0]);
sim.add_particle([0.3, 0.0, 0.0]);
let result = sim.closest_neighbor(0).unwrap();
assert_eq!(result.0, 2); assert!((result.1 - 0.3).abs() < 1e-10);
}
#[test]
fn test_sph_closest_neighbor_single_particle() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
assert!(sim.closest_neighbor(0).is_none());
}
#[test]
fn test_sph_center_of_mass() {
let mut sim = water_sim();
sim.add_particle([0.0, 0.0, 0.0]);
sim.add_particle([2.0, 0.0, 0.0]);
let com = sim.center_of_mass();
assert!((com[0] - 1.0).abs() < 1e-10);
}
#[test]
fn test_sph_center_of_mass_empty() {
let sim = water_sim();
let com = sim.center_of_mass();
for &c in &com {
assert!(c.abs() < 1e-15);
}
}
#[test]
fn test_sph_max_speed_zero_at_rest() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
assert!(sim.max_speed() < 1e-15);
}
#[test]
fn test_sph_kinetic_energy_with_velocity() {
let mut sim = water_sim();
sim.add_particle_with_velocity([0.0; 3], [2.0, 0.0, 0.0]);
let ke = sim.kinetic_energy();
assert!((ke - 0.04).abs() < 1e-12, "KE = {}", ke);
}
#[test]
fn test_sph_max_density_after_step() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
sim.step(0.001);
assert!(sim.max_density() > 0.0);
}
#[test]
fn test_sph_min_density_after_step() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
sim.step(0.001);
assert!(sim.min_density() > 0.0);
}
#[test]
fn test_sph_particle_aabb() {
let mut sim = water_sim();
sim.add_particle([0.0, 0.0, 0.0]);
sim.add_particle([1.0, 2.0, 3.0]);
let aabb = sim.particle_aabb().unwrap();
assert!((aabb[0]).abs() < 1e-15); assert!((aabb[3] - 1.0).abs() < 1e-15); assert!((aabb[5] - 3.0).abs() < 1e-15); }
#[test]
fn test_sph_particle_aabb_empty() {
let sim = water_sim();
assert!(sim.particle_aabb().is_none());
}
#[test]
fn test_sph_remove_below_y() {
let mut sim = water_sim();
sim.add_particle([0.0, 1.0, 0.0]);
sim.add_particle([0.0, -1.0, 0.0]);
sim.add_particle([0.0, 0.5, 0.0]);
let removed = sim.remove_below_y(0.0);
assert_eq!(removed, 1);
assert_eq!(sim.particle_count(), 2);
}
#[test]
fn test_sph_collect_stats() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
sim.step(0.001);
let stats = sim.collect_stats();
assert_eq!(stats.particle_count, 1);
assert!(stats.mean_density > 0.0);
}
#[test]
fn test_sph_neighbors_out_of_bounds() {
let sim = water_sim();
assert!(sim.neighbors(99, 1.0).is_none());
}
#[test]
fn test_sph_neighbor_count_out_of_bounds() {
let sim = water_sim();
assert!(sim.neighbor_count(99, 1.0).is_none());
}
}
impl PySphSimulation {
pub fn adaptive_h(&self, i: usize, k: usize) -> Option<f64> {
let pos_i = self.positions.get(i)?;
let mut dists: Vec<f64> = self
.positions
.iter()
.enumerate()
.filter(|(j, _)| *j != i)
.map(|(_, pos_j)| {
let dx = pos_i[0] - pos_j[0];
let dy = pos_i[1] - pos_j[1];
let dz = pos_i[2] - pos_j[2];
(dx * dx + dy * dy + dz * dz).sqrt()
})
.collect();
if dists.len() < k {
return None;
}
dists.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
Some(dists[k - 1])
}
pub fn global_adaptive_h(&self, k: usize) -> Option<f64> {
let n = self.positions.len();
if n < 2 {
return None;
}
let mut hs: Vec<f64> = (0..n).filter_map(|i| self.adaptive_h(i, k)).collect();
if hs.is_empty() {
return None;
}
hs.sort_by(|a, b| a.partial_cmp(b).unwrap_or(std::cmp::Ordering::Equal));
let mid = hs.len() / 2;
Some(hs[mid])
}
}
impl PySphSimulation {
pub fn detect_surface_particles(&self, r: f64, threshold: usize) -> Vec<bool> {
(0..self.positions.len())
.map(|i| {
let count = self.neighbor_count(i, r).unwrap_or(0);
count < threshold
})
.collect()
}
pub fn surface_particle_indices(&self, r: f64, threshold: usize) -> Vec<usize> {
self.detect_surface_particles(r, threshold)
.into_iter()
.enumerate()
.filter(|(_, is_surface)| *is_surface)
.map(|(i, _)| i)
.collect()
}
}
#[derive(Debug, Clone, Copy, PartialEq, Eq, serde::Serialize, serde::Deserialize)]
pub enum SphVariant {
WCSPH,
DFSPH,
}
#[allow(dead_code)]
pub struct SphSolver {
pub sim: PySphSimulation,
pub variant: SphVariant,
pub max_pressure_iters: usize,
pub pressure_tol: f64,
}
impl SphSolver {
pub fn wcsph(config: PySphConfig) -> Self {
Self {
sim: PySphSimulation::new(config),
variant: SphVariant::WCSPH,
max_pressure_iters: 50,
pressure_tol: 0.01,
}
}
pub fn dfsph(config: PySphConfig) -> Self {
Self {
sim: PySphSimulation::new(config),
variant: SphVariant::DFSPH,
max_pressure_iters: 100,
pressure_tol: 0.001,
}
}
pub fn step(&mut self, dt: f64) {
match self.variant {
SphVariant::WCSPH => self.sim.step(dt),
SphVariant::DFSPH => self.step_dfsph(dt),
}
}
fn step_dfsph(&mut self, dt: f64) {
let sub_steps = self.max_pressure_iters.clamp(1, 5);
let sub_dt = dt / sub_steps as f64;
for _ in 0..sub_steps {
self.sim.step(sub_dt);
}
}
pub fn density_error(&self) -> f64 {
let rho0 = self.sim.config.rest_density;
self.sim
.densities
.iter()
.map(|&rho| (rho - rho0).abs() / rho0)
.fold(0.0f64, f64::max)
}
pub fn particle_count(&self) -> usize {
self.sim.particle_count()
}
}
#[allow(dead_code)]
pub struct SphParticleSet {
pub name: String,
pub indices: Vec<usize>,
}
impl SphParticleSet {
pub fn new(name: impl Into<String>) -> Self {
Self {
name: name.into(),
indices: Vec::new(),
}
}
pub fn add_index(&mut self, idx: usize) {
self.indices.push(idx);
}
pub fn len(&self) -> usize {
self.indices.len()
}
pub fn is_empty(&self) -> bool {
self.indices.is_empty()
}
}
impl PySphSimulation {
pub fn add_particles(&mut self, positions: &[[f64; 3]]) -> std::ops::Range<usize> {
let start = self.positions.len();
for &pos in positions {
self.add_particle(pos);
}
start..self.positions.len()
}
pub fn set_velocity(&mut self, i: usize, vel: [f64; 3]) {
if i < self.velocities.len() {
self.velocities[i] = vel;
}
}
pub fn apply_global_impulse(&mut self, dv: [f64; 3]) {
for v in &mut self.velocities {
v[0] += dv[0];
v[1] += dv[1];
v[2] += dv[2];
}
}
pub fn clamp_velocities(&mut self, max_speed: f64) {
for v in &mut self.velocities {
let speed = (v[0] * v[0] + v[1] * v[1] + v[2] * v[2]).sqrt();
if speed > max_speed && speed > 1e-15 {
let scale = max_speed / speed;
v[0] *= scale;
v[1] *= scale;
v[2] *= scale;
}
}
}
pub fn total_momentum(&self) -> [f64; 3] {
let m = self.config.particle_mass;
let mut p = [0.0f64; 3];
for v in &self.velocities {
p[0] += m * v[0];
p[1] += m * v[1];
p[2] += m * v[2];
}
p
}
pub fn density_variance(&self) -> f64 {
if self.densities.is_empty() {
return 0.0;
}
let mean = self.mean_density();
self.densities
.iter()
.map(|&d| {
let diff = d - mean;
diff * diff
})
.sum::<f64>()
/ self.densities.len() as f64
}
}
#[cfg(test)]
mod sph_new_tests {
use crate::PySphConfig;
use crate::PySphSimulation;
use crate::sph_api::SphParticleSet;
use crate::sph_api::SphSolver;
use crate::sph_api::SphVariant;
fn water_sim() -> PySphSimulation {
PySphSimulation::new(PySphConfig::water())
}
#[test]
fn test_adaptive_h_basic() {
let mut sim = water_sim();
for i in 0..5 {
sim.add_particle([i as f64 * 0.1, 0.0, 0.0]);
}
let h = sim.adaptive_h(0, 3);
assert!(h.is_some());
assert!(h.unwrap() > 0.0);
}
#[test]
fn test_adaptive_h_not_enough_neighbors() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
sim.add_particle([1.0, 0.0, 0.0]);
let h = sim.adaptive_h(0, 5);
assert!(h.is_none());
}
#[test]
fn test_adaptive_h_out_of_bounds() {
let sim = water_sim();
assert!(sim.adaptive_h(99, 3).is_none());
}
#[test]
fn test_global_adaptive_h() {
let mut sim = water_sim();
for i in 0..10 {
sim.add_particle([i as f64 * 0.1, 0.0, 0.0]);
}
let h = sim.global_adaptive_h(3);
assert!(h.is_some());
}
#[test]
fn test_global_adaptive_h_single_particle() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
assert!(sim.global_adaptive_h(3).is_none());
}
#[test]
fn test_detect_surface_all_surface_when_isolated() {
let mut sim = water_sim();
sim.add_particle([0.0, 0.0, 0.0]);
sim.add_particle([10.0, 0.0, 0.0]); let surface = sim.detect_surface_particles(0.5, 2);
assert!(surface[0]); }
#[test]
fn test_detect_surface_not_surface_when_crowded() {
let mut sim = water_sim();
for i in 0..8 {
for j in 0..8 {
sim.add_particle([i as f64 * 0.05, j as f64 * 0.05, 0.0]);
}
}
let surface = sim.detect_surface_particles(0.3, 2);
let n_surface = surface.iter().filter(|&&s| s).count();
assert!(n_surface < sim.particle_count()); }
#[test]
fn test_surface_particle_indices_empty() {
let sim = water_sim();
let idx = sim.surface_particle_indices(1.0, 5);
assert!(idx.is_empty());
}
#[test]
fn test_sph_variant_wcsph() {
let mut solver = SphSolver::wcsph(PySphConfig::water());
solver.sim.add_particle([0.0, 1.0, 0.0]);
solver.step(0.001);
assert_eq!(solver.sim.step_count(), 1);
assert_eq!(solver.variant, SphVariant::WCSPH);
}
#[test]
fn test_sph_variant_dfsph() {
let mut solver = SphSolver::dfsph(PySphConfig::water());
solver.sim.add_particle([0.0, 1.0, 0.0]);
solver.step(0.001);
assert_eq!(solver.sim.step_count(), 5);
assert_eq!(solver.variant, SphVariant::DFSPH);
}
#[test]
fn test_density_error_at_rest() {
let mut solver = SphSolver::wcsph(PySphConfig::water());
solver.sim.add_particle([0.0; 3]);
let err = solver.density_error();
assert!(err.is_finite());
}
#[test]
fn test_sph_particle_set_empty() {
let set = SphParticleSet::new("fluid");
assert!(set.is_empty());
assert_eq!(set.len(), 0);
}
#[test]
fn test_sph_particle_set_add() {
let mut set = SphParticleSet::new("fluid");
set.add_index(0);
set.add_index(1);
assert_eq!(set.len(), 2);
}
#[test]
fn test_add_particles_bulk() {
let mut sim = water_sim();
let positions = vec![[0.0; 3], [1.0, 0.0, 0.0], [2.0, 0.0, 0.0]];
let range = sim.add_particles(&positions);
assert_eq!(range.len(), 3);
assert_eq!(sim.particle_count(), 3);
}
#[test]
fn test_set_velocity() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
sim.set_velocity(0, [1.0, 2.0, 3.0]);
let v = sim.velocity(0).unwrap();
assert!((v[0] - 1.0).abs() < 1e-15);
}
#[test]
fn test_set_velocity_out_of_bounds_no_panic() {
let mut sim = water_sim();
sim.set_velocity(99, [1.0, 0.0, 0.0]); }
#[test]
fn test_apply_global_impulse() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
sim.apply_global_impulse([10.0, 0.0, 0.0]);
let v = sim.velocity(0).unwrap();
assert!((v[0] - 10.0).abs() < 1e-15);
}
#[test]
fn test_clamp_velocities() {
let mut sim = water_sim();
sim.add_particle_with_velocity([0.0; 3], [100.0, 0.0, 0.0]);
sim.clamp_velocities(5.0);
let v = sim.velocity(0).unwrap();
let speed = (v[0] * v[0] + v[1] * v[1] + v[2] * v[2]).sqrt();
assert!((speed - 5.0).abs() < 1e-10);
}
#[test]
fn test_total_momentum_at_rest() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
let p = sim.total_momentum();
for &pi in &p {
assert!(pi.abs() < 1e-15);
}
}
#[test]
fn test_total_momentum_with_velocity() {
let mut sim = water_sim();
sim.add_particle_with_velocity([0.0; 3], [2.0, 0.0, 0.0]);
let p = sim.total_momentum();
assert!((p[0] - 0.04).abs() < 1e-12);
}
#[test]
fn test_density_variance_zero_uniform() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]);
sim.add_particle([100.0, 0.0, 0.0]); let var = sim.density_variance();
assert!((var).abs() < 1e-10);
}
#[test]
fn test_density_variance_empty() {
let sim = water_sim();
assert!((sim.density_variance()).abs() < 1e-15);
}
#[test]
fn test_sph_solver_particle_count() {
let mut solver = SphSolver::wcsph(PySphConfig::water());
solver.sim.add_particle([0.0; 3]);
assert_eq!(solver.particle_count(), 1);
}
#[test]
fn test_add_particles_returns_correct_range() {
let mut sim = water_sim();
sim.add_particle([0.0; 3]); let positions = vec![[1.0, 0.0, 0.0], [2.0, 0.0, 0.0]];
let range = sim.add_particles(&positions);
assert_eq!(range.start, 1);
assert_eq!(range.end, 3);
}
#[test]
fn test_surface_detection_length_matches_particle_count() {
let mut sim = water_sim();
sim.add_particle_block([0.0; 3], 3, 3, 3, 0.1);
let surface = sim.detect_surface_particles(0.15, 3);
assert_eq!(surface.len(), sim.particle_count());
}
#[test]
fn test_clamp_velocities_already_slow() {
let mut sim = water_sim();
sim.add_particle_with_velocity([0.0; 3], [1.0, 0.0, 0.0]);
sim.clamp_velocities(10.0); let v = sim.velocity(0).unwrap();
assert!((v[0] - 1.0).abs() < 1e-14);
}
}