1use glam::Vec3;
15use std::collections::HashMap;
16
17const PI: f32 = std::f32::consts::PI;
20const TWO_PI: f32 = 2.0 * PI;
21
22#[inline]
26pub fn cubic_kernel(r: f32, h: f32) -> f32 {
27 let q = r / h;
28 let sigma = 1.0 / (PI * h * h * h);
29 if q <= 1.0 {
30 sigma * (1.0 - 1.5 * q * q + 0.75 * q * q * q)
31 } else if q <= 2.0 {
32 sigma * 0.25 * (2.0 - q).powi(3)
33 } else {
34 0.0
35 }
36}
37
38#[inline]
40pub fn cubic_kernel_grad(r: f32, h: f32) -> f32 {
41 let q = r / h;
42 let sigma = 1.0 / (PI * h * h * h * h); if q <= 1.0 {
44 sigma * (-3.0 * q + 2.25 * q * q)
45 } else if q <= 2.0 {
46 sigma * -0.75 * (2.0 - q).powi(2)
47 } else {
48 0.0
49 }
50}
51
52#[inline]
54pub fn kernel_gradient(r_vec: Vec3, h: f32) -> Vec3 {
55 let r = r_vec.length();
56 if r < 1e-8 { return Vec3::ZERO; }
57 let dw_dr = cubic_kernel_grad(r, h);
58 r_vec / r * dw_dr
59}
60
61pub struct DensityGrid {
65 cell_size: f32,
66 cells: HashMap<(i32, i32, i32), Vec<usize>>,
67}
68
69impl DensityGrid {
70 pub fn new(cell_size: f32) -> Self {
71 Self { cell_size, cells: HashMap::new() }
72 }
73
74 fn cell_of(&self, pos: Vec3) -> (i32, i32, i32) {
75 (
76 (pos.x / self.cell_size).floor() as i32,
77 (pos.y / self.cell_size).floor() as i32,
78 (pos.z / self.cell_size).floor() as i32,
79 )
80 }
81
82 pub fn rebuild(&mut self, positions: &[Vec3]) {
84 self.cells.clear();
85 for (i, &pos) in positions.iter().enumerate() {
86 self.cells.entry(self.cell_of(pos)).or_default().push(i);
87 }
88 }
89
90 pub fn query_radius(&self, pos: Vec3, radius: f32) -> Vec<usize> {
92 let cx = (pos.x / self.cell_size).floor() as i32;
93 let cy = (pos.y / self.cell_size).floor() as i32;
94 let cz = (pos.z / self.cell_size).floor() as i32;
95 let cells_needed = (radius / self.cell_size).ceil() as i32 + 1;
96 let radius_sq = radius * radius;
97 let mut result = Vec::new();
98 for dx in -cells_needed..=cells_needed {
99 for dy in -cells_needed..=cells_needed {
100 for dz in -cells_needed..=cells_needed {
101 let key = (cx + dx, cy + dy, cz + dz);
102 if let Some(bucket) = self.cells.get(&key) {
103 for &idx in bucket {
104 result.push(idx);
106 }
107 }
108 }
109 }
110 }
111 result
112 }
113
114 pub fn particle_count(&self) -> usize {
115 self.cells.values().map(|v| v.len()).sum()
116 }
117}
118
119#[derive(Debug, Clone)]
123pub struct SphParticle {
124 pub position: Vec3,
125 pub velocity: Vec3,
126 pub rest_density: f32,
128 pub density: f32,
130 pub pressure: f32,
132 pub mass: f32,
134 pub accel: Vec3,
136 pub neighbors: Vec<usize>,
138 pub is_boundary: bool,
140 pub tag: u32,
142}
143
144impl SphParticle {
145 pub fn new(position: Vec3, mass: f32, rest_density: f32) -> Self {
146 Self {
147 position, velocity: Vec3::ZERO, rest_density, density: rest_density,
148 pressure: 0.0, mass, accel: Vec3::ZERO, neighbors: Vec::new(),
149 is_boundary: false, tag: 0,
150 }
151 }
152
153 pub fn boundary(position: Vec3, mass: f32, rest_density: f32) -> Self {
154 Self { is_boundary: true, ..Self::new(position, mass, rest_density) }
155 }
156}
157
158pub struct SphSimulation {
165 pub particles: Vec<SphParticle>,
166 pub h: f32,
168 pub rest_density: f32,
170 pub stiffness: f32,
172 pub gamma: f32,
174 pub viscosity: f32,
176 pub surface_tension: f32,
178 pub gravity: Vec3,
180 pub max_dt: f32,
182 pub min_dt: f32,
184 pub dt: f32,
186 pub cfl_factor: f32,
188 pub speed_of_sound: f32,
190 pub xsph_epsilon: f32,
192 pub boundary_repulsion: f32,
194 pub domain_min: Vec3,
196 pub domain_max: Vec3,
197 grid: DensityGrid,
198 time: f32,
199}
200
201impl SphSimulation {
202 pub fn new(h: f32, rest_density: f32) -> Self {
203 Self {
204 particles: Vec::new(),
205 h, rest_density,
206 stiffness: 100.0,
207 gamma: 7.0,
208 viscosity: 0.01,
209 surface_tension: 0.01,
210 gravity: Vec3::new(0.0, -9.81, 0.0),
211 max_dt: 0.005,
212 min_dt: 0.00001,
213 dt: 0.001,
214 cfl_factor: 0.4,
215 speed_of_sound: 88.5, xsph_epsilon: 0.5,
217 boundary_repulsion: 500.0,
218 domain_min: Vec3::new(-5.0, -5.0, -5.0),
219 domain_max: Vec3::new( 5.0, 5.0, 5.0),
220 grid: DensityGrid::new(h * 2.0),
221 time: 0.0,
222 }
223 }
224
225 pub fn add_cube(
227 &mut self, center: Vec3, half_extent: Vec3,
228 spacing: f32, mass: f32,
229 ) {
230 let nx = (2.0 * half_extent.x / spacing).ceil() as i32;
231 let ny = (2.0 * half_extent.y / spacing).ceil() as i32;
232 let nz = (2.0 * half_extent.z / spacing).ceil() as i32;
233 for ix in 0..nx {
234 for iy in 0..ny {
235 for iz in 0..nz {
236 let x = center.x - half_extent.x + ix as f32 * spacing;
237 let y = center.y - half_extent.y + iy as f32 * spacing;
238 let z = center.z - half_extent.z + iz as f32 * spacing;
239 self.particles.push(SphParticle::new(Vec3::new(x, y, z), mass, self.rest_density));
240 }
241 }
242 }
243 }
244
245 fn rebuild_grid(&mut self) {
247 let positions: Vec<Vec3> = self.particles.iter().map(|p| p.position).collect();
248 self.grid.rebuild(&positions);
249 }
250
251 fn find_neighbors(&mut self) {
253 let h2 = self.h * self.h * 4.0; let positions: Vec<Vec3> = self.particles.iter().map(|p| p.position).collect();
255 let n = self.particles.len();
256 for i in 0..n {
257 self.particles[i].neighbors.clear();
258 let candidates = self.grid.query_radius(positions[i], self.h * 2.0);
259 for j in candidates {
260 if j == i { continue; }
261 let r2 = (positions[j] - positions[i]).length_squared();
262 if r2 < h2 {
263 self.particles[i].neighbors.push(j);
264 }
265 }
266 }
267 }
268
269 fn compute_densities(&mut self) {
271 let n = self.particles.len();
272 for i in 0..n {
273 let neighbors = self.particles[i].neighbors.clone();
274 let pos_i = self.particles[i].position;
275 let mut rho = self.particles[i].mass * cubic_kernel(0.0, self.h);
276 for j in &neighbors {
277 let r = (pos_i - self.particles[*j].position).length();
278 rho += self.particles[*j].mass * cubic_kernel(r, self.h);
279 }
280 self.particles[i].density = rho.max(1e-6);
281 }
282 }
283
284 fn compute_pressures(&mut self) {
287 let b = self.stiffness * self.rest_density * self.speed_of_sound * self.speed_of_sound / self.gamma;
288 for p in &mut self.particles {
289 let ratio = p.density / p.rest_density;
290 p.pressure = b * (ratio.powf(self.gamma) - 1.0);
291 if p.pressure < 0.0 { p.pressure = 0.0; } }
293 }
294
295 fn compute_pressure_forces(&mut self) -> Vec<Vec3> {
297 let n = self.particles.len();
298 let mut forces = vec![Vec3::ZERO; n];
299 for i in 0..n {
300 if self.particles[i].is_boundary { continue; }
301 let neighbors = self.particles[i].neighbors.clone();
302 let pos_i = self.particles[i].position;
303 let rho_i = self.particles[i].density;
304 let p_i = self.particles[i].pressure;
305 let m_i = self.particles[i].mass;
306 for j in neighbors {
307 let pos_j = self.particles[j].position;
308 let r_vec = pos_i - pos_j;
309 let grad = kernel_gradient(r_vec, self.h);
310 let rho_j = self.particles[j].density;
311 let p_j = self.particles[j].pressure;
312 let m_j = self.particles[j].mass;
313 let coeff = -m_j * (p_i / (rho_i * rho_i) + p_j / (rho_j * rho_j));
315 forces[i] += grad * (coeff * m_i);
316 }
317 }
318 forces
319 }
320
321 fn compute_viscosity_forces(&mut self) -> Vec<Vec3> {
323 let n = self.particles.len();
324 let mut forces = vec![Vec3::ZERO; n];
325 for i in 0..n {
326 if self.particles[i].is_boundary { continue; }
327 let neighbors = self.particles[i].neighbors.clone();
328 let pos_i = self.particles[i].position;
329 let vel_i = self.particles[i].velocity;
330 let rho_i = self.particles[i].density;
331 let m_i = self.particles[i].mass;
332 for j in neighbors {
333 let pos_j = self.particles[j].position;
334 let vel_j = self.particles[j].velocity;
335 let r_vec = pos_i - pos_j;
336 let r = r_vec.length();
337 if r < 1e-8 { continue; }
338 let vel_diff = vel_j - vel_i;
339 let rho_j = self.particles[j].density;
340 let m_j = self.particles[j].mass;
341 let v_dot_r = vel_diff.dot(r_vec);
343 if v_dot_r < 0.0 {
344 let mu = 2.0 * self.viscosity * self.h;
345 let pi_ij = -mu * v_dot_r / (rho_i + rho_j) * 0.5 * (r * r + 0.01 * self.h * self.h).recip();
346 let grad = kernel_gradient(r_vec, self.h);
347 forces[i] += grad * (-m_j * pi_ij * m_i);
348 }
349 }
350 }
351 forces
352 }
353
354 fn compute_surface_tension_forces(&mut self) -> Vec<Vec3> {
356 let n = self.particles.len();
357 let mut forces = vec![Vec3::ZERO; n];
358 let k = self.surface_tension;
359 for i in 0..n {
360 if self.particles[i].is_boundary { continue; }
361 let neighbors = self.particles[i].neighbors.clone();
362 let pos_i = self.particles[i].position;
363 let m_i = self.particles[i].mass;
364 for j in neighbors {
365 let pos_j = self.particles[j].position;
366 let r_vec = pos_i - pos_j;
367 let r = r_vec.length();
368 if r < 1e-8 || r >= 2.0 * self.h { continue; }
369 let c_rh = cohesion_kernel(r, self.h);
371 let m_j = self.particles[j].mass;
372 let dir = r_vec / r;
373 forces[i] -= dir * (k * m_i * m_j * c_rh);
374 }
375 }
376 forces
377 }
378
379 fn compute_boundary_forces(&mut self) -> Vec<Vec3> {
381 let n = self.particles.len();
382 let mut forces = vec![Vec3::ZERO; n];
383 for i in 0..n {
384 if self.particles[i].is_boundary { continue; }
385 let pos = self.particles[i].position;
386 let m_i = self.particles[i].mass;
387 let k = self.boundary_repulsion * m_i;
389 for axis in 0..3usize {
390 let lo = match axis { 0 => self.domain_min.x, 1 => self.domain_min.y, _ => self.domain_min.z };
391 let hi = match axis { 0 => self.domain_max.x, 1 => self.domain_max.y, _ => self.domain_max.z };
392 let val = match axis { 0 => pos.x, 1 => pos.y, _ => pos.z };
393 let d_lo = val - lo;
394 let d_hi = hi - val;
395 if d_lo < self.h && d_lo > 0.0 {
396 let f = k * (self.h - d_lo) / self.h;
397 match axis { 0 => forces[i].x += f, 1 => forces[i].y += f, _ => forces[i].z += f }
398 }
399 if d_hi < self.h && d_hi > 0.0 {
400 let f = k * (self.h - d_hi) / self.h;
401 match axis { 0 => forces[i].x -= f, 1 => forces[i].y -= f, _ => forces[i].z -= f }
402 }
403 }
404 }
405 forces
406 }
407
408 fn compute_cfl_dt(&self) -> f32 {
410 let max_vel = self.particles.iter()
411 .map(|p| p.velocity.length())
412 .fold(0.0_f32, f32::max);
413 let max_speed = max_vel + self.speed_of_sound;
414 if max_speed < 1e-6 { return self.max_dt; }
415 (self.cfl_factor * self.h / max_speed).clamp(self.min_dt, self.max_dt)
416 }
417
418 fn xsph_correction(&mut self) {
420 let n = self.particles.len();
421 let mut corrections = vec![Vec3::ZERO; n];
422 for i in 0..n {
423 if self.particles[i].is_boundary { continue; }
424 let neighbors = self.particles[i].neighbors.clone();
425 let pos_i = self.particles[i].position;
426 let vel_i = self.particles[i].velocity;
427 let rho_i = self.particles[i].density;
428 for j in neighbors {
429 let r = (pos_i - self.particles[j].position).length();
430 let w = cubic_kernel(r, self.h);
431 let rho_j = self.particles[j].density;
432 let vel_j = self.particles[j].velocity;
433 let m_j = self.particles[j].mass;
434 corrections[i] += (vel_j - vel_i) * (2.0 * m_j / (rho_i + rho_j) * w);
435 }
436 }
437 for (i, p) in self.particles.iter_mut().enumerate() {
438 if !p.is_boundary {
439 p.velocity += corrections[i] * self.xsph_epsilon;
440 }
441 }
442 }
443
444 fn enforce_boundaries(&mut self) {
446 for p in &mut self.particles {
447 if p.is_boundary { continue; }
448 for axis in 0..3usize {
449 let lo = match axis { 0 => self.domain_min.x, 1 => self.domain_min.y, _ => self.domain_min.z };
450 let hi = match axis { 0 => self.domain_max.x, 1 => self.domain_max.y, _ => self.domain_max.z };
451 let val = match axis { 0 => &mut p.position.x, 1 => &mut p.position.y, _ => &mut p.position.z };
452 let vel = match axis { 0 => &mut p.velocity.x, 1 => &mut p.velocity.y, _ => &mut p.velocity.z };
453 if *val < lo { *val = lo + 1e-4; *vel = vel.abs() * 0.5; }
454 if *val > hi { *val = hi - 1e-4; *vel = -vel.abs() * 0.5; }
455 }
456 }
457 }
458
459 pub fn step(&mut self) {
461 self.dt = self.compute_cfl_dt();
462 self.step_with_dt(self.dt);
463 }
464
465 pub fn step_with_dt(&mut self, dt: f32) {
467 self.rebuild_grid();
468 self.find_neighbors();
469 self.compute_densities();
470 self.compute_pressures();
471
472 let pf = self.compute_pressure_forces();
473 let vf = self.compute_viscosity_forces();
474 let sf = self.compute_surface_tension_forces();
475 let bf = self.compute_boundary_forces();
476
477 let gravity = self.gravity;
478 let n = self.particles.len();
479 for i in 0..n {
480 if self.particles[i].is_boundary { continue; }
481 let total_force = pf[i] + vf[i] + sf[i] + bf[i];
482 let rho = self.particles[i].density.max(1e-6);
483 self.particles[i].accel = total_force / rho + gravity;
484 }
485
486 for p in &mut self.particles {
488 if p.is_boundary { continue; }
489 p.velocity += p.accel * dt;
490 p.position += p.velocity * dt;
491 }
492
493 self.xsph_correction();
494 self.enforce_boundaries();
495 self.time += dt;
496 }
497
498 pub fn advance(&mut self, total_time: f32) {
500 let mut elapsed = 0.0;
501 while elapsed < total_time {
502 let dt = self.compute_cfl_dt().min(total_time - elapsed);
503 self.step_with_dt(dt);
504 elapsed += dt;
505 }
506 }
507
508 pub fn particle_count(&self) -> usize { self.particles.len() }
509
510 pub fn total_kinetic_energy(&self) -> f32 {
512 self.particles.iter().map(|p| 0.5 * p.mass * p.velocity.length_squared()).sum()
513 }
514
515 pub fn average_density(&self) -> f32 {
517 let n = self.particles.len();
518 if n == 0 { return 0.0; }
519 self.particles.iter().map(|p| p.density).sum::<f32>() / n as f32
520 }
521
522 pub fn sample_pressure(&self, pos: Vec3) -> f32 {
524 let candidates = self.grid.query_radius(pos, self.h * 2.0);
525 let mut p_sum = 0.0;
526 let mut w_sum = 0.0;
527 for idx in candidates {
528 let r = (pos - self.particles[idx].position).length();
529 if r >= 2.0 * self.h { continue; }
530 let w = cubic_kernel(r, self.h);
531 p_sum += self.particles[idx].pressure / self.particles[idx].density.max(1e-6) * w;
532 w_sum += w;
533 }
534 if w_sum > 1e-8 { p_sum / w_sum } else { 0.0 }
535 }
536}
537
538fn cohesion_kernel(r: f32, h: f32) -> f32 {
540 let sigma = 32.0 / (PI * h.powi(9));
541 let h2 = h * h;
542 let r2 = r * r;
543 if r < h * 0.5 {
544 sigma * 2.0 * (h - r).powi(3) * r * r * r - h2 * h2 * h2 / 64.0
545 } else if r < h {
546 sigma * (h - r).powi(3) * r * r * r
547 } else {
548 0.0
549 }
550}
551
552#[derive(Debug, Clone)]
556pub struct FluidGlyph {
557 pub position: Vec3,
558 pub color: [f32; 4], pub scale: f32,
560 pub density: f32,
562}
563
564pub struct FluidRenderer {
566 pub iso_threshold: f32,
568 pub base_color: [f32; 4],
570 pub dense_color: [f32; 4],
572 pub glyph_scale: f32,
574 pub show_interior: bool,
576}
577
578impl FluidRenderer {
579 pub fn new() -> Self {
580 Self {
581 iso_threshold: 0.7,
582 base_color: [0.2, 0.5, 1.0, 0.8],
583 dense_color: [0.0, 0.1, 0.8, 1.0],
584 glyph_scale: 0.05,
585 show_interior: false,
586 }
587 }
588
589 pub fn with_color(mut self, color: [f32; 4]) -> Self { self.base_color = color; self }
590 pub fn with_iso_threshold(mut self, t: f32) -> Self { self.iso_threshold = t; self }
591 pub fn with_scale(mut self, s: f32) -> Self { self.glyph_scale = s; self }
592
593 pub fn extract(&self, sim: &SphSimulation) -> Vec<FluidGlyph> {
595 let mut glyphs = Vec::new();
596 let rho0 = sim.rest_density;
597
598 for p in &sim.particles {
599 if p.is_boundary { continue; }
600
601 let normalized_density = (p.density / rho0).clamp(0.0, 3.0);
602 let is_surface = normalized_density < self.iso_threshold + 0.3
603 && normalized_density > self.iso_threshold - 0.3;
604
605 if !is_surface && !self.show_interior { continue; }
606
607 let t = ((normalized_density - 1.0) * 0.5).clamp(0.0, 1.0);
609 let color = lerp_color(self.base_color, self.dense_color, t);
610
611 glyphs.push(FluidGlyph {
612 position: p.position,
613 color,
614 scale: self.glyph_scale * (0.8 + 0.4 * normalized_density.min(2.0)),
615 density: p.density,
616 });
617 }
618 glyphs
619 }
620
621 pub fn extract_region(&self, sim: &SphSimulation, min: Vec3, max: Vec3) -> Vec<FluidGlyph> {
623 self.extract(sim).into_iter().filter(|g| {
624 g.position.x >= min.x && g.position.x <= max.x &&
625 g.position.y >= min.y && g.position.y <= max.y &&
626 g.position.z >= min.z && g.position.z <= max.z
627 }).collect()
628 }
629}
630
631impl Default for FluidRenderer {
632 fn default() -> Self { Self::new() }
633}
634
635fn lerp_color(a: [f32; 4], b: [f32; 4], t: f32) -> [f32; 4] {
636 [
637 a[0] + (b[0] - a[0]) * t,
638 a[1] + (b[1] - a[1]) * t,
639 a[2] + (b[2] - a[2]) * t,
640 a[3] + (b[3] - a[3]) * t,
641 ]
642}
643
644pub struct WaterSurface {
654 pub width: usize,
656 pub depth: usize,
658 pub cell_size: f32,
660 pub height: Vec<f32>,
662 pub vel_x: Vec<f32>,
664 pub vel_z: Vec<f32>,
666 height_new: Vec<f32>,
668 vel_x_new: Vec<f32>,
669 vel_z_new: Vec<f32>,
670 pub rest_height: f32,
672 pub wave_speed: f32,
674 pub damping: f32,
676 pub gravity: f32,
678 time: f32,
679}
680
681impl WaterSurface {
682 pub fn new(width: usize, depth: usize, cell_size: f32, rest_height: f32) -> Self {
683 let n = width * depth;
684 Self {
685 width, depth, cell_size, rest_height,
686 height: vec![rest_height; n],
687 vel_x: vec![0.0; n],
688 vel_z: vec![0.0; n],
689 height_new: vec![rest_height; n],
690 vel_x_new: vec![0.0; n],
691 vel_z_new: vec![0.0; n],
692 wave_speed: 5.0,
693 damping: 0.005,
694 gravity: 9.81,
695 time: 0.0,
696 }
697 }
698
699 #[inline]
700 fn idx(&self, x: usize, z: usize) -> usize { z * self.width + x }
701
702 #[inline]
703 fn clamp_x(&self, x: i32) -> usize { x.clamp(0, self.width as i32 - 1) as usize }
704 #[inline]
705 fn clamp_z(&self, z: i32) -> usize { z.clamp(0, self.depth as i32 - 1) as usize }
706
707 pub fn rain_drop(&mut self, world_x: f32, world_z: f32, amplitude: f32, radius: f32) {
709 let cx = (world_x / self.cell_size) as i32;
710 let cz = (world_z / self.cell_size) as i32;
711 let r_cells = (radius / self.cell_size).ceil() as i32 + 1;
712 for dz in -r_cells..=r_cells {
713 for dx in -r_cells..=r_cells {
714 let ix = cx + dx;
715 let iz = cz + dz;
716 if ix < 0 || ix >= self.width as i32 || iz < 0 || iz >= self.depth as i32 { continue; }
717 let world_dist = ((dx as f32 * self.cell_size).powi(2) + (dz as f32 * self.cell_size).powi(2)).sqrt();
718 let gaussian = amplitude * (-world_dist * world_dist / (2.0 * radius * radius)).exp();
719 let idx = self.idx(ix as usize, iz as usize);
720 self.height[idx] += gaussian;
721 }
722 }
723 }
724
725 pub fn add_wave_source(&mut self, z_start: usize, z_end: usize, amplitude: f32, frequency: f32) {
727 let t = self.time;
728 for iz in z_start.min(self.depth - 1)..=z_end.min(self.depth - 1) {
729 let idx = self.idx(0, iz);
730 self.height[idx] = self.rest_height + amplitude * (TWO_PI * frequency * t).sin();
731 }
732 }
733
734 pub fn step(&mut self, dt: f32) {
736 let g = self.gravity;
737 let cs = self.cell_size;
738 let w = self.width;
739 let d = self.depth;
740
741 self.height_new.copy_from_slice(&self.height);
743 self.vel_x_new.copy_from_slice(&self.vel_x);
744 self.vel_z_new.copy_from_slice(&self.vel_z);
745
746 for iz in 0..d {
747 for ix in 0..w {
748 let i = self.idx(ix, iz);
749 let h = self.height[i];
750 let vx = self.vel_x[i];
751 let vz = self.vel_z[i];
752
753 let ix_p = self.clamp_x(ix as i32 + 1);
755 let ix_m = self.clamp_x(ix as i32 - 1);
756 let iz_p = self.clamp_z(iz as i32 + 1);
757 let iz_m = self.clamp_z(iz as i32 - 1);
758
759 let h_xp = self.height[self.idx(ix_p, iz)];
760 let h_xm = self.height[self.idx(ix_m, iz)];
761 let h_zp = self.height[self.idx(ix, iz_p)];
762 let h_zm = self.height[self.idx(ix, iz_m)];
763
764 let dh_dx = (h_xp - h_xm) / (2.0 * cs);
766 let dh_dz = (h_zp - h_zm) / (2.0 * cs);
767
768 let new_vx = vx - g * dh_dx * dt;
770 let new_vz = vz - g * dh_dz * dt;
771
772 let vx_p = self.vel_x[self.idx(ix_p, iz)];
774 let vx_m = self.vel_x[self.idx(ix_m, iz)];
775 let vz_p = self.vel_z[self.idx(ix, iz_p)];
776 let vz_m = self.vel_z[self.idx(ix, iz_m)];
777
778 let div_flux_x = (h * vx_p - h * vx_m) / (2.0 * cs);
779 let div_flux_z = (h * vz_p - h * vz_m) / (2.0 * cs);
780 let new_h = (h - dt * (div_flux_x + div_flux_z)).max(0.0);
781
782 self.height_new[i] = new_h * (1.0 - self.damping) + self.rest_height * self.damping;
783 self.vel_x_new[i] = new_vx * (1.0 - self.damping);
784 self.vel_z_new[i] = new_vz * (1.0 - self.damping);
785 }
786 }
787
788 std::mem::swap(&mut self.height, &mut self.height_new);
789 std::mem::swap(&mut self.vel_x, &mut self.vel_x_new);
790 std::mem::swap(&mut self.vel_z, &mut self.vel_z_new);
791
792 self.time += dt;
793 }
794
795 pub fn sample_height(&self, world_x: f32, world_z: f32) -> f32 {
797 let fx = world_x / self.cell_size;
798 let fz = world_z / self.cell_size;
799 let ix0 = (fx.floor() as i32).clamp(0, self.width as i32 - 2) as usize;
800 let iz0 = (fz.floor() as i32).clamp(0, self.depth as i32 - 2) as usize;
801 let tx = (fx - ix0 as f32).clamp(0.0, 1.0);
802 let tz = (fz - iz0 as f32).clamp(0.0, 1.0);
803 let h00 = self.height[self.idx(ix0, iz0)];
804 let h10 = self.height[self.idx(ix0 + 1, iz0)];
805 let h01 = self.height[self.idx(ix0, iz0 + 1)];
806 let h11 = self.height[self.idx(ix0 + 1, iz0 + 1)];
807 let a = h00 + tx * (h10 - h00);
808 let b = h01 + tx * (h11 - h01);
809 a + tz * (b - a)
810 }
811
812 pub fn surface_normal(&self, world_x: f32, world_z: f32) -> Vec3 {
814 let eps = self.cell_size;
815 let hx_p = self.sample_height(world_x + eps, world_z);
816 let hx_m = self.sample_height(world_x - eps, world_z);
817 let hz_p = self.sample_height(world_x, world_z + eps);
818 let hz_m = self.sample_height(world_x, world_z - eps);
819 let dx = (hx_p - hx_m) / (2.0 * eps);
820 let dz = (hz_p - hz_m) / (2.0 * eps);
821 Vec3::new(-dx, 1.0, -dz).normalize_or_zero()
822 }
823
824 pub fn to_positions(&self, y_offset: f32) -> Vec<Vec3> {
826 let mut out = Vec::with_capacity(self.width * self.depth);
827 for iz in 0..self.depth {
828 for ix in 0..self.width {
829 let h = self.height[self.idx(ix, iz)];
830 out.push(Vec3::new(
831 ix as f32 * self.cell_size,
832 h + y_offset,
833 iz as f32 * self.cell_size,
834 ));
835 }
836 }
837 out
838 }
839
840 pub fn max_amplitude(&self) -> f32 {
842 self.height.iter().map(|&h| (h - self.rest_height).abs()).fold(0.0_f32, f32::max)
843 }
844
845 pub fn total_volume(&self) -> f32 {
847 let cs2 = self.cell_size * self.cell_size;
848 self.height.iter().sum::<f32>() * cs2
849 }
850}
851
852pub struct BuoyancyForce {
859 pub fluid_density: f32,
861 pub gravity: f32,
863 pub drag: f32,
865}
866
867impl BuoyancyForce {
868 pub fn new(fluid_density: f32) -> Self {
869 Self { fluid_density, gravity: 9.81, drag: 0.5 }
870 }
871
872 pub fn compute(
875 &self,
876 position: Vec3,
877 velocity: Vec3,
878 half_extents: Vec3,
879 water_h: f32,
880 ) -> Vec3 {
881 let body_bottom = position.y - half_extents.y;
882 let body_top = position.y + half_extents.y;
883 let body_height = half_extents.y * 2.0;
884
885 if body_bottom >= water_h { return Vec3::ZERO; } let submerged_height = (water_h - body_bottom).min(body_height).max(0.0);
889 let submerged_fraction = submerged_height / body_height.max(1e-6);
890
891 let volume = 8.0 * half_extents.x * half_extents.y * half_extents.z;
893 let submerged_volume = volume * submerged_fraction;
894
895 let buoyancy = Vec3::Y * (self.fluid_density * self.gravity * submerged_volume);
897
898 let drag_force = Vec3::new(
900 -self.drag * velocity.x * submerged_fraction,
901 0.0,
902 -self.drag * velocity.z * submerged_fraction,
903 );
904
905 buoyancy + drag_force
906 }
907
908 pub fn apply_to_body(
911 &self,
912 position: Vec3,
913 velocity: Vec3,
914 half_extents: Vec3,
915 water_surface: &WaterSurface,
916 ) -> Vec3 {
917 let water_h = water_surface.sample_height(position.x, position.z);
918 self.compute(position, velocity, half_extents, water_h)
919 }
920
921 pub fn compute_with_torque(
924 &self,
925 position: Vec3,
926 velocity: Vec3,
927 half_extents: Vec3,
928 orientation: glam::Quat,
929 water_surface: &WaterSurface,
930 ) -> (Vec3, Vec3) {
931 let water_h = water_surface.sample_height(position.x, position.z);
932 let base_force = self.compute(position, velocity, half_extents, water_h);
933
934 let body_bottom = position.y - half_extents.y;
936 let submerged_height = (water_h - body_bottom).min(half_extents.y * 2.0).max(0.0);
937 let center_of_buoyancy = Vec3::new(
938 position.x,
939 body_bottom + submerged_height * 0.5,
940 position.z,
941 );
942
943 let r = center_of_buoyancy - position;
945 let torque = r.cross(base_force);
946
947 (base_force, torque)
948 }
949}
950
951#[cfg(test)]
954mod tests {
955 use super::*;
956
957 #[test]
960 fn cubic_kernel_zero_at_boundary() {
961 let w = cubic_kernel(1.0, 0.5); assert!(w.abs() < 1e-6, "kernel should be 0 at r >= h");
963 }
964
965 #[test]
966 fn cubic_kernel_positive_at_origin() {
967 let w = cubic_kernel(0.0, 0.5);
968 assert!(w > 0.0, "kernel should be positive at origin");
969 }
970
971 #[test]
972 fn kernel_gradient_zero_at_origin() {
973 let grad = kernel_gradient(Vec3::ZERO, 0.5);
974 assert!(grad.length() < 1e-6);
975 }
976
977 #[test]
980 fn density_grid_query_finds_nearby() {
981 let mut grid = DensityGrid::new(1.0);
982 let positions = vec![
983 Vec3::ZERO,
984 Vec3::new(0.5, 0.0, 0.0),
985 Vec3::new(10.0, 0.0, 0.0),
986 ];
987 grid.rebuild(&positions);
988 let candidates = grid.query_radius(Vec3::ZERO, 1.0);
989 assert!(candidates.contains(&0), "should find particle 0");
990 assert!(candidates.contains(&1), "should find particle 1");
991 let has_far = candidates.contains(&2);
992 let _ = has_far;
995 }
996
997 #[test]
1000 fn sph_particles_initialize() {
1001 let mut sim = SphSimulation::new(0.1, 1000.0);
1002 sim.add_cube(Vec3::ZERO, Vec3::splat(0.5), 0.1, 0.001);
1003 assert!(sim.particle_count() > 0);
1004 }
1005
1006 #[test]
1007 fn sph_step_does_not_nan() {
1008 let mut sim = SphSimulation::new(0.1, 1000.0);
1009 sim.domain_min = Vec3::new(-2.0, -2.0, -2.0);
1010 sim.domain_max = Vec3::new( 2.0, 2.0, 2.0);
1011 sim.add_cube(Vec3::ZERO, Vec3::new(0.3, 0.3, 0.3), 0.1, 0.001);
1012 sim.step_with_dt(0.001);
1013 for p in &sim.particles {
1014 assert!(p.position.is_finite(), "position has NaN");
1015 assert!(p.velocity.is_finite(), "velocity has NaN");
1016 }
1017 }
1018
1019 #[test]
1020 fn sph_gravity_pulls_down() {
1021 let mut sim = SphSimulation::new(0.1, 1000.0);
1022 sim.domain_min = Vec3::new(-5.0, -5.0, -5.0);
1023 sim.domain_max = Vec3::new( 5.0, 5.0, 5.0);
1024 sim.add_cube(Vec3::new(0.0, 2.0, 0.0), Vec3::new(0.2, 0.2, 0.2), 0.1, 0.001);
1025 let y_before: f32 = sim.particles.iter().map(|p| p.position.y).sum::<f32>() / sim.particle_count() as f32;
1026 for _ in 0..10 { sim.step_with_dt(0.005); }
1027 let y_after: f32 = sim.particles.iter().map(|p| p.position.y).sum::<f32>() / sim.particle_count() as f32;
1028 assert!(y_after < y_before, "particles should fall under gravity");
1029 }
1030
1031 #[test]
1032 fn sph_cfl_dt_is_bounded() {
1033 let mut sim = SphSimulation::new(0.1, 1000.0);
1034 sim.add_cube(Vec3::ZERO, Vec3::splat(0.2), 0.1, 0.001);
1035 let dt = sim.compute_cfl_dt();
1036 assert!(dt >= sim.min_dt && dt <= sim.max_dt, "CFL dt out of bounds: {dt}");
1037 }
1038
1039 #[test]
1040 fn sph_average_density_near_rest() {
1041 let mut sim = SphSimulation::new(0.1, 1000.0);
1042 sim.domain_min = Vec3::splat(-2.0);
1043 sim.domain_max = Vec3::splat( 2.0);
1044 sim.add_cube(Vec3::ZERO, Vec3::splat(0.4), 0.1, 0.001);
1045 sim.step_with_dt(0.001);
1046 let avg = sim.average_density();
1047 assert!(avg > 0.0, "average density should be positive");
1048 }
1049
1050 #[test]
1053 fn fluid_renderer_extracts_glyphs() {
1054 let mut sim = SphSimulation::new(0.1, 1000.0);
1055 sim.add_cube(Vec3::ZERO, Vec3::splat(0.3), 0.1, 0.001);
1056 sim.step_with_dt(0.001);
1057 let renderer = FluidRenderer::new().with_iso_threshold(0.5);
1058 let glyphs = renderer.extract(&sim);
1059 assert!(!glyphs.is_empty() || sim.particle_count() == 0, "renderer should produce glyphs or sim is empty");
1060 }
1061
1062 #[test]
1063 fn fluid_renderer_colors_finite() {
1064 let mut sim = SphSimulation::new(0.1, 1000.0);
1065 sim.add_cube(Vec3::ZERO, Vec3::splat(0.3), 0.1, 0.001);
1066 sim.step_with_dt(0.001);
1067 let renderer = FluidRenderer { show_interior: true, ..FluidRenderer::new() };
1068 let glyphs = renderer.extract(&sim);
1069 for g in &glyphs {
1070 for c in g.color {
1071 assert!(c.is_finite() && c >= 0.0 && c <= 1.0, "color component out of range: {c}");
1072 }
1073 }
1074 }
1075
1076 #[test]
1079 fn water_surface_rain_drop_raises_height() {
1080 let mut water = WaterSurface::new(32, 32, 0.25, 1.0);
1081 water.rain_drop(4.0, 4.0, 0.5, 0.5);
1082 let h = water.sample_height(4.0, 4.0);
1083 assert!(h > 1.0, "rain drop should raise height above rest");
1084 }
1085
1086 #[test]
1087 fn water_surface_wave_propagates() {
1088 let mut water = WaterSurface::new(64, 16, 0.1, 1.0);
1089 water.rain_drop(0.5, 0.8, 0.5, 0.3);
1090 let h0 = water.sample_height(3.0, 0.8);
1091 for _ in 0..100 { water.step(0.01); }
1092 let h1 = water.sample_height(3.0, 0.8);
1093 let changed = (h1 - h0).abs() > 0.0001;
1095 assert!(changed || water.max_amplitude() < 0.001, "wave should propagate or damp");
1096 }
1097
1098 #[test]
1099 fn water_surface_damps_to_rest() {
1100 let mut water = WaterSurface::new(16, 16, 0.25, 1.0);
1101 water.damping = 0.1;
1102 water.rain_drop(2.0, 2.0, 1.0, 0.5);
1103 for _ in 0..500 { water.step(0.01); }
1104 let amp = water.max_amplitude();
1105 assert!(amp < 0.1, "wave should damp, amplitude={amp}");
1106 }
1107
1108 #[test]
1109 fn water_surface_total_volume_stable() {
1110 let w = 32; let d = 32;
1111 let water = WaterSurface::new(w, d, 0.25, 1.0);
1112 let vol = water.total_volume();
1113 assert!(vol > 0.0, "volume should be positive");
1114 assert!(vol.is_finite());
1115 }
1116
1117 #[test]
1120 fn buoyancy_zero_above_water() {
1121 let bf = BuoyancyForce::new(1000.0);
1122 let f = bf.compute(
1123 Vec3::new(0.0, 5.0, 0.0), Vec3::ZERO,
1125 Vec3::splat(0.5),
1126 1.0, );
1128 assert_eq!(f, Vec3::ZERO, "no buoyancy above water");
1129 }
1130
1131 #[test]
1132 fn buoyancy_upward_when_submerged() {
1133 let bf = BuoyancyForce::new(1000.0);
1134 let f = bf.compute(
1135 Vec3::new(0.0, 0.0, 0.0), Vec3::ZERO,
1137 Vec3::splat(0.5),
1138 2.0, );
1140 assert!(f.y > 0.0, "buoyancy should be upward, got y={}", f.y);
1141 }
1142
1143 #[test]
1144 fn buoyancy_drag_opposes_velocity() {
1145 let bf = BuoyancyForce::new(1000.0);
1146 let f = bf.compute(
1147 Vec3::new(0.0, 0.0, 0.0),
1148 Vec3::new(5.0, 0.0, 0.0), Vec3::splat(0.5),
1150 2.0,
1151 );
1152 assert!(f.x < 0.0, "drag should oppose positive X velocity, got x={}", f.x);
1153 }
1154}