1
  2
  3
  4
  5
  6
  7
  8
  9
 10
 11
 12
 13
 14
 15
 16
 17
 18
 19
 20
 21
 22
 23
 24
 25
 26
 27
 28
 29
 30
 31
 32
 33
 34
 35
 36
 37
 38
 39
 40
 41
 42
 43
 44
 45
 46
 47
 48
 49
 50
 51
 52
 53
 54
 55
 56
 57
 58
 59
 60
 61
 62
 63
 64
 65
 66
 67
 68
 69
 70
 71
 72
 73
 74
 75
 76
 77
 78
 79
 80
 81
 82
 83
 84
 85
 86
 87
 88
 89
 90
 91
 92
 93
 94
 95
 96
 97
 98
 99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
//! # geometry-predicates
//!
//! A safe Rust port of ["Adaptive Precision Floating-Point Arithmetic and Fast Robust
//! Predicates for Computational Geometry"](https://www.cs.cmu.edu/~quake/robust.html)
//!
//! This crate provides a Rust solution to efficient exact geometry predicates
//! used widely for computational geometry.
//!
//! In addition, the building blocks of these predicates, namely the adaptive precision
//! floating-point arithmetic primitives, are also exposed in [`predicates`] to allow for extensions
//! to other predicates or exact geometric constructions.
//!
//! ## Background
//!
//! These predicates have been a staple in computational geometry for many years now
//! and are widely used in industry.   In the context of geometry algorithms, it is
//! often essential to determine the orientation of a point with respect to a line (or a
//! plane) and whether a point lies inside a circle (or a sphere) or not.  The reason
//! why these tests often need to be exact is because many geometry algorithms
//! ask questions (to determine orientation or in-circle/sphere) about point
//! configurations that require consistent answers. For instance, if `a`, `b`, and
//! `c` are three points on a 2D plane, to ask where `b` with respect to the line
//! through `a` and `c` (left-of, right-of, or coincident) is the same as asking where
//! `a` lies with respect to the line through `c` and `b`.
//! Formally this condition can be written as
//! ```ignore
//! sgn(orient2d(a,c,b)) == sgn(orient2d(c,b,a))
//! ```
//!
//! Mathematically (using MATLAB-style notation), predicates like `orient2d` are
//! defined as
//! ```ignore
//! orient2d([ax,ay], [bx,by], [cx,cy]) := det([ax ay 1; bx by 1; cx cy 1])
//! ```
//!
//! It's easy to see that these predicates solve the problem of
//! computing the determinant of small matrices with the correct sign, regardless of how
//! close the matrix is to being singular.
//!
//! For instance to compute the determinant of a matrix `[a b; c d]` with the
//! correct sign, we can invoke
//! ```ignore
//! orient2d([a,b], [c,d], [0,0])
//! ```
//!
//! For more details please refer to the [original
//! webpage](https://www.cs.cmu.edu/~quake/robust.html) for these predicates.
#![no_std]
pub mod predicates;

#[cfg(feature = "transpiled")]
mod transpiled;

// The following predicates are exposed at the top level.
// There are other alternative robust implementations (but typically slower) of these predicates.
// We use those to check the adaptive implementations.
pub use predicates::{
    incircle,
    incircle_fast,
    insphere,
    insphere_fast,
    // Adaptive robust predicates.
    orient2d,
    // Fast inexact predicates.
    orient2d_fast,
    orient3d,
    orient3d_fast,
};

#[cfg(test)]
mod tests {
    extern crate rand;
    use self::rand::{Rng, SeedableRng, StdRng};
    use super::*;

    use predicates::{
        incircle_exact, incircle_slow, insphere_exact, insphere_slow, orient2d_exact,
        orient2d_slow, orient3d_exact, orient3d_slow,
    };

    const SEED: &[usize] = &[1, 2, 3, 4];

