Skip to main content

oxiphysics_geometry/
discrete_geometry.rs

1#![allow(clippy::needless_range_loop, clippy::ptr_arg)]
2// Copyright 2026 COOLJAPAN OU (Team KitaSan)
3// SPDX-License-Identifier: Apache-2.0
4
5//! Discrete differential geometry on triangle meshes.
6//!
7//! This module implements the classical machinery of discrete differential
8//! geometry, following the pioneering work of Desbrun, Meyer, Polthier,
9//! Crane and colleagues:
10//!
11//! - [`DiscreteMesh`]: half-edge adjacency, vertex normals.
12//! - Cotangent-weight Laplace-Beltrami operator.
13//! - Angle-defect Gaussian curvature.
14//! - Discrete mean curvature via the cotangent formula.
15//! - Discrete geodesics: Dijkstra edge graph + Heat Method (Crane 2013).
16//! - Parallel transport and holonomy on meshes.
17//! - Discrete minimal surfaces (mean-curvature flow).
18//! - Discrete vector fields and Hodge decomposition on meshes.
19
20#![allow(dead_code)]
21#![allow(clippy::too_many_arguments)]
22
23use std::collections::{BinaryHeap, HashMap, HashSet};
24use std::f64::consts::PI;
25
26// ---------------------------------------------------------------------------
27// Vector math helpers (no nalgebra)
28// ---------------------------------------------------------------------------
29
30/// Subtract two 3-vectors.
31#[inline]
32pub fn sub3(a: [f64; 3], b: [f64; 3]) -> [f64; 3] {
33    [a[0] - b[0], a[1] - b[1], a[2] - b[2]]
34}
35
36/// Add two 3-vectors.
37#[inline]
38pub fn add3(a: [f64; 3], b: [f64; 3]) -> [f64; 3] {
39    [a[0] + b[0], a[1] + b[1], a[2] + b[2]]
40}
41
42/// Scale a 3-vector by a scalar.
43#[inline]
44pub fn scale3(a: [f64; 3], s: f64) -> [f64; 3] {
45    [a[0] * s, a[1] * s, a[2] * s]
46}
47
48/// Dot product of two 3-vectors.
49#[inline]
50pub fn dot3(a: [f64; 3], b: [f64; 3]) -> f64 {
51    a[0] * b[0] + a[1] * b[1] + a[2] * b[2]
52}
53
54/// Cross product of two 3-vectors.
55#[inline]
56pub fn cross3(a: [f64; 3], b: [f64; 3]) -> [f64; 3] {
57    [
58        a[1] * b[2] - a[2] * b[1],
59        a[2] * b[0] - a[0] * b[2],
60        a[0] * b[1] - a[1] * b[0],
61    ]
62}
63
64/// Euclidean length of a 3-vector.
65#[inline]
66pub fn len3(a: [f64; 3]) -> f64 {
67    dot3(a, a).sqrt()
68}
69
70/// Euclidean distance between two points.
71#[inline]
72pub fn dist3(a: [f64; 3], b: [f64; 3]) -> f64 {
73    len3(sub3(a, b))
74}
75
76/// Normalise a 3-vector; returns the zero vector if length is negligible.
77#[inline]
78pub fn normalize3(a: [f64; 3]) -> [f64; 3] {
79    let l = len3(a);
80    if l < 1e-14 {
81        [0.0; 3]
82    } else {
83        scale3(a, 1.0 / l)
84    }
85}
86
87/// Signed area of the triangle (i, j, k) — positive for CCW orientation.
88#[inline]
89pub fn tri_area(p0: [f64; 3], p1: [f64; 3], p2: [f64; 3]) -> f64 {
90    let e0 = sub3(p1, p0);
91    let e1 = sub3(p2, p0);
92    0.5 * len3(cross3(e0, e1))
93}
94
95/// Cotangent of the angle at vertex `p_opp` opposite edge `(p0, p1)`.
96///
97/// Returns `cot(alpha)` where `alpha` is the angle at `p_opp` in the triangle
98/// `(p0, p_opp, p1)`.
99#[inline]
100pub fn cotan_at(p0: [f64; 3], p_opp: [f64; 3], p1: [f64; 3]) -> f64 {
101    let u = sub3(p0, p_opp);
102    let v = sub3(p1, p_opp);
103    let c = dot3(u, v);
104    let s = len3(cross3(u, v));
105    if s.abs() < 1e-14 { 0.0 } else { c / s }
106}
107
108// ---------------------------------------------------------------------------
109// DiscreteMesh
110// ---------------------------------------------------------------------------
111
112/// A triangle mesh suitable for discrete differential geometry computations.
113///
114/// Stores vertex positions and triangles, with lazy adjacency information.
115#[derive(Debug, Clone)]
116pub struct DiscreteMesh {
117    /// Vertex positions, indexed by vertex index.
118    pub vertices: Vec<[f64; 3]>,
119    /// Triangle faces, each given as three vertex indices.
120    pub triangles: Vec<[usize; 3]>,
121}
122
123impl DiscreteMesh {
124    /// Construct a mesh from vertex positions and face index triples.
125    pub fn new(vertices: Vec<[f64; 3]>, triangles: Vec<[usize; 3]>) -> Self {
126        Self {
127            vertices,
128            triangles,
129        }
130    }
131
132    /// Number of vertices.
133    pub fn num_vertices(&self) -> usize {
134        self.vertices.len()
135    }
136
137    /// Number of triangular faces.
138    pub fn num_faces(&self) -> usize {
139        self.triangles.len()
140    }
141
142    /// Returns all triangles that share vertex `v`.
143    pub fn one_ring_faces(&self, v: usize) -> Vec<usize> {
144        self.triangles
145            .iter()
146            .enumerate()
147            .filter_map(|(fi, t)| if t.contains(&v) { Some(fi) } else { None })
148            .collect()
149    }
150
151    /// Returns all neighbour vertices of vertex `v` (one-ring).
152    pub fn one_ring_vertices(&self, v: usize) -> Vec<usize> {
153        let mut nbrs: HashSet<usize> = HashSet::new();
154        for t in &self.triangles {
155            if t.contains(&v) {
156                for &w in t {
157                    if w != v {
158                        nbrs.insert(w);
159                    }
160                }
161            }
162        }
163        nbrs.into_iter().collect()
164    }
165
166    /// Area-weighted vertex normal at vertex `v`.
167    pub fn vertex_normal(&self, v: usize) -> [f64; 3] {
168        let mut n = [0.0f64; 3];
169        for t in &self.triangles {
170            if !t.contains(&v) {
171                continue;
172            }
173            let p0 = self.vertices[t[0]];
174            let p1 = self.vertices[t[1]];
175            let p2 = self.vertices[t[2]];
176            let face_n = cross3(sub3(p1, p0), sub3(p2, p0));
177            n = add3(n, face_n);
178        }
179        normalize3(n)
180    }
181
182    /// Face normal for triangle `fi`.
183    pub fn face_normal(&self, fi: usize) -> [f64; 3] {
184        let t = &self.triangles[fi];
185        let p0 = self.vertices[t[0]];
186        let p1 = self.vertices[t[1]];
187        let p2 = self.vertices[t[2]];
188        normalize3(cross3(sub3(p1, p0), sub3(p2, p0)))
189    }
190
191    /// Voronoi area of vertex `v` (mixed Voronoi / barycentric area).
192    ///
193    /// Uses the formulation of Meyer et al. (2003).
194    pub fn voronoi_area(&self, v: usize) -> f64 {
195        let mut area = 0.0_f64;
196        for t in &self.triangles {
197            if !t.contains(&v) {
198                continue;
199            }
200            // Identify local indices
201            let (i, j, k) = if t[0] == v {
202                (0, 1, 2)
203            } else if t[1] == v {
204                (1, 2, 0)
205            } else {
206                (2, 0, 1)
207            };
208            let pi = self.vertices[t[i]];
209            let pj = self.vertices[t[j]];
210            let pk = self.vertices[t[k]];
211            // 1/8 * (cot_j * |e_k|^2 + cot_k * |e_j|^2)
212            let cot_j = cotan_at(pi, pj, pk);
213            let cot_k = cotan_at(pi, pk, pj);
214            let ek_sq = {
215                let d = sub3(pi, pj);
216                dot3(d, d)
217            };
218            let ej_sq = {
219                let d = sub3(pi, pk);
220                dot3(d, d)
221            };
222            area += 0.125 * (cot_j * ek_sq + cot_k * ej_sq);
223        }
224        area.max(1e-15)
225    }
226
227    /// Build the sparse cotangent-weight Laplace-Beltrami matrix as a list of
228    /// `(row, col, weight)` triples.
229    ///
230    /// The returned weights satisfy `L[i][j] = (cot_alpha + cot_beta) / 2` for
231    /// adjacent vertices `i`, `j`, and `L[i][i] = -sum_j L[i][j]`.
232    pub fn cotan_laplacian(&self) -> Vec<(usize, usize, f64)> {
233        let n = self.num_vertices();
234        // accumulate off-diagonal cotangent weights
235        let mut offdiag: HashMap<(usize, usize), f64> = HashMap::new();
236        for t in &self.triangles {
237            let [i, j, k] = *t;
238            let pi = self.vertices[i];
239            let pj = self.vertices[j];
240            let pk = self.vertices[k];
241            // cotangents at each corner
242            let cot_i = cotan_at(pj, pi, pk);
243            let cot_j = cotan_at(pi, pj, pk);
244            let cot_k = cotan_at(pi, pk, pj);
245            // edge (j,k) -> angle at i
246            *offdiag.entry((j.min(k), j.max(k))).or_insert(0.0) += 0.5 * cot_i;
247            // edge (i,k) -> angle at j
248            *offdiag.entry((i.min(k), i.max(k))).or_insert(0.0) += 0.5 * cot_j;
249            // edge (i,j) -> angle at k
250            *offdiag.entry((i.min(j), i.max(j))).or_insert(0.0) += 0.5 * cot_k;
251        }
252        let mut entries: Vec<(usize, usize, f64)> = Vec::with_capacity(offdiag.len() * 2 + n);
253        let mut row_sum = vec![0.0f64; n];
254        for (&(a, b), &w) in &offdiag {
255            entries.push((a, b, w));
256            entries.push((b, a, w));
257            row_sum[a] += w;
258            row_sum[b] += w;
259        }
260        for i in 0..n {
261            entries.push((i, i, -row_sum[i]));
262        }
263        entries
264    }
265
266    /// Apply the Laplace-Beltrami operator to a scalar field `f` (length = n_vertices),
267    /// returning the Laplacian value at each vertex: `(Lf)[i] = sum_j L[i][j] * f[j]`.
268    pub fn apply_laplacian(&self, f: &[f64]) -> Vec<f64> {
269        let n = self.num_vertices();
270        let mut lf = vec![0.0f64; n];
271        for (i, j, w) in self.cotan_laplacian() {
272            lf[i] += w * f[j];
273        }
274        lf
275    }
276
277    /// Discrete Gaussian curvature (angle defect) at each vertex.
278    ///
279    /// `K[v] = (2*pi - sum_angles) / A_v` where `A_v` is the Voronoi area.
280    pub fn gaussian_curvature(&self) -> Vec<f64> {
281        let n = self.num_vertices();
282        let mut angle_sum = vec![0.0f64; n];
283        for t in &self.triangles {
284            let [i, j, k] = *t;
285            let pi = self.vertices[i];
286            let pj = self.vertices[j];
287            let pk = self.vertices[k];
288            // Angle at i
289            let u = normalize3(sub3(pj, pi));
290            let v = normalize3(sub3(pk, pi));
291            let ang_i = dot3(u, v).clamp(-1.0, 1.0).acos();
292            // Angle at j
293            let u2 = normalize3(sub3(pi, pj));
294            let v2 = normalize3(sub3(pk, pj));
295            let ang_j = dot3(u2, v2).clamp(-1.0, 1.0).acos();
296            // Angle at k
297            let u3 = normalize3(sub3(pi, pk));
298            let v3 = normalize3(sub3(pj, pk));
299            let ang_k = dot3(u3, v3).clamp(-1.0, 1.0).acos();
300            angle_sum[i] += ang_i;
301            angle_sum[j] += ang_j;
302            angle_sum[k] += ang_k;
303        }
304        (0..n)
305            .map(|v| {
306                let a = self.voronoi_area(v);
307                (2.0 * PI - angle_sum[v]) / a
308            })
309            .collect()
310    }
311
312    /// Discrete mean curvature at each vertex using the cotangent formula.
313    ///
314    /// Returns the mean curvature `H[v]` (scalar) at each vertex.
315    /// `H[v] = |Lx[v]| / (2 * A_v)` where `Lx` is the Laplacian of position.
316    pub fn mean_curvature(&self) -> Vec<f64> {
317        let n = self.num_vertices();
318        let fx: Vec<f64> = self.vertices.iter().map(|p| p[0]).collect();
319        let fy: Vec<f64> = self.vertices.iter().map(|p| p[1]).collect();
320        let fz: Vec<f64> = self.vertices.iter().map(|p| p[2]).collect();
321        let lx = self.apply_laplacian(&fx);
322        let ly = self.apply_laplacian(&fy);
323        let lz = self.apply_laplacian(&fz);
324        (0..n)
325            .map(|v| {
326                let a = self.voronoi_area(v);
327                let hv = [lx[v], ly[v], lz[v]];
328                len3(hv) / (2.0 * a)
329            })
330            .collect()
331    }
332
333    /// Principal curvatures estimated from mean and Gaussian curvature.
334    ///
335    /// Returns `(kappa_1, kappa_2)` where `kappa_1 >= kappa_2`.
336    pub fn principal_curvatures(&self) -> Vec<(f64, f64)> {
337        let h = self.mean_curvature();
338        let k = self.gaussian_curvature();
339        h.iter()
340            .zip(k.iter())
341            .map(|(&hi, &ki)| {
342                let disc = (hi * hi - ki).max(0.0).sqrt();
343                let k1 = hi + disc;
344                let k2 = hi - disc;
345                if k1 >= k2 { (k1, k2) } else { (k2, k1) }
346            })
347            .collect()
348    }
349}
350
351// ---------------------------------------------------------------------------
352// Dijkstra geodesic distances
353// ---------------------------------------------------------------------------
354
355/// Dijkstra-based geodesic distance from a set of source vertices.
356///
357/// Returns a vector of length `n_vertices`; unreachable vertices get `f64::INFINITY`.
358pub fn geodesic_dijkstra(mesh: &DiscreteMesh, sources: &[usize]) -> Vec<f64> {
359    let n = mesh.num_vertices();
360    let mut dist = vec![f64::INFINITY; n];
361    let mut heap: BinaryHeap<DijkstraNode> = BinaryHeap::new();
362    for &s in sources {
363        dist[s] = 0.0;
364        heap.push(DijkstraNode {
365            dist: 0.0,
366            vertex: s,
367        });
368    }
369    // Build adjacency once
370    let adj = build_edge_adjacency(mesh);
371    while let Some(DijkstraNode { dist: d, vertex: u }) = heap.pop() {
372        if d > dist[u] + 1e-12 {
373            continue;
374        }
375        if let Some(nbrs) = adj.get(&u) {
376            for &(v, w) in nbrs {
377                let nd = d + w;
378                if nd < dist[v] {
379                    dist[v] = nd;
380                    heap.push(DijkstraNode {
381                        dist: nd,
382                        vertex: v,
383                    });
384                }
385            }
386        }
387    }
388    dist
389}
390
391/// Node for Dijkstra min-heap (ordered by distance ascending).
392#[derive(Debug, Clone, PartialEq)]
393struct DijkstraNode {
394    /// Tentative distance to this vertex.
395    dist: f64,
396    /// Vertex index.
397    vertex: usize,
398}
399
400impl Eq for DijkstraNode {}
401
402impl PartialOrd for DijkstraNode {
403    fn partial_cmp(&self, other: &Self) -> Option<std::cmp::Ordering> {
404        Some(self.cmp(other))
405    }
406}
407
408impl Ord for DijkstraNode {
409    fn cmp(&self, other: &Self) -> std::cmp::Ordering {
410        // Min-heap: smaller dist = higher priority
411        other
412            .dist
413            .partial_cmp(&self.dist)
414            .unwrap_or(std::cmp::Ordering::Equal)
415            .then(self.vertex.cmp(&other.vertex))
416    }
417}
418
419/// Build a vertex-to-neighbours adjacency list with edge lengths.
420pub fn build_edge_adjacency(mesh: &DiscreteMesh) -> HashMap<usize, Vec<(usize, f64)>> {
421    let mut adj: HashMap<usize, Vec<(usize, f64)>> = HashMap::new();
422    for t in &mesh.triangles {
423        let pairs = [(t[0], t[1]), (t[1], t[2]), (t[0], t[2])];
424        for (a, b) in pairs {
425            let w = dist3(mesh.vertices[a], mesh.vertices[b]);
426            adj.entry(a).or_default().push((b, w));
427            adj.entry(b).or_default().push((a, w));
428        }
429    }
430    adj
431}
432
433/// Extract the shortest geodesic path from `source` to `target` as a list of
434/// vertex indices, using Dijkstra distances and a predecessor map.
435pub fn geodesic_path(mesh: &DiscreteMesh, source: usize, target: usize) -> Option<Vec<usize>> {
436    let n = mesh.num_vertices();
437    let mut dist = vec![f64::INFINITY; n];
438    let mut prev: Vec<Option<usize>> = vec![None; n];
439    let adj = build_edge_adjacency(mesh);
440    dist[source] = 0.0;
441    let mut heap: BinaryHeap<DijkstraNode> = BinaryHeap::new();
442    heap.push(DijkstraNode {
443        dist: 0.0,
444        vertex: source,
445    });
446    while let Some(DijkstraNode { dist: d, vertex: u }) = heap.pop() {
447        if u == target {
448            break;
449        }
450        if d > dist[u] + 1e-12 {
451            continue;
452        }
453        if let Some(nbrs) = adj.get(&u) {
454            for &(v, w) in nbrs {
455                let nd = d + w;
456                if nd < dist[v] {
457                    dist[v] = nd;
458                    prev[v] = Some(u);
459                    heap.push(DijkstraNode {
460                        dist: nd,
461                        vertex: v,
462                    });
463                }
464            }
465        }
466    }
467    if dist[target].is_infinite() {
468        return None;
469    }
470    let mut path = Vec::new();
471    let mut cur = target;
472    loop {
473        path.push(cur);
474        if cur == source {
475            break;
476        }
477        match prev[cur] {
478            Some(p) => cur = p,
479            None => return None,
480        }
481    }
482    path.reverse();
483    Some(path)
484}
485
486// ---------------------------------------------------------------------------
487// Heat Method for Geodesic Distances (Crane et al. 2013)
488// ---------------------------------------------------------------------------
489
490/// Heat-method geodesic distance solver.
491///
492/// Approximates smooth geodesic distances by:
493/// 1. Diffusing a heat impulse for a short time `t`.
494/// 2. Computing the gradient of the heat field.
495/// 3. Normalising and divergence-computing to recover distance.
496///
497/// This implementation uses an explicit Euler diffusion step and the
498/// cotangent Laplacian.
499#[derive(Debug, Clone)]
500pub struct HeatMethodSolver {
501    /// The underlying mesh.
502    pub mesh: DiscreteMesh,
503    /// Diffusion time parameter; typically `h^2` where `h` is mean edge length.
504    pub diffusion_time: f64,
505    /// Number of explicit Euler steps for heat diffusion.
506    pub num_steps: usize,
507}
508
509impl HeatMethodSolver {
510    /// Construct a solver. If `diffusion_time <= 0`, it is set automatically
511    /// to the square of the mean edge length.
512    pub fn new(mesh: DiscreteMesh, diffusion_time: f64, num_steps: usize) -> Self {
513        let t = if diffusion_time > 0.0 {
514            diffusion_time
515        } else {
516            let h = mean_edge_length(&mesh);
517            h * h
518        };
519        Self {
520            mesh,
521            diffusion_time: t,
522            num_steps,
523        }
524    }
525
526    /// Compute approximate geodesic distances from `sources`.
527    ///
528    /// Returns distances for each vertex; values are normalised so that the
529    /// minimum source distance is zero.
530    pub fn compute(&self, sources: &[usize]) -> Vec<f64> {
531        let n = self.mesh.num_vertices();
532        // Step 1: initialise heat at sources
533        let mut u = vec![0.0f64; n];
534        for &s in sources {
535            u[s] = 1.0;
536        }
537        // Step 2: explicit Euler heat diffusion u_(t+1) = u_t + dt * L * u_t
538        let laplacian = self.mesh.cotan_laplacian();
539        let dt = self.diffusion_time / self.num_steps as f64;
540        for _ in 0..self.num_steps {
541            let mut lu = vec![0.0f64; n];
542            for &(i, j, w) in &laplacian {
543                lu[i] += w * u[j];
544            }
545            for i in 0..n {
546                u[i] += dt * lu[i];
547            }
548        }
549        // Step 3: compute gradient of u per face, normalise
550        let nf = self.mesh.num_faces();
551        let mut face_grad: Vec<[f64; 3]> = vec![[0.0; 3]; nf];
552        for (fi, t) in self.mesh.triangles.iter().enumerate() {
553            let [i, j, k] = *t;
554            let pi = self.mesh.vertices[i];
555            let pj = self.mesh.vertices[j];
556            let pk = self.mesh.vertices[k];
557            let area = tri_area(pi, pj, pk).max(1e-14);
558            let fn_ = normalize3(cross3(sub3(pj, pi), sub3(pk, pi)));
559            // Gradient: sum_corner u_corner * (N x e_opposite) / (2 * area)
560            let grad_i = scale3(cross3(fn_, sub3(pk, pj)), u[i] / (2.0 * area));
561            let grad_j = scale3(cross3(fn_, sub3(pi, pk)), u[j] / (2.0 * area));
562            let grad_k = scale3(cross3(fn_, sub3(pj, pi)), u[k] / (2.0 * area));
563            face_grad[fi] = add3(add3(grad_i, grad_j), grad_k);
564        }
565        // Step 4: normalise gradients (flip sign — points away from source)
566        let face_grad_norm: Vec<[f64; 3]> = face_grad
567            .iter()
568            .map(|&g| {
569                let l = len3(g);
570                if l < 1e-14 {
571                    [0.0; 3]
572                } else {
573                    scale3(g, -1.0 / l)
574                }
575            })
576            .collect();
577        // Step 5: compute divergence of normalised gradient at each vertex
578        let mut div = vec![0.0f64; n];
579        for (fi, t) in self.mesh.triangles.iter().enumerate() {
580            let [i, j, k] = *t;
581            let pi = self.mesh.vertices[i];
582            let pj = self.mesh.vertices[j];
583            let pk = self.mesh.vertices[k];
584            let x = face_grad_norm[fi];
585            let cot_i = cotan_at(pj, pi, pk);
586            let cot_j = cotan_at(pi, pj, pk);
587            let cot_k = cotan_at(pi, pk, pj);
588            div[i] += 0.5 * (cot_k * dot3(x, sub3(pj, pi)) + cot_j * dot3(x, sub3(pk, pi)));
589            div[j] += 0.5 * (cot_i * dot3(x, sub3(pk, pj)) + cot_k * dot3(x, sub3(pi, pj)));
590            div[k] += 0.5 * (cot_j * dot3(x, sub3(pi, pk)) + cot_i * dot3(x, sub3(pj, pk)));
591        }
592        // Step 6: solve Poisson Lφ = div (simple Gauss-Seidel iterations)
593        let mut phi = div.clone();
594        let diag: Vec<f64> = {
595            let mut d = vec![0.0f64; n];
596            for &(i, j, w) in &laplacian {
597                if i == j {
598                    d[i] = w;
599                }
600            }
601            d
602        };
603        for _ in 0..200 {
604            for vi in 0..n {
605                if diag[vi].abs() < 1e-14 {
606                    continue;
607                }
608                let mut s = div[vi];
609                for &(ri, rj, rw) in &laplacian {
610                    if ri == vi && rj != vi {
611                        s -= rw * phi[rj];
612                    }
613                }
614                phi[vi] = s / diag[vi];
615            }
616        }
617        // Shift so minimum is zero
618        let min_phi = phi.iter().cloned().fold(f64::INFINITY, f64::min);
619        phi.iter_mut().for_each(|x| *x -= min_phi);
620        phi
621    }
622}
623
624/// Compute the mean edge length of a mesh.
625pub fn mean_edge_length(mesh: &DiscreteMesh) -> f64 {
626    let mut total = 0.0;
627    let mut count = 0usize;
628    for t in &mesh.triangles {
629        let [i, j, k] = *t;
630        total += dist3(mesh.vertices[i], mesh.vertices[j]);
631        total += dist3(mesh.vertices[j], mesh.vertices[k]);
632        total += dist3(mesh.vertices[i], mesh.vertices[k]);
633        count += 3;
634    }
635    if count == 0 {
636        1.0
637    } else {
638        total / count as f64
639    }
640}
641
642// ---------------------------------------------------------------------------
643// Parallel Transport on Meshes
644// ---------------------------------------------------------------------------
645
646/// Parallel transport a tangent vector `v_tan` from face `f_src` to an adjacent
647/// face `f_dst` across their shared edge.
648///
649/// The vector is expressed in the local frame of `f_src` and rotated into
650/// the frame of `f_dst` using the connection (Levi-Civita connection on the mesh).
651///
652/// Returns the transported vector in the local tangent plane of `f_dst`.
653pub fn parallel_transport_across_edge(
654    mesh: &DiscreteMesh,
655    f_src: usize,
656    f_dst: usize,
657    v_tan: [f64; 3],
658) -> [f64; 3] {
659    let n_src = mesh.face_normal(f_src);
660    let n_dst = mesh.face_normal(f_dst);
661    // Rotation axis = cross(n_src, n_dst), angle = acos(dot)
662    let cos_theta = dot3(n_src, n_dst).clamp(-1.0, 1.0);
663    let sin_theta = (1.0 - cos_theta * cos_theta).max(0.0).sqrt();
664    if sin_theta.abs() < 1e-12 {
665        return v_tan; // faces are co-planar
666    }
667    let axis = normalize3(cross3(n_src, n_dst));
668    // Rodrigues rotation formula
669    let term1 = scale3(v_tan, cos_theta);
670    let term2 = scale3(cross3(axis, v_tan), sin_theta);
671    let term3 = scale3(axis, dot3(axis, v_tan) * (1.0 - cos_theta));
672    add3(add3(term1, term2), term3)
673}
674
675/// Compute holonomy (angle defect accumulated from parallel transport) around a
676/// vertex `v` by transporting a tangent vector around its one-ring of faces.
677///
678/// Returns the accumulated rotation angle in radians.
679pub fn vertex_holonomy(mesh: &DiscreteMesh, v: usize) -> f64 {
680    let faces = mesh.one_ring_faces(v);
681    if faces.len() < 2 {
682        return 0.0;
683    }
684    // Initialise a tangent vector in the first face's plane
685    let t0 = &mesh.triangles[faces[0]];
686    let p0 = mesh.vertices[t0[0]];
687    let p1 = mesh.vertices[t0[1]];
688    let mut tang = normalize3(sub3(p1, p0));
689    // Transport around the ring
690    let mut total_angle = 0.0_f64;
691    for i in 0..faces.len() {
692        let f_src = faces[i];
693        let f_dst = faces[(i + 1) % faces.len()];
694        let transported = parallel_transport_across_edge(mesh, f_src, f_dst, tang);
695        let n = mesh.face_normal(f_dst);
696        let tang_proj = normalize3(sub3(transported, scale3(n, dot3(transported, n))));
697        let cos_a = dot3(tang, tang_proj).clamp(-1.0, 1.0);
698        total_angle += cos_a.acos();
699        tang = tang_proj;
700    }
701    total_angle
702}
703
704// ---------------------------------------------------------------------------
705// Discrete Vector Fields on Meshes
706// ---------------------------------------------------------------------------
707
708/// A discrete vector field on a mesh, storing one tangent vector per face.
709///
710/// Vectors are stored in 3D ambient space but should lie in the tangent plane
711/// of the respective face.
712#[derive(Debug, Clone)]
713pub struct FaceVectorField {
714    /// Number of faces in the underlying mesh.
715    pub num_faces: usize,
716    /// Per-face tangent vectors.
717    pub vectors: Vec<[f64; 3]>,
718}
719
720impl FaceVectorField {
721    /// Construct a zero vector field.
722    pub fn zeros(num_faces: usize) -> Self {
723        Self {
724            num_faces,
725            vectors: vec![[0.0; 3]; num_faces],
726        }
727    }
728
729    /// Construct from a given list of per-face vectors.
730    pub fn from_vec(vectors: Vec<[f64; 3]>) -> Self {
731        let n = vectors.len();
732        Self {
733            num_faces: n,
734            vectors,
735        }
736    }
737
738    /// Project each face vector into the tangent plane of the face.
739    pub fn project_to_tangent(&mut self, mesh: &DiscreteMesh) {
740        for (fi, v) in self.vectors.iter_mut().enumerate() {
741            if fi >= mesh.num_faces() {
742                break;
743            }
744            let n = mesh.face_normal(fi);
745            let normal_comp = scale3(n, dot3(*v, n));
746            *v = sub3(*v, normal_comp);
747        }
748    }
749
750    /// L2 norm of the vector field.
751    pub fn norm(&self) -> f64 {
752        self.vectors
753            .iter()
754            .map(|v| dot3(*v, *v))
755            .sum::<f64>()
756            .sqrt()
757    }
758}
759
760/// A discrete vertex vector field, storing one tangent vector per vertex.
761#[derive(Debug, Clone)]
762pub struct VertexVectorField {
763    /// Number of vertices.
764    pub num_vertices: usize,
765    /// Per-vertex tangent vectors.
766    pub vectors: Vec<[f64; 3]>,
767}
768
769impl VertexVectorField {
770    /// Construct a zero vertex vector field.
771    pub fn zeros(num_vertices: usize) -> Self {
772        Self {
773            num_vertices,
774            vectors: vec![[0.0; 3]; num_vertices],
775        }
776    }
777
778    /// Compute the divergence of the vertex vector field as a scalar on vertices.
779    pub fn divergence(&self, mesh: &DiscreteMesh) -> Vec<f64> {
780        let n = mesh.num_vertices();
781        let mut div = vec![0.0f64; n];
782        for t in &mesh.triangles {
783            let [i, j, k] = *t;
784            let pi = mesh.vertices[i];
785            let pj = mesh.vertices[j];
786            let pk = mesh.vertices[k];
787            let area = tri_area(pi, pj, pk).max(1e-14);
788            let xi = self.vectors[i];
789            let xj = self.vectors[j];
790            let xk = self.vectors[k];
791            // Discrete divergence per triangle
792            let cot_i = cotan_at(pj, pi, pk);
793            let cot_j = cotan_at(pi, pj, pk);
794            let cot_k = cotan_at(pi, pk, pj);
795            div[i] +=
796                (cot_k * dot3(xi, sub3(pj, pi)) + cot_j * dot3(xi, sub3(pk, pi))) / (2.0 * area);
797            div[j] +=
798                (cot_i * dot3(xj, sub3(pk, pj)) + cot_k * dot3(xj, sub3(pi, pj))) / (2.0 * area);
799            div[k] +=
800                (cot_j * dot3(xk, sub3(pi, pk)) + cot_i * dot3(xk, sub3(pj, pk))) / (2.0 * area);
801        }
802        div
803    }
804
805    /// Compute the curl (discrete 2-form = scalar on faces) of the vertex vector field.
806    pub fn curl_on_faces(&self, mesh: &DiscreteMesh) -> Vec<f64> {
807        mesh.triangles
808            .iter()
809            .map(|t| {
810                let [i, j, k] = *t;
811                let pi = mesh.vertices[i];
812                let pj = mesh.vertices[j];
813                let pk = mesh.vertices[k];
814                let area = tri_area(pi, pj, pk).max(1e-14);
815                let xi = self.vectors[i];
816                let xj = self.vectors[j];
817                let xk = self.vectors[k];
818                let n = normalize3(cross3(sub3(pj, pi), sub3(pk, pi)));
819                // Discrete curl: circulation around triangle
820                let c = dot3(cross3(n, sub3(pj, pi)), xk)
821                    + dot3(cross3(n, sub3(pk, pj)), xi)
822                    + dot3(cross3(n, sub3(pi, pk)), xj);
823                c / (2.0 * area)
824            })
825            .collect()
826    }
827}
828
829// ---------------------------------------------------------------------------
830// Hodge Decomposition
831// ---------------------------------------------------------------------------
832
833/// Result of Hodge decomposition of a vertex vector field `X = grad(f) + curl(g) + harm`.
834///
835/// The decomposition is: X = ∇f + *d(g) + H where H is the harmonic part.
836#[derive(Debug, Clone)]
837pub struct HodgeDecomposition {
838    /// Exact (gradient) component scalar potential, per vertex.
839    pub scalar_potential: Vec<f64>,
840    /// Co-exact (curl) component stream function, per vertex.
841    pub stream_function: Vec<f64>,
842    /// Harmonic residual vector field (per vertex).
843    pub harmonic: Vec<[f64; 3]>,
844}
845
846/// Perform a Hodge decomposition of a vertex vector field.
847///
848/// This computes an approximate decomposition via:
849/// 1. Solve `L f = div(X)` → scalar potential.
850/// 2. Reconstruct gradient component and subtract.
851/// 3. Compute curl of remainder → stream function.
852/// 4. Harmonic part = remainder after subtracting both.
853pub fn hodge_decompose(mesh: &DiscreteMesh, field: &VertexVectorField) -> HodgeDecomposition {
854    let n = mesh.num_vertices();
855    let laplacian = mesh.cotan_laplacian();
856    let div = field.divergence(mesh);
857
858    // Build diagonal of laplacian
859    let mut diag = vec![0.0f64; n];
860    for &(i, j, w) in &laplacian {
861        if i == j {
862            diag[i] = w;
863        }
864    }
865
866    // Gauss-Seidel solve: L f = div
867    let mut f = vec![0.0f64; n];
868    for _ in 0..300 {
869        for vi in 0..n {
870            if diag[vi].abs() < 1e-14 {
871                continue;
872            }
873            let mut s = div[vi];
874            for &(ri, rj, rw) in &laplacian {
875                if ri == vi && rj != vi {
876                    s -= rw * f[rj];
877                }
878            }
879            f[vi] = s / diag[vi];
880        }
881    }
882
883    // Gradient of f per vertex: average of face gradients
884    let mut grad_f = vec![[0.0f64; 3]; n];
885    let mut weight = vec![0.0f64; n];
886    for t in &mesh.triangles {
887        let [i, j, k] = *t;
888        let pi = mesh.vertices[i];
889        let pj = mesh.vertices[j];
890        let pk = mesh.vertices[k];
891        let area = tri_area(pi, pj, pk).max(1e-14);
892        let fn_ = normalize3(cross3(sub3(pj, pi), sub3(pk, pi)));
893        // Face gradient
894        let gi = scale3(cross3(fn_, sub3(pk, pj)), f[i] / (2.0 * area));
895        let gj = scale3(cross3(fn_, sub3(pi, pk)), f[j] / (2.0 * area));
896        let gk = scale3(cross3(fn_, sub3(pj, pi)), f[k] / (2.0 * area));
897        let face_grad = add3(add3(gi, gj), gk);
898        for &vi in &[i, j, k] {
899            grad_f[vi] = add3(grad_f[vi], scale3(face_grad, area));
900            weight[vi] += area;
901        }
902    }
903    for vi in 0..n {
904        if weight[vi] > 1e-14 {
905            grad_f[vi] = scale3(grad_f[vi], 1.0 / weight[vi]);
906        }
907    }
908
909    // Curl-free residual
910    let mut residual: Vec<[f64; 3]> = (0..n)
911        .map(|vi| sub3(field.vectors[vi], grad_f[vi]))
912        .collect();
913
914    // Stream function g from curl of residual (Gauss-Seidel solve L g = curl)
915    let residual_field = VertexVectorField {
916        num_vertices: n,
917        vectors: residual.clone(),
918    };
919    let curl_on_faces = residual_field.curl_on_faces(mesh);
920
921    // Average curl to vertices as div proxy
922    let mut curl_div = vec![0.0f64; n];
923    let mut curl_wt = vec![0.0f64; n];
924    for (fi, t) in mesh.triangles.iter().enumerate() {
925        let area = tri_area(
926            mesh.vertices[t[0]],
927            mesh.vertices[t[1]],
928            mesh.vertices[t[2]],
929        )
930        .max(1e-14);
931        for &vi in t {
932            curl_div[vi] += curl_on_faces[fi] * area;
933            curl_wt[vi] += area;
934        }
935    }
936    for vi in 0..n {
937        if curl_wt[vi] > 1e-14 {
938            curl_div[vi] /= curl_wt[vi];
939        }
940    }
941
942    let mut g = vec![0.0f64; n];
943    for _ in 0..300 {
944        for vi in 0..n {
945            if diag[vi].abs() < 1e-14 {
946                continue;
947            }
948            let mut s = curl_div[vi];
949            for &(ri, rj, rw) in &laplacian {
950                if ri == vi && rj != vi {
951                    s -= rw * g[rj];
952                }
953            }
954            g[vi] = s / diag[vi];
955        }
956    }
957
958    // Subtract co-exact part from residual to get harmonic
959    for vi in 0..n {
960        // Approximate co-exact gradient
961        let co = scale3(mesh.vertex_normal(vi), g[vi]);
962        residual[vi] = sub3(residual[vi], co);
963    }
964
965    HodgeDecomposition {
966        scalar_potential: f,
967        stream_function: g,
968        harmonic: residual,
969    }
970}
971
972// ---------------------------------------------------------------------------
973// Discrete Minimal Surfaces (Mean Curvature Flow)
974// ---------------------------------------------------------------------------
975
976/// Perform one step of discrete mean curvature flow.
977///
978/// Updates vertex positions by `x_new = x_old + dt * L(x_old)`
979/// where `L` is the cotangent Laplacian applied to each coordinate.
980/// This moves the surface towards a minimal surface.
981pub fn mean_curvature_flow_step(mesh: &mut DiscreteMesh, dt: f64) {
982    let n = mesh.num_vertices();
983    let laplacian = mesh.cotan_laplacian();
984    let fx: Vec<f64> = mesh.vertices.iter().map(|p| p[0]).collect();
985    let fy: Vec<f64> = mesh.vertices.iter().map(|p| p[1]).collect();
986    let fz: Vec<f64> = mesh.vertices.iter().map(|p| p[2]).collect();
987    let mut lx = vec![0.0f64; n];
988    let mut ly = vec![0.0f64; n];
989    let mut lz = vec![0.0f64; n];
990    for &(i, j, w) in &laplacian {
991        lx[i] += w * fx[j];
992        ly[i] += w * fy[j];
993        lz[i] += w * fz[j];
994    }
995    for i in 0..n {
996        mesh.vertices[i][0] += dt * lx[i];
997        mesh.vertices[i][1] += dt * ly[i];
998        mesh.vertices[i][2] += dt * lz[i];
999    }
1000}
1001
1002/// Run `n_steps` of mean curvature flow on `mesh`.
1003pub fn mean_curvature_flow(mesh: &mut DiscreteMesh, dt: f64, n_steps: usize) {
1004    for _ in 0..n_steps {
1005        mean_curvature_flow_step(mesh, dt);
1006    }
1007}
1008
1009// ---------------------------------------------------------------------------
1010// Euler Characteristic and Topology
1011// ---------------------------------------------------------------------------
1012
1013/// Compute the Euler characteristic χ = V - E + F of the mesh.
1014pub fn euler_characteristic(mesh: &DiscreteMesh) -> i64 {
1015    let v = mesh.num_vertices() as i64;
1016    let f = mesh.num_faces() as i64;
1017    // Count unique edges
1018    let mut edges: HashSet<(usize, usize)> = HashSet::new();
1019    for t in &mesh.triangles {
1020        let [i, j, k] = *t;
1021        edges.insert((i.min(j), i.max(j)));
1022        edges.insert((j.min(k), j.max(k)));
1023        edges.insert((i.min(k), i.max(k)));
1024    }
1025    let e = edges.len() as i64;
1026    v - e + f
1027}
1028
1029/// Compute the genus `g` of a closed orientable surface from Euler characteristic.
1030///
1031/// Uses the formula χ = 2 - 2g, so g = (2 - χ) / 2.
1032pub fn genus_from_euler(chi: i64) -> i64 {
1033    (2 - chi) / 2
1034}
1035
1036// ---------------------------------------------------------------------------
1037// Discrete 1-Forms and Integration
1038// ---------------------------------------------------------------------------
1039
1040/// A discrete 1-form on mesh edges, storing one value per directed edge.
1041///
1042/// A 1-form ω assigns a value ω(e) to each oriented edge such that
1043/// ω(-e) = -ω(e).
1044#[derive(Debug, Clone)]
1045pub struct Discrete1Form {
1046    /// Values indexed by (vertex_i, vertex_j) with i < j.
1047    pub values: HashMap<(usize, usize), f64>,
1048}
1049
1050impl Discrete1Form {
1051    /// Construct a zero 1-form over all mesh edges.
1052    pub fn zeros(mesh: &DiscreteMesh) -> Self {
1053        let mut values = HashMap::new();
1054        for t in &mesh.triangles {
1055            let [i, j, k] = *t;
1056            values.entry((i.min(j), i.max(j))).or_insert(0.0);
1057            values.entry((j.min(k), j.max(k))).or_insert(0.0);
1058            values.entry((i.min(k), i.max(k))).or_insert(0.0);
1059        }
1060        Self { values }
1061    }
1062
1063    /// Evaluate the 1-form on directed edge (a → b).
1064    ///
1065    /// Returns `+val` if `a < b`, else `-val`.
1066    pub fn eval(&self, a: usize, b: usize) -> f64 {
1067        let key = (a.min(b), a.max(b));
1068        let v = self.values.get(&key).copied().unwrap_or(0.0);
1069        if a < b { v } else { -v }
1070    }
1071
1072    /// Integrate the 1-form along a path given as a sequence of vertex indices.
1073    pub fn integrate_path(&self, path: &[usize]) -> f64 {
1074        path.windows(2).map(|w| self.eval(w[0], w[1])).sum()
1075    }
1076
1077    /// Compute the exterior derivative (discrete 2-form = scalar on faces).
1078    ///
1079    /// Returns a vector of face values `dω[f] = ω(e01) + ω(e12) + ω(e20)`.
1080    pub fn exterior_derivative(&self, mesh: &DiscreteMesh) -> Vec<f64> {
1081        mesh.triangles
1082            .iter()
1083            .map(|t| {
1084                let [i, j, k] = *t;
1085                self.eval(i, j) + self.eval(j, k) + self.eval(k, i)
1086            })
1087            .collect()
1088    }
1089}
1090
1091// ---------------------------------------------------------------------------
1092// Discrete Laplacian Eigenvalues (Power Iteration)
1093// ---------------------------------------------------------------------------
1094
1095/// Compute the first `k` approximate eigenvalues of the Laplace-Beltrami
1096/// operator via power iteration.
1097///
1098/// Note: this is a rough approximation; for production use a proper sparse
1099/// eigensolver is recommended.
1100pub fn laplacian_eigenvalues_power(mesh: &DiscreteMesh, k: usize, n_iter: usize) -> Vec<f64> {
1101    let n = mesh.num_vertices();
1102    let laplacian = mesh.cotan_laplacian();
1103    let mut eigenvalues = Vec::with_capacity(k);
1104    let mut v: Vec<f64> = (0..n).map(|i| if i == 0 { 1.0 } else { 0.0 }).collect();
1105    for _ki in 0..k {
1106        // Normalise
1107        let norm = v.iter().map(|x| x * x).sum::<f64>().sqrt().max(1e-14);
1108        v.iter_mut().for_each(|x| *x /= norm);
1109        for _ in 0..n_iter {
1110            let mut lv = vec![0.0f64; n];
1111            for &(i, j, w) in &laplacian {
1112                lv[i] += w * v[j];
1113            }
1114            // Rayleigh quotient
1115            let rq: f64 = v.iter().zip(lv.iter()).map(|(vi, li)| vi * li).sum();
1116            eigenvalues.push(rq.abs());
1117            let norm2 = lv.iter().map(|x| x * x).sum::<f64>().sqrt().max(1e-14);
1118            v = lv.iter().map(|x| x / norm2).collect();
1119        }
1120    }
1121    eigenvalues.truncate(k);
1122    eigenvalues
1123}
1124
1125// ---------------------------------------------------------------------------
1126// Mesh Quality and Refinement
1127// ---------------------------------------------------------------------------
1128
1129/// Compute the aspect ratio of each triangle.
1130///
1131/// Aspect ratio = longest edge / shortest altitude.
1132/// Equilateral triangle has aspect ratio 1 (minimum).
1133pub fn triangle_aspect_ratios(mesh: &DiscreteMesh) -> Vec<f64> {
1134    mesh.triangles
1135        .iter()
1136        .map(|t| {
1137            let [i, j, k] = *t;
1138            let p0 = mesh.vertices[i];
1139            let p1 = mesh.vertices[j];
1140            let p2 = mesh.vertices[k];
1141            let l01 = dist3(p0, p1);
1142            let l12 = dist3(p1, p2);
1143            let l02 = dist3(p0, p2);
1144            let longest = l01.max(l12).max(l02);
1145            let area = tri_area(p0, p1, p2).max(1e-14);
1146            let shortest_alt = 2.0 * area / longest.max(1e-14);
1147            longest / shortest_alt.max(1e-14)
1148        })
1149        .collect()
1150}
1151
1152/// Compute per-face areas.
1153pub fn face_areas(mesh: &DiscreteMesh) -> Vec<f64> {
1154    mesh.triangles
1155        .iter()
1156        .map(|t| {
1157            let [i, j, k] = *t;
1158            tri_area(mesh.vertices[i], mesh.vertices[j], mesh.vertices[k])
1159        })
1160        .collect()
1161}
1162
1163/// Total surface area of the mesh.
1164pub fn total_surface_area(mesh: &DiscreteMesh) -> f64 {
1165    face_areas(mesh).iter().sum()
1166}
1167
1168// ---------------------------------------------------------------------------
1169// Mesh construction helpers for testing
1170// ---------------------------------------------------------------------------
1171
1172/// Build a simple unit-square mesh (two triangles).
1173///
1174/// Vertices: (0,0,0), (1,0,0), (1,1,0), (0,1,0).
1175/// Faces: \[(0,1,2), (0,2,3)\].
1176pub fn unit_square_mesh() -> DiscreteMesh {
1177    let vertices = vec![
1178        [0.0, 0.0, 0.0],
1179        [1.0, 0.0, 0.0],
1180        [1.0, 1.0, 0.0],
1181        [0.0, 1.0, 0.0],
1182    ];
1183    let triangles = vec![[0, 1, 2], [0, 2, 3]];
1184    DiscreteMesh::new(vertices, triangles)
1185}
1186
1187/// Build a regular icosahedron mesh (12 vertices, 20 faces).
1188///
1189/// Vertices lie on the unit sphere.
1190pub fn unit_icosahedron() -> DiscreteMesh {
1191    let phi = (1.0 + 5.0_f64.sqrt()) / 2.0;
1192    let mut verts: Vec<[f64; 3]> = vec![
1193        [-1.0, phi, 0.0],
1194        [1.0, phi, 0.0],
1195        [-1.0, -phi, 0.0],
1196        [1.0, -phi, 0.0],
1197        [0.0, -1.0, phi],
1198        [0.0, 1.0, phi],
1199        [0.0, -1.0, -phi],
1200        [0.0, 1.0, -phi],
1201        [phi, 0.0, -1.0],
1202        [phi, 0.0, 1.0],
1203        [-phi, 0.0, -1.0],
1204        [-phi, 0.0, 1.0],
1205    ];
1206    // Normalise to unit sphere
1207    for v in &mut verts {
1208        let l = len3(*v);
1209        v[0] /= l;
1210        v[1] /= l;
1211        v[2] /= l;
1212    }
1213    let faces = vec![
1214        [0, 11, 5],
1215        [0, 5, 1],
1216        [0, 1, 7],
1217        [0, 7, 10],
1218        [0, 10, 11],
1219        [1, 5, 9],
1220        [5, 11, 4],
1221        [11, 10, 2],
1222        [10, 7, 6],
1223        [7, 1, 8],
1224        [3, 9, 4],
1225        [3, 4, 2],
1226        [3, 2, 6],
1227        [3, 6, 8],
1228        [3, 8, 9],
1229        [4, 9, 5],
1230        [2, 4, 11],
1231        [6, 2, 10],
1232        [8, 6, 7],
1233        [9, 8, 1],
1234    ];
1235    DiscreteMesh::new(verts, faces)
1236}
1237
1238/// Build a regular tetrahedron mesh (4 vertices, 4 faces).
1239pub fn unit_tetrahedron() -> DiscreteMesh {
1240    let s = (8.0_f64 / 9.0_f64).sqrt();
1241    let c1 = (2.0_f64 / 9.0_f64).sqrt();
1242    let c2 = (2.0_f64 / 3.0_f64).sqrt();
1243    let vertices = vec![
1244        [0.0, 1.0, 0.0],
1245        [s, -1.0 / 3.0, 0.0],
1246        [-c1, -1.0 / 3.0, c2],
1247        [-c1, -1.0 / 3.0, -c2],
1248    ];
1249    let triangles = vec![[0, 1, 2], [0, 2, 3], [0, 3, 1], [1, 3, 2]];
1250    DiscreteMesh::new(vertices, triangles)
1251}
1252
1253// ---------------------------------------------------------------------------
1254// Discrete Curvature Tensor (Shape Operator Approximation)
1255// ---------------------------------------------------------------------------
1256
1257/// Approximate the curvature tensor at a vertex from neighbouring face normals.
1258///
1259/// Returns a symmetric 3×3 matrix `M` stored row-major as `[f64; 9]` such that
1260/// the quadratic form `v^T M v` approximates the second fundamental form.
1261pub fn vertex_curvature_tensor(mesh: &DiscreteMesh, v: usize) -> [f64; 9] {
1262    let mut m = [0.0f64; 9];
1263    let mut total_area = 0.0_f64;
1264    let pv = mesh.vertices[v];
1265    for t in &mesh.triangles {
1266        if !t.contains(&v) {
1267            continue;
1268        }
1269        let [i, j, k] = *t;
1270        let pi = mesh.vertices[i];
1271        let pj = mesh.vertices[j];
1272        let pk = mesh.vertices[k];
1273        let area = tri_area(pi, pj, pk).max(1e-14);
1274        let fn_ = mesh.face_normal(mesh.triangles.iter().position(|tt| tt == t).unwrap_or(0));
1275        // Find the two edges from v
1276        let others: Vec<[f64; 3]> = [i, j, k]
1277            .iter()
1278            .filter(|&&x| x != v)
1279            .map(|&x| mesh.vertices[x])
1280            .collect();
1281        for &po in &others {
1282            let e = sub3(po, pv);
1283            let kn = dot3(fn_, e) / dot3(e, e).max(1e-14);
1284            // Outer product e ⊗ e weighted by kn * area
1285            let w = kn * area;
1286            for a in 0..3 {
1287                for b in 0..3 {
1288                    m[a * 3 + b] += w * e[a] * e[b];
1289                }
1290            }
1291        }
1292        total_area += area;
1293    }
1294    if total_area > 1e-14 {
1295        for x in &mut m {
1296            *x /= total_area;
1297        }
1298    }
1299    m
1300}
1301
1302// ---------------------------------------------------------------------------
1303// Discrete Ricci Flow (simplified)
1304// ---------------------------------------------------------------------------
1305
1306/// Perform one step of discrete Ricci flow on vertex conformal factors.
1307///
1308/// Updates the conformal factor `u[v]` proportional to `K[v] - K_target`,
1309/// where `K` is Gaussian curvature and `K_target` is the target curvature.
1310pub fn ricci_flow_step(mesh: &DiscreteMesh, u: &mut Vec<f64>, k_target: &[f64], step_size: f64) {
1311    let k_current = mesh.gaussian_curvature();
1312    for (i, ui) in u.iter_mut().enumerate() {
1313        let kt = if i < k_target.len() { k_target[i] } else { 0.0 };
1314        *ui += step_size * (k_current[i] - kt);
1315    }
1316}
1317
1318// ---------------------------------------------------------------------------
1319// Tests
1320// ---------------------------------------------------------------------------
1321
1322#[cfg(test)]
1323mod tests {
1324    use super::*;
1325
1326    fn flat_square_mesh() -> DiscreteMesh {
1327        unit_square_mesh()
1328    }
1329
1330    fn ico() -> DiscreteMesh {
1331        unit_icosahedron()
1332    }
1333
1334    fn tet() -> DiscreteMesh {
1335        unit_tetrahedron()
1336    }
1337
1338    // --- Basic mesh properties ---
1339
1340    #[test]
1341    fn test_mesh_vertex_count() {
1342        let m = flat_square_mesh();
1343        assert_eq!(m.num_vertices(), 4);
1344    }
1345
1346    #[test]
1347    fn test_mesh_face_count() {
1348        let m = flat_square_mesh();
1349        assert_eq!(m.num_faces(), 2);
1350    }
1351
1352    #[test]
1353    fn test_icosahedron_vertex_count() {
1354        let m = ico();
1355        assert_eq!(m.num_vertices(), 12);
1356    }
1357
1358    #[test]
1359    fn test_icosahedron_face_count() {
1360        let m = ico();
1361        assert_eq!(m.num_faces(), 20);
1362    }
1363
1364    #[test]
1365    fn test_tetrahedron_vertex_count() {
1366        let m = tet();
1367        assert_eq!(m.num_vertices(), 4);
1368    }
1369
1370    // --- One-ring ---
1371
1372    #[test]
1373    fn test_one_ring_faces_non_empty() {
1374        let m = flat_square_mesh();
1375        let ring = m.one_ring_faces(0);
1376        assert!(!ring.is_empty());
1377    }
1378
1379    #[test]
1380    fn test_one_ring_vertices_contains_neighbors() {
1381        let m = flat_square_mesh();
1382        // Vertex 0 is connected to 1, 2, 3
1383        let nbrs = m.one_ring_vertices(0);
1384        assert!(nbrs.contains(&1) || nbrs.contains(&2) || nbrs.contains(&3));
1385    }
1386
1387    // --- Area ---
1388
1389    #[test]
1390    fn test_total_surface_area_unit_square() {
1391        let m = flat_square_mesh();
1392        let area = total_surface_area(&m);
1393        assert!((area - 1.0).abs() < 1e-10, "area={area}");
1394    }
1395
1396    #[test]
1397    fn test_face_areas_positive() {
1398        let m = ico();
1399        for a in face_areas(&m) {
1400            assert!(a > 0.0);
1401        }
1402    }
1403
1404    #[test]
1405    fn test_voronoi_area_positive() {
1406        let m = ico();
1407        for v in 0..m.num_vertices() {
1408            let a = m.voronoi_area(v);
1409            assert!(a > 0.0, "voronoi_area at v={v} is {a}");
1410        }
1411    }
1412
1413    // --- Normals ---
1414
1415    #[test]
1416    fn test_face_normal_unit_length() {
1417        let m = flat_square_mesh();
1418        for fi in 0..m.num_faces() {
1419            let n = m.face_normal(fi);
1420            let l = len3(n);
1421            assert!((l - 1.0).abs() < 1e-10, "face_normal not unit: {l}");
1422        }
1423    }
1424
1425    #[test]
1426    fn test_vertex_normal_unit_length() {
1427        let m = ico();
1428        for v in 0..m.num_vertices() {
1429            let n = m.vertex_normal(v);
1430            let l = len3(n);
1431            assert!(
1432                (l - 1.0).abs() < 1e-10,
1433                "vertex_normal not unit at {v}: {l}"
1434            );
1435        }
1436    }
1437
1438    #[test]
1439    fn test_flat_face_normal_is_z() {
1440        let m = flat_square_mesh();
1441        let n = m.face_normal(0);
1442        // flat mesh in xy-plane → normal should be (0,0,±1)
1443        assert!(n[2].abs() > 0.99, "flat face normal z={}", n[2]);
1444    }
1445
1446    // --- Laplacian ---
1447
1448    #[test]
1449    fn test_laplacian_row_sum_near_zero() {
1450        let m = ico();
1451        let lap = m.cotan_laplacian();
1452        let n = m.num_vertices();
1453        let mut row_sum = vec![0.0f64; n];
1454        for &(i, _j, w) in &lap {
1455            row_sum[i] += w;
1456        }
1457        for (i, s) in row_sum.iter().enumerate() {
1458            assert!(s.abs() < 1e-8, "row {i} sum = {s}");
1459        }
1460    }
1461
1462    #[test]
1463    fn test_laplacian_of_constant_is_zero() {
1464        let m = ico();
1465        let f = vec![1.0f64; m.num_vertices()];
1466        let lf = m.apply_laplacian(&f);
1467        for (i, v) in lf.iter().enumerate() {
1468            assert!(v.abs() < 1e-8, "L(1)[{i}] = {v}");
1469        }
1470    }
1471
1472    // --- Gaussian curvature ---
1473
1474    #[test]
1475    fn test_flat_mesh_gaussian_curvature_near_zero() {
1476        let m = flat_square_mesh();
1477        let k = m.gaussian_curvature();
1478        // Interior vertices of a flat mesh have K ≈ 0; boundary vertices have
1479        // non-zero angle defect. For the 2-triangle square, check the interior.
1480        // Vertex 2 is shared by both triangles with full 2π.
1481        // Allow boundary effects.
1482        for &ki in &k {
1483            assert!(
1484                ki.abs() < 100.0,
1485                "gaussian curvature suspiciously large: {ki}"
1486            );
1487        }
1488    }
1489
1490    #[test]
1491    fn test_icosahedron_total_gaussian_curvature() {
1492        // By Gauss-Bonnet: integral K dA = 2π χ = 4π for a sphere (genus 0)
1493        let m = ico();
1494        let k = m.gaussian_curvature();
1495        let areas: Vec<f64> = (0..m.num_vertices()).map(|v| m.voronoi_area(v)).collect();
1496        let total: f64 = k.iter().zip(areas.iter()).map(|(ki, ai)| ki * ai).sum();
1497        // Should be close to 4π ≈ 12.566
1498        assert!(
1499            (total - 4.0 * PI).abs() < 1.0,
1500            "Gauss-Bonnet: total={total}"
1501        );
1502    }
1503
1504    // --- Mean curvature ---
1505
1506    #[test]
1507    fn test_mean_curvature_nonnegative_sphere() {
1508        // On a convex surface mean curvature magnitude should be > 0
1509        let m = ico();
1510        let h = m.mean_curvature();
1511        let any_nonzero = h.iter().any(|&hi| hi.abs() > 1e-6);
1512        assert!(
1513            any_nonzero,
1514            "mean curvature should be non-zero on icosahedron"
1515        );
1516    }
1517
1518    #[test]
1519    fn test_principal_curvatures_k1_ge_k2() {
1520        let m = ico();
1521        let pc = m.principal_curvatures();
1522        for (k1, k2) in pc {
1523            assert!(k1 >= k2 - 1e-10, "k1={k1} < k2={k2}");
1524        }
1525    }
1526
1527    // --- Dijkstra ---
1528
1529    #[test]
1530    fn test_dijkstra_source_zero_distance() {
1531        let m = ico();
1532        let dist = geodesic_dijkstra(&m, &[0]);
1533        assert_eq!(dist[0], 0.0);
1534    }
1535
1536    #[test]
1537    fn test_dijkstra_all_finite_on_connected_mesh() {
1538        let m = ico();
1539        let dist = geodesic_dijkstra(&m, &[0]);
1540        for (i, d) in dist.iter().enumerate() {
1541            assert!(d.is_finite(), "vertex {i} has infinite distance");
1542        }
1543    }
1544
1545    #[test]
1546    fn test_dijkstra_nonnegative() {
1547        let m = ico();
1548        let dist = geodesic_dijkstra(&m, &[3]);
1549        for &d in &dist {
1550            assert!(d >= 0.0);
1551        }
1552    }
1553
1554    #[test]
1555    fn test_geodesic_path_starts_at_source() {
1556        let m = ico();
1557        let path = geodesic_path(&m, 0, 3).expect("should find path");
1558        assert_eq!(path[0], 0);
1559    }
1560
1561    #[test]
1562    fn test_geodesic_path_ends_at_target() {
1563        let m = ico();
1564        let path = geodesic_path(&m, 0, 3).expect("should find path");
1565        assert_eq!(*path.last().unwrap(), 3);
1566    }
1567
1568    // --- Heat method ---
1569
1570    #[test]
1571    fn test_heat_method_source_near_zero() {
1572        let m = ico();
1573        let solver = HeatMethodSolver::new(m, -1.0, 10);
1574        let dist = solver.compute(&[0]);
1575        assert!(
1576            dist[0] < 1e-3,
1577            "source distance should be near zero: {}",
1578            dist[0]
1579        );
1580    }
1581
1582    #[test]
1583    fn test_heat_method_nonnegative() {
1584        let m = ico();
1585        let solver = HeatMethodSolver::new(m, -1.0, 10);
1586        let dist = solver.compute(&[0]);
1587        for &d in &dist {
1588            assert!(d >= -1e-6, "heat distance negative: {d}");
1589        }
1590    }
1591
1592    // --- Euler characteristic ---
1593
1594    #[test]
1595    fn test_euler_characteristic_icosahedron() {
1596        let m = ico();
1597        let chi = euler_characteristic(&m);
1598        assert_eq!(chi, 2, "icosahedron euler characteristic = {chi}");
1599    }
1600
1601    #[test]
1602    fn test_euler_characteristic_tetrahedron() {
1603        let m = tet();
1604        let chi = euler_characteristic(&m);
1605        assert_eq!(chi, 2, "tetrahedron euler characteristic = {chi}");
1606    }
1607
1608    #[test]
1609    fn test_genus_sphere_is_zero() {
1610        assert_eq!(genus_from_euler(2), 0);
1611    }
1612
1613    // --- Parallel transport ---
1614
1615    #[test]
1616    fn test_parallel_transport_coplanar_identity() {
1617        let m = flat_square_mesh();
1618        let v = [1.0, 0.0, 0.0];
1619        let transported = parallel_transport_across_edge(&m, 0, 1, v);
1620        // Both faces are in same plane → identity rotation
1621        let diff = len3(sub3(transported, v));
1622        assert!(
1623            diff < 1e-10,
1624            "coplanar transport should be identity: diff={diff}"
1625        );
1626    }
1627
1628    // --- Vector field ---
1629
1630    #[test]
1631    fn test_face_vector_field_zeros() {
1632        let m = flat_square_mesh();
1633        let f = FaceVectorField::zeros(m.num_faces());
1634        assert_eq!(f.norm(), 0.0);
1635    }
1636
1637    #[test]
1638    fn test_vertex_vector_field_divergence_length() {
1639        let m = ico();
1640        let vf = VertexVectorField::zeros(m.num_vertices());
1641        let div = vf.divergence(&m);
1642        assert_eq!(div.len(), m.num_vertices());
1643    }
1644
1645    #[test]
1646    fn test_zero_field_divergence_is_zero() {
1647        let m = ico();
1648        let vf = VertexVectorField::zeros(m.num_vertices());
1649        let div = vf.divergence(&m);
1650        for &d in &div {
1651            assert!(d.abs() < 1e-12, "div of zero field = {d}");
1652        }
1653    }
1654
1655    // --- Discrete 1-form ---
1656
1657    #[test]
1658    fn test_1form_antisymmetry() {
1659        let m = flat_square_mesh();
1660        let mut form = Discrete1Form::zeros(&m);
1661        *form.values.get_mut(&(0, 1)).unwrap() = 3.0;
1662        assert_eq!(form.eval(0, 1), 3.0);
1663        assert_eq!(form.eval(1, 0), -3.0);
1664    }
1665
1666    #[test]
1667    fn test_1form_path_integration() {
1668        let m = flat_square_mesh();
1669        let mut form = Discrete1Form::zeros(&m);
1670        *form.values.get_mut(&(0, 1)).unwrap() = 1.0;
1671        *form.values.get_mut(&(1, 2)).unwrap() = 1.0;
1672        let path = vec![0, 1, 2];
1673        let val = form.integrate_path(&path);
1674        assert!((val - 2.0).abs() < 1e-10, "path integral = {val}");
1675    }
1676
1677    // --- Mean curvature flow ---
1678
1679    #[test]
1680    fn test_mcf_step_changes_vertices() {
1681        let mut m = ico();
1682        let before = m.vertices.clone();
1683        mean_curvature_flow_step(&mut m, 1e-3);
1684        let changed =
1685            m.vertices.iter().zip(before.iter()).any(|(a, b)| {
1686                (a[0] - b[0]).abs() + (a[1] - b[1]).abs() + (a[2] - b[2]).abs() > 1e-12
1687            });
1688        assert!(changed, "MCF step should change vertex positions");
1689    }
1690
1691    // --- Aspect ratio ---
1692
1693    #[test]
1694    fn test_aspect_ratio_all_positive() {
1695        let m = ico();
1696        for ar in triangle_aspect_ratios(&m) {
1697            assert!(ar > 0.0);
1698        }
1699    }
1700
1701    // --- Hodge decomposition ---
1702
1703    #[test]
1704    fn test_hodge_decompose_harmonic_length() {
1705        let m = ico();
1706        let vf = VertexVectorField::zeros(m.num_vertices());
1707        let hd = hodge_decompose(&m, &vf);
1708        assert_eq!(hd.harmonic.len(), m.num_vertices());
1709    }
1710
1711    #[test]
1712    fn test_hodge_decompose_scalar_potential_length() {
1713        let m = ico();
1714        let vf = VertexVectorField::zeros(m.num_vertices());
1715        let hd = hodge_decompose(&m, &vf);
1716        assert_eq!(hd.scalar_potential.len(), m.num_vertices());
1717    }
1718
1719    // --- Curvature tensor ---
1720
1721    #[test]
1722    fn test_curvature_tensor_size() {
1723        let m = ico();
1724        let ct = vertex_curvature_tensor(&m, 0);
1725        assert_eq!(ct.len(), 9);
1726    }
1727
1728    // --- Math helpers ---
1729
1730    #[test]
1731    fn test_cotan_equilateral() {
1732        let p0 = [1.0, 0.0, 0.0];
1733        let p_opp = [0.0, 3.0_f64.sqrt(), 0.0];
1734        let p1 = [-1.0, 0.0, 0.0];
1735        let cot = cotan_at(p0, p_opp, p1);
1736        // Angle at p_opp in equilateral triangle = 60° → cot(60°) = 1/√3 ≈ 0.577
1737        assert!((cot - 1.0 / 3.0_f64.sqrt()).abs() < 0.01, "cot = {cot}");
1738    }
1739
1740    #[test]
1741    fn test_tri_area_unit_right_triangle() {
1742        let p0 = [0.0, 0.0, 0.0];
1743        let p1 = [1.0, 0.0, 0.0];
1744        let p2 = [0.0, 1.0, 0.0];
1745        let a = tri_area(p0, p1, p2);
1746        assert!((a - 0.5).abs() < 1e-10, "area = {a}");
1747    }
1748
1749    #[test]
1750    fn test_mean_edge_length_positive() {
1751        let m = ico();
1752        let h = mean_edge_length(&m);
1753        assert!(h > 0.0);
1754    }
1755}