spade/delaunay_core/
interpolation.rs

1use core::cell::RefCell;
2
3use crate::{
4    delaunay_core::math,
5    handles::{FixedDirectedEdgeHandle, FixedVertexHandle},
6    DelaunayTriangulation, HasPosition, HintGenerator, Point2, PositionInTriangulation,
7    Triangulation,
8};
9use num_traits::{one, zero, Float};
10
11use alloc::vec::Vec;
12
13use super::VertexHandle;
14
15/// Implements methods for natural neighbor interpolation.
16///
17/// Natural neighbor interpolation is a spatial interpolation method. For a given set of 2D input points with an
18/// associated value (e.g. a height field of some terrain or temperature measurements at different locations in
19/// a country), natural neighbor interpolation allows to smoothly interpolate the associated value for every
20/// location within the convex hull of the input points.
21///
22/// Spade currently assists with 4 interpolation strategies:
23///  - **Nearest neighbor interpolation:** Fastest. Exhibits too poor quality for many tasks. Not continuous
24///    along the edges of the voronoi diagram. Use [DelaunayTriangulation::nearest_neighbor].
25///  - **Barycentric interpolation:** Fast. Not smooth on the edges of the Delaunay triangulation.
26///    See [Barycentric].
27///  - **Natural neighbor interpolation:** Slower. Smooth everywhere except the input points. The input points
28///    have no derivative. See [NaturalNeighbor::interpolate]
29///  - **Natural neighbor interpolation with gradients:** Slowest. Smooth everywhere, even at the input points.
30///    See [NaturalNeighbor::interpolate_gradient].
31///
32/// # Performance comparison
33///
34/// In general, the speed of interpolating random points in a triangulation will be determined by the speed of the
35/// lookup operation after a certain triangulation size is reached. Thus, if random sampling is required, using a
36/// [crate::HierarchyHintGenerator] may be beneficial.
37///
38/// If subsequent queries are close to each other, [crate::LastUsedVertexHintGenerator] should also work fine.
39/// In this case, interpolating a single value should result in the following (relative) run times:
40///
41/// | nearest neighbor | barycentric | natural neighbor (no gradients) | natural neighbor (gradients) |
42/// | ---------------- | ----------- | ------------------------------- | ---------------------------- |
43/// | 23ns             | 103ns       | 307ns                           | 422ns                        |
44///
45/// # Usage
46///
47/// This type is created by calling [DelaunayTriangulation::natural_neighbor]. It contains a few internal buffers
48/// that are used to prevent recurring allocations. For best performance it should be created only once per thread
49/// and then used in all interpolation activities (see example).
50///
51/// # Example
52/// ```
53/// use spade::{Point2, HasPosition, DelaunayTriangulation, InsertionError, Triangulation as _};
54///
55/// struct PointWithHeight {
56///     position: Point2<f64>,
57///     height: f64,
58/// }
59///
60/// impl HasPosition for PointWithHeight {
61///     type Scalar = f64;
62///
63///     fn position(&self) -> Point2<f64> {
64///         self.position
65///     }
66/// }
67///
68/// fn main() -> Result<(), InsertionError> {
69///     let mut t = DelaunayTriangulation::<PointWithHeight>::new();
70///     t.insert(PointWithHeight { position: Point2::new(-1.0, -1.0), height: 42.0})?;
71///     t.insert(PointWithHeight { position: Point2::new(-1.0, 1.0), height: 13.37})?;
72///     t.insert(PointWithHeight { position: Point2::new(1.0, -1.0), height: 27.18})?;
73///     t.insert(PointWithHeight { position: Point2::new(1.0, 1.0), height: 31.41})?;
74///
75///     // Set of query points (would be many more realistically):
76///     let query_points = [Point2::new(0.0, 0.1), Point2::new(0.5, 0.1)];
77///
78///     // Good: Re-use interpolation object!
79///     let nn = t.natural_neighbor();
80///     for p in &query_points {
81///         println!("Value at {:?}: {:?}", p, nn.interpolate(|v| v.data().height, *p));
82///     }
83///
84///     // Bad (slower): Don't re-use internal buffers
85///     for p in &query_points {
86///         println!("Value at {:?}: {:?}", p, t.natural_neighbor().interpolate(|v| v.data().height, *p));
87///     }
88///     Ok(())
89/// }
90/// ```
91///
92/// # Visual comparison of interpolation algorithms
93///
94/// *Note: All of these images are generated by the "interpolation" example*
95///
96/// Nearest neighbor interpolation exhibits discontinuities along the voronoi edges:
97///
98#[doc = include_str!("../../images/interpolation_nearest_neighbor.img")]
99///
100/// Barycentric interpolation, by contrast, is continuous everywhere but has no derivative on the edges of the
101/// Delaunay triangulation. These show up as sharp corners in the color gradients on the drawn edges:
102///
103#[doc = include_str!("../../images/interpolation_barycentric.img")]
104///
105/// By contrast, natural neighbor interpolation is smooth on the edges - the previously sharp angles are now rounded
106/// off. However, the vertices themselves are still not continuous and will form sharp "peaks" in the resulting
107/// surface.
108///
109#[doc = include_str!("../../images/interpolation_nn_c0.img")]
110///
111/// With a gradient, the sharp peaks are gone - the surface will smoothly approximate a linear function as defined
112/// by the gradient in the vicinity of each vertex. The gradients are estimated with the estimate_gradients function.
113/// With flatness = 0.5
114#[doc = include_str!("../../images/interpolation_nn_c1_flatness_05.img")]
115/// With flatness = 1.0
116#[doc = include_str!("../../images/interpolation_nn_c1_flatness_1.img")]
117///
118#[doc(alias = "Interpolation")]
119pub struct NaturalNeighbor<'a, T>
120where
121    T: Triangulation,
122    T::Vertex: HasPosition,
123{
124    triangulation: &'a T,
125    // Various buffers that are used for identifying natural neighbors and their weights. It's significantly faster
126    // to clear and re-use these buffers instead of allocating them anew.
127    // We also don't run the risk of mutably borrowing them twice (causing a RefCell panic) as the RefCells never leak
128    // any of the interpolation methods.
129    // Not implementing `Sync` is also fine as the workaround - creating a new `NaturalNeighbor` instance per thread -
130    // is very simple.
131    inspect_edges_buffer: RefCell<Vec<FixedDirectedEdgeHandle>>,
132    natural_neighbor_buffer: RefCell<Vec<FixedDirectedEdgeHandle>>,
133    insert_cell_buffer: RefCell<Vec<Point2<<T::Vertex as HasPosition>::Scalar>>>,
134    weight_buffer: RefCell<Vec<(FixedVertexHandle, <T::Vertex as HasPosition>::Scalar)>>,
135}
136
137/// Implements methods related to barycentric interpolation.
138///
139/// Created by calling [crate::FloatTriangulation::barycentric].
140///
141/// Refer to the documentation of [NaturalNeighbor] for an overview of different interpolation methods.
142#[doc(alias = "Interpolation")]
143pub struct Barycentric<'a, T>
144where
145    T: Triangulation,
146{
147    triangulation: &'a T,
148    weight_buffer: RefCell<Vec<(FixedVertexHandle, <T::Vertex as HasPosition>::Scalar)>>,
149}
150
151impl<'a, T> Barycentric<'a, T>
152where
153    T: Triangulation,
154    <T::Vertex as HasPosition>::Scalar: Float,
155{
156    pub(crate) fn new(triangulation: &'a T) -> Self {
157        Self {
158            triangulation,
159            weight_buffer: Default::default(),
160        }
161    }
162
163    /// Returns the barycentric coordinates and the respective vertices for a given query position.
164    ///
165    /// The resulting coordinates and vertices are stored within the given `result` `vec`` to prevent
166    /// unneeded allocations. `result` will be cleared initially.
167    ///
168    /// The number of returned elements depends on the query positions location:
169    ///  - `result` will be **empty** if the query position lies outside the triangulation's convex hull
170    ///  - `result` will contain **a single element** (with weight 1.0) if the query position lies exactly on a vertex
171    ///  - `result` will contain **two vertices** if the query point lies exactly on any edge of the triangulation.
172    ///  - `result` will contain **exactly three** elements if the query point lies on an inner face of the
173    ///    triangulation.
174    pub fn get_weights(
175        &self,
176        position: Point2<<T::Vertex as HasPosition>::Scalar>,
177        result: &mut Vec<(FixedVertexHandle, <T::Vertex as HasPosition>::Scalar)>,
178    ) {
179        result.clear();
180        match self.triangulation.locate(position) {
181            PositionInTriangulation::OnVertex(vertex) => {
182                result.push((vertex, <T::Vertex as HasPosition>::Scalar::from(1.0)))
183            }
184            PositionInTriangulation::OnEdge(edge) => {
185                let [v0, v1] = self.triangulation.directed_edge(edge).vertices();
186                let [w0, w1] = two_point_interpolation::<T>(v0, v1, position);
187                result.push((v0.fix(), w0));
188                result.push((v1.fix(), w1));
189            }
190            PositionInTriangulation::OnFace(face) => {
191                let face = self.triangulation.face(face);
192                let [v0, v1, v2] = face.vertices();
193                let [c0, c1, c2] = face.barycentric_interpolation(position);
194                result.extend([(v0.fix(), c0), (v1.fix(), c1), (v2.fix(), c2)]);
195            }
196            _ => {}
197        }
198    }
199
200    /// Performs barycentric interpolation on this triangulation at a given position.
201    ///
202    /// Returns `None` for any value outside the triangulation's convex hull.
203    /// The value to interpolate is given by the `i` parameter.
204    ///
205    /// Refer to [NaturalNeighbor] for a comparison with other interpolation methods.
206    pub fn interpolate<I>(
207        &self,
208        i: I,
209        position: Point2<<T::Vertex as HasPosition>::Scalar>,
210    ) -> Option<<T::Vertex as HasPosition>::Scalar>
211    where
212        I: Fn(
213            VertexHandle<T::Vertex, T::DirectedEdge, T::UndirectedEdge, T::Face>,
214        ) -> <T::Vertex as HasPosition>::Scalar,
215    {
216        let nns = &mut *self.weight_buffer.borrow_mut();
217        self.get_weights(position, nns);
218        if nns.is_empty() {
219            return None;
220        }
221
222        let mut total_sum = zero();
223        for (vertex, weight) in nns {
224            total_sum = total_sum + i(self.triangulation.vertex(*vertex)) * *weight;
225        }
226        Some(total_sum)
227    }
228}
229
230impl<'a, V, DE, UE, F, L> NaturalNeighbor<'a, DelaunayTriangulation<V, DE, UE, F, L>>
231where
232    V: HasPosition,
233    DE: Default,
234    UE: Default,
235    F: Default,
236    L: HintGenerator<<V as HasPosition>::Scalar>,
237    <V as HasPosition>::Scalar: Float,
238{
239    pub(crate) fn new(triangulation: &'a DelaunayTriangulation<V, DE, UE, F, L>) -> Self {
240        Self {
241            triangulation,
242            inspect_edges_buffer: Default::default(),
243            insert_cell_buffer: Default::default(),
244            natural_neighbor_buffer: Default::default(),
245            weight_buffer: Default::default(),
246        }
247    }
248
249    /// Calculates the natural neighbors and their weights (sibson coordinates) of a given query position.
250    ///
251    /// The neighbors are returned in clockwise order. The weights will add up to 1.0.
252    /// The neighbors are stored in the `result` parameter to prevent unnecessary allocations.
253    /// `result` will be cleared initially.
254    ///
255    /// The number of returned natural neighbors depends on the given query position:
256    /// - `result` will be **empty** if the query position lies outside the triangulation's convex hull
257    /// - `result` will contain **exactly one** vertex if the query position is equal to that vertex position.
258    /// - `result` will contain **exactly two** entries if the query position lies exactly *on* an edge of the
259    ///    convex hull.
260    /// - `result` will contain **at least three** `(vertex, weight)` tuples if the query point lies on an inner
261    ///    face or an inner edge.
262    ///
263    /// *Example: The natural neighbors (red vertices) of the query point (blue dot) with their weights.
264    /// The elements will be returned in clockwise order as indicated by the indices drawn within the red circles.*
265    ///
266    #[doc = include_str!("../../images/natural_neighbor_scenario.svg")]
267    pub fn get_weights(
268        &self,
269        position: Point2<<V as HasPosition>::Scalar>,
270        result: &mut Vec<(FixedVertexHandle, <V as HasPosition>::Scalar)>,
271    ) {
272        let nns = &mut *self.natural_neighbor_buffer.borrow_mut();
273        get_natural_neighbor_edges(
274            self.triangulation,
275            &mut self.inspect_edges_buffer.borrow_mut(),
276            position,
277            nns,
278        );
279        self.get_natural_neighbor_weights(position, nns, result);
280    }
281
282    /// Interpolates a value at a given position.
283    ///
284    /// Returns `None` for any point outside the triangulations convex hull.
285    /// The value to interpolate is given by the `i` parameter. The resulting interpolation will be smooth
286    /// everywhere except at the input vertices.
287    ///
288    /// Refer to [NaturalNeighbor] for an example on how to use this function.
289    pub fn interpolate<I>(
290        &self,
291        i: I,
292        position: Point2<<V as HasPosition>::Scalar>,
293    ) -> Option<<V as HasPosition>::Scalar>
294    where
295        I: Fn(VertexHandle<V, DE, UE, F>) -> <V as HasPosition>::Scalar,
296    {
297        let nns = &mut *self.weight_buffer.borrow_mut();
298        self.get_weights(position, nns);
299        if nns.is_empty() {
300            return None;
301        }
302
303        let mut total_sum = zero();
304        for (vertex, weight) in nns {
305            total_sum = total_sum + i(self.triangulation.vertex(*vertex)) * *weight;
306        }
307        Some(total_sum)
308    }
309
310    /// Interpolates a value at a given position.
311    ///
312    /// In contrast to [Self::interpolate], this method has a well-defined derivative at each vertex and will
313    /// approximate a linear function in the proximity of any vertex.
314    ///
315    /// The value to interpolate is given by the `i` parameter. The gradient that defines the derivative at
316    /// each input vertex is given by the `g` parameter.
317    ///
318    /// The `flatness` parameter (should be >= 0.0) blends between an interpolation that adheres less to the gradient (small values)
319    /// or adheres to it strongly (values larger than ~2.0) in the vicinity of any vertex. A value of 0.5 reproduces
320    /// the original Sibson C1 interpolant. A value of 1.0 is the fastest and works well in many situations.
321    /// A value of 0.0 reproduces the I1 interpolant introduced in Flötotto's thesis chapter 4.4.
322    /// If both the flatness and the gradients are 0.0 Sibsons C0 interpolant is recreated and [Self::interpolate] should be used.
323    ///
324    /// Returns `None` for any point outside the triangulation's convex hull.
325    ///
326    /// Refer to [NaturalNeighbor] for more information and a visual example.
327    ///
328    /// # Example 1
329    ///
330    /// ```
331    /// use spade::{DelaunayTriangulation, HasPosition, Point2};
332    /// # use spade::Triangulation;
333    ///
334    /// struct PointWithData {
335    ///     position: Point2<f64>,
336    ///     height: f64,
337    ///     grad: [f64; 2],
338    /// }
339    ///
340    /// impl HasPosition for PointWithData {
341    ///     type Scalar = f64;
342    ///     fn position(&self) -> Point2<f64> { self.position }
343    /// }
344    ///
345    /// let mut triangulation: DelaunayTriangulation<PointWithData> = Default::default();
346    /// // Insert some points into the triangulation
347    /// triangulation.insert(PointWithData { position: Point2::new(10.0, 10.0), height: 0.0, grad: [-1., -1.] });
348    /// triangulation.insert(PointWithData { position: Point2::new(10.0, -10.0), height: 0.0, grad: [-1., 1.] });
349    /// triangulation.insert(PointWithData { position: Point2::new(-10.0, 10.0), height: 0.0, grad: [1., -1.] });
350    /// triangulation.insert(PointWithData { position: Point2::new(-10.0, -10.0), height: 0.0, grad: [1., 1.] });
351    ///
352    /// let nn = triangulation.natural_neighbor();
353    ///
354    /// // Interpolate point at coordinates (1.0, 2.0). This example uses a known gradients at each vertex.
355    /// // The gradient is the direction of steepest ascent for the function at that point
356    /// let query_point = Point2::new(1.0, 2.0);
357    /// let value: f64 = nn.interpolate_gradient(|v| v.data().height, |v| v.data().grad, 1.0, query_point).unwrap();
358    /// ```
359    ///
360    /// /// # Example 2
361    ///
362    /// ```
363    /// use spade::{DelaunayTriangulation, HasPosition, Point2};
364    /// # use spade::Triangulation;
365    ///
366    /// struct PointWithHeight {
367    ///     position: Point2<f64>,
368    ///     height: f64,
369    /// }
370    ///
371    /// impl HasPosition for PointWithHeight {
372    ///     type Scalar = f64;
373    ///     fn position(&self) -> Point2<f64> { self.position }
374    /// }
375    ///
376    /// let mut triangulation: DelaunayTriangulation<PointWithHeight> = Default::default();
377    /// // Insert some points into the triangulation
378    /// triangulation.insert(PointWithHeight { position: Point2::new(10.0, 10.0), height: 0.0 });
379    /// triangulation.insert(PointWithHeight { position: Point2::new(10.0, -10.0), height: 0.0 });
380    /// triangulation.insert(PointWithHeight { position: Point2::new(-10.0, 10.0), height: 0.0 });
381    /// triangulation.insert(PointWithHeight { position: Point2::new(-10.0, -10.0), height: 0.0 });
382    ///
383    /// let nn = triangulation.natural_neighbor();
384    ///
385    /// // Interpolate point at coordinates (1.0, 2.0).
386    /// let query_point = Point2::new(1.0, 2.0);
387    /// // Let's interpolate the value with estimated gradients. The gradient is the direction of steepest ascent for the function at that point
388    /// let value: f64 = nn.interpolate_gradient(|v| v.data().height, |v| nn.estimate_gradient(v, |h| h.data().height), 1.0, query_point).unwrap();
389    ///
390    /// // could also have gone and calculated all grads upfront like this, notice the s in gradients
391    /// let grads = nn.estimate_gradients(|v| v.data().height);
392    /// let value: f64 = nn.interpolate_gradient(|v| v.data().height, &grads, 1.0, query_point).unwrap();
393    /// ```
394    ///
395    /// # References
396    ///
397    /// This method uses the C1 extension proposed by Sibson in
398    /// "A brief description of natural neighbor interpolation, R. Sibson, 1981"
399    ///
400    /// It is also worth looking at Julia Flötotto's thesis
401    /// A Coordinate System associated to a Point Cloud issued from a Manifold: Definition, Properties and Applications
402    /// chapter 4.1 and 4.4
403    pub fn interpolate_gradient<I, G>(
404        &self,
405        i: I,
406        g: G,
407        flatness: <V as HasPosition>::Scalar,
408        position: Point2<<V as HasPosition>::Scalar>,
409    ) -> Option<<V as HasPosition>::Scalar>
410    where
411        I: Fn(VertexHandle<V, DE, UE, F>) -> <V as HasPosition>::Scalar,
412        G: Fn(VertexHandle<V, DE, UE, F>) -> [<V as HasPosition>::Scalar; 2],
413    {
414        let nns = &mut *self.weight_buffer.borrow_mut();
415        self.get_weights(position, nns);
416        if nns.is_empty() {
417            return None;
418        }
419
420        // Variable names should make more sense after looking into the paper!
421        // Roughly speaking, this approach works by blending a smooth c1 approximation into the
422        // regular natural neighbor interpolation ("c0 contribution").
423        // The c0 / c1 contributions are stored in sum_c0 / sum_c1 and are weighted by alpha and beta
424        // respectively.
425        let mut sum_c0 = zero();
426        let mut sum_c1 = zero();
427        let mut sum_c1_weights = zero();
428        let mut alpha: <V as HasPosition>::Scalar = zero();
429        let mut beta: <V as HasPosition>::Scalar = zero();
430
431        for (handle, weight) in nns {
432            let handle = self.triangulation.vertex(*handle);
433            let pos_i = handle.position();
434            let h_i = i(handle);
435            let diff = position.sub(pos_i);
436            let r_i2 = diff.length2();
437
438            if r_i2 == zero() {
439                return Some(h_i);
440            }
441
442            let r_i = r_i2.powf(flatness);
443            let c1_weight_i = *weight / r_i;
444            let grad_i = g(handle);
445            let zeta_i = h_i + diff.dot(grad_i.into());
446
447            alpha = alpha + c1_weight_i * r_i2;
448            beta = beta + *weight * r_i2;
449            sum_c1_weights = sum_c1_weights + c1_weight_i;
450            sum_c1 = sum_c1 + zeta_i * c1_weight_i;
451            sum_c0 = sum_c0 + h_i * *weight;
452        }
453        alpha = alpha / sum_c1_weights;
454        sum_c1 = sum_c1 / sum_c1_weights;
455        let result = (alpha * sum_c0 + beta * sum_c1) / (alpha + beta);
456
457        Some(result)
458    }
459
460    /// Calculates the natural neighbor weights corresponding to a given position.
461    ///
462    /// The weight of a natural neighbor n is defined as the size of the intersection of two areas:
463    ///  - The existing voronoi cell of the natural neighbor
464    ///  - The cell that would be created if a vertex was created at the given position.
465    ///
466    /// The area of this intersection can, surprisingly, be calculated without doing the actual insertion.
467    ///
468    /// Parameters:
469    ///  - position refers to the position for which the weights should be calculated
470    ///  - nns: A list of edges that connect the natural neighbors (e.g. `nns[0].to() == nns[1].from()`).
471    ///  - result: Stores the resulting NNs (as vertex handle) and their weights.
472    ///
473    /// # Visualization
474    ///
475    /// Refer to these two .svg files for more detail on how the algorithm works, either by looking them up
476    /// directly or by running `cargo doc --document-private-items` and looking at the documentation of this
477    /// function.
478    ///
479    /// ## Insertion cell
480    ///
481    /// This .svg displays the *insertion cell* (thick orange line) which is the voronoi cell that gets
482    /// created if the query point (red dot) would be inserted.
483    /// Each point of the insertion cell lies on a circumcenter of a triangle formed by `position` and two
484    /// adjacent natural neighbors (red dots with their index shown inside).
485    #[doc = include_str!("../../images/natural_neighbor_insertion_cell.svg")]
486    ///
487    /// ## Inner loop
488    ///
489    /// This .svg illustrates the inner loop of the algorithm (see code below). The goal is to calculate
490    /// the weight of natural neighbor with index 4 which is proportional to the area of the orange polygon.
491    /// `last_edge`, `stop_edge`, `first`, `c0`, `c1` `c2` and `last` refer to variable names (see code below).
492    #[doc = include_str!("../../images/natural_neighbor_polygon.svg")]
493    fn get_natural_neighbor_weights(
494        &self,
495        position: Point2<<V as HasPosition>::Scalar>,
496        nns: &[FixedDirectedEdgeHandle],
497        result: &mut Vec<(FixedVertexHandle, <V as HasPosition>::Scalar)>,
498    ) {
499        result.clear();
500
501        if nns.is_empty() {
502            return;
503        }
504
505        if nns.len() == 1 {
506            let edge = self.triangulation.directed_edge(nns[0]);
507            result.push((edge.from().fix(), one()));
508            return;
509        }
510
511        if nns.len() == 2 {
512            let [e0, e1] = [
513                self.triangulation.directed_edge(nns[0]),
514                self.triangulation.directed_edge(nns[1]),
515            ];
516            let [v0, v1] = [e0.from(), e1.from()];
517            let [w0, w1] =
518                two_point_interpolation::<DelaunayTriangulation<V, DE, UE, F>>(v0, v1, position);
519
520            result.push((v0.fix(), w0));
521            result.push((v1.fix(), w1));
522            return;
523        }
524
525        // Get insertion cell vertices. The "insertion cell" refers to the voronoi cell that would be
526        // created if "position" would be inserted.
527        // These insertion cells happen to lie on the circumcenter of `position` and any two adjacent
528        // natural neighbors (e.g. [position, nn[2].position(), nn[3].position()]).
529        //
530        // `images/natural_neighbor_insertion_cell.svg` depicts the cell as thick orange line.
531        let mut insertion_cell = self.insert_cell_buffer.borrow_mut();
532        insertion_cell.clear();
533        for cur_nn in nns {
534            let cur_nn = self.triangulation.directed_edge(*cur_nn);
535
536            let [from, to] = cur_nn.positions();
537            insertion_cell.push(
538                math::circumcenter([
539                    to.sub(position),
540                    from.sub(position),
541                    Point2::new(zero(), zero()),
542                ])
543                .0,
544            );
545        }
546
547        let mut total_area = zero(); // Used to normalize weights at the end
548
549        let mut last_edge = self.triangulation.directed_edge(*nns.last().unwrap());
550        let mut last = *insertion_cell.last().unwrap();
551
552        for (stop_edge, first) in core::iter::zip(nns.iter(), &*insertion_cell) {
553            // Main loop
554            //
555            // Refer to images/natural_neighbor_polygon.svg for some visual aid.
556            //
557            // The outer loops calculates the weight of an individual natural neighbor.
558            // To do this, it calculates the intersection area of the insertion cell with the cell of the
559            // current natural neighbor.
560            // This intersection is a convex polygon with vertices `first, current0, current1 ... last`
561            // (ordered ccw, `currentX` refers to the variable in the inner loop)
562            //
563            // The area of a convex polygon [v0 ... vN] is given by
564            // 0.5 * ((v0.x * v1.y + ... vN.x * v0.y) - (v0.y * v1.x + ... vN.y * v0.x))
565            //        ⮤       positive_area        ⮥   ⮤         negative_area      ⮥
566            //
567            // The positive and negative contributions are calculated separately to avoid precision issues.
568            // The factor of 0.5 can be omitted as the weights are normalized anyway.
569
570            // `stop_edge` is used to know when to stop the inner loop (= once the polygon is finished)
571            let stop_edge = self.triangulation.directed_edge(*stop_edge);
572            assert!(!stop_edge.is_outer_edge());
573
574            let mut positive_area = first.x * last.y;
575            let mut negative_area = first.y * last.x;
576
577            loop {
578                // All other polygon vertices happen lie on the circumcenter of a face adjacent to an
579                // out edge of the current natural neighbor.
580                //
581                // The natural_neighbor_polygon.svg refers to this variable as `c0`, `c1`, and `c2`.
582                let vertices = last_edge
583                    .face()
584                    .as_inner()
585                    .unwrap()
586                    .positions()
587                    .map(|p| p.sub(position));
588
589                let current = math::circumcenter(vertices).0;
590
591                positive_area = positive_area + last.x * current.y;
592                negative_area = negative_area + last.y * current.x;
593
594                last_edge = last_edge.next().rev();
595                last = current;
596
597                if last_edge == stop_edge.rev() {
598                    positive_area = positive_area + current.x * first.y;
599                    negative_area = negative_area + current.y * first.x;
600                    break;
601                }
602            }
603
604            let polygon_area = positive_area - negative_area;
605
606            total_area = total_area + polygon_area;
607            result.push((stop_edge.from().fix(), polygon_area));
608
609            last = *first;
610            last_edge = stop_edge;
611        }
612
613        for tuple in result {
614            tuple.1 = tuple.1 / total_area;
615        }
616    }
617
618    /// Estimates and returns the gradient for a every vertex in this triangulation.
619    ///
620    /// This assumes that the triangulation models some kind of height field.
621    /// The gradient is calculated from the weighted and normalized average of the normals of all triangles
622    /// adjacent to the given vertex (weighted by triangle size).
623    ///
624    /// NOTE: This is not Sibson's gradient estimator
625    pub fn estimate_gradients<I>(
626        &self,
627        i: I,
628    ) -> impl Fn(VertexHandle<V, DE, UE, F>) -> [<V as HasPosition>::Scalar; 2]
629    where
630        I: Fn(VertexHandle<V, DE, UE, F>) -> <V as HasPosition>::Scalar,
631    {
632        let grads = self
633            .triangulation
634            .vertices()
635            .map(|v| self.estimate_gradient(v, &i))
636            .collect::<Vec<_>>();
637
638        move |v: VertexHandle<V, DE, UE, F>| grads[v.index()]
639    }
640
641    /// Estimates and returns the gradient for a single vertex in this triangulation.
642    ///
643    /// This assumes that the triangulation models some kind of height field.
644    /// The gradient is calculated from the weighted and normalized average of the normals of all triangles
645    /// adjacent to the given vertex (weighted by triangle size).
646    ///
647    /// NOTE: This is not Sibson's gradient estimator
648    pub fn estimate_gradient<I>(
649        &self,
650        v: VertexHandle<V, DE, UE, F>,
651        i: I,
652    ) -> [<V as HasPosition>::Scalar; 2]
653    where
654        I: Fn(VertexHandle<V, DE, UE, F>) -> <V as HasPosition>::Scalar,
655    {
656        let v_2d = v.position();
657        let v_pos = [v_2d.x, v_2d.y, i(v)];
658
659        let neighbor_positions = {
660            v.out_edges()
661                .map(|e| {
662                    let pos = e.to().position();
663                    [pos.x, pos.y, i(e.to())]
664                })
665                .collect::<Vec<_>>()
666        };
667        let mut final_normal: [<V as HasPosition>::Scalar; 3] = [zero(); 3];
668        for index in 0..neighbor_positions.len() {
669            let p0 = neighbor_positions[index];
670            let p1 = neighbor_positions[(index + 1) % neighbor_positions.len()];
671
672            let d0 = [p0[0] - v_pos[0], p0[1] - v_pos[1], p0[2] - v_pos[2]];
673            let d1 = [p1[0] - v_pos[0], p1[1] - v_pos[1], p1[2] - v_pos[2]];
674
675            let normal = [
676                d0[1] * d1[2] - d0[2] * d1[1],
677                d0[2] * d1[0] - d0[0] * d1[2],
678                d0[0] * d1[1] - d0[1] * d1[0],
679            ];
680
681            // this should be true for every normal
682            // given that the height field assumption holds
683            if normal[2] > zero() {
684                final_normal = [
685                    final_normal[0] + normal[0],
686                    final_normal[1] + normal[1],
687                    final_normal[2] + normal[2],
688                ];
689            }
690        }
691
692        // Calculate gradient from normal
693        if final_normal[2] != zero() {
694            [
695                -final_normal[0] / final_normal[2],
696                -final_normal[1] / final_normal[2],
697            ]
698        } else {
699            [zero(), zero()]
700        }
701    }
702}
703
704fn get_natural_neighbor_edges<T>(
705    triangulation: &T,
706    inspect_buffer: &mut Vec<FixedDirectedEdgeHandle>,
707    position: Point2<<T::Vertex as HasPosition>::Scalar>,
708    result: &mut Vec<FixedDirectedEdgeHandle>,
709) where
710    T: Triangulation,
711    <T::Vertex as HasPosition>::Scalar: Float,
712{
713    inspect_buffer.clear();
714    result.clear();
715    match triangulation.locate(position) {
716        PositionInTriangulation::OnFace(face) => {
717            for edge in triangulation
718                .face(face)
719                .adjacent_edges()
720                .into_iter()
721                .rev()
722                .map(|e| e.rev())
723            {
724                inspect_flips(triangulation, result, inspect_buffer, edge.fix(), position);
725            }
726        }
727        PositionInTriangulation::OnEdge(edge) => {
728            let edge = triangulation.directed_edge(edge);
729
730            if edge.is_part_of_convex_hull() {
731                result.extend([edge.fix(), edge.fix().rev()]);
732                return;
733            }
734
735            for edge in [edge, edge.rev()] {
736                inspect_flips(triangulation, result, inspect_buffer, edge.fix(), position);
737            }
738        }
739        PositionInTriangulation::OnVertex(fixed_handle) => {
740            let vertex = triangulation.vertex(fixed_handle);
741            result.push(
742                vertex
743                    .out_edge()
744                    .map(|e| e.fix())
745                    .unwrap_or(FixedDirectedEdgeHandle::new(0)),
746            )
747        }
748        _ => {}
749    }
750    result.reverse();
751}
752
753/// Identifies natural neighbors.
754///
755/// To do so, this function "simulates" the insertion of a vertex (located at `position) and keeps track of which
756/// edges would need to be flipped. A vertex is a natural neighbor if it happens to be part of an edge that would
757/// require to be flipped.
758///
759/// Similar to function `legalize_edge` (which is used for *actual* insertions).
760fn inspect_flips<T>(
761    triangulation: &T,
762    result: &mut Vec<FixedDirectedEdgeHandle>,
763    buffer: &mut Vec<FixedDirectedEdgeHandle>,
764    edge_to_validate: FixedDirectedEdgeHandle,
765    position: Point2<<T::Vertex as HasPosition>::Scalar>,
766) where
767    T: Triangulation,
768{
769    buffer.clear();
770    buffer.push(edge_to_validate);
771
772    while let Some(edge) = buffer.pop() {
773        let edge = triangulation.directed_edge(edge);
774
775        let v2 = edge.opposite_vertex();
776        let v1 = edge.from();
777
778        let mut should_flip = false;
779
780        if let Some(v2) = v2 {
781            let v0 = edge.to().position();
782            let v1 = v1.position();
783            let v3 = position;
784            debug_assert!(math::is_ordered_ccw(v2.position(), v1, v0));
785            should_flip = math::contained_in_circumference(v2.position(), v1, v0, v3);
786
787            if should_flip {
788                let e1 = edge.next().fix().rev();
789                let e2 = edge.prev().fix().rev();
790
791                buffer.push(e1);
792                buffer.push(e2);
793            }
794        }
795
796        if !should_flip {
797            result.push(edge.fix().rev());
798        }
799    }
800}
801
802fn two_point_interpolation<'a, T>(
803    v0: VertexHandle<'a, T::Vertex, T::DirectedEdge, T::UndirectedEdge, T::Face>,
804    v1: VertexHandle<'a, T::Vertex, T::DirectedEdge, T::UndirectedEdge, T::Face>,
805    position: Point2<<T::Vertex as HasPosition>::Scalar>,
806) -> [<T::Vertex as HasPosition>::Scalar; 2]
807where
808    T: Triangulation,
809    <T::Vertex as HasPosition>::Scalar: Float,
810{
811    let projection = math::project_point(v0.position(), v1.position(), position);
812    let rel = projection.relative_position();
813    let one: <T::Vertex as HasPosition>::Scalar = 1.0.into();
814    [one - rel, rel]
815}
816
817#[cfg(test)]
818mod test {
819    use approx::assert_ulps_eq;
820
821    use crate::test_utilities::{random_points_in_range, random_points_with_seed, SEED, SEED2};
822    use crate::{DelaunayTriangulation, HasPosition, InsertionError, Point2, Triangulation};
823    use alloc::vec;
824    use alloc::vec::Vec;
825
826    struct PointWithHeight {
827        position: Point2<f64>,
828        height: f64,
829    }
830
831    impl PointWithHeight {
832        fn new(position: Point2<f64>, height: f64) -> Self {
833            Self { position, height }
834        }
835    }
836
837    impl HasPosition for PointWithHeight {
838        type Scalar = f64;
839
840        fn position(&self) -> Point2<Self::Scalar> {
841            self.position
842        }
843    }
844
845    #[test]
846    fn test_natural_neighbors() -> Result<(), InsertionError> {
847        let vertices = random_points_with_seed(50, SEED);
848        let t = DelaunayTriangulation::<_>::bulk_load(vertices)?;
849
850        let mut nn_edges = Vec::new();
851
852        let mut buffer = Vec::new();
853        let query_point = Point2::new(0.5, 0.2);
854        super::get_natural_neighbor_edges(&t, &mut buffer, query_point, &mut nn_edges);
855
856        assert!(nn_edges.len() >= 3);
857
858        for edge in &nn_edges {
859            let edge = t.directed_edge(*edge);
860            assert!(edge.side_query(query_point).is_on_left_side());
861        }
862
863        let mut nns = Vec::new();
864        let natural_neighbor = &t.natural_neighbor();
865        for v in t.vertices() {
866            natural_neighbor.get_weights(v.position(), &mut nns);
867            let expected = vec![(v.fix(), 1.0f64)];
868            assert_eq!(nns, expected);
869        }
870
871        Ok(())
872    }
873
874    #[test]
875    fn test_small_interpolation() -> Result<(), InsertionError> {
876        let mut t = DelaunayTriangulation::<_>::new();
877        // Create symmetric triangle - the weights should be mirrored
878        t.insert(PointWithHeight::new(Point2::new(0.0, 2.0f64.sqrt()), 0.0))?;
879        let v1 = t.insert(PointWithHeight::new(Point2::new(-1.0, -1.0), 0.0))?;
880        let v2 = t.insert(PointWithHeight::new(Point2::new(1.0, -1.0), 0.0))?;
881
882        let mut result = Vec::new();
883        t.natural_neighbor()
884            .get_weights(Point2::new(0.0, 0.0), &mut result);
885
886        assert_eq!(result.len(), 3);
887        let mut v1_weight = None;
888        let mut v2_weight = None;
889        for (v, weight) in result {
890            if v == v1 {
891                v1_weight = Some(weight);
892            }
893            if v == v2 {
894                v2_weight = Some(weight);
895            } else {
896                assert!(weight > 0.0);
897            }
898        }
899
900        assert!(v1_weight.is_some());
901        assert!(v2_weight.is_some());
902        assert_ulps_eq!(v1_weight.unwrap(), v2_weight.unwrap());
903
904        Ok(())
905    }
906
907    #[test]
908    fn test_quad_interpolation() -> Result<(), InsertionError> {
909        let mut t = DelaunayTriangulation::<_>::new();
910        // Insert points into the corners of the unit square. The weights at the origin should then
911        // all be 0.25
912        t.insert(Point2::new(1.0, 1.0))?;
913        t.insert(Point2::new(1.0, -1.0))?;
914        t.insert(Point2::new(-1.0, 1.0))?;
915        t.insert(Point2::new(-1.0, -1.0))?;
916
917        let mut result = Vec::new();
918        t.natural_neighbor()
919            .get_weights(Point2::new(0.0, 0.0), &mut result);
920
921        assert_eq!(result.len(), 4);
922        for (_, weight) in result {
923            assert_ulps_eq!(weight, 0.25);
924        }
925
926        Ok(())
927    }
928
929    #[test]
930    fn test_constant_height_interpolation() -> Result<(), InsertionError> {
931        let mut t = DelaunayTriangulation::<_>::new();
932        let vertices = random_points_with_seed(50, SEED);
933        let fixed_height = 1.5;
934        for v in vertices {
935            t.insert(PointWithHeight::new(v, fixed_height))?;
936        }
937
938        for v in random_points_in_range(1.5, 50, SEED2) {
939            let value = t.natural_neighbor().interpolate(|p| p.data().height, v);
940            if let Some(value) = value {
941                assert_ulps_eq!(value, 1.5);
942            }
943        }
944
945        Ok(())
946    }
947
948    #[test]
949    fn test_slope_interpolation() -> Result<(), InsertionError> {
950        // Insert vertices v with
951        // v.x in [-1.0 .. 1.0]
952        // and
953        // v.height = v.x .
954        //
955        // The expected side profile of the triangulation (YZ plane through y = 0.0) looks like this:
956        //
957        //
958        //  ↑           /⇖x = 1, height = 1
959        //  height     /
960        //            /
961        //  x→       / <---- x = -1, height = -1
962        let mut t = DelaunayTriangulation::<_>::new();
963        let grid_size = 32;
964        let scale = 1.0 / grid_size as f64;
965        for x in -grid_size..=grid_size {
966            for y in -grid_size..=grid_size {
967                let coords = Point2::new(x as f64, y as f64).mul(scale);
968                t.insert(PointWithHeight::new(coords, coords.x))?;
969            }
970        }
971        let query_points = random_points_with_seed(50, SEED);
972
973        let check = |point, expected| {
974            let value = t
975                .natural_neighbor()
976                .interpolate(|v| v.data().height, point)
977                .unwrap();
978            assert_ulps_eq!(value, expected, epsilon = 0.001);
979        };
980
981        // Check inside points
982        for point in query_points {
983            check(point, point.x);
984        }
985
986        // Check positions on the boundary
987        check(Point2::new(-1.0, 0.0), -1.0);
988        check(Point2::new(-1.0, 0.2), -1.0);
989        check(Point2::new(-1.0, 0.5), -1.0);
990        check(Point2::new(-1.0, -1.0), -1.0);
991        check(Point2::new(1.0, 0.0), 1.0);
992        check(Point2::new(1.0, 0.2), 1.0);
993        check(Point2::new(1.0, 0.5), 1.0);
994        check(Point2::new(1.0, -1.0), 1.0);
995
996        for x in [-1.0, -0.8, -0.5, 0.0, 0.3, 0.7, 1.0] {
997            check(Point2::new(x, 1.0), x);
998            check(Point2::new(x, -1.0), x);
999        }
1000
1001        let expect_none = |point| {
1002            assert!(t
1003                .natural_neighbor()
1004                .interpolate(|v| v.data().height, point)
1005                .is_none());
1006        };
1007
1008        // Check positions outside the triangulation.
1009        for x in [-5.0f64, -4.0, 3.0, 2.0] {
1010            expect_none(Point2::new(x, 0.5));
1011            expect_none(Point2::new(x, -0.5));
1012            expect_none(Point2::new(x, 0.1));
1013        }
1014
1015        Ok(())
1016    }
1017
1018    #[test]
1019    fn test_parabola_interpolation() -> Result<(), InsertionError> {
1020        // Insert vertices into a regular grid and assign a height of r*r (r = distance from origin).
1021        let mut t = DelaunayTriangulation::<_>::new();
1022        let grid_size = 32;
1023        let scale = 1.0 / grid_size as f64;
1024        for x in -grid_size..=grid_size {
1025            for y in -grid_size..=grid_size {
1026                let coords = Point2::new(x as f64, y as f64).mul(scale);
1027                t.insert(PointWithHeight::new(coords, coords.length2()))?;
1028            }
1029        }
1030        let query_points = random_points_with_seed(50, SEED);
1031
1032        for point in query_points {
1033            let value = t
1034                .natural_neighbor()
1035                .interpolate(|v| v.data().height, point)
1036                .unwrap();
1037            let r2 = point.length2();
1038            assert_ulps_eq!(value, r2, epsilon = 0.001);
1039        }
1040
1041        Ok(())
1042    }
1043
1044    #[test]
1045    fn test_nan_issue() -> Result<(), InsertionError> {
1046        let vertices = vec![
1047            Point2 {
1048                x: 2273118.4147972693,
1049                y: 6205168.575915335,
1050            },
1051            Point2 {
1052                x: 2273118.2119929367,
1053                y: 6205168.6854697745,
1054            },
1055            Point2 {
1056                x: 2273118.1835989403,
1057                y: 6205168.653506873,
1058            },
1059            Point2 {
1060                x: 2273118.2647345643,
1061                y: 6205169.082690307,
1062            },
1063            Point2 {
1064                x: 2273118.105938253,
1065                y: 6205168.163882839,
1066            },
1067            Point2 {
1068                x: 2273117.7264146665,
1069                y: 6205168.718028998,
1070            },
1071            Point2 {
1072                x: 2273118.0946678743,
1073                y: 6205169.148656867,
1074            },
1075            Point2 {
1076                x: 2273118.3779900977,
1077                y: 6205168.1644559605,
1078            },
1079            Point2 {
1080                x: 2273117.71105166,
1081                y: 6205168.459413756,
1082            },
1083            Point2 {
1084                x: 2273118.30732556,
1085                y: 6205169.130634655,
1086            },
1087        ];
1088
1089        let delaunay = DelaunayTriangulation::<_>::bulk_load(vertices)?;
1090
1091        let nns = delaunay.natural_neighbor();
1092        let result: f64 = nns
1093            .interpolate(|_| 1.0, Point2::new(2273118.2, 6205168.666666667))
1094            .unwrap();
1095
1096        assert!(!result.is_nan());
1097        Ok(())
1098    }
1099
1100    #[test]
1101    fn test_gradient_estimation_planar() -> Result<(), InsertionError> {
1102        let points = vec![
1103            PointWithHeight::new(Point2::new(0., 0.), 0.),
1104            PointWithHeight::new(Point2::new(1., 0.), 0.),
1105            PointWithHeight::new(Point2::new(1., 1.), 1.),
1106            PointWithHeight::new(Point2::new(0., 1.), 1.),
1107            PointWithHeight::new(Point2::new(0.5, 0.5), 0.5),
1108        ];
1109
1110        let t = DelaunayTriangulation::<_>::bulk_load_stable(points).unwrap();
1111        let nn = t.natural_neighbor();
1112
1113        let v = nn
1114            .triangulation
1115            .locate_vertex(Point2::new(0.5, 0.5))
1116            .unwrap();
1117
1118        let grad = nn.estimate_gradient(v, |v| v.data().height);
1119
1120        assert!(grad == [0., 1.0]);
1121
1122        Ok(())
1123    }
1124
1125    #[test]
1126    fn test_gradient_estimation_flat() -> Result<(), InsertionError> {
1127        let points = vec![
1128            PointWithHeight::new(Point2::new(0., 0.), 1.),
1129            PointWithHeight::new(Point2::new(1., 0.), 1.),
1130            PointWithHeight::new(Point2::new(1., 1.), 1.),
1131            PointWithHeight::new(Point2::new(0., 1.), 1.),
1132            PointWithHeight::new(Point2::new(0.5, 0.5), 0.),
1133        ];
1134
1135        let t = DelaunayTriangulation::<_>::bulk_load_stable(points).unwrap();
1136        let nn = t.natural_neighbor();
1137
1138        let v = nn
1139            .triangulation
1140            .locate_vertex(Point2::new(0.5, 0.5))
1141            .unwrap();
1142
1143        let grad = nn.estimate_gradient(v, |v| v.data().height);
1144
1145        assert!(grad == [0., 0.]);
1146
1147        Ok(())
1148    }
1149}