    /*
     * Note on robustness testing.
     *
     * These predicates do NOT handle overflow or underflowo of the exponent.
     * Quoting Shechuk's predicates paper directly:
     *
     * > The four predicates implemented [in this library] will not overflow nor underflow if their
     * > inputs have exponents in the range [-142, 201] and IEEE 754 double precision arithmetic
     * > is used.
     * - Jonathan Shechuk 1997
     *
     * We will therefore be careful to test for inputs with exponents in that range and not beyond.
     */

    const EXP_BOUNDS: [i32; 2] = [-142, 201];

    /// Generate a tolerance value with exponent no less than -142.
    fn tol(rng: &mut StdRng) -> f64 {
        let max_exp = (EXP_BOUNDS[0] + 1) as f64;
        10.0_f64.powi((max_exp * rng.gen::<f64>()).round() as i32) * (rng.gen::<f64>() - 0.5)
    }

    /*
     * Many of the tests below ensure that the predicates produce expected results for all
     * permutations of the inputs. Thus, below we write an adhoc QuickPerm implementation to
     * generate swap based permutations.
     *
     * The following is a basic countdown implementation of QuickPerm (see quickperm.org).
     * This implementation can be refactored to use const generics when those stabilize.
     */
    macro_rules! quick_perm_impl {
        ($n:expr, $struct_name:ident) => {
            struct $struct_name<T> {
                perm: [T; $n],
                p: [usize; $n + 1],
                index: usize,
            }

            impl<T> $struct_name<T> {
                fn new(v: [T; $n]) -> Self {
                    let mut p = [0; $n + 1];
                    for i in 0..$n + 1 {
                        p[i] = i;
                    }
                    $struct_name {
                        perm: v,
                        p,
                        index: 0,
                    }
                }
            }

            impl<T: Clone> Iterator for $struct_name<T> {
                type Item = [T; $n];
                fn next(&mut self) -> Option<Self::Item> {
                    let $struct_name { perm, p, index } = self;
                    let mut i = *index;
                    let n = p.len() - 1;
                    if i == 0 {
                        *index += 1;
                        return Some(perm.clone());
                    } else if i >= n {
                        return None;
                    }
                    p[i] -= 1;
                    let j = if i % 2 == 0 { 0 } else { p[i] };
                    perm.swap(j, i);
                    i = 1;
                    while p[i] == 0 {
                        p[i] = i;
                        i += 1;
                    }
                    *index = i;
                    Some(perm.clone())
                }
            }
        };
    }

    quick_perm_impl!(3, QuickPerm3);
    quick_perm_impl!(4, QuickPerm4);
    quick_perm_impl!(5, QuickPerm5);

    /*
     * The tests below have specific purposes.
     *
     * - _robustness_ tests check that the adaptive predicates work in degenerate or close to degenerate
     *   configurations
     * - _regression_ tests verify the adpative predicates against their _slow and _exact variants.
     * - _transpiled_regression_ tests verify the predicates against the directly transpiled
     *   (unrefactored) version of the library.
     */

    #[test]
    fn orient2d_test() {
        let a = [0.0, 1.0];
        let b = [2.0, 3.0];
        let c = [4.0, 5.0];
        assert_eq!(orient2d(a, b, c), 0.0);
    }

    #[cfg(feature = "transpiled")]
    #[test]
    fn orient2d_transpiled_regression_test() {
        unsafe {
            transpiled::exactinit();
        }
        let mut rng: StdRng = SeedableRng::from_seed(SEED);

        let n = 99999;
        for _ in 0..n {
            let a = [tol(&mut rng), tol(&mut rng)];
            let b = [12.0, 12.0];
            let c = [24.0, 24.0];
            let o2d = orient2d;
            let o2d_transpiled = transpiled::orient2d;
            for p in QuickPerm3::new([a, b, c]) {
                assert_eq!(
                    o2d(p[0], p[1], p[2]),
                    o2d_transpiled(p[0], p[1], p[2]),
                    "{:?}",
                    a
                );
            }
        }
    }

