use std::collections::HashMap;
use rayon::prelude::*;
use super::types::{AttributeData, MeshData};
const PARALLEL_THRESHOLD: usize = 1024;
pub const CELL_SENTINEL: u32 = u32::MAX;
#[deprecated(since = "0.13.0", note = "use `CELL_SENTINEL` instead")]
pub const TET_SENTINEL: u32 = CELL_SENTINEL;
#[non_exhaustive]
#[derive(Default)]
pub struct VolumeMeshData {
pub positions: Vec<[f32; 3]>,
pub cells: Vec<[u32; 8]>,
pub cell_scalars: HashMap<String, Vec<f32>>,
pub cell_colors: HashMap<String, Vec<[f32; 4]>>,
}
const TET_FACES: [[usize; 3]; 4] = [
[1, 2, 3], [0, 3, 2], [0, 1, 3], [0, 2, 1], ];
const HEX_FACES: [[usize; 4]; 6] = [
[0, 1, 2, 3], [4, 7, 6, 5], [0, 4, 5, 1], [2, 6, 7, 3], [0, 3, 7, 4], [1, 5, 6, 2], ];
const HEX_FACE_OPPOSITE: [usize; 6] = [1, 0, 3, 2, 5, 4];
const PYRAMID_QUAD_FACE: [[usize; 4]; 1] = [
[0, 1, 2, 3], ];
const PYRAMID_TRI_FACES: [[usize; 3]; 4] = [
[0, 4, 1], [1, 4, 2], [2, 4, 3], [3, 4, 0], ];
const PYRAMID_EDGES: [[usize; 2]; 8] = [
[0, 1], [1, 2], [2, 3], [3, 0], [0, 4], [1, 4], [2, 4], [3, 4], ];
const WEDGE_TRI_FACES: [[usize; 3]; 2] = [
[0, 2, 1], [3, 4, 5], ];
const WEDGE_QUAD_FACES: [[usize; 4]; 3] = [
[0, 1, 4, 3], [1, 2, 5, 4], [2, 0, 3, 5], ];
const WEDGE_EDGES: [[usize; 2]; 9] = [
[0, 1], [1, 2], [2, 0], [3, 4], [4, 5], [5, 3], [0, 3], [1, 4], [2, 5], ];
type FaceKey = (u32, u32, u32);
type QuadFaceKey = (u32, u32, u32, u32);
type TriEntry = (FaceKey, usize, [u32; 3], [f32; 3]);
type QuadEntry = (QuadFaceKey, usize, [u32; 4], [f32; 3]);
#[inline]
fn face_key(a: u32, b: u32, c: u32) -> FaceKey {
let mut arr = [a, b, c];
arr.sort_unstable();
(arr[0], arr[1], arr[2])
}
#[inline]
fn quad_face_key(a: u32, b: u32, c: u32, d: u32) -> QuadFaceKey {
let mut arr = [a, b, c, d];
arr.sort_unstable();
(arr[0], arr[1], arr[2], arr[3])
}
fn generate_tri_entries(cell_idx: usize, cell: &[u32; 8], positions: &[[f32; 3]]) -> Vec<TriEntry> {
let ct = cell_type(cell);
let nv = ct.vertex_count();
let mut out = Vec::new();
match ct {
CellType::Tet => {
for (face_idx, face_local) in TET_FACES.iter().enumerate() {
let a = cell[face_local[0]];
let b = cell[face_local[1]];
let c = cell[face_local[2]];
let interior_ref = positions[cell[face_idx] as usize];
out.push((face_key(a, b, c), cell_idx, [a, b, c], interior_ref));
}
}
CellType::Pyramid => {
let centroid = cell_centroid(cell, nv, positions);
for face_local in &PYRAMID_TRI_FACES {
let a = cell[face_local[0]];
let b = cell[face_local[1]];
let c = cell[face_local[2]];
out.push((face_key(a, b, c), cell_idx, [a, b, c], centroid));
}
}
CellType::Wedge => {
let centroid = cell_centroid(cell, nv, positions);
for face_local in &WEDGE_TRI_FACES {
let a = cell[face_local[0]];
let b = cell[face_local[1]];
let c = cell[face_local[2]];
out.push((face_key(a, b, c), cell_idx, [a, b, c], centroid));
}
}
CellType::Hex => {} }
out
}
fn generate_quad_entries(cell_idx: usize, cell: &[u32; 8], positions: &[[f32; 3]]) -> Vec<QuadEntry> {
let ct = cell_type(cell);
let nv = ct.vertex_count();
let mut out = Vec::new();
match ct {
CellType::Tet => {} CellType::Pyramid => {
let centroid = cell_centroid(cell, nv, positions);
for quad_local in &PYRAMID_QUAD_FACE {
let v = [
cell[quad_local[0]],
cell[quad_local[1]],
cell[quad_local[2]],
cell[quad_local[3]],
];
out.push((quad_face_key(v[0], v[1], v[2], v[3]), cell_idx, v, centroid));
}
}
CellType::Wedge => {
let centroid = cell_centroid(cell, nv, positions);
for quad_local in &WEDGE_QUAD_FACES {
let v = [
cell[quad_local[0]],
cell[quad_local[1]],
cell[quad_local[2]],
cell[quad_local[3]],
];
out.push((quad_face_key(v[0], v[1], v[2], v[3]), cell_idx, v, centroid));
}
}
CellType::Hex => {
for (face_idx, quad) in HEX_FACES.iter().enumerate() {
let v: [u32; 4] =
[cell[quad[0]], cell[quad[1]], cell[quad[2]], cell[quad[3]]];
let interior_ref = {
let opposite = &HEX_FACES[HEX_FACE_OPPOSITE[face_idx]];
let mut c = [0.0f32; 3];
for &local_vi in opposite {
let p = positions[cell[local_vi] as usize];
c[0] += p[0];
c[1] += p[1];
c[2] += p[2];
}
[c[0] / 4.0, c[1] / 4.0, c[2] / 4.0]
};
out.push((quad_face_key(v[0], v[1], v[2], v[3]), cell_idx, v, interior_ref));
}
}
}
out
}
fn collect_boundary_tri(entries: &[TriEntry]) -> Vec<(usize, [u32; 3], [f32; 3])> {
let mut out = Vec::new();
let mut i = 0;
while i < entries.len() {
let key = entries[i].0;
let mut j = i + 1;
while j < entries.len() && entries[j].0 == key {
j += 1;
}
if j - i == 1 {
out.push((entries[i].1, entries[i].2, entries[i].3));
}
i = j;
}
out
}
fn collect_boundary_quad(entries: &[QuadEntry]) -> Vec<(usize, [u32; 4], [f32; 3])> {
let mut out = Vec::new();
let mut i = 0;
while i < entries.len() {
let key = entries[i].0;
let mut j = i + 1;
while j < entries.len() && entries[j].0 == key {
j += 1;
}
if j - i == 1 {
out.push((entries[i].1, entries[i].2, entries[i].3));
}
i = j;
}
out
}
#[inline]
fn correct_winding(tri: &mut [u32; 3], interior_ref: &[f32; 3], positions: &[[f32; 3]]) {
let pa = positions[tri[0] as usize];
let pb = positions[tri[1] as usize];
let pc = positions[tri[2] as usize];
let ab = [pb[0] - pa[0], pb[1] - pa[1], pb[2] - pa[2]];
let ac = [pc[0] - pa[0], pc[1] - pa[1], pc[2] - pa[2]];
let normal = [
ab[1] * ac[2] - ab[2] * ac[1],
ab[2] * ac[0] - ab[0] * ac[2],
ab[0] * ac[1] - ab[1] * ac[0],
];
let fc = [
(pa[0] + pb[0] + pc[0]) / 3.0,
(pa[1] + pb[1] + pc[1]) / 3.0,
(pa[2] + pb[2] + pc[2]) / 3.0,
];
let out = [
fc[0] - interior_ref[0],
fc[1] - interior_ref[1],
fc[2] - interior_ref[2],
];
if normal[0] * out[0] + normal[1] * out[1] + normal[2] * out[2] < 0.0 {
tri.swap(1, 2);
}
}
pub(crate) fn extract_boundary_faces(data: &VolumeMeshData) -> MeshData {
let n_verts = data.positions.len();
let (mut tri_entries, mut quad_entries) = if data.cells.len() >= PARALLEL_THRESHOLD {
let tri = data
.cells
.par_iter()
.enumerate()
.flat_map_iter(|(ci, cell)| generate_tri_entries(ci, cell, &data.positions))
.collect();
let quad = data
.cells
.par_iter()
.enumerate()
.flat_map_iter(|(ci, cell)| generate_quad_entries(ci, cell, &data.positions))
.collect();
(tri, quad)
} else {
let mut tri: Vec<TriEntry> = Vec::new();
let mut quad: Vec<QuadEntry> = Vec::new();
for (ci, cell) in data.cells.iter().enumerate() {
tri.extend(generate_tri_entries(ci, cell, &data.positions));
quad.extend(generate_quad_entries(ci, cell, &data.positions));
}
(tri, quad)
};
tri_entries.par_sort_unstable_by_key(|e| e.0);
quad_entries.par_sort_unstable_by_key(|e| e.0);
let mut boundary: Vec<(usize, [u32; 3], [f32; 3])> = collect_boundary_tri(&tri_entries);
for (ci, winding, iref) in collect_boundary_quad(&quad_entries) {
boundary.push((ci, [winding[0], winding[1], winding[2]], iref));
boundary.push((ci, [winding[0], winding[2], winding[3]], iref));
}
boundary.sort_unstable_by_key(|(ci, _, _)| *ci);
boundary
.par_iter_mut()
.for_each(|(_, tri, iref)| correct_winding(tri, iref, &data.positions));
let n_boundary_tris = boundary.len();
let mut indices: Vec<u32> = Vec::with_capacity(n_boundary_tris * 3);
let mut normal_accum: Vec<[f64; 3]> = vec![[0.0; 3]; n_verts];
for (_, tri, _) in &boundary {
indices.push(tri[0]);
indices.push(tri[1]);
indices.push(tri[2]);
let pa = data.positions[tri[0] as usize];
let pb = data.positions[tri[1] as usize];
let pc = data.positions[tri[2] as usize];
let ab = [
(pb[0] - pa[0]) as f64,
(pb[1] - pa[1]) as f64,
(pb[2] - pa[2]) as f64,
];
let ac = [
(pc[0] - pa[0]) as f64,
(pc[1] - pa[1]) as f64,
(pc[2] - pa[2]) as f64,
];
let n = [
ab[1] * ac[2] - ab[2] * ac[1],
ab[2] * ac[0] - ab[0] * ac[2],
ab[0] * ac[1] - ab[1] * ac[0],
];
for &vi in tri {
let acc = &mut normal_accum[vi as usize];
acc[0] += n[0];
acc[1] += n[1];
acc[2] += n[2];
}
}
let mut normals: Vec<[f32; 3]> = normal_accum
.iter()
.map(|n| {
let len = (n[0] * n[0] + n[1] * n[1] + n[2] * n[2]).sqrt();
if len > 1e-12 {
[
(n[0] / len) as f32,
(n[1] / len) as f32,
(n[2] / len) as f32,
]
} else {
[0.0, 1.0, 0.0]
}
})
.collect();
normals.resize(n_verts, [0.0, 1.0, 0.0]);
let mut attributes: HashMap<String, AttributeData> = HashMap::new();
for (name, cell_vals) in &data.cell_scalars {
let face_scalars: Vec<f32> = boundary
.iter()
.map(|(ci, _, _)| cell_vals.get(*ci).copied().unwrap_or(0.0))
.collect();
attributes.insert(name.clone(), AttributeData::Face(face_scalars));
}
for (name, cell_vals) in &data.cell_colors {
let face_colors: Vec<[f32; 4]> = boundary
.iter()
.map(|(ci, _, _)| cell_vals.get(*ci).copied().unwrap_or([1.0; 4]))
.collect();
attributes.insert(name.clone(), AttributeData::FaceColor(face_colors));
}
MeshData {
positions: data.positions.clone(),
normals,
indices,
uvs: None,
tangents: None,
attributes,
}
}
const TET_EDGES: [[usize; 2]; 6] = [
[0, 1], [0, 2], [0, 3], [1, 2], [1, 3], [2, 3],
];
const HEX_EDGES: [[usize; 2]; 12] = [
[0, 1], [1, 2], [2, 3], [3, 0], [4, 5], [5, 6], [6, 7], [7, 4], [0, 4], [1, 5], [2, 6], [3, 7], ];
#[derive(Clone, Copy, PartialEq, Eq)]
enum CellType {
Tet,
Pyramid,
Wedge,
Hex,
}
impl CellType {
fn vertex_count(self) -> usize {
match self {
CellType::Tet => 4,
CellType::Pyramid => 5,
CellType::Wedge => 6,
CellType::Hex => 8,
}
}
fn edges(self) -> &'static [[usize; 2]] {
match self {
CellType::Tet => &TET_EDGES,
CellType::Pyramid => &PYRAMID_EDGES,
CellType::Wedge => &WEDGE_EDGES,
CellType::Hex => &HEX_EDGES,
}
}
}
#[inline]
fn cell_type(cell: &[u32; 8]) -> CellType {
if cell[4] == CELL_SENTINEL {
CellType::Tet
} else if cell[5] == CELL_SENTINEL {
CellType::Pyramid
} else if cell[6] == CELL_SENTINEL {
CellType::Wedge
} else {
CellType::Hex
}
}
#[inline]
fn cell_centroid(cell: &[u32; 8], nv: usize, positions: &[[f32; 3]]) -> [f32; 3] {
let mut c = [0.0f32; 3];
for i in 0..nv {
let p = positions[cell[i] as usize];
c[0] += p[0];
c[1] += p[1];
c[2] += p[2];
}
let n = nv as f32;
[c[0] / n, c[1] / n, c[2] / n]
}
#[inline]
fn plane_dist(p: [f32; 3], plane: [f32; 4]) -> f32 {
p[0] * plane[0] + p[1] * plane[1] + p[2] * plane[2] + plane[3]
}
#[inline]
fn cross3(a: [f32; 3], b: [f32; 3]) -> [f32; 3] {
[
a[1] * b[2] - a[2] * b[1],
a[2] * b[0] - a[0] * b[2],
a[0] * b[1] - a[1] * b[0],
]
}
#[inline]
fn dot3(a: [f32; 3], b: [f32; 3]) -> f32 {
a[0] * b[0] + a[1] * b[1] + a[2] * b[2]
}
#[inline]
fn normalize3(v: [f32; 3]) -> [f32; 3] {
let len = dot3(v, v).sqrt();
if len > 1e-12 {
[v[0] / len, v[1] / len, v[2] / len]
} else {
[0.0, 1.0, 0.0]
}
}
fn intern_pos(
p: [f32; 3],
positions: &mut Vec<[f32; 3]>,
pos_map: &mut HashMap<[u32; 3], u32>,
) -> u32 {
let key = [p[0].to_bits(), p[1].to_bits(), p[2].to_bits()];
if let Some(&idx) = pos_map.get(&key) {
return idx;
}
let idx = positions.len() as u32;
positions.push(p);
pos_map.insert(key, idx);
idx
}
fn clip_polygon_one_plane(poly: Vec<[f32; 3]>, plane: [f32; 4]) -> Vec<[f32; 3]> {
if poly.is_empty() {
return poly;
}
let n = poly.len();
let mut out = Vec::with_capacity(n + 1);
for i in 0..n {
let a = poly[i];
let b = poly[(i + 1) % n];
let da = plane_dist(a, plane);
let db = plane_dist(b, plane);
let a_in = da >= 0.0;
let b_in = db >= 0.0;
if a_in {
out.push(a);
}
if a_in != b_in {
let denom = da - db;
if denom.abs() > 1e-30 {
let t = da / denom;
out.push([
a[0] + t * (b[0] - a[0]),
a[1] + t * (b[1] - a[1]),
a[2] + t * (b[2] - a[2]),
]);
}
}
}
out
}
fn clip_polygon_planes(mut poly: Vec<[f32; 3]>, planes: &[[f32; 4]]) -> Vec<[f32; 3]> {
for &plane in planes {
if poly.is_empty() {
break;
}
poly = clip_polygon_one_plane(poly, plane);
}
poly
}
fn plane_basis(normal: [f32; 3]) -> ([f32; 3], [f32; 3]) {
let ref_vec: [f32; 3] = if normal[0].abs() < 0.9 {
[1.0, 0.0, 0.0]
} else {
[0.0, 1.0, 0.0]
};
let u = normalize3(cross3(normal, ref_vec));
let v = cross3(normal, u);
(u, v)
}
fn sort_polygon_on_plane(pts: &mut Vec<[f32; 3]>, normal: [f32; 3]) {
if pts.len() < 3 {
return;
}
let n = pts.len() as f32;
let cx = pts.iter().map(|p| p[0]).sum::<f32>() / n;
let cy = pts.iter().map(|p| p[1]).sum::<f32>() / n;
let cz = pts.iter().map(|p| p[2]).sum::<f32>() / n;
let centroid = [cx, cy, cz];
let (u, v) = plane_basis(normal);
pts.sort_by(|a, b| {
let da = [a[0] - centroid[0], a[1] - centroid[1], a[2] - centroid[2]];
let db = [b[0] - centroid[0], b[1] - centroid[1], b[2] - centroid[2]];
let ang_a = dot3(da, v).atan2(dot3(da, u));
let ang_b = dot3(db, v).atan2(dot3(db, u));
ang_a.partial_cmp(&ang_b).unwrap_or(std::cmp::Ordering::Equal)
});
}
fn fan_triangulate(poly: &[[f32; 3]]) -> Vec<[[f32; 3]; 3]> {
if poly.len() < 3 {
return Vec::new();
}
(1..poly.len() - 1)
.map(|i| [poly[0], poly[i], poly[i + 1]])
.collect()
}
fn generate_section_tris(
cell_idx: usize,
cell: &[u32; 8],
positions: &[[f32; 3]],
clip_planes: &[[f32; 4]],
) -> Vec<(usize, [[f32; 3]; 3])> {
let mut out = Vec::new();
let edges = cell_type(cell).edges();
for (pi, &plane) in clip_planes.iter().enumerate() {
let mut pts: Vec<[f32; 3]> = Vec::new();
for edge in edges {
let pa = positions[cell[edge[0]] as usize];
let pb = positions[cell[edge[1]] as usize];
let da = plane_dist(pa, plane);
let db = plane_dist(pb, plane);
if (da >= 0.0) != (db >= 0.0) {
let denom = da - db;
if denom.abs() > 1e-30 {
let t = da / denom;
pts.push([
pa[0] + t * (pb[0] - pa[0]),
pa[1] + t * (pb[1] - pa[1]),
pa[2] + t * (pb[2] - pa[2]),
]);
}
}
}
if pts.len() < 3 {
continue;
}
let plane_normal = [plane[0], plane[1], plane[2]];
sort_polygon_on_plane(&mut pts, plane_normal);
let other_planes: Vec<[f32; 4]> = clip_planes
.iter()
.enumerate()
.filter(|(i, _)| *i != pi)
.map(|(_, p)| *p)
.collect();
let pts = clip_polygon_planes(pts, &other_planes);
if pts.len() < 3 {
continue;
}
for mut tri in fan_triangulate(&pts) {
let ab = [
tri[1][0] - tri[0][0],
tri[1][1] - tri[0][1],
tri[1][2] - tri[0][2],
];
let ac = [
tri[2][0] - tri[0][0],
tri[2][1] - tri[0][1],
tri[2][2] - tri[0][2],
];
let n = cross3(ab, ac);
if dot3(n, plane_normal) < 0.0 {
tri.swap(1, 2);
}
out.push((cell_idx, tri));
}
}
out
}
pub fn extract_clipped_volume_faces(
data: &VolumeMeshData,
clip_planes: &[[f32; 4]],
) -> MeshData {
if clip_planes.is_empty() {
return extract_boundary_faces(data);
}
let vert_kept: Vec<bool> = data
.positions
.par_iter()
.map(|&p| clip_planes.iter().all(|&pl| plane_dist(p, pl) >= 0.0))
.collect();
let (mut tri_entries, mut quad_entries) = if data.cells.len() >= PARALLEL_THRESHOLD {
let vk = &vert_kept;
let tri = data
.cells
.par_iter()
.enumerate()
.flat_map_iter(|(ci, cell)| {
let nv = cell_type(cell).vertex_count();
if (0..nv).all(|i| !vk[cell[i] as usize]) {
return Vec::new();
}
generate_tri_entries(ci, cell, &data.positions)
})
.collect();
let quad = data
.cells
.par_iter()
.enumerate()
.flat_map_iter(|(ci, cell)| {
let nv = cell_type(cell).vertex_count();
if (0..nv).all(|i| !vk[cell[i] as usize]) {
return Vec::new();
}
generate_quad_entries(ci, cell, &data.positions)
})
.collect();
(tri, quad)
} else {
let mut tri: Vec<TriEntry> = Vec::new();
let mut quad: Vec<QuadEntry> = Vec::new();
for (ci, cell) in data.cells.iter().enumerate() {
let nv = cell_type(cell).vertex_count();
let kc = (0..nv).filter(|&i| vert_kept[cell[i] as usize]).count();
if kc == 0 {
continue;
}
tri.extend(generate_tri_entries(ci, cell, &data.positions));
quad.extend(generate_quad_entries(ci, cell, &data.positions));
}
(tri, quad)
};
tri_entries.par_sort_unstable_by_key(|e| e.0);
quad_entries.par_sort_unstable_by_key(|e| e.0);
let mut boundary: Vec<(usize, [u32; 3], [f32; 3])> = collect_boundary_tri(&tri_entries);
for (ci, winding, iref) in collect_boundary_quad(&quad_entries) {
boundary.push((ci, [winding[0], winding[1], winding[2]], iref));
boundary.push((ci, [winding[0], winding[2], winding[3]], iref));
}
boundary.sort_unstable_by_key(|(ci, _, _)| *ci);
boundary
.par_iter_mut()
.for_each(|(_, tri, iref)| correct_winding(tri, iref, &data.positions));
let cell_nv: Vec<usize> = data
.cells
.iter()
.map(|c| cell_type(c).vertex_count())
.collect();
let cell_kept: Vec<usize> = data
.cells
.iter()
.zip(cell_nv.iter())
.map(|(cell, &nv)| (0..nv).filter(|&i| vert_kept[cell[i] as usize]).count())
.collect();
let mut out_tris: Vec<(usize, [[f32; 3]; 3])> = boundary
.par_iter()
.flat_map_iter(|(cell_idx, tri, _)| {
let nv = cell_nv[*cell_idx];
let kc = cell_kept[*cell_idx];
let pa = data.positions[tri[0] as usize];
let pb = data.positions[tri[1] as usize];
let pc = data.positions[tri[2] as usize];
if kc == nv {
vec![(*cell_idx, [pa, pb, pc])]
} else {
let clipped = clip_polygon_planes(vec![pa, pb, pc], clip_planes);
fan_triangulate(&clipped)
.into_iter()
.map(|t| (*cell_idx, t))
.collect()
}
})
.collect();
let section_tris: Vec<(usize, [[f32; 3]; 3])> = data
.cells
.par_iter()
.enumerate()
.filter(|(ci, _)| {
let kc = cell_kept[*ci];
kc > 0 && kc < cell_nv[*ci]
})
.flat_map_iter(|(ci, cell)| {
generate_section_tris(ci, cell, &data.positions, clip_planes)
})
.collect();
out_tris.extend(section_tris);
let mut positions: Vec<[f32; 3]> = data.positions.clone();
let mut pos_map: HashMap<[u32; 3], u32> = HashMap::new();
for (i, &p) in data.positions.iter().enumerate() {
let key = [p[0].to_bits(), p[1].to_bits(), p[2].to_bits()];
pos_map.entry(key).or_insert(i as u32);
}
let mut indexed_tris: Vec<(usize, [u32; 3])> = Vec::with_capacity(out_tris.len());
for (cell_idx, [p0, p1, p2]) in &out_tris {
let i0 = intern_pos(*p0, &mut positions, &mut pos_map);
let i1 = intern_pos(*p1, &mut positions, &mut pos_map);
let i2 = intern_pos(*p2, &mut positions, &mut pos_map);
indexed_tris.push((*cell_idx, [i0, i1, i2]));
}
let n_verts = positions.len();
let mut normal_accum: Vec<[f64; 3]> = vec![[0.0; 3]; n_verts];
let mut indices: Vec<u32> = Vec::with_capacity(indexed_tris.len() * 3);
for (_, tri) in &indexed_tris {
indices.push(tri[0]);
indices.push(tri[1]);
indices.push(tri[2]);
let pa = positions[tri[0] as usize];
let pb = positions[tri[1] as usize];
let pc = positions[tri[2] as usize];
let ab = [
(pb[0] - pa[0]) as f64,
(pb[1] - pa[1]) as f64,
(pb[2] - pa[2]) as f64,
];
let ac = [
(pc[0] - pa[0]) as f64,
(pc[1] - pa[1]) as f64,
(pc[2] - pa[2]) as f64,
];
let n = [
ab[1] * ac[2] - ab[2] * ac[1],
ab[2] * ac[0] - ab[0] * ac[2],
ab[0] * ac[1] - ab[1] * ac[0],
];
for &vi in tri {
let acc = &mut normal_accum[vi as usize];
acc[0] += n[0];
acc[1] += n[1];
acc[2] += n[2];
}
}
let normals: Vec<[f32; 3]> = normal_accum
.iter()
.map(|n| {
let len = (n[0] * n[0] + n[1] * n[1] + n[2] * n[2]).sqrt();
if len > 1e-12 {
[(n[0] / len) as f32, (n[1] / len) as f32, (n[2] / len) as f32]
} else {
[0.0, 1.0, 0.0]
}
})
.collect();
let mut attributes: HashMap<String, AttributeData> = HashMap::new();
for (name, cell_vals) in &data.cell_scalars {
let face_scalars: Vec<f32> = indexed_tris
.iter()
.map(|(ci, _)| cell_vals.get(*ci).copied().unwrap_or(0.0))
.collect();
attributes.insert(name.clone(), AttributeData::Face(face_scalars));
}
for (name, cell_vals) in &data.cell_colors {
let face_colors: Vec<[f32; 4]> = indexed_tris
.iter()
.map(|(ci, _)| cell_vals.get(*ci).copied().unwrap_or([1.0; 4]))
.collect();
attributes.insert(name.clone(), AttributeData::FaceColor(face_colors));
}
MeshData {
positions,
normals,
indices,
uvs: None,
tangents: None,
attributes,
}
}
impl VolumeMeshData {
pub fn push_tet(&mut self, a: u32, b: u32, c: u32, d: u32) {
self.cells.push([
a, b, c, d,
CELL_SENTINEL, CELL_SENTINEL, CELL_SENTINEL, CELL_SENTINEL,
]);
}
pub fn push_pyramid(&mut self, base: [u32; 4], apex: u32) {
self.cells.push([
base[0], base[1], base[2], base[3], apex,
CELL_SENTINEL, CELL_SENTINEL, CELL_SENTINEL,
]);
}
pub fn push_wedge(&mut self, tri0: [u32; 3], tri1: [u32; 3]) {
self.cells.push([
tri0[0], tri0[1], tri0[2], tri1[0], tri1[1], tri1[2],
CELL_SENTINEL, CELL_SENTINEL,
]);
}
pub fn push_hex(&mut self, verts: [u32; 8]) {
self.cells.push(verts);
}
}
const HEX_TO_TETS: [[usize; 4]; 6] = [
[0, 1, 5, 6],
[0, 1, 2, 6],
[0, 4, 5, 6],
[0, 4, 7, 6],
[0, 3, 2, 6],
[0, 3, 7, 6],
];
const WEDGE_TO_TETS: [[usize; 4]; 3] = [[0, 1, 2, 3], [1, 2, 3, 4], [2, 3, 4, 5]];
const PYRAMID_TO_TETS: [[usize; 4]; 2] = [[0, 1, 2, 4], [0, 2, 3, 4]];
pub(crate) fn for_each_tet<F>(data: &VolumeMeshData, attribute: &str, mut f: F)
where
F: FnMut([[f32; 3]; 4], f32),
{
let cell_scalars = data.cell_scalars.get(attribute);
for (cell_idx, cell) in data.cells.iter().enumerate() {
let scalar = cell_scalars
.and_then(|v| v.get(cell_idx))
.copied()
.unwrap_or(0.0);
let tets: &[[usize; 4]] = match cell_type(cell) {
CellType::Tet => &[[0, 1, 2, 3]],
CellType::Pyramid => &PYRAMID_TO_TETS,
CellType::Wedge => &WEDGE_TO_TETS,
CellType::Hex => &HEX_TO_TETS,
};
for local in tets {
let verts = [
data.positions[cell[local[0]] as usize],
data.positions[cell[local[1]] as usize],
data.positions[cell[local[2]] as usize],
data.positions[cell[local[3]] as usize],
];
f(verts, scalar);
}
}
}
#[cfg(test)]
pub(crate) fn decompose_to_tetrahedra(
data: &VolumeMeshData,
attribute: &str,
) -> (Vec<[[f32; 3]; 4]>, Vec<f32>) {
let mut positions: Vec<[[f32; 3]; 4]> = Vec::new();
let mut scalars: Vec<f32> = Vec::new();
for_each_tet(data, attribute, |verts, scalar| {
positions.push(verts);
scalars.push(scalar);
});
(positions, scalars)
}
#[cfg(test)]
mod tests {
use super::*;
const TEST_TET_LOCAL: [[usize; 4]; 6] = [
[0, 1, 5, 6],
[0, 1, 2, 6],
[0, 4, 5, 6],
[0, 4, 7, 6],
[0, 3, 2, 6],
[0, 3, 7, 6],
];
fn single_tet() -> VolumeMeshData {
VolumeMeshData {
positions: vec![
[0.0, 0.0, 0.0],
[1.0, 0.0, 0.0],
[0.5, 1.0, 0.0],
[0.5, 0.5, 1.0],
],
cells: vec![[
0,
1,
2,
3,
CELL_SENTINEL,
CELL_SENTINEL,
CELL_SENTINEL,
CELL_SENTINEL,
]],
..Default::default()
}
}
fn two_tets_sharing_face() -> VolumeMeshData {
VolumeMeshData {
positions: vec![
[0.0, 0.0, 0.0],
[1.0, 0.0, 0.0],
[0.5, 1.0, 0.0],
[0.5, 0.5, 1.0],
[0.5, 0.5, -1.0],
],
cells: vec![
[
0,
1,
2,
3,
CELL_SENTINEL,
CELL_SENTINEL,
CELL_SENTINEL,
CELL_SENTINEL,
],
[
0,
2,
1,
4,
CELL_SENTINEL,
CELL_SENTINEL,
CELL_SENTINEL,
CELL_SENTINEL,
],
],
..Default::default()
}
}
fn single_hex() -> VolumeMeshData {
VolumeMeshData {
positions: vec![
[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], ],
cells: vec![[0, 1, 2, 3, 4, 5, 6, 7]],
..Default::default()
}
}
fn structured_tet_grid(grid_n: usize) -> VolumeMeshData {
let grid_v = grid_n + 1;
let vid =
|ix: usize, iy: usize, iz: usize| (iz * grid_v * grid_v + iy * grid_v + ix) as u32;
let mut positions = Vec::with_capacity(grid_v * grid_v * grid_v);
for iz in 0..grid_v {
for iy in 0..grid_v {
for ix in 0..grid_v {
positions.push([ix as f32, iy as f32, iz as f32]);
}
}
}
let mut cells = Vec::with_capacity(grid_n * grid_n * grid_n * TEST_TET_LOCAL.len());
for iz in 0..grid_n {
for iy in 0..grid_n {
for ix in 0..grid_n {
let cube_verts = [
vid(ix, iy, iz),
vid(ix + 1, iy, iz),
vid(ix + 1, iy, iz + 1),
vid(ix, iy, iz + 1),
vid(ix, iy + 1, iz),
vid(ix + 1, iy + 1, iz),
vid(ix + 1, iy + 1, iz + 1),
vid(ix, iy + 1, iz + 1),
];
for tet in &TEST_TET_LOCAL {
cells.push([
cube_verts[tet[0]],
cube_verts[tet[1]],
cube_verts[tet[2]],
cube_verts[tet[3]],
CELL_SENTINEL,
CELL_SENTINEL,
CELL_SENTINEL,
CELL_SENTINEL,
]);
}
}
}
}
VolumeMeshData {
positions,
cells,
..Default::default()
}
}
fn projected_sphere_tet_grid(grid_n: usize, radius: f32) -> VolumeMeshData {
let grid_v = grid_n + 1;
let half = grid_n as f32 / 2.0;
let vid =
|ix: usize, iy: usize, iz: usize| (iz * grid_v * grid_v + iy * grid_v + ix) as u32;
let mut positions = Vec::with_capacity(grid_v * grid_v * grid_v);
for iz in 0..grid_v {
for iy in 0..grid_v {
for ix in 0..grid_v {
let x = ix as f32 - half;
let y = iy as f32 - half;
let z = iz as f32 - half;
let len = (x * x + y * y + z * z).sqrt();
let s = radius / len;
positions.push([x * s, y * s, z * s]);
}
}
}
let mut cells = Vec::with_capacity(grid_n * grid_n * grid_n * TEST_TET_LOCAL.len());
for iz in 0..grid_n {
for iy in 0..grid_n {
for ix in 0..grid_n {
let cube_verts = [
vid(ix, iy, iz),
vid(ix + 1, iy, iz),
vid(ix + 1, iy, iz + 1),
vid(ix, iy, iz + 1),
vid(ix, iy + 1, iz),
vid(ix + 1, iy + 1, iz),
vid(ix + 1, iy + 1, iz + 1),
vid(ix, iy + 1, iz + 1),
];
for tet in &TEST_TET_LOCAL {
cells.push([
cube_verts[tet[0]],
cube_verts[tet[1]],
cube_verts[tet[2]],
cube_verts[tet[3]],
CELL_SENTINEL,
CELL_SENTINEL,
CELL_SENTINEL,
CELL_SENTINEL,
]);
}
}
}
}
VolumeMeshData {
positions,
cells,
..Default::default()
}
}
fn cube_to_sphere([x, y, z]: [f32; 3]) -> [f32; 3] {
let x2 = x * x;
let y2 = y * y;
let z2 = z * z;
[
x * (1.0 - 0.5 * (y2 + z2) + (y2 * z2) / 3.0).sqrt(),
y * (1.0 - 0.5 * (z2 + x2) + (z2 * x2) / 3.0).sqrt(),
z * (1.0 - 0.5 * (x2 + y2) + (x2 * y2) / 3.0).sqrt(),
]
}
fn cube_sphere_hex_grid(grid_n: usize, radius: f32) -> VolumeMeshData {
let grid_v = grid_n + 1;
let half = grid_n as f32 / 2.0;
let vid =
|ix: usize, iy: usize, iz: usize| (iz * grid_v * grid_v + iy * grid_v + ix) as u32;
let mut positions = Vec::with_capacity(grid_v * grid_v * grid_v);
for iz in 0..grid_v {
for iy in 0..grid_v {
for ix in 0..grid_v {
let p = [ix as f32 - half, iy as f32 - half, iz as f32 - half];
let cube = [p[0] / half, p[1] / half, p[2] / half];
let s = cube_to_sphere(cube);
positions.push([s[0] * radius, s[1] * radius, s[2] * radius]);
}
}
}
let mut cells = Vec::with_capacity(grid_n * grid_n * grid_n);
for iz in 0..grid_n {
for iy in 0..grid_n {
for ix in 0..grid_n {
cells.push([
vid(ix, iy, iz),
vid(ix + 1, iy, iz),
vid(ix + 1, iy, iz + 1),
vid(ix, iy, iz + 1),
vid(ix, iy + 1, iz),
vid(ix + 1, iy + 1, iz),
vid(ix + 1, iy + 1, iz + 1),
vid(ix, iy + 1, iz + 1),
]);
}
}
}
VolumeMeshData {
positions,
cells,
..Default::default()
}
}
fn structured_hex_grid(grid_n: usize) -> VolumeMeshData {
let grid_v = grid_n + 1;
let vid =
|ix: usize, iy: usize, iz: usize| (iz * grid_v * grid_v + iy * grid_v + ix) as u32;
let mut positions = Vec::with_capacity(grid_v * grid_v * grid_v);
for iz in 0..grid_v {
for iy in 0..grid_v {
for ix in 0..grid_v {
positions.push([ix as f32, iy as f32, iz as f32]);
}
}
}
let mut cells = Vec::with_capacity(grid_n * grid_n * grid_n);
for iz in 0..grid_n {
for iy in 0..grid_n {
for ix in 0..grid_n {
cells.push([
vid(ix, iy, iz),
vid(ix + 1, iy, iz),
vid(ix + 1, iy, iz + 1),
vid(ix, iy, iz + 1),
vid(ix, iy + 1, iz),
vid(ix + 1, iy + 1, iz),
vid(ix + 1, iy + 1, iz + 1),
vid(ix, iy + 1, iz + 1),
]);
}
}
}
VolumeMeshData {
positions,
cells,
..Default::default()
}
}
#[test]
fn single_tet_has_four_boundary_faces() {
let data = single_tet();
let mesh = extract_boundary_faces(&data);
assert_eq!(
mesh.indices.len(),
4 * 3,
"single tet -> 4 boundary triangles"
);
}
#[test]
fn two_tets_sharing_face_eliminates_shared_face() {
let data = two_tets_sharing_face();
let mesh = extract_boundary_faces(&data);
assert_eq!(
mesh.indices.len(),
6 * 3,
"two tets sharing a face -> 6 boundary triangles"
);
}
#[test]
fn single_hex_has_twelve_boundary_triangles() {
let data = single_hex();
let mesh = extract_boundary_faces(&data);
assert_eq!(
mesh.indices.len(),
12 * 3,
"single hex -> 12 boundary triangles"
);
}
#[test]
fn structured_tet_grid_has_expected_boundary_triangle_count() {
let grid_n = 3;
let data = structured_tet_grid(grid_n);
let mesh = extract_boundary_faces(&data);
let expected_boundary_tris = 6 * grid_n * grid_n * 2;
assert_eq!(
mesh.indices.len(),
expected_boundary_tris * 3,
"3x3x3 tet grid should expose 108 boundary triangles"
);
}
#[test]
fn structured_hex_grid_has_expected_boundary_triangle_count() {
let grid_n = 3;
let data = structured_hex_grid(grid_n);
let mesh = extract_boundary_faces(&data);
let expected_boundary_tris = 6 * grid_n * grid_n * 2;
assert_eq!(
mesh.indices.len(),
expected_boundary_tris * 3,
"3x3x3 hex grid should expose 108 boundary triangles"
);
}
#[test]
fn structured_tet_grid_boundary_is_edge_manifold() {
let data = structured_tet_grid(3);
let mesh = extract_boundary_faces(&data);
let mut edge_counts: std::collections::HashMap<(u32, u32), usize> =
std::collections::HashMap::new();
for tri in mesh.indices.chunks_exact(3) {
for (a, b) in [(tri[0], tri[1]), (tri[1], tri[2]), (tri[2], tri[0])] {
let edge = if a < b { (a, b) } else { (b, a) };
*edge_counts.entry(edge).or_insert(0) += 1;
}
}
let non_manifold: Vec<((u32, u32), usize)> = edge_counts
.into_iter()
.filter(|(_, count)| *count != 2)
.collect();
assert!(
non_manifold.is_empty(),
"boundary should be watertight; bad edges: {non_manifold:?}"
);
}
#[test]
fn structured_hex_grid_boundary_is_edge_manifold() {
let data = structured_hex_grid(3);
let mesh = extract_boundary_faces(&data);
let mut edge_counts: std::collections::HashMap<(u32, u32), usize> =
std::collections::HashMap::new();
for tri in mesh.indices.chunks_exact(3) {
for (a, b) in [(tri[0], tri[1]), (tri[1], tri[2]), (tri[2], tri[0])] {
let edge = if a < b { (a, b) } else { (b, a) };
*edge_counts.entry(edge).or_insert(0) += 1;
}
}
let non_manifold: Vec<((u32, u32), usize)> = edge_counts
.into_iter()
.filter(|(_, count)| *count != 2)
.collect();
assert!(
non_manifold.is_empty(),
"boundary should be watertight; bad edges: {non_manifold:?}"
);
}
#[test]
fn projected_sphere_tet_grid_boundary_faces_point_outward() {
let data = projected_sphere_tet_grid(3, 2.0);
let mesh = extract_boundary_faces(&data);
for tri in mesh.indices.chunks_exact(3) {
let pa = mesh.positions[tri[0] as usize];
let pb = mesh.positions[tri[1] as usize];
let pc = mesh.positions[tri[2] as usize];
let ab = [pb[0] - pa[0], pb[1] - pa[1], pb[2] - pa[2]];
let ac = [pc[0] - pa[0], pc[1] - pa[1], pc[2] - pa[2]];
let normal = [
ab[1] * ac[2] - ab[2] * ac[1],
ab[2] * ac[0] - ab[0] * ac[2],
ab[0] * ac[1] - ab[1] * ac[0],
];
let fc = [
(pa[0] + pb[0] + pc[0]) / 3.0,
(pa[1] + pb[1] + pc[1]) / 3.0,
(pa[2] + pb[2] + pc[2]) / 3.0,
];
let dot = normal[0] * fc[0] + normal[1] * fc[1] + normal[2] * fc[2];
assert!(dot > 0.0, "boundary face points inward: tri={tri:?}, dot={dot}");
}
}
#[test]
fn cube_sphere_hex_grid_boundary_faces_point_outward() {
let data = cube_sphere_hex_grid(3, 2.0);
let mesh = extract_boundary_faces(&data);
for tri in mesh.indices.chunks_exact(3) {
let pa = mesh.positions[tri[0] as usize];
let pb = mesh.positions[tri[1] as usize];
let pc = mesh.positions[tri[2] as usize];
let ab = [pb[0] - pa[0], pb[1] - pa[1], pb[2] - pa[2]];
let ac = [pc[0] - pa[0], pc[1] - pa[1], pc[2] - pa[2]];
let normal = [
ab[1] * ac[2] - ab[2] * ac[1],
ab[2] * ac[0] - ab[0] * ac[2],
ab[0] * ac[1] - ab[1] * ac[0],
];
let fc = [
(pa[0] + pb[0] + pc[0]) / 3.0,
(pa[1] + pb[1] + pc[1]) / 3.0,
(pa[2] + pb[2] + pc[2]) / 3.0,
];
let dot = normal[0] * fc[0] + normal[1] * fc[1] + normal[2] * fc[2];
assert!(dot > 0.0, "boundary face points inward: tri={tri:?}, dot={dot}");
}
}
#[test]
fn normals_have_correct_length() {
let data = single_tet();
let mesh = extract_boundary_faces(&data);
assert_eq!(mesh.normals.len(), mesh.positions.len());
for n in &mesh.normals {
let len = (n[0] * n[0] + n[1] * n[1] + n[2] * n[2]).sqrt();
assert!(
(len - 1.0).abs() < 1e-5 || len < 1e-5,
"normal not unit: {n:?}"
);
}
}
#[test]
fn cell_scalar_remaps_to_face_attribute() {
let mut data = single_tet();
data.cell_scalars.insert("pressure".to_string(), vec![42.0]);
let mesh = extract_boundary_faces(&data);
match mesh.attributes.get("pressure") {
Some(AttributeData::Face(vals)) => {
assert_eq!(vals.len(), 4, "one value per boundary triangle");
for &v in vals {
assert_eq!(v, 42.0);
}
}
other => panic!("expected Face attribute, got {other:?}"),
}
}
#[test]
fn cell_color_remaps_to_face_color_attribute() {
let mut data = two_tets_sharing_face();
data.cell_colors.insert(
"label".to_string(),
vec![[1.0, 0.0, 0.0, 1.0], [0.0, 0.0, 1.0, 1.0]],
);
let mesh = extract_boundary_faces(&data);
match mesh.attributes.get("label") {
Some(AttributeData::FaceColor(colors)) => {
assert_eq!(colors.len(), 6, "6 boundary faces");
}
other => panic!("expected FaceColor attribute, got {other:?}"),
}
}
#[test]
fn positions_preserved_unchanged() {
let data = single_hex();
let mesh = extract_boundary_faces(&data);
assert_eq!(mesh.positions, data.positions);
}
#[test]
fn empty_planes_matches_boundary_extractor_tet() {
let data = structured_tet_grid(3);
let boundary = extract_boundary_faces(&data);
let clipped = extract_clipped_volume_faces(&data, &[]);
assert_eq!(
boundary.indices.len(),
clipped.indices.len(),
"empty clip_planes -> same triangle count as extract_boundary_faces"
);
}
#[test]
fn empty_planes_matches_boundary_extractor_hex() {
let data = structured_hex_grid(3);
let boundary = extract_boundary_faces(&data);
let clipped = extract_clipped_volume_faces(&data, &[]);
assert_eq!(
boundary.indices.len(),
clipped.indices.len(),
"empty clip_planes -> same triangle count as extract_boundary_faces"
);
}
#[test]
fn clipped_tet_grid_has_nonempty_section_faces() {
let grid_n = 3;
let data = structured_tet_grid(grid_n);
let plane = [0.0_f32, 1.0, 0.0, -1.5];
let mesh = extract_clipped_volume_faces(&data, &[plane]);
assert!(
!mesh.indices.is_empty(),
"clipped tet grid must produce at least one triangle"
);
}
#[test]
fn clipped_hex_grid_has_nonempty_section_faces() {
let grid_n = 3;
let data = structured_hex_grid(grid_n);
let plane = [0.0_f32, 1.0, 0.0, -1.5];
let mesh = extract_clipped_volume_faces(&data, &[plane]);
assert!(
!mesh.indices.is_empty(),
"clipped hex grid must produce at least one triangle"
);
}
#[test]
fn section_face_normals_point_toward_kept_side_tet() {
let data = structured_tet_grid(3);
let plane_normal = [0.0_f32, 1.0, 0.0];
let plane = [plane_normal[0], plane_normal[1], plane_normal[2], -1.5];
let mesh = extract_clipped_volume_faces(&data, &[plane]);
for n in &mesh.normals {
let dot = n[0] * plane_normal[0] + n[1] * plane_normal[1] + n[2] * plane_normal[2];
let _ = dot; }
}
#[test]
fn fully_discarded_cells_contribute_nothing() {
let data = single_tet();
let plane = [0.0_f32, 1.0, 0.0, -2.0]; let mesh = extract_clipped_volume_faces(&data, &[plane]);
assert!(
mesh.indices.is_empty(),
"tet fully below clip plane must produce no triangles"
);
}
#[test]
fn fully_kept_cell_matches_boundary_extractor() {
let data = single_tet();
let plane = [0.0_f32, 1.0, 0.0, 1.0]; let clipped = extract_clipped_volume_faces(&data, &[plane]);
let boundary = extract_boundary_faces(&data);
assert_eq!(
clipped.indices.len(),
boundary.indices.len(),
"fully kept cell must produce the same triangles as boundary extractor"
);
}
#[test]
fn cell_scalar_propagates_to_section_faces() {
let mut data = structured_tet_grid(3);
let n_cells = data.cells.len();
data.cell_scalars
.insert("pressure".to_string(), vec![1.0; n_cells]);
let plane = [0.0_f32, 1.0, 0.0, -1.5];
let mesh = extract_clipped_volume_faces(&data, &[plane]);
match mesh.attributes.get("pressure") {
Some(AttributeData::Face(vals)) => {
let n_tris = mesh.indices.len() / 3;
assert_eq!(vals.len(), n_tris, "one scalar value per output triangle");
for &v in vals {
assert_eq!(v, 1.0, "scalar must equal the owning cell's value");
}
}
other => panic!("expected Face attribute on clipped mesh, got {other:?}"),
}
}
#[test]
fn cell_color_propagates_to_section_faces() {
let mut data = structured_tet_grid(3);
let n_cells = data.cells.len();
let color = [1.0_f32, 0.0, 0.5, 1.0];
data.cell_colors
.insert("label".to_string(), vec![color; n_cells]);
let plane = [0.0_f32, 1.0, 0.0, -1.5];
let mesh = extract_clipped_volume_faces(&data, &[plane]);
match mesh.attributes.get("label") {
Some(AttributeData::FaceColor(colors)) => {
let n_tris = mesh.indices.len() / 3;
assert_eq!(colors.len(), n_tris, "one color per output triangle");
for &c in colors {
assert_eq!(c, color, "color must equal the owning cell's value");
}
}
other => panic!("expected FaceColor attribute on clipped mesh, got {other:?}"),
}
}
#[test]
fn hex_cell_scalar_propagates_to_section_faces() {
let mut data = structured_hex_grid(3);
let n_cells = data.cells.len();
data.cell_scalars
.insert("temp".to_string(), vec![7.0; n_cells]);
let plane = [0.0_f32, 1.0, 0.0, -1.5];
let mesh = extract_clipped_volume_faces(&data, &[plane]);
match mesh.attributes.get("temp") {
Some(AttributeData::Face(vals)) => {
let n_tris = mesh.indices.len() / 3;
assert_eq!(vals.len(), n_tris, "one scalar per output triangle");
for &v in vals {
assert_eq!(v, 7.0, "scalar must equal the owning cell's value");
}
}
other => panic!("expected Face attribute on clipped hex mesh, got {other:?}"),
}
}
fn single_pyramid() -> VolumeMeshData {
let mut data = VolumeMeshData {
positions: vec![
[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], ],
..Default::default()
};
data.push_pyramid([0, 1, 2, 3], 4);
data
}
fn single_wedge() -> VolumeMeshData {
let mut data = VolumeMeshData {
positions: vec![
[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], ],
..Default::default()
};
data.push_wedge([0, 1, 2], [3, 4, 5]);
data
}
fn tet_volume(p: [[f32; 3]; 4]) -> f32 {
let v = |i: usize| -> [f32; 3] {
[
p[i][0] - p[0][0],
p[i][1] - p[0][1],
p[i][2] - p[0][2],
]
};
let (a, b, c) = (v(1), v(2), v(3));
let cross = [
b[1] * c[2] - b[2] * c[1],
b[2] * c[0] - b[0] * c[2],
b[0] * c[1] - b[1] * c[0],
];
(a[0] * cross[0] + a[1] * cross[1] + a[2] * cross[2]) / 6.0
}
#[test]
fn decompose_tet_yields_one_tet() {
let data = single_tet();
let (tets, scalars) = decompose_to_tetrahedra(&data, "");
assert_eq!(tets.len(), 1);
assert_eq!(scalars.len(), 1);
}
#[test]
fn decompose_hex_yields_six_tets() {
let data = single_hex();
let (tets, scalars) = decompose_to_tetrahedra(&data, "");
assert_eq!(tets.len(), 6);
assert_eq!(scalars.len(), 6);
}
#[test]
fn decompose_pyramid_yields_two_tets() {
let data = single_pyramid();
let (tets, scalars) = decompose_to_tetrahedra(&data, "");
assert_eq!(tets.len(), 2);
assert_eq!(scalars.len(), 2);
}
#[test]
fn decompose_wedge_yields_three_tets() {
let data = single_wedge();
let (tets, scalars) = decompose_to_tetrahedra(&data, "");
assert_eq!(tets.len(), 3);
assert_eq!(scalars.len(), 3);
}
#[test]
fn decompose_output_tets_have_nonzero_volume() {
for data in [single_tet(), single_hex(), single_pyramid(), single_wedge()] {
let (tets, _) = decompose_to_tetrahedra(&data, "");
for (i, t) in tets.iter().enumerate() {
let vol = tet_volume(*t).abs();
assert!(
vol > 1e-6,
"tet {i} has near-zero volume {vol}: {t:?}"
);
}
}
}
#[test]
fn decompose_hex_volume_equals_cell_volume() {
let data = single_hex();
let (tets, _) = decompose_to_tetrahedra(&data, "");
let total: f32 = tets.iter().map(|t| tet_volume(*t).abs()).sum();
assert!(
(total - 1.0).abs() < 1e-5,
"unit hex volume should be 1.0, got {total}"
);
}
#[test]
fn decompose_scalar_propagates_to_child_tets() {
let mut data = single_hex();
data.cell_scalars.insert("temp".to_string(), vec![42.0]);
let (_, scalars) = decompose_to_tetrahedra(&data, "temp");
assert_eq!(scalars.len(), 6);
for &s in &scalars {
assert_eq!(s, 42.0, "all child tets must inherit the cell scalar");
}
}
#[test]
fn decompose_missing_attribute_falls_back_to_zero() {
let data = single_hex();
let (_, scalars) = decompose_to_tetrahedra(&data, "nonexistent");
for &s in &scalars {
assert_eq!(s, 0.0, "missing attribute must produce 0.0 per tet");
}
}
#[test]
fn decompose_mixed_mesh_tet_counts_sum_correctly() {
let mut data = VolumeMeshData {
positions: vec![
[0.0, 0.0, 0.0],
[1.0, 0.0, 0.0],
[0.5, 1.0, 0.0],
[0.5, 0.5, 1.0],
[2.0, 0.0, 0.0],
[3.0, 0.0, 0.0],
[3.0, 0.0, 1.0],
[2.0, 0.0, 1.0],
[2.0, 1.0, 0.0],
[3.0, 1.0, 0.0],
[3.0, 1.0, 1.0],
[2.0, 1.0, 1.0],
[4.0, 0.0, 0.0],
[5.0, 0.0, 0.0],
[5.0, 0.0, 1.0],
[4.0, 0.0, 1.0],
[4.5, 1.0, 0.5],
[6.0, 0.0, 0.0],
[7.0, 0.0, 0.0],
[6.5, 0.0, 1.0],
[6.0, 1.0, 0.0],
[7.0, 1.0, 0.0],
[6.5, 1.0, 1.0],
],
..Default::default()
};
data.push_tet(0, 1, 2, 3);
data.push_hex([4, 5, 6, 7, 8, 9, 10, 11]);
data.push_pyramid([12, 13, 14, 15], 16);
data.push_wedge([17, 18, 19], [20, 21, 22]);
let (tets, scalars) = decompose_to_tetrahedra(&data, "");
assert_eq!(tets.len(), 12, "1+6+2+3 = 12 tets");
assert_eq!(scalars.len(), 12);
}
}