1use std::collections::HashMap;
24
25use rayon::prelude::*;
26
27use super::types::{AttributeData, MeshData};
28
29const PARALLEL_THRESHOLD: usize = 1024;
30
31pub const CELL_SENTINEL: u32 = u32::MAX;
36
37#[deprecated(since = "0.13.0", note = "use `CELL_SENTINEL` instead")]
39pub const TET_SENTINEL: u32 = CELL_SENTINEL;
40
41#[non_exhaustive]
64#[derive(Default)]
65pub struct VolumeMeshData {
66 pub positions: Vec<[f32; 3]>,
68
69 pub cells: Vec<[u32; 8]>,
75
76 pub cell_scalars: HashMap<String, Vec<f32>>,
81
82 pub cell_colors: HashMap<String, Vec<[f32; 4]>>,
87}
88
89const TET_FACES: [[usize; 3]; 4] = [
102 [1, 2, 3], [0, 3, 2], [0, 1, 3], [0, 2, 1], ];
107
108const HEX_FACES: [[usize; 4]; 6] = [
134 [0, 1, 2, 3], [4, 7, 6, 5], [0, 4, 5, 1], [2, 6, 7, 3], [0, 3, 7, 4], [1, 5, 6, 2], ];
141
142const HEX_FACE_OPPOSITE: [usize; 6] = [1, 0, 3, 2, 5, 4];
144
145const PYRAMID_QUAD_FACE: [[usize; 4]; 1] = [
164 [0, 1, 2, 3], ];
166
167const PYRAMID_TRI_FACES: [[usize; 3]; 4] = [
169 [0, 4, 1], [1, 4, 2], [2, 4, 3], [3, 4, 0], ];
174
175const PYRAMID_EDGES: [[usize; 2]; 8] = [
177 [0, 1], [1, 2], [2, 3], [3, 0], [0, 4], [1, 4], [2, 4], [3, 4], ];
180
181const WEDGE_TRI_FACES: [[usize; 3]; 2] = [
200 [0, 2, 1], [3, 4, 5], ];
203
204const WEDGE_QUAD_FACES: [[usize; 4]; 3] = [
206 [0, 1, 4, 3], [1, 2, 5, 4], [2, 0, 3, 5], ];
210
211const WEDGE_EDGES: [[usize; 2]; 9] = [
213 [0, 1], [1, 2], [2, 0], [3, 4], [4, 5], [5, 3], [0, 3], [1, 4], [2, 5], ];
217
218type FaceKey = (u32, u32, u32);
224
225type QuadFaceKey = (u32, u32, u32, u32);
227
228type TriEntry = (FaceKey, usize, [u32; 3], [f32; 3]);
230type QuadEntry = (QuadFaceKey, usize, [u32; 4], [f32; 3]);
231
232#[inline]
234fn face_key(a: u32, b: u32, c: u32) -> FaceKey {
235 let mut arr = [a, b, c];
236 arr.sort_unstable();
237 (arr[0], arr[1], arr[2])
238}
239
240#[inline]
242fn quad_face_key(a: u32, b: u32, c: u32, d: u32) -> QuadFaceKey {
243 let mut arr = [a, b, c, d];
244 arr.sort_unstable();
245 (arr[0], arr[1], arr[2], arr[3])
246}
247
248fn generate_tri_entries(cell_idx: usize, cell: &[u32; 8], positions: &[[f32; 3]]) -> Vec<TriEntry> {
250 let ct = cell_type(cell);
251 let nv = ct.vertex_count();
252 let mut out = Vec::new();
253 match ct {
254 CellType::Tet => {
255 for (face_idx, face_local) in TET_FACES.iter().enumerate() {
256 let a = cell[face_local[0]];
257 let b = cell[face_local[1]];
258 let c = cell[face_local[2]];
259 let interior_ref = positions[cell[face_idx] as usize];
261 out.push((face_key(a, b, c), cell_idx, [a, b, c], interior_ref));
262 }
263 }
264 CellType::Pyramid => {
265 let centroid = cell_centroid(cell, nv, positions);
266 for face_local in &PYRAMID_TRI_FACES {
267 let a = cell[face_local[0]];
268 let b = cell[face_local[1]];
269 let c = cell[face_local[2]];
270 out.push((face_key(a, b, c), cell_idx, [a, b, c], centroid));
271 }
272 }
273 CellType::Wedge => {
274 let centroid = cell_centroid(cell, nv, positions);
275 for face_local in &WEDGE_TRI_FACES {
276 let a = cell[face_local[0]];
277 let b = cell[face_local[1]];
278 let c = cell[face_local[2]];
279 out.push((face_key(a, b, c), cell_idx, [a, b, c], centroid));
280 }
281 }
282 CellType::Hex => {} }
284 out
285}
286
287fn generate_quad_entries(cell_idx: usize, cell: &[u32; 8], positions: &[[f32; 3]]) -> Vec<QuadEntry> {
289 let ct = cell_type(cell);
290 let nv = ct.vertex_count();
291 let mut out = Vec::new();
292 match ct {
293 CellType::Tet => {} CellType::Pyramid => {
295 let centroid = cell_centroid(cell, nv, positions);
296 for quad_local in &PYRAMID_QUAD_FACE {
297 let v = [
298 cell[quad_local[0]],
299 cell[quad_local[1]],
300 cell[quad_local[2]],
301 cell[quad_local[3]],
302 ];
303 out.push((quad_face_key(v[0], v[1], v[2], v[3]), cell_idx, v, centroid));
304 }
305 }
306 CellType::Wedge => {
307 let centroid = cell_centroid(cell, nv, positions);
308 for quad_local in &WEDGE_QUAD_FACES {
309 let v = [
310 cell[quad_local[0]],
311 cell[quad_local[1]],
312 cell[quad_local[2]],
313 cell[quad_local[3]],
314 ];
315 out.push((quad_face_key(v[0], v[1], v[2], v[3]), cell_idx, v, centroid));
316 }
317 }
318 CellType::Hex => {
319 for (face_idx, quad) in HEX_FACES.iter().enumerate() {
320 let v: [u32; 4] =
321 [cell[quad[0]], cell[quad[1]], cell[quad[2]], cell[quad[3]]];
322 let interior_ref = {
323 let opposite = &HEX_FACES[HEX_FACE_OPPOSITE[face_idx]];
324 let mut c = [0.0f32; 3];
325 for &local_vi in opposite {
326 let p = positions[cell[local_vi] as usize];
327 c[0] += p[0];
328 c[1] += p[1];
329 c[2] += p[2];
330 }
331 [c[0] / 4.0, c[1] / 4.0, c[2] / 4.0]
332 };
333 out.push((quad_face_key(v[0], v[1], v[2], v[3]), cell_idx, v, interior_ref));
334 }
335 }
336 }
337 out
338}
339
340fn collect_boundary_tri(entries: &[TriEntry]) -> Vec<(usize, [u32; 3], [f32; 3])> {
342 let mut out = Vec::new();
343 let mut i = 0;
344 while i < entries.len() {
345 let key = entries[i].0;
346 let mut j = i + 1;
347 while j < entries.len() && entries[j].0 == key {
348 j += 1;
349 }
350 if j - i == 1 {
351 out.push((entries[i].1, entries[i].2, entries[i].3));
352 }
353 i = j;
354 }
355 out
356}
357
358fn collect_boundary_quad(entries: &[QuadEntry]) -> Vec<(usize, [u32; 4], [f32; 3])> {
360 let mut out = Vec::new();
361 let mut i = 0;
362 while i < entries.len() {
363 let key = entries[i].0;
364 let mut j = i + 1;
365 while j < entries.len() && entries[j].0 == key {
366 j += 1;
367 }
368 if j - i == 1 {
369 out.push((entries[i].1, entries[i].2, entries[i].3));
370 }
371 i = j;
372 }
373 out
374}
375
376#[inline]
379fn correct_winding(tri: &mut [u32; 3], interior_ref: &[f32; 3], positions: &[[f32; 3]]) {
380 let pa = positions[tri[0] as usize];
381 let pb = positions[tri[1] as usize];
382 let pc = positions[tri[2] as usize];
383 let ab = [pb[0] - pa[0], pb[1] - pa[1], pb[2] - pa[2]];
384 let ac = [pc[0] - pa[0], pc[1] - pa[1], pc[2] - pa[2]];
385 let normal = [
386 ab[1] * ac[2] - ab[2] * ac[1],
387 ab[2] * ac[0] - ab[0] * ac[2],
388 ab[0] * ac[1] - ab[1] * ac[0],
389 ];
390 let fc = [
391 (pa[0] + pb[0] + pc[0]) / 3.0,
392 (pa[1] + pb[1] + pc[1]) / 3.0,
393 (pa[2] + pb[2] + pc[2]) / 3.0,
394 ];
395 let out = [
396 fc[0] - interior_ref[0],
397 fc[1] - interior_ref[1],
398 fc[2] - interior_ref[2],
399 ];
400 if normal[0] * out[0] + normal[1] * out[1] + normal[2] * out[2] < 0.0 {
401 tri.swap(1, 2);
402 }
403}
404
405pub(crate) fn extract_boundary_faces(data: &VolumeMeshData) -> (MeshData, Vec<u32>) {
415 let n_verts = data.positions.len();
416
417 let (mut tri_entries, mut quad_entries) = if data.cells.len() >= PARALLEL_THRESHOLD {
419 let tri = data
420 .cells
421 .par_iter()
422 .enumerate()
423 .flat_map_iter(|(ci, cell)| generate_tri_entries(ci, cell, &data.positions))
424 .collect();
425 let quad = data
426 .cells
427 .par_iter()
428 .enumerate()
429 .flat_map_iter(|(ci, cell)| generate_quad_entries(ci, cell, &data.positions))
430 .collect();
431 (tri, quad)
432 } else {
433 let mut tri: Vec<TriEntry> = Vec::new();
434 let mut quad: Vec<QuadEntry> = Vec::new();
435 for (ci, cell) in data.cells.iter().enumerate() {
436 tri.extend(generate_tri_entries(ci, cell, &data.positions));
437 quad.extend(generate_quad_entries(ci, cell, &data.positions));
438 }
439 (tri, quad)
440 };
441
442 tri_entries.par_sort_unstable_by_key(|e| e.0);
443 quad_entries.par_sort_unstable_by_key(|e| e.0);
444
445 let mut boundary: Vec<(usize, [u32; 3], [f32; 3])> = collect_boundary_tri(&tri_entries);
447 for (ci, winding, iref) in collect_boundary_quad(&quad_entries) {
448 boundary.push((ci, [winding[0], winding[1], winding[2]], iref));
449 boundary.push((ci, [winding[0], winding[2], winding[3]], iref));
450 }
451
452 boundary.sort_unstable_by_key(|(ci, _, _)| *ci);
454
455 boundary
459 .par_iter_mut()
460 .for_each(|(_, tri, iref)| correct_winding(tri, iref, &data.positions));
461
462 let n_boundary_tris = boundary.len();
463
464 let mut indices: Vec<u32> = Vec::with_capacity(n_boundary_tris * 3);
467 let mut normal_accum: Vec<[f64; 3]> = vec![[0.0; 3]; n_verts];
468
469 for (_, tri, _) in &boundary {
470 indices.push(tri[0]);
471 indices.push(tri[1]);
472 indices.push(tri[2]);
473
474 let pa = data.positions[tri[0] as usize];
475 let pb = data.positions[tri[1] as usize];
476 let pc = data.positions[tri[2] as usize];
477 let ab = [
478 (pb[0] - pa[0]) as f64,
479 (pb[1] - pa[1]) as f64,
480 (pb[2] - pa[2]) as f64,
481 ];
482 let ac = [
483 (pc[0] - pa[0]) as f64,
484 (pc[1] - pa[1]) as f64,
485 (pc[2] - pa[2]) as f64,
486 ];
487 let n = [
488 ab[1] * ac[2] - ab[2] * ac[1],
489 ab[2] * ac[0] - ab[0] * ac[2],
490 ab[0] * ac[1] - ab[1] * ac[0],
491 ];
492 for &vi in tri {
493 let acc = &mut normal_accum[vi as usize];
494 acc[0] += n[0];
495 acc[1] += n[1];
496 acc[2] += n[2];
497 }
498 }
499
500 let mut normals: Vec<[f32; 3]> = normal_accum
501 .iter()
502 .map(|n| {
503 let len = (n[0] * n[0] + n[1] * n[1] + n[2] * n[2]).sqrt();
504 if len > 1e-12 {
505 [
506 (n[0] / len) as f32,
507 (n[1] / len) as f32,
508 (n[2] / len) as f32,
509 ]
510 } else {
511 [0.0, 1.0, 0.0]
512 }
513 })
514 .collect();
515
516 normals.resize(n_verts, [0.0, 1.0, 0.0]);
517
518 let mut attributes: HashMap<String, AttributeData> = HashMap::new();
519
520 for (name, cell_vals) in &data.cell_scalars {
521 let face_scalars: Vec<f32> = boundary
522 .iter()
523 .map(|(ci, _, _)| cell_vals.get(*ci).copied().unwrap_or(0.0))
524 .collect();
525 attributes.insert(name.clone(), AttributeData::Face(face_scalars));
526 }
527
528 for (name, cell_vals) in &data.cell_colors {
529 let face_colors: Vec<[f32; 4]> = boundary
530 .iter()
531 .map(|(ci, _, _)| cell_vals.get(*ci).copied().unwrap_or([1.0; 4]))
532 .collect();
533 attributes.insert(name.clone(), AttributeData::FaceColor(face_colors));
534 }
535
536 let face_to_cell: Vec<u32> = boundary.iter().map(|(ci, _, _)| *ci as u32).collect();
537
538 (MeshData {
539 positions: data.positions.clone(),
540 normals,
541 indices,
542 uvs: None,
543 tangents: None,
544 attributes,
545 }, face_to_cell)
546}
547
548const TET_EDGES: [[usize; 2]; 6] = [
649 [0, 1], [0, 2], [0, 3], [1, 2], [1, 3], [2, 3],
650];
651
652const HEX_EDGES: [[usize; 2]; 12] = [
663 [0, 1], [1, 2], [2, 3], [3, 0], [4, 5], [5, 6], [6, 7], [7, 4], [0, 4], [1, 5], [2, 6], [3, 7], ];
667
668#[derive(Clone, Copy, PartialEq, Eq)]
674enum CellType {
675 Tet,
676 Pyramid,
677 Wedge,
678 Hex,
679}
680
681impl CellType {
682 fn vertex_count(self) -> usize {
683 match self {
684 CellType::Tet => 4,
685 CellType::Pyramid => 5,
686 CellType::Wedge => 6,
687 CellType::Hex => 8,
688 }
689 }
690
691 fn edges(self) -> &'static [[usize; 2]] {
692 match self {
693 CellType::Tet => &TET_EDGES,
694 CellType::Pyramid => &PYRAMID_EDGES,
695 CellType::Wedge => &WEDGE_EDGES,
696 CellType::Hex => &HEX_EDGES,
697 }
698 }
699}
700
701#[inline]
703fn cell_type(cell: &[u32; 8]) -> CellType {
704 if cell[4] == CELL_SENTINEL {
705 CellType::Tet
706 } else if cell[5] == CELL_SENTINEL {
707 CellType::Pyramid
708 } else if cell[6] == CELL_SENTINEL {
709 CellType::Wedge
710 } else {
711 CellType::Hex
712 }
713}
714
715#[inline]
717fn cell_centroid(cell: &[u32; 8], nv: usize, positions: &[[f32; 3]]) -> [f32; 3] {
718 let mut c = [0.0f32; 3];
719 for i in 0..nv {
720 let p = positions[cell[i] as usize];
721 c[0] += p[0];
722 c[1] += p[1];
723 c[2] += p[2];
724 }
725 let n = nv as f32;
726 [c[0] / n, c[1] / n, c[2] / n]
727}
728
729#[inline]
732fn plane_dist(p: [f32; 3], plane: [f32; 4]) -> f32 {
733 p[0] * plane[0] + p[1] * plane[1] + p[2] * plane[2] + plane[3]
734}
735
736#[inline]
738fn cross3(a: [f32; 3], b: [f32; 3]) -> [f32; 3] {
739 [
740 a[1] * b[2] - a[2] * b[1],
741 a[2] * b[0] - a[0] * b[2],
742 a[0] * b[1] - a[1] * b[0],
743 ]
744}
745
746#[inline]
748fn dot3(a: [f32; 3], b: [f32; 3]) -> f32 {
749 a[0] * b[0] + a[1] * b[1] + a[2] * b[2]
750}
751
752#[inline]
754fn normalize3(v: [f32; 3]) -> [f32; 3] {
755 let len = dot3(v, v).sqrt();
756 if len > 1e-12 {
757 [v[0] / len, v[1] / len, v[2] / len]
758 } else {
759 [0.0, 1.0, 0.0]
760 }
761}
762
763fn intern_pos(
767 p: [f32; 3],
768 positions: &mut Vec<[f32; 3]>,
769 pos_map: &mut HashMap<[u32; 3], u32>,
770) -> u32 {
771 let key = [p[0].to_bits(), p[1].to_bits(), p[2].to_bits()];
772 if let Some(&idx) = pos_map.get(&key) {
773 return idx;
774 }
775 let idx = positions.len() as u32;
776 positions.push(p);
777 pos_map.insert(key, idx);
778 idx
779}
780
781fn clip_polygon_one_plane(poly: Vec<[f32; 3]>, plane: [f32; 4]) -> Vec<[f32; 3]> {
784 if poly.is_empty() {
785 return poly;
786 }
787 let n = poly.len();
788 let mut out = Vec::with_capacity(n + 1);
789 for i in 0..n {
790 let a = poly[i];
791 let b = poly[(i + 1) % n];
792 let da = plane_dist(a, plane);
793 let db = plane_dist(b, plane);
794 let a_in = da >= 0.0;
795 let b_in = db >= 0.0;
796 if a_in {
797 out.push(a);
798 }
799 if a_in != b_in {
800 let denom = da - db;
801 if denom.abs() > 1e-30 {
802 let t = da / denom;
803 out.push([
804 a[0] + t * (b[0] - a[0]),
805 a[1] + t * (b[1] - a[1]),
806 a[2] + t * (b[2] - a[2]),
807 ]);
808 }
809 }
810 }
811 out
812}
813
814fn clip_polygon_planes(mut poly: Vec<[f32; 3]>, planes: &[[f32; 4]]) -> Vec<[f32; 3]> {
816 for &plane in planes {
817 if poly.is_empty() {
818 break;
819 }
820 poly = clip_polygon_one_plane(poly, plane);
821 }
822 poly
823}
824
825fn plane_basis(normal: [f32; 3]) -> ([f32; 3], [f32; 3]) {
827 let ref_vec: [f32; 3] = if normal[0].abs() < 0.9 {
828 [1.0, 0.0, 0.0]
829 } else {
830 [0.0, 1.0, 0.0]
831 };
832 let u = normalize3(cross3(normal, ref_vec));
833 let v = cross3(normal, u);
834 (u, v)
835}
836
837fn sort_polygon_on_plane(pts: &mut Vec<[f32; 3]>, normal: [f32; 3]) {
842 if pts.len() < 3 {
843 return;
844 }
845 let n = pts.len() as f32;
846 let cx = pts.iter().map(|p| p[0]).sum::<f32>() / n;
847 let cy = pts.iter().map(|p| p[1]).sum::<f32>() / n;
848 let cz = pts.iter().map(|p| p[2]).sum::<f32>() / n;
849 let centroid = [cx, cy, cz];
850 let (u, v) = plane_basis(normal);
851 pts.sort_by(|a, b| {
852 let da = [a[0] - centroid[0], a[1] - centroid[1], a[2] - centroid[2]];
853 let db = [b[0] - centroid[0], b[1] - centroid[1], b[2] - centroid[2]];
854 let ang_a = dot3(da, v).atan2(dot3(da, u));
855 let ang_b = dot3(db, v).atan2(dot3(db, u));
856 ang_a.partial_cmp(&ang_b).unwrap_or(std::cmp::Ordering::Equal)
857 });
858}
859
860fn fan_triangulate(poly: &[[f32; 3]]) -> Vec<[[f32; 3]; 3]> {
862 if poly.len() < 3 {
863 return Vec::new();
864 }
865 (1..poly.len() - 1)
866 .map(|i| [poly[0], poly[i], poly[i + 1]])
867 .collect()
868}
869
870fn generate_section_tris(
872 cell_idx: usize,
873 cell: &[u32; 8],
874 positions: &[[f32; 3]],
875 clip_planes: &[[f32; 4]],
876) -> Vec<(usize, [[f32; 3]; 3])> {
877 let mut out = Vec::new();
878 let edges = cell_type(cell).edges();
879
880 for (pi, &plane) in clip_planes.iter().enumerate() {
881 let mut pts: Vec<[f32; 3]> = Vec::new();
882 for edge in edges {
883 let pa = positions[cell[edge[0]] as usize];
884 let pb = positions[cell[edge[1]] as usize];
885 let da = plane_dist(pa, plane);
886 let db = plane_dist(pb, plane);
887 if (da >= 0.0) != (db >= 0.0) {
888 let denom = da - db;
889 if denom.abs() > 1e-30 {
890 let t = da / denom;
891 pts.push([
892 pa[0] + t * (pb[0] - pa[0]),
893 pa[1] + t * (pb[1] - pa[1]),
894 pa[2] + t * (pb[2] - pa[2]),
895 ]);
896 }
897 }
898 }
899 if pts.len() < 3 {
900 continue;
901 }
902 let plane_normal = [plane[0], plane[1], plane[2]];
903 sort_polygon_on_plane(&mut pts, plane_normal);
904 let other_planes: Vec<[f32; 4]> = clip_planes
905 .iter()
906 .enumerate()
907 .filter(|(i, _)| *i != pi)
908 .map(|(_, p)| *p)
909 .collect();
910 let pts = clip_polygon_planes(pts, &other_planes);
911 if pts.len() < 3 {
912 continue;
913 }
914 for mut tri in fan_triangulate(&pts) {
915 let ab = [
916 tri[1][0] - tri[0][0],
917 tri[1][1] - tri[0][1],
918 tri[1][2] - tri[0][2],
919 ];
920 let ac = [
921 tri[2][0] - tri[0][0],
922 tri[2][1] - tri[0][1],
923 tri[2][2] - tri[0][2],
924 ];
925 let n = cross3(ab, ac);
926 if dot3(n, plane_normal) < 0.0 {
927 tri.swap(1, 2);
928 }
929 out.push((cell_idx, tri));
930 }
931 }
932 out
933}
934
935pub fn extract_clipped_volume_faces(
945 data: &VolumeMeshData,
946 clip_planes: &[[f32; 4]],
947) -> (MeshData, Vec<u32>) {
948 if clip_planes.is_empty() {
949 return extract_boundary_faces(data);
950 }
951
952 let vert_kept: Vec<bool> = data
954 .positions
955 .par_iter()
956 .map(|&p| clip_planes.iter().all(|&pl| plane_dist(p, pl) >= 0.0))
957 .collect();
958
959 let (mut tri_entries, mut quad_entries) = if data.cells.len() >= PARALLEL_THRESHOLD {
961 let vk = &vert_kept;
962 let tri = data
963 .cells
964 .par_iter()
965 .enumerate()
966 .flat_map_iter(|(ci, cell)| {
967 let nv = cell_type(cell).vertex_count();
968 if (0..nv).all(|i| !vk[cell[i] as usize]) {
969 return Vec::new();
970 }
971 generate_tri_entries(ci, cell, &data.positions)
972 })
973 .collect();
974 let quad = data
975 .cells
976 .par_iter()
977 .enumerate()
978 .flat_map_iter(|(ci, cell)| {
979 let nv = cell_type(cell).vertex_count();
980 if (0..nv).all(|i| !vk[cell[i] as usize]) {
981 return Vec::new();
982 }
983 generate_quad_entries(ci, cell, &data.positions)
984 })
985 .collect();
986 (tri, quad)
987 } else {
988 let mut tri: Vec<TriEntry> = Vec::new();
989 let mut quad: Vec<QuadEntry> = Vec::new();
990 for (ci, cell) in data.cells.iter().enumerate() {
991 let nv = cell_type(cell).vertex_count();
992 let kc = (0..nv).filter(|&i| vert_kept[cell[i] as usize]).count();
993 if kc == 0 {
994 continue;
995 }
996 tri.extend(generate_tri_entries(ci, cell, &data.positions));
997 quad.extend(generate_quad_entries(ci, cell, &data.positions));
998 }
999 (tri, quad)
1000 };
1001
1002 tri_entries.par_sort_unstable_by_key(|e| e.0);
1003 quad_entries.par_sort_unstable_by_key(|e| e.0);
1004
1005 let mut boundary: Vec<(usize, [u32; 3], [f32; 3])> = collect_boundary_tri(&tri_entries);
1006 for (ci, winding, iref) in collect_boundary_quad(&quad_entries) {
1007 boundary.push((ci, [winding[0], winding[1], winding[2]], iref));
1008 boundary.push((ci, [winding[0], winding[2], winding[3]], iref));
1009 }
1010 boundary.sort_unstable_by_key(|(ci, _, _)| *ci);
1011
1012 boundary
1013 .par_iter_mut()
1014 .for_each(|(_, tri, iref)| correct_winding(tri, iref, &data.positions));
1015
1016 let cell_nv: Vec<usize> = data
1018 .cells
1019 .iter()
1020 .map(|c| cell_type(c).vertex_count())
1021 .collect();
1022 let cell_kept: Vec<usize> = data
1023 .cells
1024 .iter()
1025 .zip(cell_nv.iter())
1026 .map(|(cell, &nv)| (0..nv).filter(|&i| vert_kept[cell[i] as usize]).count())
1027 .collect();
1028
1029 let mut out_tris: Vec<(usize, [[f32; 3]; 3])> = boundary
1031 .par_iter()
1032 .flat_map_iter(|(cell_idx, tri, _)| {
1033 let nv = cell_nv[*cell_idx];
1034 let kc = cell_kept[*cell_idx];
1035 let pa = data.positions[tri[0] as usize];
1036 let pb = data.positions[tri[1] as usize];
1037 let pc = data.positions[tri[2] as usize];
1038 if kc == nv {
1039 vec![(*cell_idx, [pa, pb, pc])]
1040 } else {
1041 let clipped = clip_polygon_planes(vec![pa, pb, pc], clip_planes);
1042 fan_triangulate(&clipped)
1043 .into_iter()
1044 .map(|t| (*cell_idx, t))
1045 .collect()
1046 }
1047 })
1048 .collect();
1049
1050 let section_tris: Vec<(usize, [[f32; 3]; 3])> = data
1052 .cells
1053 .par_iter()
1054 .enumerate()
1055 .filter(|(ci, _)| {
1056 let kc = cell_kept[*ci];
1057 kc > 0 && kc < cell_nv[*ci]
1058 })
1059 .flat_map_iter(|(ci, cell)| {
1060 generate_section_tris(ci, cell, &data.positions, clip_planes)
1061 })
1062 .collect();
1063 out_tris.extend(section_tris);
1064
1065 let mut positions: Vec<[f32; 3]> = data.positions.clone();
1067 let mut pos_map: HashMap<[u32; 3], u32> = HashMap::new();
1068 for (i, &p) in data.positions.iter().enumerate() {
1069 let key = [p[0].to_bits(), p[1].to_bits(), p[2].to_bits()];
1070 pos_map.entry(key).or_insert(i as u32);
1071 }
1072
1073 let mut indexed_tris: Vec<(usize, [u32; 3])> = Vec::with_capacity(out_tris.len());
1074 for (cell_idx, [p0, p1, p2]) in &out_tris {
1075 let i0 = intern_pos(*p0, &mut positions, &mut pos_map);
1076 let i1 = intern_pos(*p1, &mut positions, &mut pos_map);
1077 let i2 = intern_pos(*p2, &mut positions, &mut pos_map);
1078 indexed_tris.push((*cell_idx, [i0, i1, i2]));
1079 }
1080
1081 let n_verts = positions.len();
1082 let mut normal_accum: Vec<[f64; 3]> = vec![[0.0; 3]; n_verts];
1083 let mut indices: Vec<u32> = Vec::with_capacity(indexed_tris.len() * 3);
1084
1085 for (_, tri) in &indexed_tris {
1086 indices.push(tri[0]);
1087 indices.push(tri[1]);
1088 indices.push(tri[2]);
1089
1090 let pa = positions[tri[0] as usize];
1091 let pb = positions[tri[1] as usize];
1092 let pc = positions[tri[2] as usize];
1093 let ab = [
1094 (pb[0] - pa[0]) as f64,
1095 (pb[1] - pa[1]) as f64,
1096 (pb[2] - pa[2]) as f64,
1097 ];
1098 let ac = [
1099 (pc[0] - pa[0]) as f64,
1100 (pc[1] - pa[1]) as f64,
1101 (pc[2] - pa[2]) as f64,
1102 ];
1103 let n = [
1104 ab[1] * ac[2] - ab[2] * ac[1],
1105 ab[2] * ac[0] - ab[0] * ac[2],
1106 ab[0] * ac[1] - ab[1] * ac[0],
1107 ];
1108 for &vi in tri {
1109 let acc = &mut normal_accum[vi as usize];
1110 acc[0] += n[0];
1111 acc[1] += n[1];
1112 acc[2] += n[2];
1113 }
1114 }
1115
1116 let normals: Vec<[f32; 3]> = normal_accum
1117 .iter()
1118 .map(|n| {
1119 let len = (n[0] * n[0] + n[1] * n[1] + n[2] * n[2]).sqrt();
1120 if len > 1e-12 {
1121 [(n[0] / len) as f32, (n[1] / len) as f32, (n[2] / len) as f32]
1122 } else {
1123 [0.0, 1.0, 0.0]
1124 }
1125 })
1126 .collect();
1127
1128 let mut attributes: HashMap<String, AttributeData> = HashMap::new();
1129 for (name, cell_vals) in &data.cell_scalars {
1130 let face_scalars: Vec<f32> = indexed_tris
1131 .iter()
1132 .map(|(ci, _)| cell_vals.get(*ci).copied().unwrap_or(0.0))
1133 .collect();
1134 attributes.insert(name.clone(), AttributeData::Face(face_scalars));
1135 }
1136 for (name, cell_vals) in &data.cell_colors {
1137 let face_colors: Vec<[f32; 4]> = indexed_tris
1138 .iter()
1139 .map(|(ci, _)| cell_vals.get(*ci).copied().unwrap_or([1.0; 4]))
1140 .collect();
1141 attributes.insert(name.clone(), AttributeData::FaceColor(face_colors));
1142 }
1143
1144 let face_to_cell: Vec<u32> = indexed_tris.iter().map(|(ci, _)| *ci as u32).collect();
1145
1146 (MeshData {
1147 positions,
1148 normals,
1149 indices,
1150 uvs: None,
1151 tangents: None,
1152 attributes,
1153 }, face_to_cell)
1154}
1155
1156impl VolumeMeshData {
1161 pub fn push_tet(&mut self, a: u32, b: u32, c: u32, d: u32) {
1165 self.cells.push([
1166 a, b, c, d,
1167 CELL_SENTINEL, CELL_SENTINEL, CELL_SENTINEL, CELL_SENTINEL,
1168 ]);
1169 }
1170
1171 pub fn push_pyramid(&mut self, base: [u32; 4], apex: u32) {
1177 self.cells.push([
1178 base[0], base[1], base[2], base[3], apex,
1179 CELL_SENTINEL, CELL_SENTINEL, CELL_SENTINEL,
1180 ]);
1181 }
1182
1183 pub fn push_wedge(&mut self, tri0: [u32; 3], tri1: [u32; 3]) {
1189 self.cells.push([
1190 tri0[0], tri0[1], tri0[2], tri1[0], tri1[1], tri1[2],
1191 CELL_SENTINEL, CELL_SENTINEL,
1192 ]);
1193 }
1194
1195 pub fn push_hex(&mut self, verts: [u32; 8]) {
1197 self.cells.push(verts);
1198 }
1199}
1200
1201const HEX_TO_TETS: [[usize; 4]; 6] = [
1209 [0, 1, 5, 6],
1210 [0, 1, 2, 6],
1211 [0, 4, 5, 6],
1212 [0, 4, 7, 6],
1213 [0, 3, 2, 6],
1214 [0, 3, 7, 6],
1215];
1216
1217const WEDGE_TO_TETS: [[usize; 4]; 3] = [[0, 1, 2, 3], [1, 2, 3, 4], [2, 3, 4, 5]];
1221
1222const PYRAMID_TO_TETS: [[usize; 4]; 2] = [[0, 1, 2, 4], [0, 2, 3, 4]];
1226
1227pub(crate) fn for_each_tet<F>(data: &VolumeMeshData, attribute: &str, mut f: F)
1239where
1240 F: FnMut([[f32; 3]; 4], f32),
1241{
1242 let cell_scalars = data.cell_scalars.get(attribute);
1243 for (cell_idx, cell) in data.cells.iter().enumerate() {
1244 let scalar = cell_scalars
1245 .and_then(|v| v.get(cell_idx))
1246 .copied()
1247 .unwrap_or(0.0);
1248 let tets: &[[usize; 4]] = match cell_type(cell) {
1249 CellType::Tet => &[[0, 1, 2, 3]],
1250 CellType::Pyramid => &PYRAMID_TO_TETS,
1251 CellType::Wedge => &WEDGE_TO_TETS,
1252 CellType::Hex => &HEX_TO_TETS,
1253 };
1254 for local in tets {
1255 let verts = [
1256 data.positions[cell[local[0]] as usize],
1257 data.positions[cell[local[1]] as usize],
1258 data.positions[cell[local[2]] as usize],
1259 data.positions[cell[local[3]] as usize],
1260 ];
1261 f(verts, scalar);
1262 }
1263 }
1264}
1265
1266#[cfg(test)]
1276pub(crate) fn decompose_to_tetrahedra(
1277 data: &VolumeMeshData,
1278 attribute: &str,
1279) -> (Vec<[[f32; 3]; 4]>, Vec<f32>) {
1280 let mut positions: Vec<[[f32; 3]; 4]> = Vec::new();
1281 let mut scalars: Vec<f32> = Vec::new();
1282 for_each_tet(data, attribute, |verts, scalar| {
1283 positions.push(verts);
1284 scalars.push(scalar);
1285 });
1286 (positions, scalars)
1287}
1288
1289#[cfg(test)]
1294mod tests {
1295 use super::*;
1296
1297 const TEST_TET_LOCAL: [[usize; 4]; 6] = [
1298 [0, 1, 5, 6],
1299 [0, 1, 2, 6],
1300 [0, 4, 5, 6],
1301 [0, 4, 7, 6],
1302 [0, 3, 2, 6],
1303 [0, 3, 7, 6],
1304 ];
1305
1306 fn single_tet() -> VolumeMeshData {
1307 VolumeMeshData {
1308 positions: vec![
1309 [0.0, 0.0, 0.0],
1310 [1.0, 0.0, 0.0],
1311 [0.5, 1.0, 0.0],
1312 [0.5, 0.5, 1.0],
1313 ],
1314 cells: vec![[
1315 0,
1316 1,
1317 2,
1318 3,
1319 CELL_SENTINEL,
1320 CELL_SENTINEL,
1321 CELL_SENTINEL,
1322 CELL_SENTINEL,
1323 ]],
1324 ..Default::default()
1325 }
1326 }
1327
1328 fn two_tets_sharing_face() -> VolumeMeshData {
1329 VolumeMeshData {
1332 positions: vec![
1333 [0.0, 0.0, 0.0],
1334 [1.0, 0.0, 0.0],
1335 [0.5, 1.0, 0.0],
1336 [0.5, 0.5, 1.0],
1337 [0.5, 0.5, -1.0],
1338 ],
1339 cells: vec![
1340 [
1341 0,
1342 1,
1343 2,
1344 3,
1345 CELL_SENTINEL,
1346 CELL_SENTINEL,
1347 CELL_SENTINEL,
1348 CELL_SENTINEL,
1349 ],
1350 [
1351 0,
1352 2,
1353 1,
1354 4,
1355 CELL_SENTINEL,
1356 CELL_SENTINEL,
1357 CELL_SENTINEL,
1358 CELL_SENTINEL,
1359 ],
1360 ],
1361 ..Default::default()
1362 }
1363 }
1364
1365 fn single_hex() -> VolumeMeshData {
1366 VolumeMeshData {
1367 positions: vec![
1368 [0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [1.0, 0.0, 1.0], [0.0, 0.0, 1.0], [0.0, 1.0, 0.0], [1.0, 1.0, 0.0], [1.0, 1.0, 1.0], [0.0, 1.0, 1.0], ],
1377 cells: vec![[0, 1, 2, 3, 4, 5, 6, 7]],
1378 ..Default::default()
1379 }
1380 }
1381
1382 fn structured_tet_grid(grid_n: usize) -> VolumeMeshData {
1383 let grid_v = grid_n + 1;
1384 let vid =
1385 |ix: usize, iy: usize, iz: usize| (iz * grid_v * grid_v + iy * grid_v + ix) as u32;
1386
1387 let mut positions = Vec::with_capacity(grid_v * grid_v * grid_v);
1388 for iz in 0..grid_v {
1389 for iy in 0..grid_v {
1390 for ix in 0..grid_v {
1391 positions.push([ix as f32, iy as f32, iz as f32]);
1392 }
1393 }
1394 }
1395
1396 let mut cells = Vec::with_capacity(grid_n * grid_n * grid_n * TEST_TET_LOCAL.len());
1397 for iz in 0..grid_n {
1398 for iy in 0..grid_n {
1399 for ix in 0..grid_n {
1400 let cube_verts = [
1401 vid(ix, iy, iz),
1402 vid(ix + 1, iy, iz),
1403 vid(ix + 1, iy, iz + 1),
1404 vid(ix, iy, iz + 1),
1405 vid(ix, iy + 1, iz),
1406 vid(ix + 1, iy + 1, iz),
1407 vid(ix + 1, iy + 1, iz + 1),
1408 vid(ix, iy + 1, iz + 1),
1409 ];
1410 for tet in &TEST_TET_LOCAL {
1411 cells.push([
1412 cube_verts[tet[0]],
1413 cube_verts[tet[1]],
1414 cube_verts[tet[2]],
1415 cube_verts[tet[3]],
1416 CELL_SENTINEL,
1417 CELL_SENTINEL,
1418 CELL_SENTINEL,
1419 CELL_SENTINEL,
1420 ]);
1421 }
1422 }
1423 }
1424 }
1425
1426 VolumeMeshData {
1427 positions,
1428 cells,
1429 ..Default::default()
1430 }
1431 }
1432
1433 fn projected_sphere_tet_grid(grid_n: usize, radius: f32) -> VolumeMeshData {
1434 let grid_v = grid_n + 1;
1435 let half = grid_n as f32 / 2.0;
1436 let vid =
1437 |ix: usize, iy: usize, iz: usize| (iz * grid_v * grid_v + iy * grid_v + ix) as u32;
1438
1439 let mut positions = Vec::with_capacity(grid_v * grid_v * grid_v);
1440 for iz in 0..grid_v {
1441 for iy in 0..grid_v {
1442 for ix in 0..grid_v {
1443 let x = ix as f32 - half;
1444 let y = iy as f32 - half;
1445 let z = iz as f32 - half;
1446 let len = (x * x + y * y + z * z).sqrt();
1447 let s = radius / len;
1448 positions.push([x * s, y * s, z * s]);
1449 }
1450 }
1451 }
1452
1453 let mut cells = Vec::with_capacity(grid_n * grid_n * grid_n * TEST_TET_LOCAL.len());
1454 for iz in 0..grid_n {
1455 for iy in 0..grid_n {
1456 for ix in 0..grid_n {
1457 let cube_verts = [
1458 vid(ix, iy, iz),
1459 vid(ix + 1, iy, iz),
1460 vid(ix + 1, iy, iz + 1),
1461 vid(ix, iy, iz + 1),
1462 vid(ix, iy + 1, iz),
1463 vid(ix + 1, iy + 1, iz),
1464 vid(ix + 1, iy + 1, iz + 1),
1465 vid(ix, iy + 1, iz + 1),
1466 ];
1467 for tet in &TEST_TET_LOCAL {
1468 cells.push([
1469 cube_verts[tet[0]],
1470 cube_verts[tet[1]],
1471 cube_verts[tet[2]],
1472 cube_verts[tet[3]],
1473 CELL_SENTINEL,
1474 CELL_SENTINEL,
1475 CELL_SENTINEL,
1476 CELL_SENTINEL,
1477 ]);
1478 }
1479 }
1480 }
1481 }
1482
1483 VolumeMeshData {
1484 positions,
1485 cells,
1486 ..Default::default()
1487 }
1488 }
1489
1490 fn cube_to_sphere([x, y, z]: [f32; 3]) -> [f32; 3] {
1491 let x2 = x * x;
1492 let y2 = y * y;
1493 let z2 = z * z;
1494 [
1495 x * (1.0 - 0.5 * (y2 + z2) + (y2 * z2) / 3.0).sqrt(),
1496 y * (1.0 - 0.5 * (z2 + x2) + (z2 * x2) / 3.0).sqrt(),
1497 z * (1.0 - 0.5 * (x2 + y2) + (x2 * y2) / 3.0).sqrt(),
1498 ]
1499 }
1500
1501 fn cube_sphere_hex_grid(grid_n: usize, radius: f32) -> VolumeMeshData {
1502 let grid_v = grid_n + 1;
1503 let half = grid_n as f32 / 2.0;
1504 let vid =
1505 |ix: usize, iy: usize, iz: usize| (iz * grid_v * grid_v + iy * grid_v + ix) as u32;
1506
1507 let mut positions = Vec::with_capacity(grid_v * grid_v * grid_v);
1508 for iz in 0..grid_v {
1509 for iy in 0..grid_v {
1510 for ix in 0..grid_v {
1511 let p = [ix as f32 - half, iy as f32 - half, iz as f32 - half];
1512 let cube = [p[0] / half, p[1] / half, p[2] / half];
1513 let s = cube_to_sphere(cube);
1514 positions.push([s[0] * radius, s[1] * radius, s[2] * radius]);
1515 }
1516 }
1517 }
1518
1519 let mut cells = Vec::with_capacity(grid_n * grid_n * grid_n);
1520 for iz in 0..grid_n {
1521 for iy in 0..grid_n {
1522 for ix in 0..grid_n {
1523 cells.push([
1524 vid(ix, iy, iz),
1525 vid(ix + 1, iy, iz),
1526 vid(ix + 1, iy, iz + 1),
1527 vid(ix, iy, iz + 1),
1528 vid(ix, iy + 1, iz),
1529 vid(ix + 1, iy + 1, iz),
1530 vid(ix + 1, iy + 1, iz + 1),
1531 vid(ix, iy + 1, iz + 1),
1532 ]);
1533 }
1534 }
1535 }
1536
1537 VolumeMeshData {
1538 positions,
1539 cells,
1540 ..Default::default()
1541 }
1542 }
1543
1544 fn structured_hex_grid(grid_n: usize) -> VolumeMeshData {
1545 let grid_v = grid_n + 1;
1546 let vid =
1547 |ix: usize, iy: usize, iz: usize| (iz * grid_v * grid_v + iy * grid_v + ix) as u32;
1548
1549 let mut positions = Vec::with_capacity(grid_v * grid_v * grid_v);
1550 for iz in 0..grid_v {
1551 for iy in 0..grid_v {
1552 for ix in 0..grid_v {
1553 positions.push([ix as f32, iy as f32, iz as f32]);
1554 }
1555 }
1556 }
1557
1558 let mut cells = Vec::with_capacity(grid_n * grid_n * grid_n);
1559 for iz in 0..grid_n {
1560 for iy in 0..grid_n {
1561 for ix in 0..grid_n {
1562 cells.push([
1563 vid(ix, iy, iz),
1564 vid(ix + 1, iy, iz),
1565 vid(ix + 1, iy, iz + 1),
1566 vid(ix, iy, iz + 1),
1567 vid(ix, iy + 1, iz),
1568 vid(ix + 1, iy + 1, iz),
1569 vid(ix + 1, iy + 1, iz + 1),
1570 vid(ix, iy + 1, iz + 1),
1571 ]);
1572 }
1573 }
1574 }
1575
1576 VolumeMeshData {
1577 positions,
1578 cells,
1579 ..Default::default()
1580 }
1581 }
1582
1583 #[test]
1584 fn single_tet_has_four_boundary_faces() {
1585 let data = single_tet();
1586 let (mesh, _) = extract_boundary_faces(&data);
1587 assert_eq!(
1588 mesh.indices.len(),
1589 4 * 3,
1590 "single tet -> 4 boundary triangles"
1591 );
1592 }
1593
1594 #[test]
1595 fn two_tets_sharing_face_eliminates_shared_face() {
1596 let data = two_tets_sharing_face();
1597 let (mesh, _) = extract_boundary_faces(&data);
1598 assert_eq!(
1601 mesh.indices.len(),
1602 6 * 3,
1603 "two tets sharing a face -> 6 boundary triangles"
1604 );
1605 }
1606
1607 #[test]
1608 fn single_hex_has_twelve_boundary_triangles() {
1609 let data = single_hex();
1610 let (mesh, _) = extract_boundary_faces(&data);
1611 assert_eq!(
1613 mesh.indices.len(),
1614 12 * 3,
1615 "single hex -> 12 boundary triangles"
1616 );
1617 }
1618
1619 #[test]
1620 fn structured_tet_grid_has_expected_boundary_triangle_count() {
1621 let grid_n = 3;
1622 let data = structured_tet_grid(grid_n);
1623 let (mesh, _) = extract_boundary_faces(&data);
1624 let expected_boundary_tris = 6 * grid_n * grid_n * 2;
1625 assert_eq!(
1626 mesh.indices.len(),
1627 expected_boundary_tris * 3,
1628 "3x3x3 tet grid should expose 108 boundary triangles"
1629 );
1630 }
1631
1632 #[test]
1633 fn structured_hex_grid_has_expected_boundary_triangle_count() {
1634 let grid_n = 3;
1635 let data = structured_hex_grid(grid_n);
1636 let (mesh, _) = extract_boundary_faces(&data);
1637 let expected_boundary_tris = 6 * grid_n * grid_n * 2;
1638 assert_eq!(
1639 mesh.indices.len(),
1640 expected_boundary_tris * 3,
1641 "3x3x3 hex grid should expose 108 boundary triangles"
1642 );
1643 }
1644
1645 #[test]
1646 fn structured_tet_grid_boundary_is_edge_manifold() {
1647 let data = structured_tet_grid(3);
1648 let (mesh, _) = extract_boundary_faces(&data);
1649
1650 let mut edge_counts: std::collections::HashMap<(u32, u32), usize> =
1651 std::collections::HashMap::new();
1652 for tri in mesh.indices.chunks_exact(3) {
1653 for (a, b) in [(tri[0], tri[1]), (tri[1], tri[2]), (tri[2], tri[0])] {
1654 let edge = if a < b { (a, b) } else { (b, a) };
1655 *edge_counts.entry(edge).or_insert(0) += 1;
1656 }
1657 }
1658
1659 let non_manifold: Vec<((u32, u32), usize)> = edge_counts
1660 .into_iter()
1661 .filter(|(_, count)| *count != 2)
1662 .collect();
1663
1664 assert!(
1665 non_manifold.is_empty(),
1666 "boundary should be watertight; bad edges: {non_manifold:?}"
1667 );
1668 }
1669
1670 #[test]
1671 fn structured_hex_grid_boundary_is_edge_manifold() {
1672 let data = structured_hex_grid(3);
1673 let (mesh, _) = extract_boundary_faces(&data);
1674
1675 let mut edge_counts: std::collections::HashMap<(u32, u32), usize> =
1676 std::collections::HashMap::new();
1677 for tri in mesh.indices.chunks_exact(3) {
1678 for (a, b) in [(tri[0], tri[1]), (tri[1], tri[2]), (tri[2], tri[0])] {
1679 let edge = if a < b { (a, b) } else { (b, a) };
1680 *edge_counts.entry(edge).or_insert(0) += 1;
1681 }
1682 }
1683
1684 let non_manifold: Vec<((u32, u32), usize)> = edge_counts
1685 .into_iter()
1686 .filter(|(_, count)| *count != 2)
1687 .collect();
1688
1689 assert!(
1690 non_manifold.is_empty(),
1691 "boundary should be watertight; bad edges: {non_manifold:?}"
1692 );
1693 }
1694
1695 #[test]
1696 fn projected_sphere_tet_grid_boundary_faces_point_outward() {
1697 let data = projected_sphere_tet_grid(3, 2.0);
1698 let (mesh, _) = extract_boundary_faces(&data);
1699
1700 for tri in mesh.indices.chunks_exact(3) {
1701 let pa = mesh.positions[tri[0] as usize];
1702 let pb = mesh.positions[tri[1] as usize];
1703 let pc = mesh.positions[tri[2] as usize];
1704
1705 let ab = [pb[0] - pa[0], pb[1] - pa[1], pb[2] - pa[2]];
1706 let ac = [pc[0] - pa[0], pc[1] - pa[1], pc[2] - pa[2]];
1707 let normal = [
1708 ab[1] * ac[2] - ab[2] * ac[1],
1709 ab[2] * ac[0] - ab[0] * ac[2],
1710 ab[0] * ac[1] - ab[1] * ac[0],
1711 ];
1712 let fc = [
1713 (pa[0] + pb[0] + pc[0]) / 3.0,
1714 (pa[1] + pb[1] + pc[1]) / 3.0,
1715 (pa[2] + pb[2] + pc[2]) / 3.0,
1716 ];
1717 let dot = normal[0] * fc[0] + normal[1] * fc[1] + normal[2] * fc[2];
1718 assert!(dot > 0.0, "boundary face points inward: tri={tri:?}, dot={dot}");
1719 }
1720 }
1721
1722 #[test]
1723 fn cube_sphere_hex_grid_boundary_faces_point_outward() {
1724 let data = cube_sphere_hex_grid(3, 2.0);
1725 let (mesh, _) = extract_boundary_faces(&data);
1726
1727 for tri in mesh.indices.chunks_exact(3) {
1728 let pa = mesh.positions[tri[0] as usize];
1729 let pb = mesh.positions[tri[1] as usize];
1730 let pc = mesh.positions[tri[2] as usize];
1731
1732 let ab = [pb[0] - pa[0], pb[1] - pa[1], pb[2] - pa[2]];
1733 let ac = [pc[0] - pa[0], pc[1] - pa[1], pc[2] - pa[2]];
1734 let normal = [
1735 ab[1] * ac[2] - ab[2] * ac[1],
1736 ab[2] * ac[0] - ab[0] * ac[2],
1737 ab[0] * ac[1] - ab[1] * ac[0],
1738 ];
1739 let fc = [
1740 (pa[0] + pb[0] + pc[0]) / 3.0,
1741 (pa[1] + pb[1] + pc[1]) / 3.0,
1742 (pa[2] + pb[2] + pc[2]) / 3.0,
1743 ];
1744 let dot = normal[0] * fc[0] + normal[1] * fc[1] + normal[2] * fc[2];
1745 assert!(dot > 0.0, "boundary face points inward: tri={tri:?}, dot={dot}");
1746 }
1747 }
1748
1749 #[test]
1750 fn normals_have_correct_length() {
1751 let data = single_tet();
1752 let (mesh, _) = extract_boundary_faces(&data);
1753 assert_eq!(mesh.normals.len(), mesh.positions.len());
1754 for n in &mesh.normals {
1755 let len = (n[0] * n[0] + n[1] * n[1] + n[2] * n[2]).sqrt();
1756 assert!(
1757 (len - 1.0).abs() < 1e-5 || len < 1e-5,
1758 "normal not unit: {n:?}"
1759 );
1760 }
1761 }
1762
1763 #[test]
1764 fn cell_scalar_remaps_to_face_attribute() {
1765 let mut data = single_tet();
1766 data.cell_scalars.insert("pressure".to_string(), vec![42.0]);
1767 let (mesh, _) = extract_boundary_faces(&data);
1768 match mesh.attributes.get("pressure") {
1769 Some(AttributeData::Face(vals)) => {
1770 assert_eq!(vals.len(), 4, "one value per boundary triangle");
1771 for &v in vals {
1772 assert_eq!(v, 42.0);
1773 }
1774 }
1775 other => panic!("expected Face attribute, got {other:?}"),
1776 }
1777 }
1778
1779 #[test]
1780 fn cell_color_remaps_to_face_color_attribute() {
1781 let mut data = two_tets_sharing_face();
1782 data.cell_colors.insert(
1783 "label".to_string(),
1784 vec![[1.0, 0.0, 0.0, 1.0], [0.0, 0.0, 1.0, 1.0]],
1785 );
1786 let (mesh, _) = extract_boundary_faces(&data);
1787 match mesh.attributes.get("label") {
1788 Some(AttributeData::FaceColor(colors)) => {
1789 assert_eq!(colors.len(), 6, "6 boundary faces");
1790 }
1791 other => panic!("expected FaceColor attribute, got {other:?}"),
1792 }
1793 }
1794
1795 #[test]
1796 fn positions_preserved_unchanged() {
1797 let data = single_hex();
1798 let (mesh, _) = extract_boundary_faces(&data);
1799 assert_eq!(mesh.positions, data.positions);
1800 }
1801
1802 #[test]
1812
1813 fn empty_planes_matches_boundary_extractor_tet() {
1814 let data = structured_tet_grid(3);
1815 let (boundary, _) = extract_boundary_faces(&data);
1816 let (clipped, _) = extract_clipped_volume_faces(&data, &[]);
1817 assert_eq!(
1818 boundary.indices.len(),
1819 clipped.indices.len(),
1820 "empty clip_planes -> same triangle count as extract_boundary_faces"
1821 );
1822 }
1823
1824 #[test]
1827
1828 fn empty_planes_matches_boundary_extractor_hex() {
1829 let data = structured_hex_grid(3);
1830 let (boundary, _) = extract_boundary_faces(&data);
1831 let (clipped, _) = extract_clipped_volume_faces(&data, &[]);
1832 assert_eq!(
1833 boundary.indices.len(),
1834 clipped.indices.len(),
1835 "empty clip_planes -> same triangle count as extract_boundary_faces"
1836 );
1837 }
1838
1839 #[test]
1842
1843 fn clipped_tet_grid_has_nonempty_section_faces() {
1844 let grid_n = 3;
1845 let data = structured_tet_grid(grid_n);
1846 let plane = [0.0_f32, 1.0, 0.0, -1.5];
1849 let (mesh, _) = extract_clipped_volume_faces(&data, &[plane]);
1850 assert!(
1852 !mesh.indices.is_empty(),
1853 "clipped tet grid must produce at least one triangle"
1854 );
1855 }
1856
1857 #[test]
1859
1860 fn clipped_hex_grid_has_nonempty_section_faces() {
1861 let grid_n = 3;
1862 let data = structured_hex_grid(grid_n);
1863 let plane = [0.0_f32, 1.0, 0.0, -1.5];
1864 let (mesh, _) = extract_clipped_volume_faces(&data, &[plane]);
1865 assert!(
1866 !mesh.indices.is_empty(),
1867 "clipped hex grid must produce at least one triangle"
1868 );
1869 }
1870
1871 #[test]
1874
1875 fn section_face_normals_point_toward_kept_side_tet() {
1876 let data = structured_tet_grid(3);
1877 let plane_normal = [0.0_f32, 1.0, 0.0];
1878 let plane = [plane_normal[0], plane_normal[1], plane_normal[2], -1.5];
1879 let (mesh, _) = extract_clipped_volume_faces(&data, &[plane]);
1880
1881 for n in &mesh.normals {
1882 let dot = n[0] * plane_normal[0] + n[1] * plane_normal[1] + n[2] * plane_normal[2];
1883 let _ = dot; }
1889 }
1890
1891 #[test]
1893
1894 fn fully_discarded_cells_contribute_nothing() {
1895 let data = single_tet();
1897 let plane = [0.0_f32, 1.0, 0.0, -2.0]; let (mesh, _) = extract_clipped_volume_faces(&data, &[plane]);
1899 assert!(
1900 mesh.indices.is_empty(),
1901 "tet fully below clip plane must produce no triangles"
1902 );
1903 }
1904
1905 #[test]
1908
1909 fn fully_kept_cell_matches_boundary_extractor() {
1910 let data = single_tet();
1912 let plane = [0.0_f32, 1.0, 0.0, 1.0]; let (clipped, _) = extract_clipped_volume_faces(&data, &[plane]);
1914 let (boundary, _) = extract_boundary_faces(&data);
1915 assert_eq!(
1916 clipped.indices.len(),
1917 boundary.indices.len(),
1918 "fully kept cell must produce the same triangles as boundary extractor"
1919 );
1920 }
1921
1922 #[test]
1925 fn cell_scalar_propagates_to_section_faces() {
1926 let mut data = structured_tet_grid(3);
1927 let n_cells = data.cells.len();
1928 data.cell_scalars
1929 .insert("pressure".to_string(), vec![1.0; n_cells]);
1930 let plane = [0.0_f32, 1.0, 0.0, -1.5];
1931 let (mesh, _) = extract_clipped_volume_faces(&data, &[plane]);
1932 match mesh.attributes.get("pressure") {
1933 Some(AttributeData::Face(vals)) => {
1934 let n_tris = mesh.indices.len() / 3;
1935 assert_eq!(vals.len(), n_tris, "one scalar value per output triangle");
1936 for &v in vals {
1937 assert_eq!(v, 1.0, "scalar must equal the owning cell's value");
1938 }
1939 }
1940 other => panic!("expected Face attribute on clipped mesh, got {other:?}"),
1941 }
1942 }
1943
1944 #[test]
1947 fn cell_color_propagates_to_section_faces() {
1948 let mut data = structured_tet_grid(3);
1949 let n_cells = data.cells.len();
1950 let color = [1.0_f32, 0.0, 0.5, 1.0];
1951 data.cell_colors
1952 .insert("label".to_string(), vec![color; n_cells]);
1953 let plane = [0.0_f32, 1.0, 0.0, -1.5];
1954 let (mesh, _) = extract_clipped_volume_faces(&data, &[plane]);
1955 match mesh.attributes.get("label") {
1956 Some(AttributeData::FaceColor(colors)) => {
1957 let n_tris = mesh.indices.len() / 3;
1958 assert_eq!(colors.len(), n_tris, "one color per output triangle");
1959 for &c in colors {
1960 assert_eq!(c, color, "color must equal the owning cell's value");
1961 }
1962 }
1963 other => panic!("expected FaceColor attribute on clipped mesh, got {other:?}"),
1964 }
1965 }
1966
1967 #[test]
1969 fn hex_cell_scalar_propagates_to_section_faces() {
1970 let mut data = structured_hex_grid(3);
1971 let n_cells = data.cells.len();
1972 data.cell_scalars
1973 .insert("temp".to_string(), vec![7.0; n_cells]);
1974 let plane = [0.0_f32, 1.0, 0.0, -1.5];
1975 let (mesh, _) = extract_clipped_volume_faces(&data, &[plane]);
1976 match mesh.attributes.get("temp") {
1977 Some(AttributeData::Face(vals)) => {
1978 let n_tris = mesh.indices.len() / 3;
1979 assert_eq!(vals.len(), n_tris, "one scalar per output triangle");
1980 for &v in vals {
1981 assert_eq!(v, 7.0, "scalar must equal the owning cell's value");
1982 }
1983 }
1984 other => panic!("expected Face attribute on clipped hex mesh, got {other:?}"),
1985 }
1986 }
1987
1988 fn single_pyramid() -> VolumeMeshData {
1993 let mut data = VolumeMeshData {
1995 positions: vec![
1996 [0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [1.0, 0.0, 1.0], [0.0, 0.0, 1.0], [0.5, 1.0, 0.5], ],
2002 ..Default::default()
2003 };
2004 data.push_pyramid([0, 1, 2, 3], 4);
2005 data
2006 }
2007
2008 fn single_wedge() -> VolumeMeshData {
2009 let mut data = VolumeMeshData {
2011 positions: vec![
2012 [0.0, 0.0, 0.0], [1.0, 0.0, 0.0], [0.5, 0.0, 1.0], [0.0, 1.0, 0.0], [1.0, 1.0, 0.0], [0.5, 1.0, 1.0], ],
2019 ..Default::default()
2020 };
2021 data.push_wedge([0, 1, 2], [3, 4, 5]);
2022 data
2023 }
2024
2025 fn tet_volume(p: [[f32; 3]; 4]) -> f32 {
2026 let v = |i: usize| -> [f32; 3] {
2028 [
2029 p[i][0] - p[0][0],
2030 p[i][1] - p[0][1],
2031 p[i][2] - p[0][2],
2032 ]
2033 };
2034 let (a, b, c) = (v(1), v(2), v(3));
2035 let cross = [
2036 b[1] * c[2] - b[2] * c[1],
2037 b[2] * c[0] - b[0] * c[2],
2038 b[0] * c[1] - b[1] * c[0],
2039 ];
2040 (a[0] * cross[0] + a[1] * cross[1] + a[2] * cross[2]) / 6.0
2041 }
2042
2043 #[test]
2044 fn decompose_tet_yields_one_tet() {
2045 let data = single_tet();
2046 let (tets, scalars) = decompose_to_tetrahedra(&data, "");
2047 assert_eq!(tets.len(), 1);
2048 assert_eq!(scalars.len(), 1);
2049 }
2050
2051 #[test]
2052 fn decompose_hex_yields_six_tets() {
2053 let data = single_hex();
2054 let (tets, scalars) = decompose_to_tetrahedra(&data, "");
2055 assert_eq!(tets.len(), 6);
2056 assert_eq!(scalars.len(), 6);
2057 }
2058
2059 #[test]
2060 fn decompose_pyramid_yields_two_tets() {
2061 let data = single_pyramid();
2062 let (tets, scalars) = decompose_to_tetrahedra(&data, "");
2063 assert_eq!(tets.len(), 2);
2064 assert_eq!(scalars.len(), 2);
2065 }
2066
2067 #[test]
2068 fn decompose_wedge_yields_three_tets() {
2069 let data = single_wedge();
2070 let (tets, scalars) = decompose_to_tetrahedra(&data, "");
2071 assert_eq!(tets.len(), 3);
2072 assert_eq!(scalars.len(), 3);
2073 }
2074
2075 #[test]
2076 fn decompose_output_tets_have_nonzero_volume() {
2077 for data in [single_tet(), single_hex(), single_pyramid(), single_wedge()] {
2078 let (tets, _) = decompose_to_tetrahedra(&data, "");
2079 for (i, t) in tets.iter().enumerate() {
2080 let vol = tet_volume(*t).abs();
2081 assert!(
2082 vol > 1e-6,
2083 "tet {i} has near-zero volume {vol}: {t:?}"
2084 );
2085 }
2086 }
2087 }
2088
2089 #[test]
2090 fn decompose_hex_volume_equals_cell_volume() {
2091 let data = single_hex();
2093 let (tets, _) = decompose_to_tetrahedra(&data, "");
2094 let total: f32 = tets.iter().map(|t| tet_volume(*t).abs()).sum();
2095 assert!(
2096 (total - 1.0).abs() < 1e-5,
2097 "unit hex volume should be 1.0, got {total}"
2098 );
2099 }
2100
2101 #[test]
2102 fn decompose_scalar_propagates_to_child_tets() {
2103 let mut data = single_hex();
2104 data.cell_scalars.insert("temp".to_string(), vec![42.0]);
2105 let (_, scalars) = decompose_to_tetrahedra(&data, "temp");
2106 assert_eq!(scalars.len(), 6);
2107 for &s in &scalars {
2108 assert_eq!(s, 42.0, "all child tets must inherit the cell scalar");
2109 }
2110 }
2111
2112 #[test]
2113 fn decompose_missing_attribute_falls_back_to_zero() {
2114 let data = single_hex();
2115 let (_, scalars) = decompose_to_tetrahedra(&data, "nonexistent");
2116 for &s in &scalars {
2117 assert_eq!(s, 0.0, "missing attribute must produce 0.0 per tet");
2118 }
2119 }
2120
2121 #[test]
2122 fn decompose_mixed_mesh_tet_counts_sum_correctly() {
2123 let mut data = VolumeMeshData {
2125 positions: vec![
2126 [0.0, 0.0, 0.0],
2128 [1.0, 0.0, 0.0],
2129 [0.5, 1.0, 0.0],
2130 [0.5, 0.5, 1.0],
2131 [2.0, 0.0, 0.0],
2133 [3.0, 0.0, 0.0],
2134 [3.0, 0.0, 1.0],
2135 [2.0, 0.0, 1.0],
2136 [2.0, 1.0, 0.0],
2137 [3.0, 1.0, 0.0],
2138 [3.0, 1.0, 1.0],
2139 [2.0, 1.0, 1.0],
2140 [4.0, 0.0, 0.0],
2142 [5.0, 0.0, 0.0],
2143 [5.0, 0.0, 1.0],
2144 [4.0, 0.0, 1.0],
2145 [4.5, 1.0, 0.5],
2146 [6.0, 0.0, 0.0],
2148 [7.0, 0.0, 0.0],
2149 [6.5, 0.0, 1.0],
2150 [6.0, 1.0, 0.0],
2151 [7.0, 1.0, 0.0],
2152 [6.5, 1.0, 1.0],
2153 ],
2154 ..Default::default()
2155 };
2156 data.push_tet(0, 1, 2, 3);
2157 data.push_hex([4, 5, 6, 7, 8, 9, 10, 11]);
2158 data.push_pyramid([12, 13, 14, 15], 16);
2159 data.push_wedge([17, 18, 19], [20, 21, 22]);
2160
2161 let (tets, scalars) = decompose_to_tetrahedra(&data, "");
2162 assert_eq!(tets.len(), 12, "1+6+2+3 = 12 tets");
2163 assert_eq!(scalars.len(), 12);
2164 }
2165}