    #[test]
    fn orient2d_robustness_test() {
        let mut rng: StdRng = SeedableRng::from_seed(SEED);
        let n = 99999;

        for o2d in &[orient2d, orient2d_exact, orient2d_slow] {
            for _ in 0..n {
                let pa = [tol(&mut rng), tol(&mut rng)];
                let pb = [12.0, 12.0];
                let pc = [24.0, 24.0];

                let main = o2d(pa, pb, pc);
                if main == 0.0 {
                    for [a, b, c] in QuickPerm3::new([pa, pb, pc]) {
                        assert_eq!(o2d(a, b, c), 0.0);
                    }
                }

                let pred = main > 0.0;
                for (i, [a, b, c]) in QuickPerm3::new([pa, pb, pc]).enumerate() {
                    let t = o2d(a, b, c);
                    assert_eq!(pred, if i % 2 == 0 { t > 0.0 } else { t < 0.0 });
                }
            }
        }
    }

    // The following test verifies equivalence of all of the robust orient2d variants (including
    // exact and slow variants).
    #[test]
    fn orient2d_regression_test() {
        let mut rng: StdRng = SeedableRng::from_seed(SEED);

        let n = 99999;
        for _ in 0..n {
            let pa = [tol(&mut rng), tol(&mut rng)];
            let pb = [12.0, 12.0];
            let pc = [24.0, 24.0];

            let o2d = predicates::orient2d;
            let o2de = predicates::orient2d_exact;
            let o2ds = predicates::orient2d_slow;

            // Test all permutations.

            for [a, b, c] in QuickPerm3::new([pa, pb, pc]) {
                assert_eq!(
                    o2d(a, b, c).partial_cmp(&0.0),
                    o2de(a, b, c).partial_cmp(&0.0),
                    "{:?}",
                    pa
                );
                assert_eq!(
                    o2d(a, b, c).partial_cmp(&0.0),
                    o2ds(a, b, c).partial_cmp(&0.0),
                    "{:?}",
                    pa
                );
            }
        }
    }

    #[test]
    fn orient2d_fast_test() {
        let a = [0.0, 1.0];
        let b = [2.0, 3.0];
        let c = [4.0, 5.0];
        assert_eq!(orient2d_fast(a, b, c), 0.0);

        // The fast orientation test should also work when the given points are sufficiently
        // non-colinear.
        let mut rng: StdRng = SeedableRng::from_seed(SEED);
        let tol = 5.0e-10; // will not work with 5.0e-14

        for _ in 0..999 {
            let a = [0.5 + tol * rng.gen::<f64>(), 0.5 + tol * rng.gen::<f64>()];
            let b = [12.0, 12.0];
            let c = [24.0, 24.0];
            assert_eq!(orient2d_fast(a, b, c) > 0.0, orient2d_fast(b, c, a) > 0.0);
            assert_eq!(orient2d_fast(b, c, a) > 0.0, orient2d_fast(c, a, b) > 0.0);
            assert_eq!(orient2d_fast(a, b, c) > 0.0, orient2d_fast(b, a, c) < 0.0);
            assert_eq!(orient2d_fast(a, b, c) > 0.0, orient2d_fast(a, c, b) < 0.0);
        }
    }

    #[test]
    fn orient3d_test() {
        let a = [0.0, 1.0, 6.0];
        let b = [2.0, 3.0, 4.0];
        let c = [4.0, 5.0, 1.0];
        let d = [6.0, 2.0, 5.3];
        assert_eq!(orient3d(a, b, c, d), 10.0);
    }

