nd_triangulation/
lib.rs

1//! nd-triangulation provides an interface to the
2//! [dD-Triangulation](https://doc.cgal.org/latest/Triangulation/index.html)
3//! Library of [CGAL](https://cgal.org). It allows to create
4//! Triangulations of points in arbitrary dimensions and to traverse
5//! the whole triangulation or just its convex hull.
6//!
7//! The main entry point to the crate is the [`Triangulation`](struct.Triangulation.html)  struct.
8
9#[macro_use]
10extern crate cpp;
11
12#[cfg(not(feature = "docs-rs"))]
13cpp! {{
14    #include <iterator>
15    #include <CGAL/Epick_d.h>
16    #include <CGAL/Triangulation.h>
17
18    using DynDimension = CGAL::Dynamic_dimension_tag;
19    using K = CGAL::Epick_d<DynDimension>;
20
21    using Vertex = CGAL::Triangulation_vertex<K, size_t>;
22    using FullCell = CGAL::Triangulation_full_cell<K, size_t>;
23    using TDS = CGAL::Triangulation_data_structure<DynDimension, Vertex, FullCell>;
24    using Triangulation = CGAL::Triangulation<K, TDS>;
25
26}}
27
28#[cfg(not(feature = "docs-rs"))]
29cpp! {{
30    using Point = Triangulation::Point;
31    using Facet_iterator = Triangulation::Facet_iterator;
32    using Facet = Triangulation::Facet;
33
34    using Full_cell_handle = Triangulation::Full_cell_handle;
35    using Vertex_handle = Triangulation::Vertex_handle;
36    using Full_cells = std::vector<Full_cell_handle>;
37}}
38
39use std::fmt;
40
41mod cell;
42mod vertex;
43
44pub use cell::*;
45pub use vertex::*;
46
47#[non_exhaustive]
48#[derive(Debug, PartialEq, Eq)]
49/// Error Type for Triangulation Errors
50pub enum TriangulationError {
51    /// Returned if a vertex of the wrong dimension is added.
52    WrongDimension {
53        actual_dim: usize,
54        expected_dim: usize,
55    },
56}
57
58impl fmt::Display for TriangulationError {
59    fn fmt(&self, f: &mut fmt::Formatter) -> fmt::Result {
60        use TriangulationError::*;
61        match self {
62            WrongDimension {
63                actual_dim,
64                expected_dim,
65            } => write!(
66                f,
67                "The vertex could not be added because its dimension was {} instead of {} ",
68                actual_dim, expected_dim
69            ),
70        }
71    }
72}
73
74impl std::error::Error for TriangulationError {}
75
76/// Triangulation
77///
78/// Uses the dD triangulation package from CGAL internally.
79#[derive(Debug, PartialEq, Eq)]
80pub struct Triangulation {
81    /// Pointer to CGAL triangulation
82    ptr: *mut u8, //c++ type: Triangulation*
83    /// Dimension of the triangulation
84    dim: usize,
85    next_vertex_id: usize,
86    next_cell_id: usize,
87}
88
89impl Triangulation {
90    /// Create new triangulation for vertices of size/dimension `dim`
91    pub fn new(dim: usize) -> Triangulation {
92        let ptr = unsafe { Self::init_triangulation_ptr(dim) };
93        Triangulation {
94            ptr,
95            dim,
96            next_vertex_id: 0,
97            next_cell_id: 1,
98        }
99    }
100
101    unsafe fn init_triangulation_ptr(dim: usize) -> *mut u8 {
102        #[cfg(not(feature = "docs-rs"))]
103        return cpp!([dim as "size_t"] -> *mut u8 as "Triangulation*"{
104            return new Triangulation(dim);
105        });
106
107        #[cfg(feature = "docs-rs")]
108        std::ptr::null_mut()
109    }
110
111    /// Add vertex to the triangulation.
112    ///
113    /// The operation fails if `coords` has the wrong dimension.
114    pub fn add_vertex(&mut self, coords: &[f64]) -> Result<usize, TriangulationError> {
115        if coords.len() != self.dim {
116            return Err(TriangulationError::WrongDimension {
117                actual_dim: coords.len(),
118                expected_dim: self.dim,
119            });
120        }
121        let id = unsafe { self.add_vertex_unchecked(coords) };
122        Ok(id)
123    }
124
125    /// Add vertex to triangulation without veryfing the dimension
126    ///
127    /// # Safety
128    /// If the dimension of `coords` is too small undefined behavior might be triggered on the c++ side.
129    /// If the dimension of `coords` is too large only the first `dim` values will be considered.
130    pub unsafe fn add_vertex_unchecked(&mut self, coords: &[f64]) -> usize {
131        let tri = self.ptr;
132        let dim = self.dim;
133        let coords = coords.as_ptr();
134        let vertex_id = self.next_vertex_id;
135
136        #[cfg(not(feature = "docs-rs"))]
137        cpp!([tri as "Triangulation*", dim as "size_t", coords as "double*", vertex_id as "size_t"] {
138            auto p = Point(dim, &coords[0], &coords[dim]);
139            auto vertex = tri->insert(p);
140            auto& id = vertex->data();
141            id = vertex_id;
142        });
143        self.next_vertex_id += 1;
144        vertex_id
145    }
146
147    /// Returns a iterator over all convex hull cells/facets.
148    ///
149    /// This allocates a vector with cell handles internally and is
150    /// not implemented in the typical streaming fashion of rust iterators.
151    pub fn convex_hull_cells(&mut self) -> CellIter {
152        let cells = unsafe { self.gather_ch_cells() };
153        CellIter::new(self, cells)
154    }
155
156    #[rustfmt::skip]
157    unsafe fn gather_ch_cells(&mut self) -> *mut u8 {
158        let tri = self.ptr;
159        let cell_id = &mut self.next_cell_id;
160
161        #[cfg(not(feature = "docs-rs"))]
162	return cpp!([tri as "Triangulation*", cell_id as "size_t*"] -> *mut u8 as "Full_cells*" {
163	    auto infinite_full_cells = new Full_cells();
164	    std::back_insert_iterator<Full_cells> out(*infinite_full_cells);
165	    tri->incident_full_cells(tri->infinite_vertex(), out);
166	    for (auto& cell: *infinite_full_cells){
167		auto& id = cell->data();
168		if(id == 0){
169		    id = *cell_id;
170		    (*cell_id)++;
171		}
172	    }
173	    return infinite_full_cells;
174        });
175
176	#[cfg(feature = "docs-rs")]
177	std::ptr::null_mut()
178    }
179
180    /// Returns a iterator over all cells/facets of the triangulation.
181    ///
182    /// This allocates a vector with cell handles internally and is
183    /// not implemented in the typical streaming fashion of rust iterators.
184    pub fn cells(&mut self) -> CellIter {
185        let cells = unsafe { self.gather_all_cells() };
186        CellIter::new(self, cells)
187    }
188
189    #[rustfmt::skip]
190    unsafe fn gather_all_cells(&mut self) -> *mut u8 {
191        let tri = self.ptr;
192        let cell_id = &mut self.next_cell_id;
193
194        #[cfg(not(feature = "docs-rs"))]
195	return cpp!([tri as "Triangulation*", cell_id as "size_t*"] -> *mut u8 as "Full_cells*" {
196	    auto full_cells = new Full_cells();
197	    std::back_insert_iterator<Full_cells> out(*full_cells);
198	    for (auto cit = tri->full_cells_begin(); cit != tri->full_cells_end(); ++cit){
199		auto cell = cit;
200		auto& id = cell->data();
201		if(id == 0){
202		    id = *cell_id;
203		    (*cell_id)++;
204		}
205		full_cells->push_back(cell);
206	    }
207	    return full_cells;
208        });
209
210	#[cfg(feature = "docs-rs")]
211	std::ptr::null_mut()
212    }
213}
214
215impl Drop for Triangulation {
216    fn drop(&mut self) {
217        let ptr = self.ptr;
218        unsafe {
219            #[cfg(not(feature = "docs-rs"))]
220            cpp!([ptr as "Triangulation*"] {
221                delete ptr;
222            })
223        }
224    }
225}
226
227#[test]
228fn test_triangulation_can_be_created_and_dropped_safely() {
229    let tri = Triangulation::new(3);
230    assert_eq!(3, tri.dim);
231}
232
233#[test]
234fn test_vertices_have_to_be_of_right_dimension() {
235    let mut tri = Triangulation::new(3);
236    assert!(tri.add_vertex(&[1.0]).is_err());
237    assert!(tri.add_vertex(&[1.0, 2.0]).is_err());
238    assert!(tri.add_vertex(&[1.0, 2.0, 3.0]).is_ok());
239    assert!(tri.add_vertex(&[4.0, 5.0, 6.0]).is_ok());
240    assert_eq!(
241        tri.add_vertex(&[1.0, 2.0, 3.0, 4.0]),
242        Err(TriangulationError::WrongDimension {
243            actual_dim: 4,
244            expected_dim: 3
245        })
246    );
247}
248
249#[test]
250fn test_empty_triangulation_has_pseudo_cell() {
251    let mut tri = Triangulation::new(3);
252    let ch_cells = tri.convex_hull_cells();
253
254    assert_eq!(1, ch_cells.count());
255}
256
257#[test]
258fn test_convex_hull_has_right_size() {
259    let mut tri = Triangulation::new(2);
260
261    tri.add_vertex(&[1.0, 1.0]).unwrap();
262    tri.add_vertex(&[2.0, 1.0]).unwrap();
263    tri.add_vertex(&[1.5, 1.5]).unwrap();
264
265    let ch_cells = tri.convex_hull_cells();
266    assert_eq!(3, ch_cells.count());
267}
268
269#[test]
270fn test_convex_hull_has_right_cells() {
271    let mut tri = Triangulation::new(2);
272
273    let p1 = &[1.0, 1.0];
274    let p2 = &[2.0, 1.0];
275    let p3 = &[1.5, 1.5];
276
277    let id1 = tri.add_vertex(p1).unwrap();
278    let id2 = tri.add_vertex(p2).unwrap();
279    let id3 = tri.add_vertex(p3).unwrap();
280
281    let ch_cells = tri.convex_hull_cells();
282
283    for cell in ch_cells {
284        let mut all_vertices: Vec<_> = cell.vertices().collect();
285
286        all_vertices.dedup_by_key(|p| p.id());
287
288        assert_eq!(2, all_vertices.len());
289
290        let only_input_vertices = all_vertices
291            .iter()
292            .map(Vertex::id)
293            .all(|id| id == id1 || id == id2 || id == id3);
294        assert!(only_input_vertices);
295    }
296}
297
298#[test]
299fn test_triangulation_has_right_size() {
300    let mut tri = Triangulation::new(2);
301
302    tri.add_vertex(&[1.0, 1.0]).unwrap();
303    tri.add_vertex(&[2.0, 1.0]).unwrap();
304    tri.add_vertex(&[1.5, 1.5]).unwrap();
305
306    let cells = tri.cells();
307    assert_eq!(4, cells.count());
308}