Skip to main content

oxiphysics_geometry/
discrete_geometry.rs

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