qhull_enhanced/
builder.rs

1use std::{cell::{RefCell, UnsafeCell}, ffi::{CString}, marker::PhantomData, ptr};
2
3use crate::{helpers::{collect_coords, CollectedCoords}, io_buffers::IOBuffers, log, sys, Qh, QhError};
4
5type QhConfigurator = Box<dyn for<'b> Fn(&'b mut Qh) -> Result<(), QhError<'b>> + 'static>;
6
7#[derive(Debug)]
8#[derive(thiserror::Error)]
9pub enum InvalidStringError {
10    #[error(transparent)]
11    NulError(#[from] std::ffi::NulError),
12    #[error("string is too long ({actual_length}), max length is {max_length} for {field}")]
13    TooLong {
14        field: &'static str,
15        max_length: usize,
16        actual_length: usize,
17    },
18}
19
20/// Builder for a Qhull instance
21///
22/// # Example
23/// ```
24/// # use qhull_enhanced::*;
25/// let mut points = [
26///    0.0, 0.0,
27///    1.0, 0.0,
28///    0.0, 1.0,
29///    0.25, 0.25,
30/// ];
31///
32/// let qh = QhBuilder::default()
33///     .build(2, &mut points)
34///     .unwrap();
35///
36/// assert_eq!(qh.num_facets(), 3);
37/// ```
38#[must_use]
39pub struct QhBuilder {
40    dim: Option<usize>,
41    capture_stdout: bool,
42    capture_stderr: bool,
43    compute: bool,
44    check_output: bool,
45    check_points: bool,
46    args: Vec<CString>,
47    configs: Vec<QhConfigurator>,
48}
49
50/// Default settings:
51/// * No [dimension hint](QhBuilder::dim)
52/// * [stdout](QhBuilder::capture_stdout) is not captured
53/// * [stderr](QhBuilder::capture_stderr) is captured
54/// * [compute](QhBuilder::compute) is `true`
55impl Default for QhBuilder {
56    fn default() -> Self {
57        let args = CString::new("qhull").unwrap();
58
59        Self {
60            dim: None,
61            capture_stdout: false,
62            capture_stderr: true,
63            compute: true,
64            check_output: false,
65            check_points: false,
66            args: vec![args],
67            configs: Vec::new(),
68        }
69    }
70}
71
72impl QhBuilder {
73    /// Add arguments for the qhull library
74    ///
75    /// These arguments are used by [`qhull_sys_enhanced::qh_init_A`] before setting any other option.
76    pub fn qhull_args<S: AsRef<str>>(mut self, args: impl IntoIterator<Item = S>) -> Result<Self, InvalidStringError> {
77        self.args.extend(
78            args
79                .into_iter()
80                .map(|arg| CString::new(arg.as_ref()))
81                .collect::<Result<Vec<_>, _>>()?
82        );
83        Ok(self)
84    }
85
86    /// Sets a dimension hint for the data
87    ///
88    /// This is useful when you want to **assert** that the data
89    /// has the correct dimensionality before building the Qhull instance.
90    ///
91    /// This build will panic if the dimensionality of the data does not match the hint.
92    pub fn dim(mut self, dim: usize) -> Self {
93        assert!(dim > 0, "dim must be > 0");
94        self.dim = Some(dim);
95        self
96    }
97
98    /// Capture stdout
99    ///
100    /// When enabled, the output of the qhull library will be captured instead of printed to the console.
101    pub fn capture_stdout(mut self, capture: bool) -> Self {
102        self.capture_stdout = capture;
103        self
104    }
105
106    /// Capture stderr
107    ///
108    /// When enabled, the error output of the qhull library will be captured instead of printed to the console.
109    pub fn capture_stderr(mut self, capture: bool) -> Self {
110        self.capture_stderr = capture;
111        self
112    }
113
114    /// Set whether to compute the hull when building the Qhull instance
115    ///
116    /// When enabled, [`Qh::compute`] and [`Qh::prepare_output`] will be called.
117    /// When disabled, you will have to call this method manually.
118    ///
119    /// # Example
120    /// ```
121    /// # use qhull_enhanced::*;
122    /// let points = [
123    ///     [0.0, 0.0],
124    ///     [1.0, 0.0],
125    ///     [0.0, 1.0],
126    ///     [0.25, 0.25],
127    /// ];
128    /// let qh = QhBuilder::default()
129    ///     .compute(true) // this is the default
130    ///     .build_from_iter(points)
131    ///     .unwrap();
132    /// assert_eq!(qh.num_facets(), 3);
133    ///
134    /// let mut qh = QhBuilder::default()
135    ///     .compute(false)
136    ///     .build_from_iter(points)
137    ///     .unwrap();
138    /// assert_eq!(qh.num_facets(), 0);
139    /// qh.compute().unwrap();
140    /// assert_eq!(qh.num_facets(), 3);
141    /// ```
142    pub fn compute(mut self, compute: bool) -> Self {
143        self.compute = compute;
144        self
145    }
146
147    /// Set whether to check the output when building the Qhull instance
148    ///
149    /// When enabled, [`Qh::check_output`] will be called after computing the hull.  
150    /// If [`compute`](QhBuilder::compute) is disabled, this setting will have no effect.
151    pub fn check_output(mut self, check: bool) -> Self {
152        self.check_output = check;
153        self
154    }
155
156    /// Set whether to check the points when building the Qhull instance
157    ///
158    /// When enabled, [`Qh::check_points`] will be called after computing the hull.
159    /// If [`compute`](QhBuilder::compute) is disabled, this setting will have no effect.
160    pub fn check_points(mut self, check: bool) -> Self {
161        self.check_points = check;
162        self
163    }
164
165    /// Build a Qhull instance
166    ///
167    /// # Example
168    /// ```
169    /// # use qhull_enhanced::*;
170    /// let mut points = [
171    ///    0.0, 0.0,
172    ///    1.0, 0.0,
173    ///    0.0, 1.0,
174    ///    0.25, 0.25,
175    /// ];
176    ///
177    /// let qh = QhBuilder::default()
178    ///     .build(2, &mut points).unwrap();
179    ///
180    /// assert_eq!(qh.num_facets(), 3);
181    /// ```
182    ///
183    /// # Panics
184    /// * If the number of points is not divisible by the dimension
185    /// * If the dimensionality of the points does not match the hint
186    /// * Cannot create a temporary file for capturing stdout or stderr
187    pub fn build(self, dim: usize, points: &mut [f64]) -> Result<Qh, QhError> {
188        if let Some(dim_hint) = self.dim {
189            assert_eq!(
190                dim, dim_hint,
191                "data dimensionality does not match hint that was given with QhBuilder::dim"
192            );
193        }
194
195        assert_eq!(points.len() % dim, 0, "points.len() % dim != 0");
196        let num_points = points.len() / dim;
197
198        unsafe {
199            let mut qh: sys::qhT = std::mem::zeroed();
200            let buffers = IOBuffers::new(self.capture_stdout, self.capture_stderr);
201
202            let argv = self
203                .args
204                .iter()
205                .map(|arg| arg.as_ptr())
206                .chain(std::iter::once(ptr::null()))
207                .collect::<Vec<_>>();
208            log!("argv {:?}", argv);
209            let argc = argv.len() as i32 - 1;
210            let argv = argv.as_ptr() as *mut _;
211
212            // Note: this function cannot be called
213            // inside of a try
214            sys::qh_init_A(
215                &mut qh,
216                buffers.in_file(),
217                buffers.out_file(),
218                buffers.err_file(),
219                argc,
220                argv,
221            );
222
223            let mut qh = Qh {
224                qh: UnsafeCell::new(qh),
225                coords_holder: None,
226                dim,
227                buffers: RefCell::new(buffers),
228                owned_values: Default::default(),
229                phantom: PhantomData,
230            };
231
232            // from `qconvex_r.c`
233            let hidden_options = CString::new(" d v H Qbb Qf Qg Qm Qr Qu Qv Qx Qz TR E V Fp Gt Q0 Q1 Q2 Q3 Q4 Q5 Q6 Q7 Q8 Q9 Q10 Q11 Q15 ").unwrap();
234
235            let qh_ptr = qh.qh.get();
236
237            QhError::try_3(
238                qh_ptr,
239                &mut qh.buffers.borrow_mut().err_file,
240                sys::qh_checkflags,
241                (
242                    qh_ptr,
243                    (*qh_ptr).qhull_command.as_ptr() as *mut _,
244                    hidden_options.as_ptr() as *mut _,
245                ),
246            )
247            .map_err(|e| e.into_static())?;
248            QhError::try_2(
249                qh_ptr,
250                &mut qh.buffers.borrow_mut().err_file,
251                sys::qh_initflags,
252                (
253                    qh_ptr,
254                    (*qh_ptr).qhull_command.as_ptr() as *mut _,
255                ),
256            ).map_err(|e| e.into_static())?;
257
258            for config in self.configs {
259                config(&mut qh).map_err(|e| e.into_static())?;
260            }
261
262            QhError::try_5(
263                qh_ptr,
264                &mut qh.buffers.borrow_mut().err_file,
265                sys::qh_init_B,
266                (
267                    qh_ptr,
268                    points.as_ptr() as *mut f64,
269                    num_points as _,
270                    dim as _,
271                    false as _,
272                ),
273            ).map_err(|e| e.into_static())?;
274
275            if self.compute {
276                qh.compute().map_err(|e| e.into_static())?;
277                if self.check_output {
278                    qh.check_output().map_err(|e| e.into_static())?;
279                }
280                qh.prepare_output().map_err(|e| e.into_static())?;
281                if self.check_points {
282                    qh.check_points().map_err(|e| e.into_static())?;
283                }
284            }
285
286            Ok(qh)
287        }
288    }
289
290    /// Build a Qhull instance with managed points
291    ///
292    /// This is useful when you want to keep the points alive
293    /// for the lifetime of the Qhull instance.
294    ///
295    /// # Example
296    /// ```
297    /// # use qhull_enhanced::*;
298    ///
299    /// let qh = QhBuilder::default()
300    ///     .build_managed(2, vec![
301    ///         0.0, 0.0,
302    ///         1.0, 0.0,
303    ///         0.0, 1.0,
304    ///         0.25, 0.25,
305    ///     ]).unwrap();
306    ///
307    /// assert_eq!(qh.num_facets(), 3);
308    /// ```
309    pub fn build_managed(
310        self,
311        dim: usize,
312        points: impl ToOwned<Owned = Vec<f64>>,
313    ) -> Result<Qh<'static>, QhError<'static>> {
314        let mut points = points.to_owned();
315        let points_ptr = points.as_mut_ptr();
316        let mut qh: Qh<'static> = self.build(dim, unsafe {
317            std::slice::from_raw_parts_mut(points_ptr, points.len())
318        })?;
319        assert!(qh.coords_holder.is_none());
320        qh.coords_holder = Some(points);
321        Ok(qh)
322    }
323
324    /// Build a Qhull instance from an iterator of points
325    ///
326    /// # Example
327    /// ```
328    /// # use qhull_enhanced::*;
329    /// let qh = QhBuilder::default()
330    ///     .build_from_iter([
331    ///         [0.0, 0.0],
332    ///         [1.0, 0.0],
333    ///         [0.0, 1.0],
334    ///         [0.25, 0.25],
335    ///     ]).unwrap();
336    ///
337    /// assert_eq!(qh.num_facets(), 3);
338    /// ```
339    pub fn build_from_iter<I>(
340        self,
341        points: impl IntoIterator<Item = I>,
342    ) -> Result<Qh<'static>, QhError<'static>>
343    where
344        I: IntoIterator<Item = f64>,
345    {
346        let CollectedCoords {
347            coords,
348            count: _,
349            dim,
350        } = collect_coords(points);
351        self.build_managed(dim, coords)
352    }
353
354    /// Configure the qhull instance with a closure
355    ///
356    /// # Safety
357    /// * closure must not panic
358    /// * closure shall not invalidate the qhull instance
359    /// * closure shall not keep references to the qhull instance
360    /// * closure shall not initialize the qhull instance
361    /// * closure shall not modify the error handling state of the qhull instance
362    ///
363    /// # Example
364    /// ```
365    /// # use qhull_enhanced::*;
366    /// let builder = unsafe {
367    ///     QhBuilder::default()
368    ///         .with_configure(|qh| {
369    ///             (*Qh::raw_ptr(qh)).DELAUNAY = true as _;
370    ///             Ok(())
371    ///         });
372    /// };
373    /// ```
374    pub unsafe fn with_configure(
375        mut self,
376        configurator: impl for<'a> Fn(&'a mut Qh) -> Result<(), QhError<'a>> + 'static,
377    ) -> Self {
378        self.configs.push(Box::new(configurator));
379        self
380    }
381
382    // TODO args and checkflags
383}
384
385// https://doc.rust-lang.org/book/ch03-02-data-types.html
386
387macro_rules! add_setting {
388    ($(
389        $(#[$meta:meta])*
390        $class:ident $($unsafe:ident)? ($($ty:tt)*) $setter:ident => $qhull_name:ident $($orig_doc:literal)?,
391    )*) => {
392        /// Raw settings
393        ///
394        /// Methods to set options that have a direct mapping to the qhull library.
395        ///
396        /// Some methods may invalidate the qhull instance if used incorrectly, they are marked as `unsafe`.
397        impl QhBuilder {
398            $(
399                add_setting! {
400                    $(#[$meta])*
401                    $class $($unsafe)? ($($ty)*) $setter => $qhull_name $($orig_doc)?
402                }
403            )*
404        }
405    };
406    (
407        basic documentation: $setter:ident => $qhull_name:ident: $(#[$meta:meta])* $($orig_doc:literal)?
408    ) => {
409        concat!(
410            "setter for [`", stringify!($qhull_name), "`](crate::sys::qhT::", stringify!($qhull_name), ")",
411            "\n\n",
412            $(
413                concat!("Original documentation:\n> <em>", $orig_doc, "</em>\n"),
414            )?
415            "\n\n",
416        )
417    };
418    (
419        safety documentation: unsafe
420    ) => {
421        "# Safety
422* This setting in unsafe because it can lead to undefined behavior if used incorrectly.
423
424"
425    };
426    (
427        safety documentation:
428    ) => { "\n\n" };
429    (
430        $(#[$meta:meta])*
431        scalar $($unsafe:ident)? ($ty:ident) $setter:ident => $qhull_name:ident $($orig_doc:literal)?
432    ) => {
433        #[doc = add_setting!(basic documentation: $setter => $qhull_name: $($orig_doc)?)]
434        $(#[$meta])*
435        #[doc = add_setting!(safety documentation: $($unsafe)?)]
436        pub $($unsafe)? fn $setter(mut self, $setter: type_mapping::$ty) -> Self {
437            self = unsafe {
438                self.with_configure(move |qh| {
439                    (qh.qh.get_mut()).$qhull_name = $setter as _;
440                    Ok(())
441                })
442            };
443            self
444        }
445    };
446    (
447        $(#[$meta:meta])*
448        array $($unsafe:ident)? ($ty:ident[$N:expr]) $setter:ident => $qhull_name:ident $($orig_doc:literal)?
449    ) => {
450        #[doc = add_setting!(basic documentation: $setter => $qhull_name: $($orig_doc)?)]
451        $(#[$meta])*
452        #[doc = add_setting!(safety documentation: unsafe)]
453        pub unsafe fn $setter(mut self, $setter: [type_mapping::$ty; $N]) -> Self {
454            self = unsafe {
455                self.with_configure(move |qh| {
456                    (qh.qh.get_mut()).$qhull_name = $setter;
457                    Ok(())
458                })
459            };
460            self
461        }
462    };
463    (
464        $(#[$meta:meta])*
465        array_string $($unsafe:ident)? (char[$N:expr]) $setter:ident => $qhull_name:ident $($orig_doc:literal)?
466    ) => {
467        #[doc = add_setting!(basic documentation: $setter => $qhull_name: $($orig_doc)?)]
468        $(#[$meta])*
469        #[doc = add_setting!(safety documentation: unsafe)]
470        pub fn $setter(mut self, $setter: &str) -> Result<Self, InvalidStringError> { // or maybe QhError or something else?
471            let bytes = std::ffi::CString::new($setter)?
472                .as_bytes_with_nul()
473                .iter()
474                .map(|&b| b as core::ffi::c_char)
475                .collect::<Vec<_>>();
476            if bytes.len() > $N {
477                return Err(InvalidStringError::TooLong {
478                    field: stringify!($setter),
479                    max_length: $N,
480                    actual_length: bytes.len(),
481                });
482            }
483            self = unsafe {
484                self.with_configure(move |qh| {
485                    (qh.qh.get_mut()).$qhull_name = [0; $N];
486                    (qh.qh.get_mut()).$qhull_name[..bytes.len()].copy_from_slice(&bytes[..]);
487                    Ok(())
488                })
489            };
490            Ok(self)
491        }
492    };
493    (
494        $(#[$meta:meta])*
495        dyn_string $($unsafe:ident)? (char*) $setter:ident => $qhull_name:ident $($orig_doc:literal)?
496    ) => {
497        #[doc = add_setting!(basic documentation: $setter => $qhull_name: $($orig_doc)?)]
498        $(#[$meta])*
499        #[doc = add_setting!(safety documentation: unsafe)]
500        pub fn $setter(mut self, $setter: &str) -> Result<Self, InvalidStringError> { // or maybe QhError or something else?
501            let $setter = std::ffi::CString::new($setter)?;
502            self = unsafe {
503                self.with_configure(move |qh| {
504                    let ptr = $setter.as_ptr();
505                    qh.owned_values.$setter = Some($setter.clone());
506                    (qh.qh.get_mut()).$qhull_name = ptr as *mut _;
507                    Ok(())
508                })
509            };
510            Ok(self)
511        }
512    };
513    (
514        $(#[$meta:meta])*
515        point $($unsafe:ident)? ($ty:ident*) $setter:ident => $qhull_name:ident $($orig_doc:literal)?
516    ) => {
517        #[doc = add_setting!(basic documentation: $setter => $qhull_name: $($orig_doc)?)]
518        $(#[$meta])*
519        #[doc = add_setting!(safety documentation: $($unsafe)?)]
520        pub $($unsafe)? fn $setter(mut self, $setter: impl IntoIterator<Item = type_mapping::$ty>) -> Self {
521            let dim = self.dim.expect(concat!("dimension hint is required for ", stringify!($setter), " setter"));
522            let $setter = $setter.into_iter().collect::<Vec<_>>();
523            assert_eq!($setter.len() % dim, 0, concat!("number of elements in ", stringify!($setter), " must be divisible by dim"));
524            self = unsafe {
525                self.with_configure(move |qh| {
526                    let ptr = $setter.as_ptr();
527                    qh.owned_values.$setter = Some($setter.clone());
528                    (qh.qh.get_mut()).$qhull_name = ptr as *mut _;
529                    Ok(())
530                })
531            };
532            self
533        }
534    };
535}
536
537add_setting!(
538    scalar unsafe(boolT)  all_points => ALLpoints "true 'Qs' if search all points for initial simplex",
539    scalar(boolT)  allow_short => ALLOWshort "true 'Qa' allow input with fewer or more points than coordinates",
540    scalar(boolT)  allow_warning => ALLOWwarning "true 'Qw' if allow option warnings",
541    scalar(boolT)  allow_wide => ALLOWwide "true 'Q12' if allow wide facets and wide dupridges, c.f. qh_WIDEmaxoutside",
542    scalar(boolT)  angle_merge => ANGLEmerge "true 'Q1' if sort potential merges by type/angle instead of type/distance ",
543    scalar(boolT)  approx_hull => APPROXhull "true 'Wn' if MINoutside set",
544    scalar(realT)  min_outside => MINoutside "  Minimum distance for an outside point ('Wn' or 2*qh.MINvisible)",
545    scalar(boolT)  annotate_output => ANNOTATEoutput "true 'Ta' if annotate output with message codes",
546    scalar(boolT)  at_infinity => ATinfinity "true 'Qz' if point num_points-1 is \"at-infinity\" for improving precision in Delaunay triangulations",
547    scalar(boolT)  avoid_old => AVOIDold "true 'Q4' if avoid old->new merges",
548    scalar(boolT)  best_outside => BESToutside "true 'Qf' if partition points into best outsideset",
549    scalar(boolT)  cdd_input => CDDinput "true 'Pc' if input uses CDD format (1.0/offset first)",
550    scalar(boolT)  cdd_output => CDDoutput "true 'PC' if print normals in CDD format (offset first)",
551    scalar(boolT)  check_duplicates => CHECKduplicates "true 'Q15' if qh_maybe_duplicateridges after each qh_mergefacet",
552    scalar(boolT)  check_frequently => CHECKfrequently "true 'Tc' if checking frequently",
553    scalar(realT)  premerge_cos => premerge_cos "'A-n'   cos_max when pre merging",
554    scalar(realT)  postmerge_cos => postmerge_cos "'An'    cos_max when post merging",
555    scalar(boolT)  delaunay => DELAUNAY "true 'd' or 'v' if computing DELAUNAY triangulation",
556    scalar(boolT)  do_intersections => DOintersections "true 'Gh' if print hyperplane intersections",
557    scalar(int)    drop_dim => DROPdim "drops dim 'GDn' for 4-d -> 3-d output",
558    scalar(boolT)  flush_print => FLUSHprint "true 'Tf' if flush after qh_fprintf for segfaults",
559    scalar(boolT)  force_output => FORCEoutput "true 'Po' if forcing output despite degeneracies",
560    scalar(int)    good_point => GOODpoint "'QGn' or 'QG-n' (n+1, n-1), good facet if visible from point n (or not)",
561    point(pointT*) good_point_coords => GOODpointp "the actual point",
562    scalar(boolT)  good_threshold => GOODthreshold "true 'Pd/PD' if qh.lower_threshold/upper_threshold defined set if qh.UPPERdelaunay (qh_initbuild) false if qh.SPLITthreshold",
563    scalar(int)    good_vertex => GOODvertex "'QVn' or 'QV-n' (n+1, n-1), good facet if vertex for point n (or not)",
564    point(pointT*) good_vertex_coords => GOODvertexp "the actual point",
565    scalar(boolT) half_space => HALFspace "true 'Hn,n,n' if halfspace intersection",
566    // scalar(boolT) is_qhull_qh => ISqhullQh "Set by Qhull.cpp on initialization",
567    scalar(int)  is_tracing => IStracing "'Tn' trace execution, 0=none, 1=least, 4=most, -1=events",
568    scalar(int)  keep_area => KEEParea "'PAn' number of largest facets to keep",
569    scalar(boolT) keep_coplanar => KEEPcoplanar "true 'Qc' if keeping nearest facet for coplanar points",
570    scalar(boolT) keep_inside => KEEPinside "true 'Qi' if keeping nearest facet for inside points set automatically if 'd Qc'",
571    scalar(int)   keep_merge => KEEPmerge "'PMn' number of facets to keep with most merges ",
572    scalar(realT) keep_min_area => KEEPminArea "'PFn' minimum facet area to keep ",
573    scalar(realT) max_coplanar => MAXcoplanar "'Un' max distance below a facet to be coplanar",
574    scalar(int)   max_wide => MAXwide "'QWn' max ratio for wide facet, otherwise error unless Q12-allow-wide ",
575    scalar(boolT) merge_exact => MERGEexact "true 'Qx' if exact merges (concave, degen, dupridge, flipped) tested by qh_checkzero and qh_test_*_merge",
576    scalar(boolT) merge_independent => MERGEindependent "true if merging independent sets of coplanar facets. 'Q2' disables ",
577    scalar(boolT) merging => MERGING "true if exact-, pre- or post-merging, with angle and centrum tests ",
578    scalar(realT) premerge_centrum => premerge_centrum "  'C-n' centrum_radius when pre merging.  Default is round-off ",
579    scalar(realT) postmerge_centrum => postmerge_centrum "  'Cn' centrum_radius when post merging.  Default is round-off ",
580    scalar(boolT) merge_pinched => MERGEpinched "true 'Q14' if merging pinched vertices due to dupridge ",
581    scalar(boolT) merge_vertices => MERGEvertices "true if merging redundant vertices, 'Q3' disables or qh.hull_dim > qh_DIMmergeVertex ",
582    scalar(realT) min_visible => MINvisible "'Vn' min. distance for a facet to be visible ",
583    scalar(boolT) no_narrow => NOnarrow "true 'Q10' if no special processing for narrow distributions ",
584    scalar(boolT) no_near_inside => NOnearinside "true 'Q8' if ignore near-inside points when partitioning, qh_check_points may fail ",
585    scalar(boolT) no_premerge => NOpremerge "true 'Q0' if no defaults for C-0 or Qx ",
586    scalar(boolT) only_good => ONLYgood "true 'Qg' if process points with good visible or horizon facets ",
587    scalar(boolT) only_max => ONLYmax "true 'Qm' if only process points that increase max_outside ",
588    scalar(boolT) pick_furthest => PICKfurthest "true 'Q9' if process furthest of furthest points",
589    scalar(boolT) post_merge => POSTmerge "true if merging after buildhull ('Cn' or 'An') ",
590    scalar(boolT) pre_merge => PREmerge "true if merging during buildhull ('C-n' or 'A-n') ",
591    scalar(boolT) print_centrums => PRINTcentrums "true 'Gc' if printing centrums",
592    scalar(boolT) print_coplanar => PRINTcoplanar "true 'Gp' if printing coplanar points",
593    scalar(int)   print_dim => PRINTdim "print dimension for Geomview output",
594    scalar(boolT) print_dots => PRINTdots "true 'Ga' if printing all points as dots",
595    scalar(boolT) print_good => PRINTgood "true 'Pg' if printing good facets PGood set if 'd', 'PAn', 'PFn', 'PMn', 'QGn', 'QG-n', 'QVn', or 'QV-n'",
596    scalar(boolT) print_inner => PRINTinner "true 'Gi' if printing inner planes",
597    scalar(boolT) print_neighbors => PRINTneighbors "true 'PG' if printing neighbors of good facets",
598    scalar(boolT) print_no_planes => PRINTnoplanes "true 'Gn' if printing no planes",
599    scalar(boolT) print_options_1st => PRINToptions1st "true 'FO' if printing options to stderr",
600    scalar(boolT) print_outer => PRINTouter "true 'Go' if printing outer planes",
601    scalar(boolT) print_precision => PRINTprecision "false 'Pp' if not reporting precision problems",
602    array(qh_PRINT[sys::qh_PRINT_qh_PRINTEND as usize]) print_out => PRINTout "list of output formats to print",
603    scalar(boolT) print_ridges => PRINTridges "true 'Gr' if print ridges",
604    scalar(boolT) print_spheres => PRINTspheres "true 'Gv' if print vertices as spheres",
605    scalar(boolT) print_statistics => PRINTstatistics "true 'Ts' if printing statistics to stderr",
606    scalar(boolT) print_summary => PRINTsummary "true 's' if printing summary to stderr",
607    scalar(boolT) print_transparent => PRINTtransparent "true 'Gt' if print transparent outer ridges",
608    scalar(boolT) project_delaunay => PROJECTdelaunay "true if DELAUNAY, no readpoints() and need projectinput() for Delaunay in qh_init_B",
609    scalar(int)   project_input => PROJECTinput "number of projected dimensions 'bn:0Bn:0'",
610    scalar(boolT) random_dist => RANDOMdist "true 'Rn' if randomly change distplane and setfacetplane",
611    scalar(realT) random_factor => RANDOMfactor "maximum random perturbation",
612    scalar(realT) random_a => RANDOMa "qh_randomfactor is randr * RANDOMa + RANDOMb",
613    scalar(realT) random_b => RANDOMb "boolT) RANDOMoutside \"true 'Qr' if select a random outside point",
614    scalar(int)   report_freq => REPORTfreq "TFn' buildtracing reports every n facets",
615    scalar(int)   report_freq_2 => REPORTfreq2 "tracemerging reports every REPORTfreq/2 facets",
616    scalar(int)   rerun => RERUN "TRn' rerun qhull n times (qh.build_cnt)",
617    scalar(int)   rotate_random => ROTATErandom "QRn' n<-1 random seed, n==-1 time is seed, n==0 random rotation by time, n>0 rotate input",
618    scalar(boolT) scale_input => SCALEinput "true 'Qbk' if scaling input",
619    scalar(boolT) scale_last => SCALElast "true 'Qbb' if scale last coord to max prev coord",
620    scalar(boolT) set_roundoff => SETroundoff "true 'En' if qh.DISTround is predefined",
621    scalar(boolT) skip_check_max => SKIPcheckmax "true 'Q5' if skip qh_check_maxout, qh_check_points may fail",
622    scalar(boolT) skip_convex => SKIPconvex "true 'Q6' if skip convexity testing during pre-merge",
623    scalar(boolT) split_thresholds => SPLITthresholds "true 'Pd/PD' if upper_/lower_threshold defines a region else qh.GOODthresholds set if qh.DELAUNAY (qh_initbuild) used  only for printing (!for qh.ONLYgood)",
624    scalar(int)   stop_add => STOPadd "'TAn' 1+n for stop after adding n vertices",
625    scalar(int)   stop_clone => STOPcone "'TCn' 1+n for stopping after cone for point n also used by qh_build_withresart for err exi",
626    scalar(int)   stop_point => STOPpoint "'TVn' 'TV-n' 1+n for stopping after/before(-) adding point n",
627    scalar(int)   test_point => TESTpoints "'QTn' num of test points after qh.num_points.  Test points always coplanar.",
628    scalar(boolT) test_v_neighbors => TESTvneighbors " true 'Qv' if test vertex neighbors at end",
629    scalar(int)   trace_level => TRACElevel "'Tn' conditional IStracing level",
630    scalar(int)   trace_last_run => TRACElastrun " qh.TRACElevel applies to last qh.RERUN",
631    scalar(int)   trace_point => TRACEpoint "'TPn' start tracing when point n is a vertex, use qh_IDunknown (-1) after qh_buildhull and qh_postmerge",
632    scalar(realT) trace_dist => TRACEdist "'TWn' start tracing when merge distance too big",
633    scalar(int)   trace_merge => TRACEmerge "'TMn' start tracing before this merge",
634    scalar(boolT) triangulate => TRIangulate "true 'Qt' if triangulate non-simplicial facets",
635    scalar(boolT) tri_normals => TRInormals "true 'Q11' if triangulate duplicates ->normal and ->center (sets Qt)",
636    scalar(boolT) upper_delaunay => UPPERdelaunay "true 'Qu' if computing furthest-site Delaunay",
637    scalar(boolT) use_stdout => USEstdout "true 'Tz' if using stdout instead of stderr",
638    scalar(boolT) verify_output => VERIFYoutput "true 'Tv' if verify output at end of qhull",
639    scalar(boolT) virtual_memory => VIRTUALmemory "true 'Q7' if depth-first processing in buildhull",
640    scalar(boolT) voronoi => VORONOI "true 'v' if computing Voronoi diagram, also sets qh.DELAUNAY",
641
642  /*--------input constants ---------*/
643    scalar(realT) area_factor => AREAfactor "1/(hull_dim-1)! for converting det's to area",
644    scalar(boolT) do_check_max => DOcheckmax "true if calling qh_check_maxout (!qh.SKIPcheckmax && qh.MERGING)",
645    dyn_string(char*) feasible_string => feasible_string "feasible point 'Hn,n,n' for halfspace intersection",
646    point(coordT*) feasible_point => feasible_point "   as coordinates, both malloc'd",
647    scalar(boolT) get_area => GETarea "true 'Fa', 'FA', 'FS', 'PAn', 'PFn' if compute facet area/Voronoi volume in io_r.c",
648    scalar(boolT) keep_near_inside => KEEPnearinside "true if near-inside points in coplanarset",
649    //scalar(int)   hull_dim => hull_dim "dimension of hull, set by initbuffers",
650    //scalar(int)   input_dim "dimension of input, set by initbuffers",
651    scalar(int)   num_points => num_points "number of input points",
652    // point(pointT*) first_point => first_point "array of input points, see POINTSmalloc",
653    scalar(boolT) points_malloc => POINTSmalloc "  true if qh.first_point/num_points allocated",
654    // TODO ???(pointT*) input_points => input_points "copy of original qh.first_point for input points for qh_joggleinput",
655    scalar(boolT) input_malloc => input_malloc "true if qh.input_points malloc'd",
656    array_string(char[256])  qhull_command => qhull_command "command line that invoked this program",
657    scalar(int)   qhull_command_size_2 => qhull_commandsiz2 "size of qhull_command at qh_clear_outputflags",
658    array_string(char[256]) rbox_command => rbox_command "command line that produced the input points",
659    array_string(char[512]) qhull_options => qhull_options "descriptive list of options",
660    scalar(int)  qhull_option_len => qhull_optionlen "length of last line",
661    scalar(int)  qhull_option_size => qhull_optionsiz "size of qhull_options at qh_build_withrestart",
662    scalar(int)  qhull_option_size_2 => qhull_optionsiz2 "size of qhull_options at qh_clear_outputflags",
663    scalar(int)   run_id => run_id "non-zero, random identifier for this instance of qhull",
664    scalar(boolT) vertex_neighbors => VERTEXneighbors "true if maintaining vertex neighbors",
665    scalar(boolT) zero_centrum => ZEROcentrum "true if 'C-0' or 'C-0 Qx' and not post-merging or 'A-n'.  Sets ZEROall_ok",
666    point(realT*) upper_threshold => upper_threshold "don't print if facet->normal\\[k\\]>=upper_threshold\\[k\\] must set either GOODthreshold or SPLITthreshold if qh.DELAUNAY, default is 0.0 for upper envelope (qh_initbuild)",
667    point(realT*) lower_threshold => lower_threshold "don't print if facet->normal\\[k\\] <=lower_threshold\\[k\\]",
668    point(realT*) upper_bound => upper_bound "scale point\\[k\\] to new upper bound",
669    point(realT*) lower_bound => lower_bound "scale point\\[k\\] to new lower bound project if both upper_ and lower_bound == 0",
670
671    /* precision constants */
672    scalar(realT) angle_round => ANGLEround "max round off error for angles",
673    scalar(realT) centrum_radius => centrum_radius "max centrum radius for convexity ('Cn' + 2*qh.DISTround)",
674    scalar(realT) cos_max => cos_max "max cosine for convexity (roundoff added)",
675    scalar(realT) dist_round => DISTround "max round off error for distances, qh.SETroundoff ('En') overrides qh_distround",
676    scalar(realT) max_abs_coors => MAXabs_coord "max absolute coordinate",
677    scalar(realT) max_last_coord => MAXlastcoord "max last coordinate for qh_scalelast",
678    scalar(realT) max_outside => MAXoutside "max target for qh.max_outside/f.maxoutside, base for qh_RATIO... recomputed at qh_addpoint, unrelated to qh_MAXoutside",
679    scalar(realT) max_sum_coord => MAXsumcoord "max sum of coordinates",
680    scalar(realT) max_width => MAXwidth "max rectilinear width of point coordinates",
681    scalar(realT) min_denom_1 => MINdenom_1 "min. abs. value for 1/x",
682    scalar(realT) min_denom => MINdenom "   use divzero if denominator < MINdenom",
683    scalar(realT) min_denom_1_2 => MINdenom_1_2 "min. abs. val for 1/x that allows normalization",
684    scalar(realT) min_denom_2 => MINdenom_2 "   use divzero if denominator < MINdenom_2",
685    scalar(realT) min_last_coord => MINlastcoord "min. last coordinate for qh_scalelast",
686    point(realT*) near_zero => NEARzero "hull_dim array for near zero in gausselim",
687    scalar(realT) near_inside => NEARinside "keep points for qh_check_maxout if close to facet",
688    scalar(realT) one_merge => ONEmerge "max distance for merging simplicial facets",
689    scalar(realT) outside_err => outside_err "application's epsilon for coplanar points qh_check_bestdist() qh_check_points() reports error if point outside",
690    scalar(realT) wide_face => WIDEfacet "size of wide facet for skipping ridge in area computation and locking centrum",
691    scalar(boolT) narrow_hull => NARROWhull "set in qh_initialhull if angle < qh_MAXnarrow", /**/
692);
693
694#[allow(non_camel_case_types)]
695mod type_mapping {
696    use crate::sys;
697    pub type boolT = bool;
698    pub type realT = f64;
699    pub type int = i32;
700    pub type pointT = f64;
701    pub type coordT = f64;
702    pub type qh_PRINT = sys::qh_PRINT;
703}
704
705#[cfg(test)]
706mod tests {
707    use crate::Qh;
708
709    // from https://github.com/LucaCiucci/qhull-rs/issues/15
710    const POINTS_WITH_COPLANAR_SUBSET: &[[f64; 3]] = &[
711        [10.0, -1000.0, 0.0],
712        [-838.3645082073471, -1000.0, 0.0],
713        [0.0, -1000.0, 0.5],
714        [0.0, -1000.0, -1.0],
715        [0.0, 0.0, 0.0],
716    ];
717
718    #[test]
719    fn triangulate() {
720        let qh = Qh::builder()
721            .compute(true)
722            .triangulate(false)
723            .build_from_iter(POINTS_WITH_COPLANAR_SUBSET.iter().cloned())
724            .expect("Failed to compute convex hull");
725
726        let faces = qh
727            .facets()
728            .map(|face| face
729                .vertices()
730                .unwrap()
731                .iter()
732                .map(|v| v.point_id(&qh).unwrap())
733                .collect::<Vec<_>>()
734            )
735            .collect::<Vec<_>>();
736
737        assert_eq!(
738            faces,
739            vec![
740                vec![3, 4, 1],
741                vec![3, 4, 0],
742                vec![2, 4, 1],
743                vec![2, 4, 0],
744                vec![2, 3, 0, 1],
745            ]
746        );
747
748        let qh = Qh::builder()
749            .compute(true)
750            .triangulate(true)
751            .build_from_iter(POINTS_WITH_COPLANAR_SUBSET.iter().cloned())
752            .expect("Failed to compute convex hull");
753
754        let faces = qh
755            .facets()
756            .map(|face| face
757                .vertices()
758                .unwrap()
759                .iter()
760                .map(|v| v.point_id(&qh).unwrap())
761                .collect::<Vec<_>>()
762            )
763            .collect::<Vec<_>>();
764
765        assert_eq!(
766            faces,
767            vec![
768                vec![3, 4, 1],
769                vec![3, 4, 0],
770                vec![2, 4, 1],
771                vec![2, 4, 0],
772                vec![2, 3, 0],
773                vec![2, 3, 1],
774            ]
775        );
776    }
777
778    #[test]
779    fn triangulate_with_args() {
780        let qh = Qh::builder()
781            .compute(true)
782            .build_from_iter(POINTS_WITH_COPLANAR_SUBSET.iter().cloned())
783            .expect("Failed to compute convex hull");
784
785        let faces = qh
786            .facets()
787            .map(|face| face
788                .vertices()
789                .unwrap()
790                .iter()
791                .map(|v| v.point_id(&qh).unwrap())
792                .collect::<Vec<_>>()
793            )
794            .collect::<Vec<_>>();
795
796        assert_eq!(
797            faces,
798            vec![
799                vec![3, 4, 1],
800                vec![3, 4, 0],
801                vec![2, 4, 1],
802                vec![2, 4, 0],
803                vec![2, 3, 0, 1],
804            ]
805        );
806
807        let qh = Qh::builder()
808            .compute(true)
809            .qhull_args(&["Qt"]).unwrap()
810            .build_from_iter(POINTS_WITH_COPLANAR_SUBSET.iter().cloned())
811            .expect("Failed to compute convex hull");
812
813        let faces = qh
814            .facets()
815            .map(|face| face
816                .vertices()
817                .unwrap()
818                .iter()
819                .map(|v| v.point_id(&qh).unwrap())
820                .collect::<Vec<_>>()
821            )
822            .collect::<Vec<_>>();
823
824        assert_eq!(
825            faces,
826            vec![
827                vec![3, 4, 1],
828                vec![3, 4, 0],
829                vec![2, 4, 1],
830                vec![2, 4, 0],
831                vec![2, 3, 0],
832                vec![2, 3, 1],
833            ]
834        );
835    }
836}