sfcgal/
geometry.rs

1use std::{
2    ffi::{c_void, CString},
3    mem::MaybeUninit,
4    os::raw::c_char,
5    ptr::NonNull,
6};
7
8use num_traits::{FromPrimitive, ToPrimitive};
9use sfcgal_sys::{
10    initialize, sfcgal_alloc_handler_t, sfcgal_error_handler_t, sfcgal_free_handler_t,
11    sfcgal_geometry_alpha_shapes, sfcgal_geometry_approximate_medial_axis, sfcgal_geometry_area,
12    sfcgal_geometry_area_3d, sfcgal_geometry_as_text, sfcgal_geometry_as_text_decim,
13    sfcgal_geometry_clone, sfcgal_geometry_collection_add_geometry,
14    sfcgal_geometry_collection_create, sfcgal_geometry_collection_geometry_n,
15    sfcgal_geometry_collection_num_geometries, sfcgal_geometry_convexhull,
16    sfcgal_geometry_convexhull_3d, sfcgal_geometry_covers, sfcgal_geometry_covers_3d,
17    sfcgal_geometry_delete, sfcgal_geometry_difference, sfcgal_geometry_difference_3d,
18    sfcgal_geometry_distance, sfcgal_geometry_distance_3d, sfcgal_geometry_extrude,
19    sfcgal_geometry_extrude_polygon_straight_skeleton, sfcgal_geometry_extrude_straight_skeleton,
20    sfcgal_geometry_intersection, sfcgal_geometry_intersection_3d, sfcgal_geometry_intersects,
21    sfcgal_geometry_intersects_3d, sfcgal_geometry_is_3d, sfcgal_geometry_is_empty,
22    sfcgal_geometry_is_measured, sfcgal_geometry_is_planar, sfcgal_geometry_is_valid,
23    sfcgal_geometry_is_valid_detail, sfcgal_geometry_line_sub_string,
24    sfcgal_geometry_minkowski_sum, sfcgal_geometry_offset_polygon,
25    sfcgal_geometry_optimal_alpha_shapes, sfcgal_geometry_orientation,
26    sfcgal_geometry_straight_skeleton, sfcgal_geometry_straight_skeleton_distance_in_m,
27    sfcgal_geometry_t, sfcgal_geometry_tesselate, sfcgal_geometry_triangulate_2dz,
28    sfcgal_geometry_type_id, sfcgal_geometry_union, sfcgal_geometry_union_3d,
29    sfcgal_geometry_volume, sfcgal_io_read_wkt, sfcgal_multi_linestring_create,
30    sfcgal_multi_point_create, sfcgal_multi_polygon_create, sfcgal_prepared_geometry_t, srid_t,
31};
32use sfcgal_sys::{
33    /* sfcgal_solid_set_exterior_shell, */
34    sfcgal_approx_convex_partition_2, sfcgal_full_version, sfcgal_geometry_as_hexwkb,
35    sfcgal_geometry_as_obj, sfcgal_geometry_as_obj_file, sfcgal_geometry_as_vtk,
36    sfcgal_geometry_as_vtk_file, sfcgal_geometry_as_wkb, sfcgal_geometry_buffer3d,
37    sfcgal_geometry_force_lhr, sfcgal_geometry_force_rhr, sfcgal_geometry_force_valid,
38    sfcgal_geometry_has_validity_flag, sfcgal_geometry_make_solid, sfcgal_geometry_rotate,
39    sfcgal_geometry_rotate_2d, sfcgal_geometry_rotate_3d, sfcgal_geometry_rotate_3d_around_center,
40    sfcgal_geometry_rotate_x, sfcgal_geometry_rotate_y, sfcgal_geometry_rotate_z,
41    sfcgal_geometry_round, sfcgal_geometry_scale, sfcgal_geometry_scale_3d,
42    sfcgal_geometry_scale_3d_around_center, sfcgal_geometry_straight_skeleton_partition,
43    sfcgal_geometry_translate_2d, sfcgal_geometry_translate_3d, sfcgal_geometry_visibility_point,
44    sfcgal_geometry_visibility_segment, sfcgal_greene_approx_convex_partition_2, sfcgal_init,
45    sfcgal_io_read_binary_prepared, sfcgal_io_read_ewkt, sfcgal_io_read_wkb,
46    sfcgal_io_write_binary_prepared, sfcgal_linestring_add_point, sfcgal_linestring_create,
47    sfcgal_linestring_num_points, sfcgal_linestring_point_n, sfcgal_multi_solid_create,
48    sfcgal_optimal_convex_partition_2, sfcgal_point_create, sfcgal_point_create_from_xy,
49    sfcgal_point_create_from_xym, sfcgal_point_create_from_xyz, sfcgal_point_create_from_xyzm,
50    sfcgal_point_m, sfcgal_point_x, sfcgal_point_y, sfcgal_point_z,
51    sfcgal_polygon_add_interior_ring, sfcgal_polygon_create,
52    sfcgal_polygon_create_from_exterior_ring, sfcgal_polygon_exterior_ring,
53    sfcgal_polygon_interior_ring_n, sfcgal_polygon_num_interior_rings,
54    sfcgal_polyhedral_surface_add_polygon, sfcgal_polyhedral_surface_create,
55    sfcgal_polyhedral_surface_num_polygons, sfcgal_polyhedral_surface_polygon_n,
56    sfcgal_prepared_geometry_as_ewkt, sfcgal_prepared_geometry_create,
57    sfcgal_prepared_geometry_create_from_geometry, sfcgal_prepared_geometry_delete,
58    sfcgal_prepared_geometry_geometry, sfcgal_prepared_geometry_set_geometry,
59    sfcgal_prepared_geometry_set_srid, sfcgal_prepared_geometry_srid, sfcgal_set_alloc_handlers,
60    sfcgal_set_error_handlers, sfcgal_set_geometry_validation, sfcgal_solid_add_interior_shell,
61    sfcgal_solid_create, sfcgal_solid_create_from_exterior_shell, sfcgal_solid_num_shells,
62    sfcgal_solid_shell_n, sfcgal_triangle_create, sfcgal_triangle_create_from_points,
63    sfcgal_triangle_set_vertex, sfcgal_triangle_set_vertex_from_xy,
64    sfcgal_triangle_set_vertex_from_xyz, sfcgal_triangle_vertex,
65    sfcgal_triangulated_surface_add_triangle, sfcgal_triangulated_surface_create,
66    sfcgal_triangulated_surface_num_triangles, sfcgal_triangulated_surface_triangle_n,
67    sfcgal_version, sfcgal_y_monotone_partition_2,
68};
69
70use crate::{
71    conversion::{CoordSeq, CoordType, ToSFCGALGeom},
72    errors::get_last_error,
73    utils::{
74        _c_string_with_size, _string, check_computed_value, check_nan_value,
75        check_null_prepared_geom, check_predicate, get_raw_bytes,
76    },
77    Result, ToSFCGAL,
78};
79
80/// Represent the buffer types usable with the `buffer3d` method.
81#[repr(C)]
82#[derive(PartialEq, Eq, PartialOrd, Ord, Debug, Primitive)]
83pub enum BufferType {
84    Round = 0,
85    CylSphere = 1,
86    Flat = 2,
87}
88
89/// SFCGAL Geometry types.
90///
91/// Indicates the type of shape represented by a `SFCGeometry`.
92#[repr(C)]
93#[derive(PartialEq, Eq, PartialOrd, Ord, Debug, Primitive)]
94pub enum GeomType {
95    Point = 1,
96    Linestring = 2,
97    Polygon = 3,
98    Multipoint = 4,
99    Multilinestring = 5,
100    Multipolygon = 6,
101    Geometrycollection = 7,
102    Polyhedralsurface = 15,
103    Triangulatedsurface = 16,
104    Triangle = 17,
105    Solid = 101,
106    Multisolid = 102,
107}
108
109impl GeomType {
110    fn is_collection_type(&self) -> bool {
111        matches!(
112            &self,
113            GeomType::Multipoint
114                | GeomType::Multilinestring
115                | GeomType::Multipolygon
116                | GeomType::Multisolid
117                | GeomType::Geometrycollection
118        )
119    }
120}
121
122/// Represents the orientation of a `SFCGeometry`.
123#[derive(PartialEq, Eq, Debug, Primitive)]
124pub enum Orientation {
125    CounterClockWise = -1isize,
126    ClockWise = 1isize,
127    Undetermined = 0isize,
128}
129
130macro_rules! precondition_match_type {
131    ($value : expr, $input_type: expr) => {{
132        if $value._type()? != $input_type {
133            bail!("Wrong input Geometry type")
134        }
135    }};
136}
137
138macro_rules! precondition_match_type_other {
139    ($value : expr, $input_type: expr) => {{
140        if $value._type()? != $input_type {
141            bail!("Wrong input Geometry type for other")
142        }
143    }};
144}
145
146macro_rules! precondition_match_validity {
147    ($value : expr) => {{
148        if !$value.is_valid()? {
149            bail!("Input geometry is not valid")
150        }
151    }};
152}
153
154macro_rules! precondition_match_validity_other {
155    ($value : expr) => {{
156        if !$value.is_valid()? {
157            bail!("Input geometry is not valid for other")
158        }
159    }};
160}
161
162macro_rules! precondition_match_not_empty {
163    ($value : expr) => {{
164        if $value.is_empty()? {
165            bail!("Input geometry is empty")
166        }
167    }};
168}
169
170macro_rules! precondition_index_in_range {
171    ($index: expr, $range: expr) => {{
172        if !$range.contains(&$index) {
173            bail!("Index not in the expected range")
174        }
175    }};
176}
177
178macro_rules! precondition_index_in_result_value {
179    ($index: expr, $range_result: expr) => {{
180        match $range_result {
181            Ok(max_value) => {
182                if !(0..max_value).contains(&$index) {
183                    bail!("Index not in the expected range")
184                }
185            }
186            Err(err) => {
187                bail!(err)
188            }
189        }
190    }};
191}
192
193/// Object representing a SFCGAL Geometry.
194#[repr(C)]
195pub struct SFCGeometry {
196    pub(crate) c_geom: NonNull<sfcgal_geometry_t>,
197    pub(crate) owned: bool,
198}
199
200impl Drop for SFCGeometry {
201    fn drop(&mut self) {
202        if self.owned {
203            unsafe { sfcgal_geometry_delete(self.c_geom.as_mut()) }
204        }
205    }
206}
207
208impl Clone for SFCGeometry {
209    fn clone(&self) -> SFCGeometry {
210        SFCGeometry {
211            c_geom: NonNull::new(unsafe { sfcgal_geometry_clone(self.c_geom.as_ref()) }).unwrap(),
212            owned: true,
213        }
214    }
215}
216
217impl std::fmt::Debug for SFCGeometry {
218    fn fmt(&self, f: &mut std::fmt::Formatter) -> std::fmt::Result {
219        write!(f, "{}", self.to_wkt_decim(8).unwrap())
220    }
221}
222
223impl SFCGeometry {
224    // TODO: reactivate when this binding is added in sfcgal-sys
225    /*pub fn solid_set_exterior_shell(&self, shell: &SFCGeometry) -> Result<()> {
226        precondition_match_type!(self, GeomType::Solid);
227        precondition_match_type_other!(shell, GeomType::Polyhedralsurface);
228
229        unsafe { sfcgal_solid_set_exterior_shell(self.c_geom.as_ptr(), shell.c_geom.as_ptr()) };
230
231        Ok(())
232    }*/
233
234    /// Create a geometry by parsing a [WKT](https://en.wikipedia.org/wiki/Well-known_text) string.
235    pub fn new(wkt: &str) -> Result<SFCGeometry> {
236        initialize();
237
238        let c_str = CString::new(wkt)?;
239
240        let obj = unsafe { sfcgal_io_read_wkt(c_str.as_ptr(), wkt.len()) };
241
242        unsafe { SFCGeometry::new_from_raw(obj, true) }
243    }
244
245    pub(crate) fn _init() {
246        // Right now it is a no-op
247        // Don't bother
248        unsafe {
249            sfcgal_init();
250        }
251    }
252
253    pub(crate) unsafe fn new_from_raw(
254        g: *mut sfcgal_geometry_t,
255        owned: bool,
256    ) -> Result<SFCGeometry> {
257        Ok(SFCGeometry {
258            owned,
259            c_geom: NonNull::new(g).ok_or_else(|| {
260                format_err!(
261                    "Obtained null pointer when creating geometry: {}",
262                    get_last_error()
263                )
264            })?,
265        })
266    }
267
268    pub fn new_from_coordinates<T>(coords: &CoordSeq<T>) -> Result<SFCGeometry>
269    where
270        T: ToSFCGALGeom + CoordType,
271    {
272        coords.to_sfcgal()
273    }
274
275    /// Returns a WKT representation of the given `SFCGeometry` using CGAL
276    /// exact integer fractions as coordinate values.
277    pub fn to_wkt(&self) -> Result<String> {
278        let mut ptr = MaybeUninit::<*mut c_char>::uninit();
279
280        let mut length: usize = 0;
281
282        unsafe {
283            sfcgal_geometry_as_text(self.c_geom.as_ref(), ptr.as_mut_ptr(), &mut length);
284
285            Ok(_c_string_with_size(ptr.assume_init(), length))
286        }
287    }
288
289    /// Returns a WKT representation of the given `SFCGeometry` using
290    /// floating point coordinate values with the desired number of
291    /// decimals.
292    pub fn to_wkt_decim(&self, nb_decim: i32) -> Result<String> {
293        let mut ptr = MaybeUninit::<*mut c_char>::uninit();
294
295        let mut length: usize = 0;
296
297        unsafe {
298            sfcgal_geometry_as_text_decim(
299                self.c_geom.as_ref(),
300                nb_decim,
301                ptr.as_mut_ptr(),
302                &mut length,
303            );
304
305            Ok(_c_string_with_size(ptr.assume_init(), length))
306        }
307    }
308
309    /// Creates a Wkb string of the given geometry. In memory version.
310    pub fn to_wkb_in_memory(&self) -> Result<Vec<u8>> {
311        let mut ptr = MaybeUninit::<*mut c_char>::uninit();
312
313        let mut length: usize = 0;
314
315        unsafe {
316            sfcgal_geometry_as_wkb(self.c_geom.as_ref(), ptr.as_mut_ptr(), &mut length);
317
318            Ok(get_raw_bytes(ptr.assume_init(), length))
319        }
320    }
321
322    /// # Safety
323    /// Returns an EWKT representation of the given PreparedGeometry
324    pub unsafe fn to_ewkt_in_memory(
325        prepared: *const sfcgal_prepared_geometry_t,
326        num_decimals: i32,
327    ) -> Result<Vec<u8>> {
328        let mut ptr = MaybeUninit::<*mut c_char>::uninit();
329
330        let mut length: usize = 0;
331
332        unsafe {
333            sfcgal_prepared_geometry_as_ewkt(prepared, num_decimals, ptr.as_mut_ptr(), &mut length);
334
335            Ok(get_raw_bytes(ptr.assume_init(), length))
336        }
337    }
338
339    /// Creates a Hexwkb string of the given geometry. In memory version.
340    pub fn to_hexwkb_in_memory(&self) -> Result<Vec<u8>> {
341        let mut ptr = MaybeUninit::<*mut c_char>::uninit();
342
343        let mut length: usize = 0;
344
345        unsafe {
346            sfcgal_geometry_as_hexwkb(self.c_geom.as_ref(), ptr.as_mut_ptr(), &mut length);
347
348            Ok(get_raw_bytes(ptr.assume_init(), length))
349        }
350    }
351
352    /// Read the binary prepared geometry
353    pub fn io_read_binary_prepared(data: &[u8]) -> Result<*mut sfcgal_prepared_geometry_t> {
354        unsafe {
355            let ptr = data.as_ptr() as *const i8;
356
357            let length = data.len();
358
359            let result = sfcgal_io_read_binary_prepared(ptr, length);
360
361            check_null_prepared_geom(result)?;
362
363            Ok(result)
364        }
365    }
366
367    /// Read the binary prepared geometry
368    pub fn io_read_wkb(data: &[u8]) -> Result<*mut sfcgal_prepared_geometry_t> {
369        unsafe {
370            let ptr = data.as_ptr() as *const i8;
371
372            let length = data.len();
373
374            let result = sfcgal_io_read_wkb(ptr, length);
375
376            check_null_prepared_geom(result)?;
377
378            Ok(result)
379        }
380    }
381
382    /// Read the binary prepared geometry
383    pub fn io_read_ewkt(data: &[u8]) -> Result<*mut sfcgal_prepared_geometry_t> {
384        unsafe {
385            let ptr = data.as_ptr() as *const i8;
386
387            let length = data.len();
388
389            let result = sfcgal_io_read_ewkt(ptr, length);
390
391            check_null_prepared_geom(result)?;
392
393            Ok(result)
394        }
395    }
396
397    /// # Safety
398    /// Read the binary prepared geometry in memory
399    pub unsafe fn io_write_binary_prepared(
400        prepared_geometry: *mut sfcgal_prepared_geometry_t,
401    ) -> Result<Vec<u8>> {
402        unsafe {
403            let mut ptr = MaybeUninit::<*mut c_char>::uninit();
404
405            let mut length: usize = 0;
406
407            sfcgal_io_write_binary_prepared(prepared_geometry, ptr.as_mut_ptr(), &mut length);
408
409            Ok(get_raw_bytes(ptr.assume_init(), length))
410        }
411    }
412
413    /// Creates a OBJ file of the given geometry
414    pub fn to_obj_file(&self, filename: &str) -> Result<()> {
415        unsafe {
416            let c_string = CString::new(filename)?;
417
418            let raw: *mut c_char = c_string.into_raw();
419
420            sfcgal_geometry_as_obj_file(self.c_geom.as_ptr(), raw);
421        };
422
423        Ok(())
424    }
425
426    /// Creates a OBJ string of the given geometry. In memory version.
427    pub fn to_obj_in_memory(&self) -> Result<Vec<u8>> {
428        let mut ptr = MaybeUninit::<*mut c_char>::uninit();
429
430        let mut length: usize = 0;
431
432        unsafe {
433            sfcgal_geometry_as_obj(self.c_geom.as_ref(), ptr.as_mut_ptr(), &mut length);
434
435            Ok(get_raw_bytes(ptr.assume_init(), length))
436        }
437    }
438
439    /// Creates a VTK file of the given geometry
440    pub fn to_vtk_file(&self, filename: &str) -> Result<()> {
441        unsafe {
442            let c_string = CString::new(filename)?;
443
444            let raw: *mut c_char = c_string.into_raw();
445
446            sfcgal_geometry_as_vtk_file(self.c_geom.as_ptr(), raw);
447        };
448
449        Ok(())
450    }
451
452    /// Creates a VTK string of the given geometry. In memory version.
453    pub fn to_vtk_in_memory(&self) -> Result<Vec<u8>> {
454        let mut ptr = MaybeUninit::<*mut c_char>::uninit();
455
456        let mut length: usize = 0;
457
458        unsafe {
459            sfcgal_geometry_as_vtk(self.c_geom.as_ref(), ptr.as_mut_ptr(), &mut length);
460
461            Ok(get_raw_bytes(ptr.assume_init(), length))
462        }
463    }
464
465    /// Test if the given `SFCGeometry` is empty or not.
466    pub fn is_empty(&self) -> Result<bool> {
467        let rv = unsafe { sfcgal_geometry_is_empty(self.c_geom.as_ptr()) };
468
469        check_predicate(rv)
470    }
471
472    /// Test if the given `SFCGeometry` is valid or not.
473    pub fn is_valid(&self) -> Result<bool> {
474        let rv = unsafe { sfcgal_geometry_is_valid(self.c_geom.as_ptr()) };
475
476        check_predicate(rv)
477    }
478
479    /// Test if the given `SFCGeometry` is measured (has an 'm' coordinates)
480    pub fn is_measured(&self) -> Result<bool> {
481        let rv = unsafe { sfcgal_geometry_is_measured(self.c_geom.as_ptr()) };
482
483        check_predicate(rv)
484    }
485
486    /// Test if the given `SFCGeometry` is planar or not.
487    pub fn is_planar(&self) -> Result<bool> {
488        precondition_match_validity!(self);
489
490        let rv = unsafe { sfcgal_geometry_is_planar(self.c_geom.as_ptr()) };
491
492        check_predicate(rv)
493    }
494
495    /// Test if the given `SFCGeometry` is a 3d geometry or not.
496    pub fn is_3d(&self) -> Result<bool> {
497        let rv = unsafe { sfcgal_geometry_is_3d(self.c_geom.as_ptr()) };
498
499        check_predicate(rv)
500    }
501
502    /// Returns reason for the invalidity or None in case of validity.
503    pub fn validity_detail(&self) -> Result<Option<String>> {
504        let mut ptr = MaybeUninit::<*mut c_char>::uninit();
505
506        unsafe {
507            let rv = sfcgal_geometry_is_valid_detail(
508                self.c_geom.as_ptr(),
509                ptr.as_mut_ptr(),
510                std::ptr::null::<sfcgal_geometry_t>() as *mut *mut sfcgal_geometry_t,
511            );
512
513            match rv {
514                1 => Ok(None),
515                0 => Ok(Some(_string(ptr.assume_init()))),
516                _ => Err(format_err!("SFCGAL error: {}", get_last_error())),
517            }
518        }
519    }
520
521    /// Returns the SFCGAL type of the given `SFCGeometry`.
522    pub fn _type(&self) -> Result<GeomType> {
523        let type_geom = unsafe { sfcgal_geometry_type_id(self.c_geom.as_ptr()) };
524
525        GeomType::from_u32(type_geom)
526            .ok_or_else(|| format_err!("Unknown geometry type (val={})", type_geom))
527    }
528
529    /// Computes the distance to an other `SFCGeometry`.
530    pub fn distance(&self, other: &SFCGeometry) -> Result<f64> {
531        precondition_match_validity!(self);
532
533        precondition_match_validity_other!(other);
534
535        let distance =
536            unsafe { sfcgal_geometry_distance(self.c_geom.as_ptr(), other.c_geom.as_ptr()) };
537
538        check_computed_value(distance)
539    }
540
541    /// Computes the 3d distance to an other `SFCGeometry`.
542    pub fn distance_3d(&self, other: &SFCGeometry) -> Result<f64> {
543        precondition_match_validity!(self);
544
545        precondition_match_validity_other!(other);
546
547        let distance =
548            unsafe { sfcgal_geometry_distance_3d(self.c_geom.as_ptr(), other.c_geom.as_ptr()) };
549
550        check_computed_value(distance)
551    }
552
553    /// Computes the area of the given `SFCGeometry`.
554    pub fn area(&self) -> Result<f64> {
555        precondition_match_validity!(self);
556
557        let area = unsafe { sfcgal_geometry_area(self.c_geom.as_ptr()) };
558
559        check_computed_value(area)
560    }
561
562    /// Computes the 3d area of the given `SFCGeometry`.
563    pub fn area_3d(&self) -> Result<f64> {
564        precondition_match_validity!(self);
565
566        let area = unsafe { sfcgal_geometry_area_3d(self.c_geom.as_ptr()) };
567
568        check_computed_value(area)
569    }
570
571    /// Computes the volume of the given `SFCGeometry` (must be a volume).
572    pub fn volume(&self) -> Result<f64> {
573        precondition_match_validity!(self);
574
575        let volume = unsafe { sfcgal_geometry_volume(self.c_geom.as_ptr()) };
576
577        check_computed_value(volume)
578    }
579
580    /// Computes the orientation of the given `SFCGeometry` (must be a
581    /// Polygon)
582    pub fn orientation(&self) -> Result<Orientation> {
583        precondition_match_type!(self, GeomType::Polygon);
584
585        precondition_match_validity!(self);
586
587        let orientation = unsafe { sfcgal_geometry_orientation(self.c_geom.as_ptr()) };
588
589        Orientation::from_i32(orientation)
590            .ok_or_else(|| format_err!("Error while retrieving orientation (val={})", orientation))
591    }
592
593    /// Test the intersection with an other `SFCGeometry`.
594    pub fn intersects(&self, other: &SFCGeometry) -> Result<bool> {
595        precondition_match_validity!(self);
596
597        precondition_match_validity_other!(other);
598
599        let rv = unsafe { sfcgal_geometry_intersects(self.c_geom.as_ptr(), other.c_geom.as_ptr()) };
600
601        check_predicate(rv)
602    }
603
604    /// Test the 3d intersection with an other `SFCGeometry`.
605    pub fn intersects_3d(&self, other: &SFCGeometry) -> Result<bool> {
606        precondition_match_validity!(self);
607
608        precondition_match_validity_other!(other);
609
610        let rv =
611            unsafe { sfcgal_geometry_intersects_3d(self.c_geom.as_ptr(), other.c_geom.as_ptr()) };
612
613        check_predicate(rv)
614    }
615
616    /// Returns the intersection of the given `SFCGeometry` to an other
617    /// `SFCGeometry`.
618    pub fn intersection(&self, other: &SFCGeometry) -> Result<SFCGeometry> {
619        precondition_match_validity!(self);
620
621        precondition_match_validity_other!(other);
622
623        let result =
624            unsafe { sfcgal_geometry_intersection(self.c_geom.as_ptr(), other.c_geom.as_ptr()) };
625
626        unsafe { SFCGeometry::new_from_raw(result, true) }
627    }
628
629    /// Returns the 3d intersection of the given `SFCGeometry` to an other
630    /// `SFCGeometry`.
631    pub fn intersection_3d(&self, other: &SFCGeometry) -> Result<SFCGeometry> {
632        precondition_match_validity!(self);
633
634        precondition_match_validity_other!(other);
635
636        let result =
637            unsafe { sfcgal_geometry_intersection_3d(self.c_geom.as_ptr(), other.c_geom.as_ptr()) };
638
639        unsafe { SFCGeometry::new_from_raw(result, true) }
640    }
641
642    /// Test if the geometry covers an other `SFCGeometry`.
643    pub fn covers(&self, other: &SFCGeometry) -> Result<bool> {
644        precondition_match_validity!(self);
645
646        precondition_match_validity_other!(other);
647
648        let rv = unsafe { sfcgal_geometry_covers(self.c_geom.as_ptr(), other.c_geom.as_ptr()) };
649
650        check_predicate(rv)
651    }
652
653    /// Test if the geometry covers an other `SFCGeometry` in 3D.
654    pub fn covers_3d(&self, other: &SFCGeometry) -> Result<bool> {
655        precondition_match_validity!(self);
656
657        precondition_match_validity_other!(other);
658
659        let rv = unsafe { sfcgal_geometry_covers_3d(self.c_geom.as_ptr(), other.c_geom.as_ptr()) };
660
661        check_predicate(rv)
662    }
663
664    /// Returns the difference of the given `SFCGeometry` to an other
665    /// `SFCGeometry`.
666    pub fn difference(&self, other: &SFCGeometry) -> Result<SFCGeometry> {
667        precondition_match_validity!(self);
668
669        precondition_match_validity_other!(other);
670
671        let result =
672            unsafe { sfcgal_geometry_difference(self.c_geom.as_ptr(), other.c_geom.as_ptr()) };
673
674        unsafe { SFCGeometry::new_from_raw(result, true) }
675    }
676
677    /// Returns the 3d difference of the given `SFCGeometry` to an other
678    /// `SFCGeometry`.
679    pub fn difference_3d(&self, other: &SFCGeometry) -> Result<SFCGeometry> {
680        precondition_match_validity!(self);
681
682        precondition_match_validity_other!(other);
683
684        let result =
685            unsafe { sfcgal_geometry_difference_3d(self.c_geom.as_ptr(), other.c_geom.as_ptr()) };
686
687        unsafe { SFCGeometry::new_from_raw(result, true) }
688    }
689
690    /// Returns the union of the given `SFCGeometry` to an other
691    /// `SFCGeometry`.
692    pub fn union(&self, other: &SFCGeometry) -> Result<SFCGeometry> {
693        precondition_match_validity!(self);
694
695        precondition_match_validity_other!(other);
696
697        let result = unsafe { sfcgal_geometry_union(self.c_geom.as_ptr(), other.c_geom.as_ptr()) };
698
699        unsafe { SFCGeometry::new_from_raw(result, true) }
700    }
701
702    /// Returns the 3d union of the given `SFCGeometry` to an other
703    /// `SFCGeometry`.
704    pub fn union_3d(&self, other: &SFCGeometry) -> Result<SFCGeometry> {
705        precondition_match_validity!(self);
706
707        precondition_match_validity_other!(other);
708
709        let result =
710            unsafe { sfcgal_geometry_union_3d(self.c_geom.as_ptr(), other.c_geom.as_ptr()) };
711
712        unsafe { SFCGeometry::new_from_raw(result, true) }
713    }
714
715    /// Returns the minkowski sum of the given `SFCGeometry` and an other
716    /// `SFCGEOMETRY`.
717    pub fn minkowski_sum(&self, other: &SFCGeometry) -> Result<SFCGeometry> {
718        precondition_match_validity!(self);
719
720        precondition_match_validity_other!(other);
721
722        let result =
723            unsafe { sfcgal_geometry_minkowski_sum(self.c_geom.as_ptr(), other.c_geom.as_ptr()) };
724
725        unsafe { SFCGeometry::new_from_raw(result, true) }
726    }
727
728    /// Returns the straight skeleton of the given `SFCGeometry`.
729    pub fn straight_skeleton(&self) -> Result<SFCGeometry> {
730        precondition_match_validity!(self);
731
732        let result = unsafe { sfcgal_geometry_straight_skeleton(self.c_geom.as_ptr()) };
733
734        unsafe { SFCGeometry::new_from_raw(result, true) }
735    }
736
737    /// Returns the straight skeleton of the given `SFCGeometry` with the
738    /// distance to the border as M coordinate.
739    pub fn straight_skeleton_distance_in_m(&self) -> Result<SFCGeometry> {
740        precondition_match_validity!(self);
741
742        let result =
743            unsafe { sfcgal_geometry_straight_skeleton_distance_in_m(self.c_geom.as_ptr()) };
744
745        unsafe { SFCGeometry::new_from_raw(result, true) }
746    }
747
748    /// Returns the extrude straight skeleton of the given Polygon.
749    pub fn extrude_straight_skeleton(&self, height: f64) -> Result<SFCGeometry> {
750        precondition_match_type!(self, GeomType::Polygon);
751
752        precondition_match_validity!(self);
753
754        if height == 0. {
755            bail!("Height cannot be 0.");
756        }
757
758        let result =
759            unsafe { sfcgal_geometry_extrude_straight_skeleton(self.c_geom.as_ptr(), height) };
760
761        unsafe { SFCGeometry::new_from_raw(result, true) }
762    }
763
764    /// Returns the union of the polygon z-extrusion (with respect to
765    /// building_height) and the extrude straight skeleton (with
766    /// respect to roof_height) of the given Polygon.
767    pub fn extrude_polygon_straight_skeleton(
768        &self,
769        building_height: f64,
770        roof_height: f64,
771    ) -> Result<SFCGeometry> {
772        precondition_match_type!(self, GeomType::Polygon);
773
774        precondition_match_validity!(self);
775
776        let result = unsafe {
777            sfcgal_geometry_extrude_polygon_straight_skeleton(
778                self.c_geom.as_ptr(),
779                building_height,
780                roof_height,
781            )
782        };
783
784        unsafe { SFCGeometry::new_from_raw(result, true) }
785    }
786
787    /// Returns the approximate medial axis for the given `SFCGeometry`
788    /// Polygon.
789    pub fn approximate_medial_axis(&self) -> Result<SFCGeometry> {
790        precondition_match_validity!(self);
791
792        let result = unsafe { sfcgal_geometry_approximate_medial_axis(self.c_geom.as_ptr()) };
793
794        unsafe { SFCGeometry::new_from_raw(result, true) }
795    }
796
797    /// Returns the offset polygon of the given `SFCGeometry`.
798    pub fn offset_polygon(&self, radius: f64) -> Result<SFCGeometry> {
799        precondition_match_validity!(self);
800
801        let result = unsafe { sfcgal_geometry_offset_polygon(self.c_geom.as_ptr(), radius) };
802
803        unsafe { SFCGeometry::new_from_raw(result, true) }
804    }
805
806    /// Returns the extrusion of the given `SFCGeometry` (not supported on
807    /// Solid and Multisolid).
808    pub fn extrude(&self, ex: f64, ey: f64, ez: f64) -> Result<SFCGeometry> {
809        precondition_match_validity!(self);
810
811        let result = unsafe { sfcgal_geometry_extrude(self.c_geom.as_ptr(), ex, ey, ez) };
812
813        unsafe { SFCGeometry::new_from_raw(result, true) }
814    }
815
816    /// Returns a tesselation of the given `SFCGeometry`.
817    pub fn tesselate(&self) -> Result<SFCGeometry> {
818        precondition_match_validity!(self);
819
820        let result = unsafe { sfcgal_geometry_tesselate(self.c_geom.as_ptr()) };
821
822        unsafe { SFCGeometry::new_from_raw(result, true) }
823    }
824
825    /// Returns a triangulation of the given `SFCGeometry`.
826    pub fn triangulate_2dz(&self) -> Result<SFCGeometry> {
827        precondition_match_validity!(self);
828
829        let result = unsafe { sfcgal_geometry_triangulate_2dz(self.c_geom.as_ptr()) };
830
831        unsafe { SFCGeometry::new_from_raw(result, true) }
832    }
833
834    /// Returns the convex hull of the given `SFCGeometry`.
835    pub fn convexhull(&self) -> Result<SFCGeometry> {
836        precondition_match_validity!(self);
837
838        let result = unsafe { sfcgal_geometry_convexhull(self.c_geom.as_ptr()) };
839
840        unsafe { SFCGeometry::new_from_raw(result, true) }
841    }
842
843    /// Returns the 3d convex hull of the given `SFCGeometry`.
844    pub fn convexhull_3d(&self) -> Result<SFCGeometry> {
845        precondition_match_validity!(self);
846
847        let result = unsafe { sfcgal_geometry_convexhull_3d(self.c_geom.as_ptr()) };
848
849        unsafe { SFCGeometry::new_from_raw(result, true) }
850    }
851
852    /// Returns the substring of the given `SFCGeometry` LineString between
853    /// fractional distances.
854    pub fn line_substring(&self, start: f64, end: f64) -> Result<SFCGeometry> {
855        precondition_match_type!(self, GeomType::Linestring);
856
857        precondition_match_validity!(self);
858
859        precondition_index_in_range!(start, (-1.0..=1.0));
860
861        precondition_index_in_range!(end, (-1.0..=1.0));
862
863        let result = unsafe { sfcgal_geometry_line_sub_string(self.c_geom.as_ptr(), start, end) };
864
865        unsafe { SFCGeometry::new_from_raw(result, true) }
866    }
867
868    /// Returns the alpha shape of the given `SFCGeometry` Point set.
869    pub fn alpha_shapes(&self, alpha: f64, allow_holes: bool) -> Result<SFCGeometry> {
870        if !self.is_valid().unwrap() {
871            return Err(format_err!(
872                "Error: alpha shapes can only be computed on valid geometries"
873            ));
874        }
875
876        if alpha < 0.0 || !alpha.is_finite() {
877            return Err(format_err!(
878                "Error: alpha parameter must be positive or equal to 0.0, got {}",
879                alpha,
880            ));
881        }
882
883        let result =
884            unsafe { sfcgal_geometry_alpha_shapes(self.c_geom.as_ptr(), alpha, allow_holes) };
885
886        unsafe { SFCGeometry::new_from_raw(result, true) }
887    }
888
889    /// Return the optimal alpha shape of the given `SFCGeometry` Point set.
890    pub fn optimal_alpha_shapes(
891        &self,
892        allow_holes: bool,
893        nb_components: usize,
894    ) -> Result<SFCGeometry> {
895        precondition_match_validity!(self);
896
897        let result = unsafe {
898            sfcgal_geometry_optimal_alpha_shapes(self.c_geom.as_ptr(), allow_holes, nb_components)
899        };
900
901        unsafe { SFCGeometry::new_from_raw(result, true) }
902    }
903
904    /// Create a SFCGeometry collection type (MultiPoint, MultiLineString,
905    /// MultiPolygon, MultiSolid or GeometryCollection) given a
906    /// mutable slice of `SFCGeometry`'s (this is a destructive
907    /// operation).
908    /// ```rust
909    /// use sfcgal::SFCGeometry;
910    /// let a = SFCGeometry::new("POINT (1.0 1.0)").unwrap();
911    /// let b = SFCGeometry::new("POINT (2.0 2.0)").unwrap();
912    /// let g = SFCGeometry::create_collection(&mut[a, b]).unwrap();
913    /// assert_eq!(
914    ///     g.to_wkt_decim(1).unwrap(),
915    ///     "MULTIPOINT ((1.0 1.0),(2.0 2.0))",
916    /// );
917    /// ```
918    pub fn create_collection(geoms: &mut [SFCGeometry]) -> Result<SFCGeometry> {
919        if geoms.is_empty() {
920            let res_geom = unsafe { sfcgal_geometry_collection_create() };
921
922            return unsafe { SFCGeometry::new_from_raw(res_geom, true) };
923        }
924
925        let types = geoms
926            .iter()
927            .map(|g| g._type().unwrap())
928            .collect::<Vec<GeomType>>();
929
930        let multis = types
931            .iter()
932            .map(|gt| gt.is_collection_type())
933            .collect::<Vec<bool>>();
934
935        if !is_all_same(&types) || multis.iter().any(|&x| x) {
936            let res_geom = unsafe { sfcgal_geometry_collection_create() };
937
938            make_multi_geom(res_geom, geoms)
939        } else if types[0] == GeomType::Point {
940            let res_geom = unsafe { sfcgal_multi_point_create() };
941
942            make_multi_geom(res_geom, geoms)
943        } else if types[0] == GeomType::Linestring {
944            let res_geom = unsafe { sfcgal_multi_linestring_create() };
945
946            make_multi_geom(res_geom, geoms)
947        } else if types[0] == GeomType::Polygon {
948            let res_geom = unsafe { sfcgal_multi_polygon_create() };
949
950            make_multi_geom(res_geom, geoms)
951        } else if types[0] == GeomType::Solid {
952            let mut res_geom = SFCGeometry::new("MULTISOLID EMPTY")?;
953
954            res_geom.owned = false;
955
956            make_multi_geom(res_geom.c_geom.as_ptr(), geoms)
957        } else {
958            unreachable!();
959        }
960    }
961
962    /// Get the members of a SFCGeometry.
963    /// Returns Err if the SFCGeometry if not a collection (i.e. if it's
964    /// type is not in { MultiPoint, MultiLineString, MultiPolygon,
965    /// MultiSolid, GeometryCollection }). The original geometry
966    /// stay untouched.
967    /// ```rust
968    /// use sfcgal::SFCGeometry;
969    /// let g = SFCGeometry::new("MULTIPOINT ((1.0 1.0),(2.0
970    /// 2.0))").unwrap(); let members =
971    /// g.get_collection_members().unwrap(); assert_eq!(
972    ///     members[0].to_wkt_decim(1).unwrap(),
973    ///     "POINT (1.0 1.0)",
974    /// );
975    /// assert_eq!(
976    ///     members[1].to_wkt_decim(1).unwrap(),
977    ///     "POINT (2.0 2.0)",
978    /// );
979    /// ```
980    pub fn get_collection_members(self) -> Result<Vec<SFCGeometry>> {
981        let _type = self._type()?;
982
983        if !_type.is_collection_type() {
984            return Err(format_err!(
985                "Error: the given geometry doesn't have any member ({:?} is not a collection type)",
986                _type,
987            ));
988        }
989
990        unsafe {
991            let ptr = self.c_geom.as_ptr();
992
993            let n_geom = sfcgal_geometry_collection_num_geometries(ptr);
994
995            let mut result = Vec::new();
996
997            for n in 0..n_geom {
998                let _original_c_geom = sfcgal_geometry_collection_geometry_n(ptr, n);
999
1000                let clone_c_geom = sfcgal_geometry_clone(_original_c_geom);
1001
1002                result.push(SFCGeometry::new_from_raw(clone_c_geom, true)?);
1003            }
1004
1005            Ok(result)
1006        }
1007    }
1008
1009    /// Creates an empty point
1010    pub fn point_create() -> Result<SFCGeometry> {
1011        let result = unsafe { sfcgal_point_create() };
1012
1013        unsafe { SFCGeometry::new_from_raw(result, true) }
1014    }
1015
1016    /// Creates a point from two X and Y coordinates
1017    pub fn point_create_from_xy(x: f64, y: f64) -> Result<SFCGeometry> {
1018        let result = unsafe { sfcgal_point_create_from_xy(x, y) };
1019
1020        unsafe { SFCGeometry::new_from_raw(result, true) }
1021    }
1022
1023    /// Creates a point from three X, Y and M coordinates
1024    pub fn point_create_from_xym(x: f64, y: f64, m: f64) -> Result<SFCGeometry> {
1025        let result = unsafe { sfcgal_point_create_from_xym(x, y, m) };
1026
1027        unsafe { SFCGeometry::new_from_raw(result, true) }
1028    }
1029
1030    /// Creates a point from three X, Y and Z coordinates
1031    pub fn point_create_from_xyz(x: f64, y: f64, z: f64) -> Result<SFCGeometry> {
1032        let result = unsafe { sfcgal_point_create_from_xyz(x, y, z) };
1033
1034        unsafe { SFCGeometry::new_from_raw(result, true) }
1035    }
1036
1037    /// Create a point from x, y, z, m components.
1038    pub fn point_create_from_xyzm(&self, x: f64, y: f64, z: f64, m: f64) -> Result<SFCGeometry> {
1039        let result = unsafe { sfcgal_point_create_from_xyzm(x, y, z, m) };
1040
1041        unsafe { SFCGeometry::new_from_raw(result, true) }
1042    }
1043
1044    /// Returns the X coordinate of the given Point
1045    pub fn point_x(&self) -> Result<f64> {
1046        precondition_match_type!(self, GeomType::Point);
1047
1048        precondition_match_not_empty!(self);
1049
1050        unsafe { Ok(sfcgal_point_x(self.c_geom.as_ptr())) }
1051    }
1052
1053    /// Returns the Y coordinate of the given Point
1054    pub fn point_y(&self) -> Result<f64> {
1055        precondition_match_type!(self, GeomType::Point);
1056
1057        precondition_match_not_empty!(self);
1058
1059        unsafe { Ok(sfcgal_point_y(self.c_geom.as_ptr())) }
1060    }
1061
1062    /// Returns the Z coordinate of the given Point
1063    pub fn point_z(&self) -> Result<f64> {
1064        precondition_match_type!(self, GeomType::Point);
1065
1066        precondition_match_not_empty!(self);
1067
1068        unsafe {
1069            let result = sfcgal_point_z(self.c_geom.as_ptr());
1070
1071            check_nan_value(result)
1072        }
1073    }
1074
1075    /// Returns the M coordinate of the given Point
1076    pub fn point_m(&self) -> Result<f64> {
1077        precondition_match_type!(self, GeomType::Point);
1078
1079        precondition_match_not_empty!(self);
1080
1081        unsafe {
1082            let result = sfcgal_point_m(self.c_geom.as_ptr());
1083
1084            check_nan_value(result)
1085        }
1086    }
1087
1088    /// Sets one vertex of a Triangle
1089    pub fn triangle_set_vertex(&self, index: i32, vertex: &SFCGeometry) -> Result<()> {
1090        precondition_match_type!(self, GeomType::Triangle);
1091
1092        precondition_index_in_range!(index, (0..3));
1093
1094        unsafe {
1095            let explicit_converted_int: ::std::os::raw::c_int = index;
1096
1097            sfcgal_triangle_set_vertex(
1098                self.c_geom.as_ptr(),
1099                explicit_converted_int,
1100                vertex.c_geom.as_ptr(),
1101            );
1102
1103            Ok(())
1104        }
1105    }
1106
1107    /// Sets one vertex of a Triangle from two coordinates
1108    pub fn triangle_set_vertex_from_xy(&self, index: i32, x: f64, y: f64) -> Result<()> {
1109        precondition_match_type!(self, GeomType::Triangle);
1110
1111        precondition_index_in_range!(index, (0..3));
1112
1113        unsafe {
1114            let explicit_converted_int: ::std::os::raw::c_int = index;
1115
1116            sfcgal_triangle_set_vertex_from_xy(self.c_geom.as_ptr(), explicit_converted_int, x, y);
1117
1118            Ok(())
1119        }
1120    }
1121
1122    /// Sets one vertex of a Triangle from three coordinates
1123    pub fn triangle_set_vertex_from_xyz(&self, index: i32, x: f64, y: f64, z: f64) -> Result<()> {
1124        precondition_match_type!(self, GeomType::Triangle);
1125
1126        precondition_index_in_range!(index, (0..3));
1127
1128        unsafe {
1129            let explicit_converted_int: ::std::os::raw::c_int = index;
1130
1131            sfcgal_triangle_set_vertex_from_xyz(
1132                self.c_geom.as_ptr(),
1133                explicit_converted_int,
1134                x,
1135                y,
1136                z,
1137            );
1138
1139            Ok(())
1140        }
1141    }
1142
1143    /// Returns the straight skeleton partition for the given Polygon
1144    pub fn straight_skeleton_partition(&self, auto_orientation: bool) -> Result<SFCGeometry> {
1145        match self._type()? {
1146            GeomType::Polygon | GeomType::Multipolygon | GeomType::Triangle => (),
1147            _ => {
1148                bail!("Wrong input Geometry type");
1149            }
1150        }
1151
1152        precondition_match_validity!(self);
1153
1154        let result = unsafe {
1155            sfcgal_geometry_straight_skeleton_partition(self.c_geom.as_ptr(), auto_orientation)
1156        };
1157
1158        unsafe { SFCGeometry::new_from_raw(result, true) }
1159    }
1160
1161    /// Returns the visibility polygon of a Point inside a Polygon
1162    pub fn visibility_point(&self, point: &SFCGeometry) -> Result<SFCGeometry> {
1163        precondition_match_type!(self, GeomType::Polygon);
1164
1165        precondition_match_type_other!(point, GeomType::Point);
1166
1167        precondition_match_validity!(self);
1168
1169        if !self.covers(point)? {
1170            bail!("Point must be inside the Polygon");
1171        }
1172
1173        let result = unsafe {
1174            sfcgal_geometry_visibility_point(self.c_geom.as_ptr(), point.c_geom.as_ptr())
1175        };
1176
1177        unsafe { SFCGeometry::new_from_raw(result, true) }
1178    }
1179
1180    /// Rotates a geometry around the origin (0,0,0) by a given angle
1181    pub fn rotate(&self, angle: f64) -> Result<SFCGeometry> {
1182        let result = unsafe { sfcgal_geometry_rotate(self.c_geom.as_ptr(), angle) };
1183
1184        unsafe { SFCGeometry::new_from_raw(result, true) }
1185    }
1186
1187    /// Rotates a geometry around the X axis by a given angle
1188    pub fn rotate_x(&self, angle: f64) -> Result<SFCGeometry> {
1189        let result = unsafe { sfcgal_geometry_rotate_x(self.c_geom.as_ptr(), angle) };
1190
1191        unsafe { SFCGeometry::new_from_raw(result, true) }
1192    }
1193
1194    /// Rotates a geometry around the Y axis by a given angle
1195    pub fn rotate_y(&self, angle: f64) -> Result<SFCGeometry> {
1196        let result = unsafe { sfcgal_geometry_rotate_y(self.c_geom.as_ptr(), angle) };
1197
1198        unsafe { SFCGeometry::new_from_raw(result, true) }
1199    }
1200
1201    /// Rotates a geometry around the Z axis by a given angle
1202    pub fn rotate_z(&self, angle: f64) -> Result<SFCGeometry> {
1203        let result = unsafe { sfcgal_geometry_rotate_z(self.c_geom.as_ptr(), angle) };
1204
1205        unsafe { SFCGeometry::new_from_raw(result, true) }
1206    }
1207
1208    /// Rotates a geometry around a specified point by a given angle
1209    pub fn rotate_2d(&self, angle: f64, origin_x: f64, origin_y: f64) -> Result<SFCGeometry> {
1210        let result =
1211            unsafe { sfcgal_geometry_rotate_2d(self.c_geom.as_ptr(), angle, origin_x, origin_y) };
1212
1213        unsafe { SFCGeometry::new_from_raw(result, true) }
1214    }
1215
1216    /// Rotates a 3D geometry around a specified axis by a given angle
1217    pub fn rotate_3d(
1218        &self,
1219        angle: f64,
1220        axis_x_angle: f64,
1221        axis_y_angle: f64,
1222        axis_z_angle: f64,
1223    ) -> Result<SFCGeometry> {
1224        let result = unsafe {
1225            sfcgal_geometry_rotate_3d(
1226                self.c_geom.as_ptr(),
1227                angle,
1228                axis_x_angle,
1229                axis_y_angle,
1230                axis_z_angle,
1231            )
1232        };
1233
1234        unsafe { SFCGeometry::new_from_raw(result, true) }
1235    }
1236
1237    /// Rotates a 3D geometry around a specified axis and center point by a
1238    /// given
1239    #[allow(clippy::too_many_arguments)]
1240    pub fn rotate_3d_around_center(
1241        &self,
1242        angle: f64,
1243        axis_x_angle: f64,
1244        axis_y_angle: f64,
1245        axis_z_angle: f64,
1246        center_x: f64,
1247        center_y: f64,
1248        center_z: f64,
1249    ) -> Result<SFCGeometry> {
1250        let result = unsafe {
1251            sfcgal_geometry_rotate_3d_around_center(
1252                self.c_geom.as_ptr(),
1253                angle,
1254                axis_x_angle,
1255                axis_y_angle,
1256                axis_z_angle,
1257                center_x,
1258                center_y,
1259                center_z,
1260            )
1261        };
1262
1263        unsafe { SFCGeometry::new_from_raw(result, true) }
1264    }
1265
1266    /// Force a Right Handed Rule on the given Geometry
1267    pub fn force_rhr(&self) -> Result<SFCGeometry> {
1268        precondition_match_validity!(self);
1269
1270        let result = unsafe { sfcgal_geometry_force_rhr(self.c_geom.as_ptr()) };
1271
1272        unsafe { SFCGeometry::new_from_raw(result, true) }
1273    }
1274
1275    /// Force a Left Handed Rule on the given Geometry
1276    pub fn force_lhr(&self) -> Result<SFCGeometry> {
1277        precondition_match_validity!(self);
1278
1279        let result = unsafe { sfcgal_geometry_force_lhr(self.c_geom.as_ptr()) };
1280
1281        unsafe { SFCGeometry::new_from_raw(result, true) }
1282    }
1283
1284    /// Scale a geometry by a given factor
1285    pub fn scale(&self, scale: f64) -> Result<SFCGeometry> {
1286        let result = unsafe { sfcgal_geometry_scale(self.c_geom.as_ptr(), scale) };
1287
1288        unsafe { SFCGeometry::new_from_raw(result, true) }
1289    }
1290
1291    /// Scale a geometry by different factors for each dimension
1292    pub fn scale_3d(&self, scale_x: f64, scale_y: f64, scale_z: f64) -> Result<SFCGeometry> {
1293        let result =
1294            unsafe { sfcgal_geometry_scale_3d(self.c_geom.as_ptr(), scale_x, scale_y, scale_z) };
1295
1296        unsafe { SFCGeometry::new_from_raw(result, true) }
1297    }
1298
1299    /// Scale a geometry by different factors for each dimension around a
1300    /// center
1301    pub fn scale_3d_around_center(
1302        &self,
1303        factor_x: f64,
1304        factor_y: f64,
1305        factor_z: f64,
1306        center_x: f64,
1307        center_y: f64,
1308        center_z: f64,
1309    ) -> Result<SFCGeometry> {
1310        let result = unsafe {
1311            sfcgal_geometry_scale_3d_around_center(
1312                self.c_geom.as_ptr(),
1313                factor_x,
1314                factor_y,
1315                factor_z,
1316                center_x,
1317                center_y,
1318                center_z,
1319            )
1320        };
1321
1322        unsafe { SFCGeometry::new_from_raw(result, true) }
1323    }
1324
1325    /// Round coordinates of the given Geometry
1326    pub fn round(&self, value: i32) -> Result<SFCGeometry> {
1327        precondition_match_validity!(self);
1328
1329        let explicit_conversion: ::std::os::raw::c_int = value;
1330
1331        let result = unsafe { sfcgal_geometry_round(self.c_geom.as_ptr(), explicit_conversion) };
1332
1333        unsafe { SFCGeometry::new_from_raw(result, true) }
1334    }
1335
1336    /// Computes a 3D buffer around a geometry
1337    pub fn buffer3d(
1338        &self,
1339        radius: f64,
1340        segments: i32,
1341        buffer_type: BufferType,
1342    ) -> Result<SFCGeometry> {
1343        match self._type()? {
1344            GeomType::Point | GeomType::Linestring => (),
1345            _ => bail!("Geometry must be a Point or a Linestring"),
1346        }
1347
1348        precondition_match_validity!(self);
1349
1350        if radius <= 0. {
1351            bail!("Radius must be greater than 0.0");
1352        }
1353
1354        if segments <= 3 {
1355            bail!("The algorithm needs at least 3 segments");
1356        }
1357
1358        unsafe {
1359            let result = sfcgal_geometry_buffer3d(
1360                self.c_geom.as_ptr(),
1361                radius,
1362                segments,
1363                buffer_type.to_u32().unwrap(),
1364            );
1365
1366            SFCGeometry::new_from_raw(result, true)
1367        }
1368    }
1369
1370    /// Gets the validity flag of the geometry
1371    pub fn has_validity_flag(&self) -> i32 {
1372        unsafe { sfcgal_geometry_has_validity_flag(self.c_geom.as_ptr()) }
1373    }
1374
1375    /// Build the visibility polygon of the segment [pointA ; pointB] on a
1376    /// Polygon
1377    pub fn visibility_segment(
1378        &self,
1379        point_a: &SFCGeometry,
1380        point_b: &SFCGeometry,
1381    ) -> Result<SFCGeometry> {
1382        precondition_match_validity!(self);
1383
1384        let result = unsafe {
1385            sfcgal_geometry_visibility_segment(
1386                self.c_geom.as_ptr(),
1387                point_a.c_geom.as_ptr(),
1388                point_b.c_geom.as_ptr(),
1389            )
1390        };
1391
1392        unsafe { SFCGeometry::new_from_raw(result, true) }
1393    }
1394
1395    /// Convert a PolyhedralSurface to a Solid
1396    pub fn make_solid(&self) -> Result<SFCGeometry> {
1397        precondition_match_validity!(self);
1398
1399        let result = unsafe { sfcgal_geometry_make_solid(self.c_geom.as_ptr()) };
1400
1401        unsafe { SFCGeometry::new_from_raw(result, true) }
1402    }
1403
1404    /// Translate a geometry by a 2D vector
1405    pub fn translate_2d(&self, dx: f64, dy: f64) -> Result<SFCGeometry> {
1406        let result = unsafe { sfcgal_geometry_translate_2d(self.c_geom.as_ptr(), dx, dy) };
1407
1408        unsafe { SFCGeometry::new_from_raw(result, true) }
1409    }
1410
1411    /// Translate a geometry by a 3D vector
1412    pub fn translate_3d(
1413        &self,
1414        translation_x: f64,
1415        translation_y: f64,
1416        translation_z: f64,
1417    ) -> Result<SFCGeometry> {
1418        let result = unsafe {
1419            sfcgal_geometry_translate_3d(
1420                self.c_geom.as_ptr(),
1421                translation_x,
1422                translation_y,
1423                translation_z,
1424            )
1425        };
1426
1427        unsafe { SFCGeometry::new_from_raw(result, true) }
1428    }
1429
1430    /// Sets the validity flag of the geometry.
1431    pub fn force_valid(&self, validity: i32) {
1432        unsafe { sfcgal_geometry_force_valid(self.c_geom.as_ptr(), validity) };
1433    }
1434
1435    /// Set the geometry validation mode
1436    pub fn set_geometry_validation(enabled: i32) {
1437        unsafe { sfcgal_set_geometry_validation(enabled) };
1438    }
1439
1440    /// Creates an empty Triangle
1441    pub fn triangle_create() -> Result<SFCGeometry> {
1442        let result = unsafe { sfcgal_triangle_create() };
1443
1444        unsafe { SFCGeometry::new_from_raw(result, true) }
1445    }
1446
1447    /// Returns one the Triangle's vertex as a Point
1448    pub fn triangle_vertex(&self, index: i32) -> Result<SFCGeometry> {
1449        precondition_match_type!(self, GeomType::Triangle);
1450
1451        precondition_index_in_range!(index, (0..3));
1452
1453        let result = unsafe { sfcgal_triangle_vertex(self.c_geom.as_ptr(), index) };
1454
1455        let convert_mutability = result as *mut c_void;
1456
1457        unsafe { SFCGeometry::new_from_raw(convert_mutability, true) }
1458    }
1459
1460    /// Creates a Triangle from three given Point
1461    pub fn triangle_create_from_points(
1462        point_a: &SFCGeometry,
1463        point_b: &SFCGeometry,
1464        point_c: &SFCGeometry,
1465    ) -> Result<SFCGeometry> {
1466        precondition_match_type!(point_a, GeomType::Point);
1467
1468        precondition_match_type!(point_a, GeomType::Point);
1469
1470        precondition_match_type!(point_a, GeomType::Point);
1471
1472        let result = unsafe {
1473            sfcgal_triangle_create_from_points(
1474                point_a.c_geom.as_ptr(),
1475                point_b.c_geom.as_ptr(),
1476                point_c.c_geom.as_ptr(),
1477            )
1478        };
1479
1480        unsafe { SFCGeometry::new_from_raw(result, true) }
1481    }
1482
1483    /// Creates an empty TriangulatedSurface
1484    pub fn triangulated_surface_create() -> Result<SFCGeometry> {
1485        let result = unsafe { sfcgal_triangulated_surface_create() };
1486
1487        unsafe { SFCGeometry::new_from_raw(result, true) }
1488    }
1489
1490    /// Returns the ith Triangle of a given TriangulatedSurface
1491    pub fn triangulated_surface_triangle_n(&self, index: usize) -> Result<SFCGeometry> {
1492        precondition_match_type!(self, GeomType::Triangulatedsurface);
1493
1494        precondition_index_in_result_value!(index, self.triangulated_surface_num_triangles());
1495
1496        unsafe {
1497            let result = sfcgal_triangulated_surface_triangle_n(self.c_geom.as_ptr(), index);
1498
1499            let convert_mutability = result as *mut c_void;
1500
1501            SFCGeometry::new_from_raw(convert_mutability, true)
1502        }
1503    }
1504
1505    /// Adds a Triangle to a given TriangulatedSurface
1506    pub fn triangulated_surface_add_triangle(&self, triangle: &SFCGeometry) -> Result<()> {
1507        precondition_match_type!(self, GeomType::Triangulatedsurface);
1508
1509        precondition_match_type_other!(triangle, GeomType::Triangle);
1510
1511        unsafe {
1512            sfcgal_triangulated_surface_add_triangle(self.c_geom.as_ptr(), triangle.c_geom.as_ptr())
1513        };
1514
1515        Ok(())
1516    }
1517
1518    /// Returns the number of triangles of a given TriangulatedSurface
1519    pub fn triangulated_surface_num_triangles(&self) -> Result<usize> {
1520        precondition_match_type!(self, GeomType::Triangulatedsurface);
1521
1522        unsafe {
1523            Ok(sfcgal_triangulated_surface_num_triangles(
1524                self.c_geom.as_ptr(),
1525            ))
1526        }
1527    }
1528
1529    /// Returns the optimal convex partition of a geometry (polygon without
1530    /// hole)
1531    pub fn optimal_convex_partition_2(&self) -> Result<SFCGeometry> {
1532        precondition_match_validity!(self);
1533
1534        let result = unsafe { sfcgal_optimal_convex_partition_2(self.c_geom.as_ptr()) };
1535
1536        unsafe { SFCGeometry::new_from_raw(result, true) }
1537    }
1538
1539    /// Returns the approximal convex partition of a geometry (polygon
1540    /// without hole)
1541    pub fn approx_convex_partition_2(&self) -> Result<SFCGeometry> {
1542        precondition_match_validity!(self);
1543
1544        let result = unsafe { sfcgal_approx_convex_partition_2(self.c_geom.as_ptr()) };
1545
1546        unsafe { SFCGeometry::new_from_raw(result, true) }
1547    }
1548
1549    /// Returns the greene approximal convex partition of a geometry
1550    /// (polygon without hole)
1551    pub fn greene_approx_convex_partition_2(&self) -> Result<SFCGeometry> {
1552        precondition_match_validity!(self);
1553
1554        let result = unsafe { sfcgal_greene_approx_convex_partition_2(self.c_geom.as_ptr()) };
1555
1556        unsafe { SFCGeometry::new_from_raw(result, true) }
1557    }
1558
1559    /// Returns the y monotone partition of a geometry (polygon without
1560    /// hole)
1561    pub fn y_monotone_partition_2(&self) -> Result<SFCGeometry> {
1562        precondition_match_validity!(self);
1563
1564        let result = unsafe { sfcgal_y_monotone_partition_2(self.c_geom.as_ptr()) };
1565
1566        unsafe { SFCGeometry::new_from_raw(result, true) }
1567    }
1568
1569    /// Returns the number of points of the given LineString
1570    pub fn linestring_num_points(&self) -> Result<usize> {
1571        precondition_match_type!(self, GeomType::Linestring);
1572
1573        unsafe { Ok(sfcgal_linestring_num_points(self.c_geom.as_ptr())) }
1574    }
1575
1576    /// Adds a point to a LineString
1577    pub fn linestring_add_point(&self, point: &SFCGeometry) -> Result<()> {
1578        precondition_match_type!(self, GeomType::Linestring);
1579
1580        precondition_match_type!(point, GeomType::Point);
1581
1582        unsafe { sfcgal_linestring_add_point(self.c_geom.as_ptr(), point.c_geom.as_ptr()) };
1583
1584        Ok(())
1585    }
1586
1587    /// Returns the ith point of a given LineString
1588    pub fn linestring_point_n(&self, index: usize) -> Result<SFCGeometry> {
1589        precondition_match_type!(self, GeomType::Point);
1590
1591        match self.linestring_num_points() {
1592            Ok(num) => {
1593                if index >= num {
1594                    bail!("Overflowing index")
1595                }
1596            }
1597            Err(e) => bail!(e),
1598        }
1599
1600        unsafe {
1601            let result = sfcgal_linestring_point_n(self.c_geom.as_ptr(), index);
1602
1603            let convert_mutability = result as *mut c_void;
1604
1605            SFCGeometry::new_from_raw(convert_mutability, true)
1606        }
1607    }
1608
1609    /// Creates an empty LineString
1610    pub fn linestring_create() -> Result<SFCGeometry> {
1611        let result = unsafe { sfcgal_linestring_create() };
1612
1613        unsafe { SFCGeometry::new_from_raw(result, true) }
1614    }
1615
1616    /// Creates an empty PolyhedralSurface
1617    pub fn polyhedral_surface_create() -> Result<SFCGeometry> {
1618        let result = unsafe { sfcgal_polyhedral_surface_create() };
1619
1620        unsafe { SFCGeometry::new_from_raw(result, true) }
1621    }
1622
1623    /// Returns the number of polygons of a given PolyhedralSurface
1624    pub fn polyhedral_surface_num_polygons(&self) -> Result<usize> {
1625        precondition_match_type!(self, GeomType::Polyhedralsurface);
1626
1627        unsafe { Ok(sfcgal_polyhedral_surface_num_polygons(self.c_geom.as_ptr())) }
1628    }
1629
1630    /// Adds a Polygon to a given PolyhedralSurface
1631    pub fn polyhedral_surface_add_polygon(&self, other: &SFCGeometry) -> Result<()> {
1632        precondition_match_type!(self, GeomType::Polyhedralsurface);
1633
1634        precondition_match_type_other!(other, GeomType::Polygon);
1635
1636        unsafe {
1637            sfcgal_polyhedral_surface_add_polygon(self.c_geom.as_ptr(), other.c_geom.as_ptr())
1638        };
1639
1640        Ok(())
1641    }
1642
1643    /// Returns the ith polygon of a given PolyhedralSurface
1644    pub fn polyhedral_surface_polygon_n(&self, index: usize) -> Result<SFCGeometry> {
1645        precondition_match_type!(self, GeomType::Polyhedralsurface);
1646
1647        precondition_index_in_result_value!(index, self.polyhedral_surface_num_polygons());
1648
1649        let result = unsafe { sfcgal_polyhedral_surface_polygon_n(self.c_geom.as_ptr(), index) };
1650
1651        let convert_mutability = result as *mut c_void;
1652
1653        unsafe { SFCGeometry::new_from_raw(convert_mutability, true) }
1654    }
1655
1656    /// Sets the error handlers. These callbacks are called on warning or
1657    /// error
1658    pub fn set_error_handlers(
1659        warning_handler: sfcgal_error_handler_t,
1660        error_handler: sfcgal_error_handler_t,
1661    ) {
1662        unsafe { sfcgal_set_error_handlers(warning_handler, error_handler) }
1663    }
1664
1665    /// Sets the error handlers. These callbacks are called on warning or
1666    /// error
1667    pub fn set_alloc_handlers(
1668        alloc_handler: sfcgal_alloc_handler_t,
1669        free_handler: sfcgal_free_handler_t,
1670    ) {
1671        unsafe { sfcgal_set_alloc_handlers(alloc_handler, free_handler) };
1672    }
1673
1674    /// Creates an empty Solid
1675    pub fn solid_create() -> Result<SFCGeometry> {
1676        let result = unsafe { sfcgal_solid_create() };
1677
1678        unsafe { SFCGeometry::new_from_raw(result, true) }
1679    }
1680
1681    /// Returns the number of shells of a given Solid
1682    pub fn solid_num_shells(&self) -> Result<usize> {
1683        precondition_match_type!(self, GeomType::Solid);
1684
1685        unsafe { Ok(sfcgal_solid_num_shells(self.c_geom.as_ptr())) }
1686    }
1687
1688    /// Creates a Solid from an exterior shell
1689    pub fn solid_create_from_exterior_shell(exterior_shell: &SFCGeometry) -> Result<SFCGeometry> {
1690        precondition_match_type_other!(exterior_shell, GeomType::Polyhedralsurface);
1691
1692        unsafe {
1693            let result = sfcgal_solid_create_from_exterior_shell(exterior_shell.c_geom.as_ptr());
1694
1695            SFCGeometry::new_from_raw(result, true)
1696        }
1697    }
1698
1699    /// Returns the ith shell of a given Solid
1700    pub fn solid_shell_n(&self, index: usize) -> Result<SFCGeometry> {
1701        precondition_match_type!(self, GeomType::Solid);
1702
1703        precondition_index_in_result_value!(index, self.solid_num_shells());
1704
1705        let result = unsafe { sfcgal_solid_shell_n(self.c_geom.as_ptr(), index) };
1706
1707        let convert_mutability = result as *mut c_void;
1708
1709        unsafe { SFCGeometry::new_from_raw(convert_mutability, true) }
1710    }
1711
1712    /// Adds a shell to a given Solid
1713    pub fn solid_add_interior_shell(&self, shell: &SFCGeometry) -> Result<()> {
1714        precondition_match_type!(self, GeomType::Solid);
1715
1716        precondition_match_type_other!(shell, GeomType::Polyhedralsurface);
1717
1718        unsafe { sfcgal_solid_add_interior_shell(self.c_geom.as_ptr(), shell.c_geom.as_ptr()) };
1719
1720        Ok(())
1721    }
1722
1723    /// Creates an empty Polygon
1724    pub fn polygon_create() -> Result<SFCGeometry> {
1725        let result = unsafe { sfcgal_polygon_create() };
1726
1727        unsafe { SFCGeometry::new_from_raw(result, true) }
1728    }
1729
1730    /// Returns the number of interior rings of a given Polygon
1731    pub fn polygon_num_interior_rings(&self) -> Result<usize> {
1732        precondition_match_type!(self, GeomType::Polygon);
1733
1734        unsafe { Ok(sfcgal_polygon_num_interior_rings(self.c_geom.as_ptr())) }
1735    }
1736
1737    /// Returns the exterior ring of a given Polygon
1738    pub fn polygon_exterior_ring(&self) -> Result<SFCGeometry> {
1739        let result = unsafe { sfcgal_polygon_exterior_ring(self.c_geom.as_ptr()) };
1740
1741        let convert_mutability = result as *mut c_void;
1742
1743        unsafe { SFCGeometry::new_from_raw(convert_mutability, true) }
1744    }
1745
1746    /// Creates an empty Polygon from an extrior ring
1747    pub fn polygon_create_from_exterior_ring(&self) -> Result<SFCGeometry> {
1748        precondition_match_type!(self, GeomType::Linestring);
1749
1750        let result = unsafe { sfcgal_polygon_create_from_exterior_ring(self.c_geom.as_ptr()) };
1751
1752        unsafe { SFCGeometry::new_from_raw(result, true) }
1753    }
1754
1755    /// Adds an interior ring to a given Polygon
1756    pub fn polygon_add_interior_ring(&self, ring: &SFCGeometry) -> Result<()> {
1757        precondition_match_type!(self, GeomType::Polygon);
1758
1759        precondition_match_type!(ring, GeomType::Linestring);
1760
1761        unsafe { sfcgal_polygon_add_interior_ring(self.c_geom.as_ptr(), ring.c_geom.as_ptr()) };
1762
1763        Ok(())
1764    }
1765
1766    /// Returns the ith interior ring of a given Polygon
1767    pub fn polygon_interior_ring_n(&self, index: usize) -> Result<SFCGeometry> {
1768        precondition_match_type!(self, GeomType::Polygon);
1769
1770        precondition_index_in_result_value!(index, self.polygon_num_interior_rings());
1771
1772        let result = unsafe { sfcgal_polygon_interior_ring_n(self.c_geom.as_ptr(), index) };
1773
1774        let convert_mutability = result as *mut c_void;
1775
1776        unsafe { SFCGeometry::new_from_raw(convert_mutability, true) }
1777    }
1778
1779    /// Creates an empty MultiSolid
1780    pub fn multi_solid_create() -> Result<SFCGeometry> {
1781        let result = unsafe { sfcgal_multi_solid_create() };
1782
1783        unsafe { SFCGeometry::new_from_raw(result, true) }
1784    }
1785
1786    // ------------------------------
1787    // Prepared geometries section
1788    // ------------------------------
1789
1790    /// # Safety
1791    /// Sets SRID associated with a given PreparedGeometry
1792    pub unsafe fn prepared_geometry_set_srid(
1793        prepared_geometry: *mut sfcgal_prepared_geometry_t,
1794        srid: srid_t,
1795    ) {
1796        sfcgal_prepared_geometry_set_srid(prepared_geometry, srid);
1797    }
1798
1799    /// # Safety
1800    /// Deletes a given PreparedGeometry
1801    pub unsafe fn prepared_geometry_delete(prepared_geometry: *mut sfcgal_prepared_geometry_t) {
1802        sfcgal_prepared_geometry_delete(prepared_geometry);
1803    }
1804
1805    /// # Safety
1806    /// Creates a PreparedGeometry from a Geometry and an SRID
1807    pub fn prepared_geometry_create_from_geometry(
1808        &self,
1809        srid: srid_t,
1810    ) -> Result<*mut sfcgal_prepared_geometry_t> {
1811        let result =
1812            unsafe { sfcgal_prepared_geometry_create_from_geometry(self.c_geom.as_ptr(), srid) };
1813
1814        check_null_prepared_geom(result)?;
1815
1816        Ok(result)
1817    }
1818
1819    /// # Safety
1820    /// Sets the Geometry associated with the given PreparedGeometry
1821    pub unsafe fn prepared_geometry_set_geometry(&self, prepared: *mut sfcgal_prepared_geometry_t) {
1822        sfcgal_prepared_geometry_set_geometry(prepared, self.c_geom.as_ptr());
1823    }
1824
1825    /// # Safety
1826    /// Returns SRID associated with a given PreparedGeometry
1827    pub unsafe fn prepared_geometry_srid(prepared: *mut sfcgal_prepared_geometry_t) -> u32 {
1828        sfcgal_prepared_geometry_srid(prepared)
1829    }
1830
1831    /// Creates an empty PreparedGeometry
1832    pub fn prepared_geometry_create() -> Result<*mut sfcgal_prepared_geometry_t> {
1833        let result = unsafe { sfcgal_prepared_geometry_create() };
1834
1835        check_null_prepared_geom(result)?;
1836
1837        Ok(result)
1838    }
1839
1840    /// # Safety
1841    /// Returns the Geometry associated with a given PreparedGeometry
1842    pub unsafe fn prepared_geometry_geometry(
1843        prepared_geometry: *mut sfcgal_prepared_geometry_t,
1844    ) -> Result<SFCGeometry> {
1845        let result = unsafe { sfcgal_prepared_geometry_geometry(prepared_geometry) };
1846
1847        // NOTE: Not proud..
1848        let converted = result as *mut c_void;
1849
1850        SFCGeometry::new_from_raw(converted, true)
1851    }
1852}
1853fn is_all_same<T>(arr: &[T]) -> bool
1854where
1855    T: Ord + Eq,
1856{
1857    arr.iter().min() == arr.iter().max()
1858}
1859fn make_multi_geom(
1860    out_multi: *mut sfcgal_geometry_t,
1861    geoms: &mut [SFCGeometry],
1862) -> Result<SFCGeometry> {
1863    for sfcgal_geom in geoms.iter_mut() {
1864        unsafe {
1865            sfcgal_geom.owned = false;
1866
1867            sfcgal_geometry_collection_add_geometry(
1868                out_multi,
1869                sfcgal_geom.c_geom.as_ptr() as *mut sfcgal_geometry_t,
1870            )
1871        };
1872    }
1873
1874    unsafe { SFCGeometry::new_from_raw(out_multi, true) }
1875}
1876
1877pub fn _sfcgal_get_full_version() -> String {
1878    let result = unsafe { sfcgal_full_version() };
1879
1880    _string(result)
1881}
1882
1883pub fn _sfcgal_get_version() -> String {
1884    let result = unsafe { sfcgal_version() };
1885
1886    _string(result)
1887}
1888
1889#[cfg(test)]
1890
1891mod tests {
1892
1893    use std::{env, f64::consts::PI};
1894
1895    use num_traits::abs;
1896
1897    use super::*;
1898    use crate::{Point2d, Point3d, ToCoordinates};
1899
1900    #[test]
1901
1902    fn creation_point_from_wkt() {
1903        let geom = SFCGeometry::new("POINT (1.0 1.0)");
1904
1905        assert!(geom.is_ok());
1906    }
1907
1908    #[test]
1909
1910    fn creation_polygon_from_wkt() {
1911        let geom = SFCGeometry::new("POLYGON((0.0 0.0, 1.0 0.0, 1.0 1.0, 0.0 0.0))");
1912
1913        assert!(geom.is_ok());
1914
1915        let geom = geom.unwrap();
1916
1917        assert!(geom.is_valid().unwrap());
1918
1919        let geom1 = SFCGeometry::new("POINT (1.0 1.0)").unwrap();
1920
1921        assert!(geom.intersects(&geom1).unwrap());
1922    }
1923
1924    #[test]
1925
1926    fn writing_to_wkt() {
1927        let geom = SFCGeometry::new("POINT (1.0 1.0)");
1928
1929        assert!(geom.is_ok());
1930
1931        let wkt = geom.unwrap().to_wkt();
1932
1933        assert!(wkt.is_ok());
1934
1935        assert_eq!(wkt.unwrap(), String::from("POINT (1/1 1/1)"));
1936    }
1937
1938    #[test]
1939
1940    fn writing_to_wkt_with_decimals() {
1941        let geom = SFCGeometry::new("POINT (1.0 1.0)");
1942
1943        assert!(geom.is_ok());
1944
1945        let wkt = geom.unwrap().to_wkt_decim(1);
1946
1947        assert!(wkt.is_ok());
1948
1949        assert_eq!(wkt.unwrap(), String::from("POINT (1.0 1.0)"));
1950    }
1951
1952    #[test]
1953
1954    fn creation_failed_with_error_message() {
1955        let geom = SFCGeometry::new("POINT (1, 1)");
1956
1957        assert!(geom.is_err());
1958
1959        assert_eq!(geom.err().unwrap().to_string(), "Obtained null pointer when creating geometry: WKT parse error, Coordinate dimension < 2 (, 1))",)
1960    }
1961
1962    #[test]
1963
1964    fn distance_to_other() {
1965        let pt1 = SFCGeometry::new("POINT (1.0 1.0)").unwrap();
1966
1967        let pt2 = SFCGeometry::new("POINT (10.0 1.0)").unwrap();
1968
1969        let distance = pt1.distance(&pt2).unwrap();
1970
1971        assert_eq!(distance, 9.0);
1972    }
1973
1974    #[test]
1975
1976    fn distance_3d_to_other() {
1977        let pt1 = SFCGeometry::new("POINT (1.0 1.0 2.0)").unwrap();
1978
1979        let pt2 = SFCGeometry::new("POINT (10.0 1.0 2.0)").unwrap();
1980
1981        let distance = pt1.distance_3d(&pt2).unwrap();
1982
1983        assert_eq!(distance, 9.0);
1984    }
1985
1986    #[test]
1987
1988    fn measured_geometry() {
1989        let pt1 = SFCGeometry::new("POINT (1.0 1.0)").unwrap();
1990
1991        let pt2 = SFCGeometry::new("POINTM(1.0 1.0 2.0)").unwrap();
1992
1993        assert!(!pt1.is_measured().unwrap());
1994
1995        assert!(pt2.is_measured().unwrap());
1996    }
1997
1998    #[test]
1999
2000    fn area() {
2001        let polygon = SFCGeometry::new("POLYGON((1 1, 3 1, 4 4, 1 3, 1 1))").unwrap();
2002
2003        assert_eq!(polygon.area().unwrap(), 6.0);
2004    }
2005
2006    #[test]
2007
2008    fn area_3d() {
2009        let polygon = SFCGeometry::new("POLYGON((1 1 1, 3 1 1, 4 4 1, 1 3 1, 1 1 1))").unwrap();
2010
2011        assert_ulps_eq!(polygon.area_3d().unwrap(), 6.0);
2012    }
2013
2014    #[test]
2015
2016    fn volume() {
2017        let cube = SFCGeometry::new(
2018            "SOLID((((0 0 0,0 0 1,0 1 1,0 1 0,0 0 0)),\
2019             ((0 0 0,0 1 0,1 1 0,1 0 0,0 0 0)),\
2020             ((0 0 0,1 0 0,1 0 1,0 0 1,0 0 0)),\
2021             ((1 0 0,1 1 0,1 1 1,1 0 1,1 0 0)),\
2022             ((0 0 1,1 0 1,1 1 1,0 1 1,0 0 1)),\
2023             ((0 1 0,0 1 1,1 1 1,1 1 0,0 1 0))))",
2024        )
2025        .unwrap();
2026
2027        assert_eq!(cube.volume().unwrap(), 1.);
2028    }
2029
2030    #[test]
2031
2032    fn volume_on_not_volume_geometry() {
2033        let surface = SFCGeometry::new(
2034            "POLYHEDRALSURFACE Z \
2035             (((0 0 0, 0 1 0, 1 1 0, 1 0 0, 0 0 0)),\
2036             ((0 0 0, 0 1 0, 0 1 1, 0 0 1, 0 0 0)),\
2037             ((0 0 0, 1 0 0, 1 0 1, 0 0 1, 0 0 0)),\
2038             ((1 1 1, 1 0 1, 0 0 1, 0 1 1, 1 1 1)),\
2039             ((1 1 1, 1 0 1, 1 0 0, 1 1 0, 1 1 1)),\
2040             ((1 1 1, 1 1 0, 0 1 0, 0 1 1, 1 1 1)))",
2041        )
2042        .unwrap();
2043
2044        assert!(surface.volume().is_err());
2045    }
2046
2047    #[test]
2048
2049    fn predicates() {
2050        let pt = SFCGeometry::new("POINT (1.0 1.0)").unwrap();
2051
2052        assert!(pt.is_valid().unwrap());
2053
2054        assert!(!pt.is_3d().unwrap());
2055
2056        assert!(!pt.is_empty().unwrap());
2057
2058        assert_eq!(
2059            pt.is_planar().err().unwrap().to_string(),
2060            "SFCGAL error: is_planar() only applies to polygons",
2061        );
2062
2063        let linestring_3d = SFCGeometry::new("LINESTRING (10.0 1.0 2.0, 1.0 2.0 1.7)").unwrap();
2064
2065        assert!(linestring_3d.is_valid().unwrap());
2066
2067        assert!(linestring_3d.is_3d().unwrap());
2068
2069        assert!(!linestring_3d.is_empty().unwrap());
2070
2071        assert_eq!(
2072            linestring_3d.is_planar().err().unwrap().to_string(),
2073            "SFCGAL error: is_planar() only applies to polygons",
2074        );
2075
2076        let empty_geom = SFCGeometry::new("LINESTRING EMPTY").unwrap();
2077
2078        assert!(empty_geom.is_valid().unwrap());
2079
2080        assert!(!empty_geom.is_3d().unwrap());
2081
2082        assert!(empty_geom.is_empty().unwrap());
2083
2084        assert_eq!(
2085            linestring_3d.is_planar().err().unwrap().to_string(),
2086            "SFCGAL error: is_planar() only applies to polygons",
2087        );
2088
2089        let polyg = SFCGeometry::new("POLYGON ((1 1, 3 1, 4 4, 1 3, 1 1))").unwrap();
2090
2091        assert!(polyg.is_valid().unwrap());
2092
2093        assert!(!polyg.is_3d().unwrap());
2094
2095        assert!(!polyg.is_empty().unwrap());
2096
2097        assert!(polyg.is_planar().unwrap());
2098
2099        assert!(pt.intersects(&polyg).unwrap());
2100
2101        assert!(!pt.intersects_3d(&linestring_3d).unwrap());
2102    }
2103
2104    #[test]
2105
2106    fn validity_detail_on_valid_geom() {
2107        let line = SFCGeometry::new("LINESTRING (10.0 1.0 2.0, 1.0 2.0 1.7)").unwrap();
2108
2109        assert!(line.is_valid().unwrap());
2110
2111        assert_eq!(line.validity_detail().unwrap(), None);
2112    }
2113
2114    #[test]
2115
2116    fn validity_detail_on_invalid_geom_1() {
2117        let surface = SFCGeometry::new(
2118            "POLYHEDRALSURFACE Z \
2119             (((0 0 0, 0 1 0, 1 1 0, 1 0 0, 0 0 0)),\
2120             ((0 0 0, 0 1 0, 0 1 1, 0 0 1, 0 0 0)),\
2121             ((0 0 0, 1 0 0, 1 0 1, 0 0 1, 0 0 0)),\
2122             ((1 1 1, 1 0 1, 0 0 1, 0 1 1, 1 1 1)),\
2123             ((1 1 1, 1 0 1, 1 0 0, 1 1 0, 1 1 1)),\
2124             ((1 1 1, 1 1 0, 0 1 0, 0 1 1, 1 1 1)))",
2125        )
2126        .unwrap();
2127
2128        assert!(!surface.is_valid().unwrap());
2129
2130        assert_eq!(surface.validity_detail().unwrap(), Some(String::from("inconsistent orientation of PolyhedralSurface detected at edge 3 (4-7) of polygon 5")),);
2131    }
2132
2133    #[test]
2134
2135    fn validity_detail_on_invalid_geom_2() {
2136        let surface = SFCGeometry::new("POLYGON ((1 2,1 2,1 2,1 2))").unwrap();
2137
2138        assert!(!surface.is_valid().unwrap());
2139
2140        assert_eq!(
2141            surface.validity_detail().unwrap(),
2142            Some(String::from("ring 0 degenerated to a point")),
2143        );
2144    }
2145
2146    #[test]
2147
2148    fn validity_detail_on_invalid_geom_3() {
2149        let surface = SFCGeometry::new("LINESTRING (1 2, 1 2, 1 2)").unwrap();
2150
2151        assert!(!surface.is_valid().unwrap());
2152
2153        assert_eq!(
2154            surface.validity_detail().unwrap(),
2155            Some(String::from("no length")),
2156        );
2157    }
2158
2159    #[test]
2160
2161    fn straight_skeleton() {
2162        let geom = SFCGeometry::new("POLYGON ((0 0,1 0,1 1,0 1,0 0))").unwrap();
2163
2164        let result = geom.straight_skeleton().unwrap();
2165
2166        let wkt = result.to_wkt_decim(1).unwrap();
2167
2168        assert_eq!(wkt, "MULTILINESTRING ((0.0 0.0,0.5 0.5),(1.0 0.0,0.5 0.5),(1.0 1.0,0.5 0.5),(0.0 1.0,0.5 0.5))",);
2169    }
2170
2171    #[test]
2172
2173    fn straight_skeleton_distance_in_m() {
2174        let geom = SFCGeometry::new("POLYGON ((0 0,1 0,1 1,0 1,0 0))").unwrap();
2175
2176        let result = geom.straight_skeleton_distance_in_m().unwrap();
2177
2178        let wkt = result.to_wkt_decim(1).unwrap();
2179
2180        assert_eq!(
2181            wkt,
2182            "MULTILINESTRING M (\
2183             (0.0 0.0 0.0,0.5 0.5 0.5),\
2184             (1.0 0.0 0.0,0.5 0.5 0.5),\
2185             (1.0 1.0 0.0,0.5 0.5 0.5),\
2186             (0.0 1.0 0.0,0.5 0.5 0.5))",
2187        );
2188    }
2189
2190    #[test]
2191
2192    fn extrude_straight_skeleton() {
2193        // Test adapted from https://gitlab.com/sfcgal/SFCGAL/-/blob/master/test/unit/SFCGAL/algorithm/StraightSkeletonTest.cpp#L275
2194        let geom = SFCGeometry::new("POLYGON ((0 0, 5 0, 5 5, 4 5, 4 4, 0 4, 0 0))").unwrap();
2195
2196        let result = geom.extrude_straight_skeleton(2.).unwrap();
2197
2198        let wkt = result.to_wkt_decim(2).unwrap();
2199
2200        assert_eq!(
2201            wkt,
2202            "POLYHEDRALSURFACE Z (((4.00 5.00 0.00,5.00 5.00 0.00,4.00 4.00 0.00,4.00 \
2203              5.00 0.00)),((0.00 4.00 0.00,4.00 4.00 0.00,0.00 0.00 0.00,0.00 4.00 \
2204              0.00)),((4.00 4.00 0.00,5.00 0.00 0.00,0.00 0.00 0.00,4.00 4.00 \
2205              0.00)),((5.00 5.00 0.00,5.00 0.00 0.00,4.00 4.00 0.00,5.00 5.00 \
2206              0.00)),((0.00 4.00 0.00,0.00 0.00 0.00,2.00 2.00 2.00,0.00 4.00 \
2207              0.00)),((0.00 0.00 0.00,5.00 0.00 0.00,3.00 2.00 2.00,0.00 0.00 \
2208              0.00)),((2.00 2.00 2.00,0.00 0.00 0.00,3.00 2.00 2.00,2.00 2.00 \
2209              2.00)),((4.50 3.50 0.50,5.00 5.00 0.00,4.50 4.50 0.50,4.50 3.50 \
2210              0.50)),((3.00 2.00 2.00,5.00 0.00 0.00,4.50 3.50 0.50,3.00 2.00 \
2211              2.00)),((4.50 3.50 0.50,5.00 0.00 0.00,5.00 5.00 0.00,4.50 3.50 \
2212              0.50)),((5.00 5.00 0.00,4.00 5.00 0.00,4.50 4.50 0.50,5.00 5.00 \
2213              0.00)),((4.50 4.50 0.50,4.00 4.00 0.00,4.50 3.50 0.50,4.50 4.50 \
2214              0.50)),((4.50 4.50 0.50,4.00 5.00 0.00,4.00 4.00 0.00,4.50 4.50 \
2215              0.50)),((4.00 4.00 0.00,0.00 4.00 0.00,2.00 2.00 2.00,4.00 4.00 \
2216              0.00)),((4.50 3.50 0.50,4.00 4.00 0.00,3.00 2.00 2.00,4.50 3.50 \
2217              0.50)),((3.00 2.00 2.00,4.00 4.00 0.00,2.00 2.00 2.00,3.00 2.00 \
2218              2.00)))"
2219        );
2220    }
2221
2222    #[test]
2223
2224    fn extrude_polygon_straight_skeleton() {
2225        // Test adapted from https://gitlab.com/sfcgal/SFCGAL/-/blob/master/test/unit/SFCGAL/algorithm/StraightSkeletonTest.cpp#L354
2226        let geom = SFCGeometry::new(
2227            "POLYGON (( 0 0, 5 0, 5 5, 4 5, 4 4, 0 4, 0 0 ), (1 1, 1 2, 2 2, 2 1, 1 1))",
2228        )
2229        .unwrap();
2230
2231        let result = geom.extrude_polygon_straight_skeleton(9., 2.).unwrap();
2232
2233        let wkt = result.to_wkt_decim(1).unwrap();
2234
2235        assert_eq!(
2236            wkt,
2237            "POLYHEDRALSURFACE Z (((0.0 0.0 0.0,0.0 4.0 0.0,4.0 4.0 \
2238             0.0,4.0 5.0 0.0,5.0 5.0 0.0,5.0 0.0 0.0,0.0 0.0 0.0),\
2239             (1.0 1.0 0.0,2.0 1.0 0.0,2.0 2.0 0.0,1.0 2.0 0.0,1.0 1.0 0.0)),\
2240             ((0.0 0.0 0.0,0.0 0.0 9.0,0.0 4.0 9.0,0.0 4.0 0.0,0.0 0.0 0.0)),\
2241             ((0.0 4.0 0.0,0.0 4.0 9.0,4.0 4.0 9.0,4.0 4.0 0.0,0.0 4.0 0.0)),\
2242             ((4.0 4.0 0.0,4.0 4.0 9.0,4.0 5.0 9.0,4.0 5.0 0.0,4.0 4.0 0.0)),\
2243             ((4.0 5.0 0.0,4.0 5.0 9.0,5.0 5.0 9.0,5.0 5.0 0.0,4.0 5.0 0.0)),\
2244             ((5.0 5.0 0.0,5.0 5.0 9.0,5.0 0.0 9.0,5.0 0.0 0.0,5.0 5.0 0.0)),\
2245             ((5.0 0.0 0.0,5.0 0.0 9.0,0.0 0.0 9.0,0.0 0.0 0.0,5.0 0.0 0.0)),\
2246             ((1.0 1.0 0.0,1.0 1.0 9.0,2.0 1.0 9.0,2.0 1.0 0.0,1.0 1.0 0.0)),\
2247             ((2.0 1.0 0.0,2.0 1.0 9.0,2.0 2.0 9.0,2.0 2.0 0.0,2.0 1.0 0.0)),\
2248             ((2.0 2.0 0.0,2.0 2.0 9.0,1.0 2.0 9.0,1.0 2.0 0.0,2.0 2.0 0.0)),\
2249             ((1.0 2.0 0.0,1.0 2.0 9.0,1.0 1.0 9.0,1.0 1.0 0.0,1.0 2.0 0.0)),\
2250             ((4.0 5.0 9.0,5.0 5.0 9.0,4.0 4.0 9.0,4.0 5.0 9.0)),\
2251             ((2.0 1.0 9.0,5.0 0.0 9.0,0.0 0.0 9.0,2.0 1.0 9.0)),\
2252             ((5.0 5.0 9.0,5.0 0.0 9.0,4.0 4.0 9.0,5.0 5.0 9.0)),\
2253             ((2.0 1.0 9.0,0.0 0.0 9.0,1.0 1.0 9.0,2.0 1.0 9.0)),\
2254             ((1.0 2.0 9.0,1.0 1.0 9.0,0.0 0.0 9.0,1.0 2.0 9.0)),\
2255             ((0.0 4.0 9.0,2.0 2.0 9.0,1.0 2.0 9.0,0.0 4.0 9.0)),\
2256             ((0.0 4.0 9.0,1.0 2.0 9.0,0.0 0.0 9.0,0.0 4.0 9.0)),\
2257             ((4.0 4.0 9.0,5.0 0.0 9.0,2.0 2.0 9.0,4.0 4.0 9.0)),\
2258             ((4.0 4.0 9.0,2.0 2.0 9.0,0.0 4.0 9.0,4.0 4.0 9.0)),\
2259             ((2.0 2.0 9.0,5.0 0.0 9.0,2.0 1.0 9.0,2.0 2.0 9.0)),\
2260             ((0.5 2.5 9.5,0.0 0.0 9.0,0.5 0.5 9.5,0.5 2.5 9.5)),\
2261             ((1.0 3.0 10.0,0.0 4.0 9.0,0.5 2.5 9.5,1.0 3.0 10.0)),\
2262             ((0.5 2.5 9.5,0.0 4.0 9.0,0.0 0.0 9.0,0.5 2.5 9.5)),\
2263             ((2.5 0.5 9.5,5.0 0.0 9.0,3.5 1.5 10.5,2.5 0.5 9.5)),\
2264             ((0.0 0.0 9.0,5.0 0.0 9.0,2.5 0.5 9.5,0.0 0.0 9.0)),\
2265             ((0.5 0.5 9.5,0.0 0.0 9.0,2.5 0.5 9.5,0.5 0.5 9.5)),\
2266             ((4.5 3.5 9.5,5.0 5.0 9.0,4.5 4.5 9.5,4.5 3.5 9.5)),\
2267             ((3.5 2.5 10.5,3.5 1.5 10.5,4.5 3.5 9.5,3.5 2.5 10.5)),\
2268             ((4.5 3.5 9.5,5.0 0.0 9.0,5.0 5.0 9.0,4.5 3.5 9.5)),\
2269             ((3.5 1.5 10.5,5.0 0.0 9.0,4.5 3.5 9.5,3.5 1.5 10.5)),\
2270             ((5.0 5.0 9.0,4.0 5.0 9.0,4.5 4.5 9.5,5.0 5.0 9.0)),\
2271             ((4.5 4.5 9.5,4.0 4.0 9.0,4.5 3.5 9.5,4.5 4.5 9.5)),\
2272             ((4.5 4.5 9.5,4.0 5.0 9.0,4.0 4.0 9.0,4.5 4.5 9.5)),\
2273             ((3.0 3.0 10.0,0.0 4.0 9.0,1.0 3.0 10.0,3.0 3.0 10.0)),\
2274             ((3.5 2.5 10.5,4.5 3.5 9.5,3.0 3.0 10.0,3.5 2.5 10.5)),\
2275             ((3.0 3.0 10.0,4.0 4.0 9.0,0.0 4.0 9.0,3.0 3.0 10.0)),\
2276             ((4.5 3.5 9.5,4.0 4.0 9.0,3.0 3.0 10.0,4.5 3.5 9.5)),\
2277             ((2.0 1.0 9.0,1.0 1.0 9.0,0.5 0.5 9.5,2.0 1.0 9.0)),\
2278             ((2.5 0.5 9.5,2.0 1.0 9.0,0.5 0.5 9.5,2.5 0.5 9.5)),\
2279             ((1.0 1.0 9.0,1.0 2.0 9.0,0.5 2.5 9.5,1.0 1.0 9.0)),\
2280             ((0.5 0.5 9.5,1.0 1.0 9.0,0.5 2.5 9.5,0.5 0.5 9.5)),\
2281             ((1.0 3.0 10.0,2.0 2.0 9.0,3.0 3.0 10.0,1.0 3.0 10.0)),\
2282             ((0.5 2.5 9.5,1.0 2.0 9.0,1.0 3.0 10.0,0.5 2.5 9.5)),\
2283             ((1.0 3.0 10.0,1.0 2.0 9.0,2.0 2.0 9.0,1.0 3.0 10.0)),\
2284             ((2.0 2.0 9.0,2.0 1.0 9.0,2.5 0.5 9.5,2.0 2.0 9.0)),\
2285             ((3.5 2.5 10.5,3.0 3.0 10.0,3.5 1.5 10.5,3.5 2.5 10.5)),\
2286             ((3.5 1.5 10.5,2.0 2.0 9.0,2.5 0.5 9.5,3.5 1.5 10.5)),\
2287             ((3.0 3.0 10.0,2.0 2.0 9.0,3.5 1.5 10.5,3.0 3.0 10.0)))"
2288        );
2289    }
2290
2291    #[test]
2292
2293    fn tesselate() {
2294        let geom = SFCGeometry::new("POLYGON ((0.0 0.0,1.0 0.0,1.0 1.0,0.0 1.0,0.0 0.0))").unwrap();
2295
2296        let result = geom.tesselate().unwrap();
2297
2298        let output_wkt = result.to_wkt_decim(1).unwrap();
2299
2300        assert_eq!(
2301            output_wkt,
2302            "TIN (((0.0 1.0,1.0 0.0,1.0 1.0,0.0 1.0)),((0.0 1.0,0.0 0.0,1.0 0.0,0.0 1.0)))",
2303        );
2304    }
2305
2306    #[test]
2307
2308    fn offset_polygon() {
2309        let geom = SFCGeometry::new("POLYGON ((0.0 0.0,1.0 0.0,1.0 1.0,0.0 1.0,0.0 0.0))").unwrap();
2310
2311        let buff = geom.offset_polygon(1.).unwrap();
2312
2313        assert!(buff.is_valid().unwrap());
2314
2315        assert!(!buff.is_empty().unwrap());
2316    }
2317
2318    #[test]
2319
2320    fn extrude_polygon() {
2321        let geom = SFCGeometry::new("POLYGON ((0.0 0.0,1.0 0.0,1.0 1.0,0.0 1.0,0.0 0.0))").unwrap();
2322
2323        let extr = geom.extrude(0., 0., 1.).unwrap();
2324
2325        assert!(extr.is_valid().unwrap());
2326
2327        assert!(!extr.is_empty().unwrap());
2328
2329        assert_eq!(extr._type().unwrap(), GeomType::Solid);
2330    }
2331
2332    #[test]
2333
2334    fn tesselate_invariant_geom() {
2335        let input_wkt = String::from("POINT (1.0 1.0)");
2336
2337        let pt = SFCGeometry::new(&input_wkt).unwrap();
2338
2339        let result = pt.tesselate().unwrap();
2340
2341        let output_wkt = result.to_wkt_decim(1).unwrap();
2342
2343        assert_eq!(input_wkt, output_wkt);
2344    }
2345
2346    #[test]
2347
2348    fn line_substring() {
2349        let g = SFCGeometry::new("LINESTRING Z (10.0 1.0 2.0, 1.0 2.0 1.7)").unwrap();
2350
2351        let result = g.line_substring(-0.2, 0.2).unwrap();
2352
2353        assert_eq!(
2354            result.to_wkt_decim(1).unwrap(),
2355            "LINESTRING Z (2.8 1.8 1.8,8.2 1.2 1.9)"
2356        );
2357
2358        // With "start" or "end" point not in [-1; 1]
2359        assert_eq!(
2360            g.line_substring(-2., 0.2).err().unwrap().to_string(),
2361            "Index not in the expected range"
2362        );
2363    }
2364
2365    #[test]
2366
2367    fn difference_3d() {
2368        let cube1 = SFCGeometry::new(
2369            "
2370            SOLID((((0 0 0, 0 1 0, 1 1 0, 1 0 0, 0 0 0)),\
2371            ((0 0 0, 0 0 1, 0 1 1, 0 1 0, 0 0 0)),\
2372            ((0 0 0, 1 0 0, 1 0 1, 0 0 1, 0 0 0)),\
2373            ((1 1 1, 0 1 1, 0 0 1, 1 0 1, 1 1 1)),\
2374            ((1 1 1, 1 0 1, 1 0 0, 1 1 0, 1 1 1)),\
2375            ((1 1 1, 1 1 0, 0 1 0, 0 1 1, 1 1 1))))",
2376        )
2377        .unwrap();
2378
2379        let cube2 = SFCGeometry::new(
2380            "
2381            SOLID((((0 0 0.5, 0 1 0.5, 1 1 0.5, 1 0 0.5, 0 0 0.5)),\
2382            ((0 0 0.5, 0 0 1, 0 1 1, 0 1 0.5, 0 0 0.5)),\
2383            ((0 0 0.5, 1 0 0.5, 1 0 1, 0 0 1, 0 0 0.5)),\
2384            ((1 1 1, 0 1 1, 0 0 1, 1 0 1, 1 1 1)),\
2385            ((1 1 1, 1 0 1, 1 0 0.5, 1 1 0.5, 1 1 1)),\
2386            ((1 1 1, 1 1 0.5, 0 1 0.5, 0 1 1, 1 1 1))))",
2387        )
2388        .unwrap();
2389
2390        let diff = cube1.difference_3d(&cube2).unwrap();
2391
2392        assert!(diff.is_valid().unwrap());
2393
2394        assert_ulps_eq!(diff.volume().unwrap(), 0.5);
2395    }
2396
2397    #[test]
2398
2399    fn intersection_3d() {
2400        let cube1 = SFCGeometry::new(
2401            "
2402            SOLID((((0 0 0, 0 1 0, 1 1 0, 1 0 0, 0 0 0)),\
2403            ((0 0 0, 0 0 1, 0 1 1, 0 1 0, 0 0 0)),\
2404            ((0 0 0, 1 0 0, 1 0 1, 0 0 1, 0 0 0)),\
2405            ((1 1 1, 0 1 1, 0 0 1, 1 0 1, 1 1 1)),\
2406            ((1 1 1, 1 0 1, 1 0 0, 1 1 0, 1 1 1)),\
2407            ((1 1 1, 1 1 0, 0 1 0, 0 1 1, 1 1 1))))",
2408        )
2409        .unwrap();
2410
2411        let cube2 = SFCGeometry::new(
2412            "
2413            SOLID((((0 0 0.5, 0 1 0.5, 1 1 0.5, 1 0 0.5, 0 0 0.5)),\
2414            ((0 0 0.5, 0 0 1, 0 1 1, 0 1 0.5, 0 0 0.5)),\
2415            ((0 0 0.5, 1 0 0.5, 1 0 1, 0 0 1, 0 0 0.5)),\
2416            ((1 1 1, 0 1 1, 0 0 1, 1 0 1, 1 1 1)),\
2417            ((1 1 1, 1 0 1, 1 0 0.5, 1 1 0.5, 1 1 1)),\
2418            ((1 1 1, 1 1 0.5, 0 1 0.5, 0 1 1, 1 1 1))))",
2419        )
2420        .unwrap();
2421
2422        let diff = cube1.intersection_3d(&cube2).unwrap();
2423
2424        assert!(diff.is_valid().unwrap());
2425
2426        assert_ulps_eq!(diff.volume().unwrap(), 0.5);
2427    }
2428
2429    #[test]
2430
2431    fn alpha_shapes_on_point() {
2432        let multipoint = SFCGeometry::new("POINT (1 3)").unwrap();
2433
2434        let res = multipoint.alpha_shapes(20.0, false).unwrap();
2435
2436        assert_eq!(res.to_wkt_decim(1).unwrap(), "GEOMETRYCOLLECTION EMPTY");
2437    }
2438
2439    #[test]
2440
2441    fn alpha_shapes_on_multipoint() {
2442        let multipoint = SFCGeometry::new("MULTIPOINT ((1 2),(2 2),(3 0),(1 3))").unwrap();
2443
2444        let res = multipoint.alpha_shapes(20.0, false).unwrap();
2445
2446        assert_eq!(
2447            res.to_wkt_decim(1).unwrap(),
2448            "POLYGON ((1.0 2.0,1.0 3.0,2.0 2.0,3.0 0.0,1.0 2.0))"
2449        );
2450    }
2451
2452    #[test]
2453
2454    fn alpha_shapes_on_linestring() {
2455        let multipoint = SFCGeometry::new("LINESTRING (1 2,2 2,3 0,1 3)").unwrap();
2456
2457        let res = multipoint.alpha_shapes(20.0, false).unwrap();
2458
2459        assert_eq!(
2460            res.to_wkt_decim(1).unwrap(),
2461            "POLYGON ((1.0 2.0,1.0 3.0,2.0 2.0,3.0 0.0,1.0 2.0))"
2462        );
2463    }
2464
2465    #[test]
2466
2467    fn alpha_shapes_on_multilinestring() {
2468        let multipoint =
2469            SFCGeometry::new("MULTILINESTRING ((1 2,2 2,3 0,1 3), (2 6,3 5,4 2))").unwrap();
2470
2471        let res = multipoint.alpha_shapes(20.0, false).unwrap();
2472
2473        assert_eq!(
2474            res.to_wkt_decim(1).unwrap(),
2475            "POLYGON ((1.0 2.0,1.0 3.0,2.0 6.0,3.0 5.0,4.0 2.0,3.0 0.0,1.0 2.0))"
2476        );
2477    }
2478
2479    #[test]
2480
2481    fn alpha_shapes_on_polygon() {
2482        let pol1 = SFCGeometry::new("POLYGON ((0 0, 0 4, 4 4, 4 0, 0 0))").unwrap();
2483
2484        let res = pol1.alpha_shapes(10.0, false).unwrap();
2485
2486        assert_eq!(
2487            res.to_wkt_decim(1).unwrap(),
2488            "POLYGON ((0.0 0.0,0.0 4.0,4.0 4.0,4.0 0.0,0.0 0.0))"
2489        );
2490    }
2491
2492    #[test]
2493
2494    fn alpha_shapes_on_invalid() {
2495        let pol1 = SFCGeometry::new("LINESTRING (1 2, 1 2, 1 2, 1 2)").unwrap();
2496
2497        let res = pol1.alpha_shapes(10.0, false);
2498
2499        assert!(res.is_err());
2500    }
2501
2502    #[test]
2503
2504    fn create_collection_empty() {
2505        let g = SFCGeometry::create_collection(&mut []).unwrap();
2506
2507        assert_eq!(g.to_wkt_decim(1).unwrap(), "GEOMETRYCOLLECTION EMPTY",);
2508    }
2509
2510    #[test]
2511
2512    fn create_collection_heterogenous() {
2513        let a = SFCGeometry::new("POINT (1.0 1.0)").unwrap();
2514
2515        let b = SFCGeometry::new("LINESTRING Z (10.0 1.0 2.0, 1.0 2.0 1.7)").unwrap();
2516
2517        let g = SFCGeometry::create_collection(&mut [a, b]).unwrap();
2518
2519        assert_eq!(
2520            g.to_wkt_decim(1).unwrap(),
2521            "GEOMETRYCOLLECTION (POINT (1.0 1.0),LINESTRING Z (10.0 1.0 2.0,1.0 2.0 1.7))",
2522        );
2523    }
2524
2525    #[test]
2526
2527    fn create_collection_multipoint_from_points() {
2528        let a = SFCGeometry::new("POINT (1.0 1.0)").unwrap();
2529
2530        let b = SFCGeometry::new("POINT (2.0 2.0)").unwrap();
2531
2532        let g = SFCGeometry::create_collection(&mut [a, b]).unwrap();
2533
2534        assert_eq!(
2535            g.to_wkt_decim(1).unwrap(),
2536            "MULTIPOINT ((1.0 1.0),(2.0 2.0))",
2537        );
2538    }
2539
2540    #[test]
2541
2542    fn create_collection_multilinestring_from_linestrings() {
2543        let a = SFCGeometry::new("LINESTRING (10.0 1.0 2.0, 1.0 2.0 1.7)").unwrap();
2544
2545        let b = SFCGeometry::new("LINESTRING (10.0 1.0 2.0, 1.0 2.0 1.7)").unwrap();
2546
2547        let g = SFCGeometry::create_collection(&mut [a, b]).unwrap();
2548
2549        assert_eq!(
2550            g.to_wkt_decim(1).unwrap(),
2551            "MULTILINESTRING Z ((10.0 1.0 2.0,1.0 2.0 1.7),(10.0 1.0 2.0,1.0 2.0 1.7))",
2552        );
2553    }
2554
2555    #[test]
2556
2557    fn create_collection_multisolid_from_solids() {
2558        let a = SFCGeometry::new(
2559            "SOLID((((0 0 0,0 0 1,0 1 1,0 1 0,0 0 0)),\
2560             ((0 0 0,0 1 0,1 1 0,1 0 0,0 0 0)),\
2561             ((0 0 0,1 0 0,1 0 1,0 0 1,0 0 0)),\
2562             ((1 0 0,1 1 0,1 1 1,1 0 1,1 0 0)),\
2563             ((0 0 1,1 0 1,1 1 1,0 1 1,0 0 1)),\
2564             ((0 1 0,0 1 1,1 1 1,1 1 0,0 1 0))))",
2565        )
2566        .unwrap();
2567
2568        let b = SFCGeometry::new(
2569            "SOLID((((0 0 0,0 0 1,0 1 1,0 1 0,0 0 0)),\
2570             ((0 0 0,0 1 0,1 1 0,1 0 0,0 0 0)),\
2571             ((0 0 0,1 0 0,1 0 1,0 0 1,0 0 0)),\
2572             ((1 0 0,1 1 0,1 1 1,1 0 1,1 0 0)),\
2573             ((0 0 1,1 0 1,1 1 1,0 1 1,0 0 1)),\
2574             ((0 1 0,0 1 1,1 1 1,1 1 0,0 1 0))))",
2575        )
2576        .unwrap();
2577
2578        let g = SFCGeometry::create_collection(&mut [a, b]).unwrap();
2579
2580        assert_eq!(
2581            g.to_wkt_decim(1).unwrap(),
2582            "MULTISOLID Z (\
2583             ((\
2584             ((0.0 0.0 0.0,0.0 0.0 1.0,0.0 1.0 1.0,0.0 1.0 0.0,0.0 0.0 0.0)),\
2585             ((0.0 0.0 0.0,0.0 1.0 0.0,1.0 1.0 0.0,1.0 0.0 0.0,0.0 0.0 0.0)),\
2586             ((0.0 0.0 0.0,1.0 0.0 0.0,1.0 0.0 1.0,0.0 0.0 1.0,0.0 0.0 0.0)),\
2587             ((1.0 0.0 0.0,1.0 1.0 0.0,1.0 1.0 1.0,1.0 0.0 1.0,1.0 0.0 0.0)),\
2588             ((0.0 0.0 1.0,1.0 0.0 1.0,1.0 1.0 1.0,0.0 1.0 1.0,0.0 0.0 1.0)),\
2589             ((0.0 1.0 0.0,0.0 1.0 1.0,1.0 1.0 1.0,1.0 1.0 0.0,0.0 1.0 0.0))\
2590             )),\
2591             ((\
2592             ((0.0 0.0 0.0,0.0 0.0 1.0,0.0 1.0 1.0,0.0 1.0 0.0,0.0 0.0 0.0)),\
2593             ((0.0 0.0 0.0,0.0 1.0 0.0,1.0 1.0 0.0,1.0 0.0 0.0,0.0 0.0 0.0)),\
2594             ((0.0 0.0 0.0,1.0 0.0 0.0,1.0 0.0 1.0,0.0 0.0 1.0,0.0 0.0 0.0)),\
2595             ((1.0 0.0 0.0,1.0 1.0 0.0,1.0 1.0 1.0,1.0 0.0 1.0,1.0 0.0 0.0)),\
2596             ((0.0 0.0 1.0,1.0 0.0 1.0,1.0 1.0 1.0,0.0 1.0 1.0,0.0 0.0 1.0)),\
2597             ((0.0 1.0 0.0,0.0 1.0 1.0,1.0 1.0 1.0,1.0 1.0 0.0,0.0 1.0 0.0))\
2598             ))\
2599             )",
2600        );
2601    }
2602
2603    fn points_are_close2d(p1: (f64, f64), p2: (f64, f64)) -> bool {
2604        let tolerance = 1e-10;
2605
2606        abs(p1.0 - p2.0) < tolerance && abs(p1.1 - p2.1) < tolerance
2607    }
2608
2609    fn points_are_close3d(p1: (f64, f64, f64), p2: (f64, f64, f64)) -> bool {
2610        let tolerance = 1e-10;
2611
2612        abs(p1.0 - p2.0) < tolerance && abs(p1.1 - p2.1) < tolerance && abs(p1.2 - p2.2) < tolerance
2613    }
2614
2615    #[test]
2616
2617    fn test_polygon_translation_2d() {
2618        let line = CoordSeq::<Point2d>::Linestring(vec![(0., 0.), (1.0, 1.0)])
2619            .to_sfcgal()
2620            .unwrap();
2621
2622        let line = line.translate_2d(2., 3.).unwrap();
2623
2624        let coords: CoordSeq<Point2d> = line.to_coordinates().unwrap();
2625
2626        match coords {
2627            CoordSeq::Linestring(vec) => {
2628                assert!(points_are_close2d(vec[0], (2.0, 3.0)));
2629
2630                assert!(points_are_close2d(vec[1], (3.0, 4.0)));
2631            }
2632            _ => panic!("Bad coordinates variant"),
2633        }
2634    }
2635
2636    #[test]
2637
2638    fn test_polygon_translation_3d() {
2639        let poly = CoordSeq::<Point3d>::Polygon(vec![vec![
2640            (0., 0., 0.),
2641            (1., 0., 0.),
2642            (1., 1., 0.),
2643            (0., 1., 0.),
2644            (0., 0., 0.),
2645        ]])
2646        .to_sfcgal()
2647        .unwrap();
2648
2649        let poly = poly.translate_3d(1.0, 2.0, 3.0).unwrap();
2650
2651        let coords: CoordSeq<Point3d> = poly.to_coordinates().unwrap();
2652
2653        match coords {
2654            CoordSeq::Polygon(vec) => {
2655                let poly_test = &vec[0];
2656
2657                assert!(points_are_close3d(poly_test[0], (1., 2., 3.)));
2658
2659                assert!(points_are_close3d(poly_test[2], (2., 3., 3.)));
2660            }
2661            _ => panic!("Bad coordinate variant"),
2662        }
2663    }
2664
2665    #[test]
2666
2667    fn test_rotate_point() {
2668        let point = CoordSeq::Point((1., 0., 0.)).to_sfcgal().unwrap();
2669
2670        let point = point.rotate_z(5. * PI / 2.).unwrap();
2671
2672        let coords: CoordSeq<Point3d> = point.to_coordinates().unwrap();
2673
2674        match coords {
2675            CoordSeq::Point(pt) => {
2676                assert!(points_are_close3d(pt, (0., 1., 0.)))
2677            }
2678            _ => panic!("Bad coordinate variant"),
2679        }
2680    }
2681
2682    #[test]
2683
2684    fn test_obj_export() {
2685        // ----------------------------------------
2686        // Open with Blender or any other 3d soft
2687        // to visually check
2688        // ----------------------------------------
2689        let temp_dir = env::temp_dir();
2690
2691        let final_path = format!("{}/sfcgal_test.obj", temp_dir.to_str().unwrap());
2692
2693        println!("Writing to {:?}", temp_dir);
2694
2695        let input = "MULTISOLID Z (((((0 0 0,0 1 0,1 1 0,1 0 0,0 0 0)),((0 0 1,1 0 1,1 1 1,0 1 1,0 0 1)),((0 0 0,1 0 0,1 0 1,0 0 1,0 0 0)),((1 1 0,0 1 0,0 1 1,1 1 1,1 1 0)),((1 0 0,1 1 0,1 1 1,1 0 1,1 0 0)),((0 0 0,0 0 1,0 1 1,0 1 0,0 0 0)))),((((2 4 6,2 5 6,3 5 6,3 4 6,2 4 6)),((2 4 7,3 4 7,3 5 7,2 5 7,2 4 7)),((2 4 6,3 4 6,3 4 7,2 4 7,2 4 6)),((3 5 6,2 5 6,2 5 7,3 5 7,3 5 6)),((3 4 6,3 5 6,3 5 7,3 4 7,3 4 6)),((2 4 6,2 4 7,2 5 7,2 5 6,2 4 6)))))";
2696
2697        let shape = SFCGeometry::new(input).unwrap();
2698
2699        assert!(shape.to_obj_file(final_path.as_str()).is_ok());
2700    }
2701
2702    #[test]
2703
2704    fn show_full_version() {
2705        println!("{:?}", _sfcgal_get_full_version());
2706    }
2707}