    #[cfg(feature = "transpiled")]
    #[test]
    fn orient3d_transpiled_regression_test() {
        unsafe {
            transpiled::exactinit();
        }
        let mut rng: StdRng = SeedableRng::from_seed(SEED);

        let n = 9999;
        for _ in 0..n {
            let pa = [tol(&mut rng), tol(&mut rng), tol(&mut rng)];
            let pb = [12.0, 12.0, 12.0];
            let pc = [24.0, 24.0, 24.0];
            let pd = [48.0, 48.0, 48.0];

            let o3d = predicates::orient3d;
            let o3d_transpiled = transpiled::orient3d;

            // Test all permutations.
            for [a, b, c, d] in QuickPerm4::new([pa, pb, pc, pd]) {
                assert_eq!(o3d(a, c, d, b), o3d_transpiled(a, c, d, b), "{:?}", pa);
            }
        }
    }

    // The following test verifies equivalence of all of the robust orient3d variants.
    #[test]
    fn orient3d_regression_test() {
        let mut rng: StdRng = SeedableRng::from_seed(SEED);

        let n = 5000;
        for _ in 0..n {
            let pa = [tol(&mut rng), tol(&mut rng), tol(&mut rng)];
            let pb = [12.0, 12.0, 12.0];
            let pc = [24.0, 24.0, 24.0];
            let pd = [48.0, 48.0, 48.0];

            let o3d = predicates::orient3d;
            let o3de = predicates::orient3d_exact;
            let o3ds = predicates::orient3d_slow;

            // Test all permutations.
            // Actually these don't need to be exactly equal, they just need to compare equally to
            // 0.0. It just so happens that they are also equal.
            for [a, b, c, d] in QuickPerm4::new([pa, pb, pc, pd]) {
                assert_eq!(o3d(a, c, d, b), o3de(a, c, d, b), "{:?}", pa);
                assert_eq!(o3d(a, c, d, b), o3ds(a, c, d, b), "{:?}", pa);
            }
        }
    }

    #[test]
    fn orient3d_robustness_test() {
        let mut rng: StdRng = SeedableRng::from_seed(SEED);

        let n = 999;

        for o3d in &[orient3d, orient3d_exact, orient3d_slow] {
            for _ in 0..n {
                let pa = [tol(&mut rng), tol(&mut rng), tol(&mut rng)];
                let pb = [12.0, 12.0, 12.0];
                let pc = [24.0, 24.0, 24.0];
                let pd = [48.0, 48.0, 48.0];

                // Test all permutations.

                let main = o3d(pa, pb, pc, pd);

                if main == 0.0 {
                    for [a, b, c, d] in QuickPerm4::new([pa, pb, pc, pd]) {
                        assert_eq!(o3d(a, b, c, d), 0.0);
                    }
                }

                let pred = main > 0.0;
                for (i, [a, b, c, d]) in QuickPerm4::new([pa, pb, pc, pd]).enumerate() {
                    let t = o3d(a, b, c, d);
                    assert_eq!(
                        pred,
                        if i % 2 == 0 { t > 0.0 } else { t < 0.0 },
                        "{} vs. {} at {:?}",
                        t,
                        -main,
                        pa
                    );
                }
            }
        }
    }

    #[test]
    fn orient3d_fast_test() {
        let a = [0.0, 1.0, 6.0];
        let b = [2.0, 3.0, 4.0];
        let c = [4.0, 5.0, 1.0];
        let d = [6.0, 2.0, 5.3];
        assert!(f64::abs(orient3d_fast(a, b, c, d) - 10.0) < 1.0e-4);

        // The fast orientation test should also work when the given points are sufficiently
        // non-colinear.
        let mut rng: StdRng = SeedableRng::from_seed(SEED);
        let tol = 5.0;

        let a = [
            0.5 + tol * rng.gen::<f64>(),
            0.5 + tol * rng.gen::<f64>(),
            0.5 + tol * rng.gen::<f64>(),
        ];
        let b = [12.0, 12.0, 12.0];
        let c = [24.0, 24.0, 24.0];
        let d = [48.0, 48.0, 48.0];
        assert_eq!(
            orient3d_fast(a, b, c, d) > 0.0,
            orient3d_fast(b, c, a, d) > 0.0
        );
        assert_eq!(
            orient3d_fast(b, a, c, d) < 0.0,
            orient3d_fast(c, b, d, a) < 0.0
        );
        // The following orientations are expected to be inconsistent
        // (this is why we need exact predicates)
        assert_ne!(
            orient3d_fast(b, c, d, a) > 0.0,
            orient3d_fast(c, d, a, b) > 0.0
        );
        assert_ne!(
            orient3d_fast(b, d, c, a) < 0.0,
            orient3d_fast(c, d, b, a) < 0.0
        );
    }

    #[test]
    fn incircle_test() {
        let a = [0.0, 1.0];
        let b = [1.0, 0.0];
        let c = [1.0, 1.0];
        let d = [0.0, 0.0];
        assert_eq!(incircle(a, b, c, d), 0.0);
        let d = [0.1, 0.1];
        assert!(incircle(a, b, c, d) > 0.0);
        let d = [-0.1, -0.1];
        assert!(incircle(a, b, c, d) < 0.0);
    }

    #[test]
    fn incircle_robustness_test() {
        let mut rng: StdRng = SeedableRng::from_seed(SEED);

        let n = 999;
        for ic in &[incircle, incircle_exact, incircle_slow] {
            for _ in 0..n {
                let pa = [0.0, 1.0];
                let pb = [1.0, 0.0];
                let pc = [1.0, 1.0];
                let pd = [tol(&mut rng), tol(&mut rng)];

                let main = ic(pa, pb, pc, pd);

                if main == 0.0 {
                    for [a, b, c, d] in QuickPerm4::new([pa, pb, pc, pd]) {
                        assert_eq!(ic(a, b, c, d), 0.0);
                    }
                }

                let pred = main > 0.0;
                for (i, [a, b, c, d]) in QuickPerm4::new([pa, pb, pc, pd]).enumerate() {
                    let t = ic(a, b, c, d);
                    assert_eq!(
                        pred,
                        if i % 2 == 0 { t > 0.0 } else { t < 0.0 },
                        "{} vs. {} at {:?}",
                        t,
                        -main,
                        pd
                    );
                }
            }
        }
    }

    #[test]
    fn insphere_test() {
        let a = [0.0, 1.0, 0.0];
        let b = [1.0, 0.0, 0.0];
        let c = [0.0, 0.0, 1.0];
        let d = [1.0, 1.0, 1.0];
        let e = [0.0, 0.0, 0.0];
        assert_eq!(insphere(a, b, c, d, e), 0.0);
        let e = [0.1, 0.1, 0.1];
        assert!(insphere(a, b, c, d, e) > 0.0);
        let e = [-0.1, -0.1, -0.1];
        assert!(insphere(a, b, c, d, e) < 0.0);
    }

    #[test]
    fn insphere_robustness_test() {
        let mut rng: StdRng = SeedableRng::from_seed(SEED);

        let n = 99;
        for is in &[insphere, insphere_exact, insphere_slow] {
            for _ in 0..n {
                let pa = [0.0, 1.0, 0.0];
                let pb = [1.0, 0.0, 0.0];
                let pc = [0.0, 0.0, 1.0];
                let pd = [1.0, 1.0, 1.0];
                let pe = [tol(&mut rng), tol(&mut rng), tol(&mut rng)];

                let main = is(pa, pb, pc, pd, pe);

                if main == 0.0 {
                    for [a, b, c, d, e] in QuickPerm5::new([pa, pb, pc, pd, pe]) {
                        assert_eq!(is(a, b, c, d, e), 0.0);
                    }
                }

                let pred = main > 0.0;
                for (i, [a, b, c, d, e]) in QuickPerm5::new([pa, pb, pc, pd, pe]).enumerate() {
                    let t = is(a, b, c, d, e);
                    assert_eq!(
                        pred,
                        if i % 2 == 0 { t > 0.0 } else { t < 0.0 },
                        "{} vs. {} at {:?}",
                        t,
                        -main,
                        pe
                    );
                }
            }
        }
    }
}