-
Notifications
You must be signed in to change notification settings - Fork 6
/
predicates.c
4273 lines (3892 loc) · 168 KB
/
predicates.c
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
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/*****************************************************************************/
/* */
/* Routines for Arbitrary Precision Floating-point Arithmetic */
/* and Fast Robust Geometric Predicates */
/* (predicates.c) */
/* */
/* May 18, 1996 */
/* */
/* Placed in the public domain by */
/* Jonathan Richard Shewchuk */
/* School of Computer Science */
/* Carnegie Mellon University */
/* 5000 Forbes Avenue */
/* Pittsburgh, Pennsylvania 15213-3891 */
/* [email protected] */
/* */
/* This file contains C implementation of algorithms for exact addition */
/* and multiplication of floating-point numbers, and predicates for */
/* robustly performing the orientation and incircle tests used in */
/* computational geometry. The algorithms and underlying theory are */
/* described in Jonathan Richard Shewchuk. "Adaptive Precision Floating- */
/* Point Arithmetic and Fast Robust Geometric Predicates." Technical */
/* Report CMU-CS-96-140, School of Computer Science, Carnegie Mellon */
/* University, Pittsburgh, Pennsylvania, May 1996. (Submitted to */
/* Discrete & Computational Geometry.) */
/* */
/* This file, the paper listed above, and other information are available */
/* from the Web page http://www.cs.cmu.edu/~quake/robust.html . */
/* */
/*****************************************************************************/
/*****************************************************************************/
/* */
/* Using this code: */
/* */
/* First, read the short or long version of the paper (from the Web page */
/* above). */
/* */
/* Be sure to call exactinit() once, before calling any of the arithmetic */
/* functions or geometric predicates. Also be sure to turn on the */
/* optimizer when compiling this file. */
/* */
/* */
/* Several geometric predicates are defined. Their parameters are all */
/* points. Each point is an array of two or three floating-point */
/* numbers. The geometric predicates, described in the papers, are */
/* */
/* orient2d(pa, pb, pc) */
/* orient2dfast(pa, pb, pc) */
/* orient3d(pa, pb, pc, pd) */
/* orient3dfast(pa, pb, pc, pd) */
/* incircle(pa, pb, pc, pd) */
/* incirclefast(pa, pb, pc, pd) */
/* insphere(pa, pb, pc, pd, pe) */
/* inspherefast(pa, pb, pc, pd, pe) */
/* */
/* Those with suffix "fast" are approximate, non-robust versions. Those */
/* without the suffix are adaptive precision, robust versions. There */
/* are also versions with the suffices "exact" and "slow", which are */
/* non-adaptive, exact arithmetic versions, which I use only for timings */
/* in my arithmetic papers. */
/* */
/* */
/* An expansion is represented by an array of floating-point numbers, */
/* sorted from smallest to largest magnitude (possibly with interspersed */
/* zeros). The length of each expansion is stored as a separate integer, */
/* and each arithmetic function returns an integer which is the length */
/* of the expansion it created. */
/* */
/* Several arithmetic functions are defined. Their parameters are */
/* */
/* e, f Input expansions */
/* elen, flen Lengths of input expansions (must be >= 1) */
/* h Output expansion */
/* b Input scalar */
/* */
/* The arithmetic functions are */
/* */
/* grow_expansion(elen, e, b, h) */
/* grow_expansion_zeroelim(elen, e, b, h) */
/* expansion_sum(elen, e, flen, f, h) */
/* expansion_sum_zeroelim1(elen, e, flen, f, h) */
/* expansion_sum_zeroelim2(elen, e, flen, f, h) */
/* fast_expansion_sum(elen, e, flen, f, h) */
/* fast_expansion_sum_zeroelim(elen, e, flen, f, h) */
/* linear_expansion_sum(elen, e, flen, f, h) */
/* linear_expansion_sum_zeroelim(elen, e, flen, f, h) */
/* scale_expansion(elen, e, b, h) */
/* scale_expansion_zeroelim(elen, e, b, h) */
/* compress(elen, e, h) */
/* */
/* All of these are described in the long version of the paper; some are */
/* described in the short version. All return an integer that is the */
/* length of h. Those with suffix _zeroelim perform zero elimination, */
/* and are recommended over their counterparts. The procedure */
/* fast_expansion_sum_zeroelim() (or linear_expansion_sum_zeroelim() on */
/* processors that do not use the round-to-even tiebreaking rule) is */
/* recommended over expansion_sum_zeroelim(). Each procedure has a */
/* little note next to it (in the code below) that tells you whether or */
/* not the output expansion may be the same array as one of the input */
/* expansions. */
/* */
/* */
/* If you look around below, you'll also find macros for a bunch of */
/* simple unrolled arithmetic operations, and procedures for printing */
/* expansions (commented out because they don't work with all C */
/* compilers) and for generating random floating-point numbers whose */
/* significand bits are all random. Most of the macros have undocumented */
/* requirements that certain of their parameters should not be the same */
/* variable; for safety, better to make sure all the parameters are */
/* distinct variables. Feel free to send email to [email protected] if you */
/* have questions. */
/* */
/*****************************************************************************/
#include <stdio.h>
#include <stdlib.h>
#include <math.h>
//#include <sys/time.h>
/* On some machines, the exact arithmetic routines might be defeated by the */
/* use of internal extended precision floating-point registers. Sometimes */
/* this problem can be fixed by defining certain values to be volatile, */
/* thus forcing them to be stored to memory and rounded off. This isn't */
/* a great solution, though, as it slows the arithmetic down. */
/* */
/* To try this out, write "#define INEXACT volatile" below. Normally, */
/* however, INEXACT should be defined to be nothing. ("#define INEXACT".) */
#define INEXACT /* Nothing */
/* #define INEXACT volatile */
#ifdef LIBIGL_PREDICATES_USE_FLOAT
#define REAL float // float
#define REALPRINT floatprint
#define REALRAND floatrand
#define NARROWRAND narrowfloatrand
#define UNIFORMRAND uniformfloatrand
#else
#define REAL double // double
#define REALPRINT doubleprint
#define REALRAND doublerand
#define NARROWRAND narrowdoublerand
#define UNIFORMRAND uniformdoublerand
#endif
/* Which of the following two methods of finding the absolute values is */
/* fastest is compiler-dependent. A few compilers can inline and optimize */
/* the fabs() call; but most will incur the overhead of a function call, */
/* which is disastrously slow. A faster way on IEEE machines might be to */
/* mask the appropriate bit, but that's difficult to do in C. */
#define Absolute(a) ((a) >= 0.0 ? (a) : -(a))
/* #define Absolute(a) fabs(a) */
/* Many of the operations are broken up into two pieces, a main part that */
/* performs an approximate operation, and a "tail" that computes the */
/* roundoff error of that operation. */
/* */
/* The operations Fast_Two_Sum(), Fast_Two_Diff(), Two_Sum(), Two_Diff(), */
/* Split(), and Two_Product() are all implemented as described in the */
/* reference. Each of these macros requires certain variables to be */
/* defined in the calling routine. The variables `bvirt', `c', `abig', */
/* `_i', `_j', `_k', `_l', `_m', and `_n' are declared `INEXACT' because */
/* they store the result of an operation that may incur roundoff error. */
/* The input parameter `x' (or the highest numbered `x_' parameter) must */
/* also be declared `INEXACT'. */
#define Fast_Two_Sum_Tail(a, b, x, y) \
bvirt = x - a; \
y = b - bvirt
#define Fast_Two_Sum(a, b, x, y) \
x = (REAL) (a + b); \
Fast_Two_Sum_Tail(a, b, x, y)
#define Fast_Two_Diff_Tail(a, b, x, y) \
bvirt = a - x; \
y = bvirt - b
#define Fast_Two_Diff(a, b, x, y) \
x = (REAL) (a - b); \
Fast_Two_Diff_Tail(a, b, x, y)
#define Two_Sum_Tail(a, b, x, y) \
bvirt = (REAL) (x - a); \
avirt = x - bvirt; \
bround = b - bvirt; \
around = a - avirt; \
y = around + bround
#define Two_Sum(a, b, x, y) \
x = (REAL) (a + b); \
Two_Sum_Tail(a, b, x, y)
#define Two_Diff_Tail(a, b, x, y) \
bvirt = (REAL) (a - x); \
avirt = x + bvirt; \
bround = bvirt - b; \
around = a - avirt; \
y = around + bround
#define Two_Diff(a, b, x, y) \
x = (REAL) (a - b); \
Two_Diff_Tail(a, b, x, y)
#define Split(a, ahi, alo) \
c = (REAL) (splitter * a); \
abig = (REAL) (c - a); \
ahi = c - abig; \
alo = a - ahi
#define Two_Product_Tail(a, b, x, y) \
Split(a, ahi, alo); \
Split(b, bhi, blo); \
err1 = x - (ahi * bhi); \
err2 = err1 - (alo * bhi); \
err3 = err2 - (ahi * blo); \
y = (alo * blo) - err3
#define Two_Product(a, b, x, y) \
x = (REAL) (a * b); \
Two_Product_Tail(a, b, x, y)
/* Two_Product_Presplit() is Two_Product() where one of the inputs has */
/* already been split. Avoids redundant splitting. */
#define Two_Product_Presplit(a, b, bhi, blo, x, y) \
x = (REAL) (a * b); \
Split(a, ahi, alo); \
err1 = x - (ahi * bhi); \
err2 = err1 - (alo * bhi); \
err3 = err2 - (ahi * blo); \
y = (alo * blo) - err3
/* Two_Product_2Presplit() is Two_Product() where both of the inputs have */
/* already been split. Avoids redundant splitting. */
#define Two_Product_2Presplit(a, ahi, alo, b, bhi, blo, x, y) \
x = (REAL) (a * b); \
err1 = x - (ahi * bhi); \
err2 = err1 - (alo * bhi); \
err3 = err2 - (ahi * blo); \
y = (alo * blo) - err3
/* Square() can be done more quickly than Two_Product(). */
#define Square_Tail(a, x, y) \
Split(a, ahi, alo); \
err1 = x - (ahi * ahi); \
err3 = err1 - ((ahi + ahi) * alo); \
y = (alo * alo) - err3
#define Square(a, x, y) \
x = (REAL) (a * a); \
Square_Tail(a, x, y)
/* Macros for summing expansions of various fixed lengths. These are all */
/* unrolled versions of Expansion_Sum(). */
#define Two_One_Sum(a1, a0, b, x2, x1, x0) \
Two_Sum(a0, b , _i, x0); \
Two_Sum(a1, _i, x2, x1)
#define Two_One_Diff(a1, a0, b, x2, x1, x0) \
Two_Diff(a0, b , _i, x0); \
Two_Sum( a1, _i, x2, x1)
#define Two_Two_Sum(a1, a0, b1, b0, x3, x2, x1, x0) \
Two_One_Sum(a1, a0, b0, _j, _0, x0); \
Two_One_Sum(_j, _0, b1, x3, x2, x1)
#define Two_Two_Diff(a1, a0, b1, b0, x3, x2, x1, x0) \
Two_One_Diff(a1, a0, b0, _j, _0, x0); \
Two_One_Diff(_j, _0, b1, x3, x2, x1)
#define Four_One_Sum(a3, a2, a1, a0, b, x4, x3, x2, x1, x0) \
Two_One_Sum(a1, a0, b , _j, x1, x0); \
Two_One_Sum(a3, a2, _j, x4, x3, x2)
#define Four_Two_Sum(a3, a2, a1, a0, b1, b0, x5, x4, x3, x2, x1, x0) \
Four_One_Sum(a3, a2, a1, a0, b0, _k, _2, _1, _0, x0); \
Four_One_Sum(_k, _2, _1, _0, b1, x5, x4, x3, x2, x1)
#define Four_Four_Sum(a3, a2, a1, a0, b4, b3, b1, b0, x7, x6, x5, x4, x3, x2, \
x1, x0) \
Four_Two_Sum(a3, a2, a1, a0, b1, b0, _l, _2, _1, _0, x1, x0); \
Four_Two_Sum(_l, _2, _1, _0, b4, b3, x7, x6, x5, x4, x3, x2)
#define Eight_One_Sum(a7, a6, a5, a4, a3, a2, a1, a0, b, x8, x7, x6, x5, x4, \
x3, x2, x1, x0) \
Four_One_Sum(a3, a2, a1, a0, b , _j, x3, x2, x1, x0); \
Four_One_Sum(a7, a6, a5, a4, _j, x8, x7, x6, x5, x4)
#define Eight_Two_Sum(a7, a6, a5, a4, a3, a2, a1, a0, b1, b0, x9, x8, x7, \
x6, x5, x4, x3, x2, x1, x0) \
Eight_One_Sum(a7, a6, a5, a4, a3, a2, a1, a0, b0, _k, _6, _5, _4, _3, _2, \
_1, _0, x0); \
Eight_One_Sum(_k, _6, _5, _4, _3, _2, _1, _0, b1, x9, x8, x7, x6, x5, x4, \
x3, x2, x1)
#define Eight_Four_Sum(a7, a6, a5, a4, a3, a2, a1, a0, b4, b3, b1, b0, x11, \
x10, x9, x8, x7, x6, x5, x4, x3, x2, x1, x0) \
Eight_Two_Sum(a7, a6, a5, a4, a3, a2, a1, a0, b1, b0, _l, _6, _5, _4, _3, \
_2, _1, _0, x1, x0); \
Eight_Two_Sum(_l, _6, _5, _4, _3, _2, _1, _0, b4, b3, x11, x10, x9, x8, \
x7, x6, x5, x4, x3, x2)
/* Macros for multiplying expansions of various fixed lengths. */
#define Two_One_Product(a1, a0, b, x3, x2, x1, x0) \
Split(b, bhi, blo); \
Two_Product_Presplit(a0, b, bhi, blo, _i, x0); \
Two_Product_Presplit(a1, b, bhi, blo, _j, _0); \
Two_Sum(_i, _0, _k, x1); \
Fast_Two_Sum(_j, _k, x3, x2)
#define Four_One_Product(a3, a2, a1, a0, b, x7, x6, x5, x4, x3, x2, x1, x0) \
Split(b, bhi, blo); \
Two_Product_Presplit(a0, b, bhi, blo, _i, x0); \
Two_Product_Presplit(a1, b, bhi, blo, _j, _0); \
Two_Sum(_i, _0, _k, x1); \
Fast_Two_Sum(_j, _k, _i, x2); \
Two_Product_Presplit(a2, b, bhi, blo, _j, _0); \
Two_Sum(_i, _0, _k, x3); \
Fast_Two_Sum(_j, _k, _i, x4); \
Two_Product_Presplit(a3, b, bhi, blo, _j, _0); \
Two_Sum(_i, _0, _k, x5); \
Fast_Two_Sum(_j, _k, x7, x6)
#define Two_Two_Product(a1, a0, b1, b0, x7, x6, x5, x4, x3, x2, x1, x0) \
Split(a0, a0hi, a0lo); \
Split(b0, bhi, blo); \
Two_Product_2Presplit(a0, a0hi, a0lo, b0, bhi, blo, _i, x0); \
Split(a1, a1hi, a1lo); \
Two_Product_2Presplit(a1, a1hi, a1lo, b0, bhi, blo, _j, _0); \
Two_Sum(_i, _0, _k, _1); \
Fast_Two_Sum(_j, _k, _l, _2); \
Split(b1, bhi, blo); \
Two_Product_2Presplit(a0, a0hi, a0lo, b1, bhi, blo, _i, _0); \
Two_Sum(_1, _0, _k, x1); \
Two_Sum(_2, _k, _j, _1); \
Two_Sum(_l, _j, _m, _2); \
Two_Product_2Presplit(a1, a1hi, a1lo, b1, bhi, blo, _j, _0); \
Two_Sum(_i, _0, _n, _0); \
Two_Sum(_1, _0, _i, x2); \
Two_Sum(_2, _i, _k, _1); \
Two_Sum(_m, _k, _l, _2); \
Two_Sum(_j, _n, _k, _0); \
Two_Sum(_1, _0, _j, x3); \
Two_Sum(_2, _j, _i, _1); \
Two_Sum(_l, _i, _m, _2); \
Two_Sum(_1, _k, _i, x4); \
Two_Sum(_2, _i, _k, x5); \
Two_Sum(_m, _k, x7, x6)
/* An expansion of length two can be squared more quickly than finding the */
/* product of two different expansions of length two, and the result is */
/* guaranteed to have no more than six (rather than eight) components. */
#define Two_Square(a1, a0, x5, x4, x3, x2, x1, x0) \
Square(a0, _j, x0); \
_0 = a0 + a0; \
Two_Product(a1, _0, _k, _1); \
Two_One_Sum(_k, _1, _j, _l, _2, x1); \
Square(a1, _j, _1); \
Two_Two_Sum(_j, _1, _l, _2, x5, x4, x3, x2)
static REAL splitter; /* = 2^ceiling(p / 2) + 1. Used to split floats in half. */
static REAL epsilon; /* = 2^(-p). Used to estimate roundoff errors. */
/* A set of coefficients used to calculate maximum roundoff errors. */
static REAL resulterrbound;
static REAL ccwerrboundA, ccwerrboundB, ccwerrboundC;
static REAL o3derrboundA, o3derrboundB, o3derrboundC;
static REAL iccerrboundA, iccerrboundB, iccerrboundC;
static REAL isperrboundA, isperrboundB, isperrboundC;
/*****************************************************************************/
/* */
/* doubleprint() Print the bit representation of a double. */
/* */
/* Useful for debugging exact arithmetic routines. */
/* */
/*****************************************************************************/
/*
void doubleprint(number)
double number;
{
unsigned long long no;
unsigned long long sign, expo;
int exponent;
int i, bottomi;
no = *(unsigned long long *) &number;
sign = no & 0x8000000000000000ll;
expo = (no >> 52) & 0x7ffll;
exponent = (int) expo;
exponent = exponent - 1023;
if (sign) {
printf("-");
} else {
printf(" ");
}
if (exponent == -1023) {
printf(
"0.0000000000000000000000000000000000000000000000000000_ ( )");
} else {
printf("1.");
bottomi = -1;
for (i = 0; i < 52; i++) {
if (no & 0x0008000000000000ll) {
printf("1");
bottomi = i;
} else {
printf("0");
}
no <<= 1;
}
printf("_%d (%d)", exponent, exponent - 1 - bottomi);
}
}
*/
/*****************************************************************************/
/* */
/* floatprint() Print the bit representation of a float. */
/* */
/* Useful for debugging exact arithmetic routines. */
/* */
/*****************************************************************************/
/*
void floatprint(number)
float number;
{
unsigned no;
unsigned sign, expo;
int exponent;
int i, bottomi;
no = *(unsigned *) &number;
sign = no & 0x80000000;
expo = (no >> 23) & 0xff;
exponent = (int) expo;
exponent = exponent - 127;
if (sign) {
printf("-");
} else {
printf(" ");
}
if (exponent == -127) {
printf("0.00000000000000000000000_ ( )");
} else {
printf("1.");
bottomi = -1;
for (i = 0; i < 23; i++) {
if (no & 0x00400000) {
printf("1");
bottomi = i;
} else {
printf("0");
}
no <<= 1;
}
printf("_%3d (%3d)", exponent, exponent - 1 - bottomi);
}
}
*/
/*****************************************************************************/
/* */
/* expansion_print() Print the bit representation of an expansion. */
/* */
/* Useful for debugging exact arithmetic routines. */
/* */
/*****************************************************************************/
/*
void expansion_print(elen, e)
int elen;
REAL *e;
{
int i;
for (i = elen - 1; i >= 0; i--) {
REALPRINT(e[i]);
if (i > 0) {
printf(" +\n");
} else {
printf("\n");
}
}
}
*/
/*****************************************************************************/
/* */
/* doublerand() Generate a double with random 53-bit significand and a */
/* random exponent in [0, 511]. */
/* */
/*****************************************************************************/
static double doublerand()
{
double result;
double expo;
long a, b, c;
long i;
a = rand();
b = rand();
c = rand();
result = (double) (a - 1073741824) * 8388608.0 + (double) (b >> 8);
for (i = 512, expo = 2; i <= 131072; i *= 2, expo = expo * expo) {
if (c & i) {
result *= expo;
}
}
return result;
}
/*****************************************************************************/
/* */
/* narrowdoublerand() Generate a double with random 53-bit significand */
/* and a random exponent in [0, 7]. */
/* */
/*****************************************************************************/
static double narrowdoublerand()
{
double result;
double expo;
long a, b, c;
long i;
a = rand();
b = rand();
c = rand();
result = (double) (a - 1073741824) * 8388608.0 + (double) (b >> 8);
for (i = 512, expo = 2; i <= 2048; i *= 2, expo = expo * expo) {
if (c & i) {
result *= expo;
}
}
return result;
}
/*****************************************************************************/
/* */
/* uniformdoublerand() Generate a double with random 53-bit significand. */
/* */
/*****************************************************************************/
static double uniformdoublerand()
{
double result;
long a, b;
a = rand();
b = rand();
result = (double) (a - 1073741824) * 8388608.0 + (double) (b >> 8);
return result;
}
/*****************************************************************************/
/* */
/* floatrand() Generate a float with random 24-bit significand and a */
/* random exponent in [0, 63]. */
/* */
/*****************************************************************************/
static float floatrand()
{
float result;
float expo;
long a, c;
long i;
a = rand();
c = rand();
result = (float) ((a - 1073741824) >> 6);
for (i = 512, expo = 2; i <= 16384; i *= 2, expo = expo * expo) {
if (c & i) {
result *= expo;
}
}
return result;
}
/*****************************************************************************/
/* */
/* narrowfloatrand() Generate a float with random 24-bit significand and */
/* a random exponent in [0, 7]. */
/* */
/*****************************************************************************/
static float narrowfloatrand()
{
float result;
float expo;
long a, c;
long i;
a = rand();
c = rand();
result = (float) ((a - 1073741824) >> 6);
for (i = 512, expo = 2; i <= 2048; i *= 2, expo = expo * expo) {
if (c & i) {
result *= expo;
}
}
return result;
}
/*****************************************************************************/
/* */
/* uniformfloatrand() Generate a float with random 24-bit significand. */
/* */
/*****************************************************************************/
static float uniformfloatrand()
{
float result;
long a;
a = rand();
result = (float) ((a - 1073741824) >> 6);
return result;
}
/*****************************************************************************/
/* */
/* exactinit() Initialize the variables used for exact arithmetic. */
/* */
/* `epsilon' is the largest power of two such that 1.0 + epsilon = 1.0 in */
/* floating-point arithmetic. `epsilon' bounds the relative roundoff */
/* error. It is used for floating-point error analysis. */
/* */
/* `splitter' is used to split floating-point numbers into two half- */
/* length significands for exact multiplication. */
/* */
/* I imagine that a highly optimizing compiler might be too smart for its */
/* own good, and somehow cause this routine to fail, if it pretends that */
/* floating-point arithmetic is too much like real arithmetic. */
/* */
/* Don't change this routine unless you fully understand it. */
/* */
/*****************************************************************************/
void exactinit()
{
REAL half;
REAL check, lastcheck;
int every_other;
every_other = 1;
half = 0.5;
epsilon = 1.0;
splitter = 1.0;
check = 1.0;
/* Repeatedly divide `epsilon' by two until it is too small to add to */
/* one without causing roundoff. (Also check if the sum is equal to */
/* the previous sum, for machines that round up instead of using exact */
/* rounding. Not that this library will work on such machines anyway. */
do {
lastcheck = check;
epsilon *= half;
if (every_other) {
splitter *= 2.0;
}
every_other = !every_other;
check = 1.0 + epsilon;
} while ((check != 1.0) && (check != lastcheck));
splitter += 1.0;
/* Error bounds for orientation and incircle tests. */
resulterrbound = (3.0 + 8.0 * epsilon) * epsilon;
ccwerrboundA = (3.0 + 16.0 * epsilon) * epsilon;
ccwerrboundB = (2.0 + 12.0 * epsilon) * epsilon;
ccwerrboundC = (9.0 + 64.0 * epsilon) * epsilon * epsilon;
o3derrboundA = (7.0 + 56.0 * epsilon) * epsilon;
o3derrboundB = (3.0 + 28.0 * epsilon) * epsilon;
o3derrboundC = (26.0 + 288.0 * epsilon) * epsilon * epsilon;
iccerrboundA = (10.0 + 96.0 * epsilon) * epsilon;
iccerrboundB = (4.0 + 48.0 * epsilon) * epsilon;
iccerrboundC = (44.0 + 576.0 * epsilon) * epsilon * epsilon;
isperrboundA = (16.0 + 224.0 * epsilon) * epsilon;
isperrboundB = (5.0 + 72.0 * epsilon) * epsilon;
isperrboundC = (71.0 + 1408.0 * epsilon) * epsilon * epsilon;
}
/*****************************************************************************/
/* */
/* grow_expansion() Add a scalar to an expansion. */
/* */
/* Sets h = e + b. See the long version of my paper for details. */
/* */
/* Maintains the nonoverlapping property. If round-to-even is used (as */
/* with IEEE 754), maintains the strongly nonoverlapping and nonadjacent */
/* properties as well. (That is, if e has one of these properties, so */
/* will h.) */
/* */
/*****************************************************************************/
static int grow_expansion(elen, e, b, h) /* e and h can be the same. */
int elen;
REAL *e;
REAL b;
REAL *h;
{
REAL Q;
INEXACT REAL Qnew;
int eindex;
REAL enow;
INEXACT REAL bvirt;
REAL avirt, bround, around;
Q = b;
for (eindex = 0; eindex < elen; eindex++) {
enow = e[eindex];
Two_Sum(Q, enow, Qnew, h[eindex]);
Q = Qnew;
}
h[eindex] = Q;
return eindex + 1;
}
/*****************************************************************************/
/* */
/* grow_expansion_zeroelim() Add a scalar to an expansion, eliminating */
/* zero components from the output expansion. */
/* */
/* Sets h = e + b. See the long version of my paper for details. */
/* */
/* Maintains the nonoverlapping property. If round-to-even is used (as */
/* with IEEE 754), maintains the strongly nonoverlapping and nonadjacent */
/* properties as well. (That is, if e has one of these properties, so */
/* will h.) */
/* */
/*****************************************************************************/
static int grow_expansion_zeroelim(elen, e, b, h) /* e and h can be the same. */
int elen;
REAL *e;
REAL b;
REAL *h;
{
REAL Q, hh;
INEXACT REAL Qnew;
int eindex, hindex;
REAL enow;
INEXACT REAL bvirt;
REAL avirt, bround, around;
hindex = 0;
Q = b;
for (eindex = 0; eindex < elen; eindex++) {
enow = e[eindex];
Two_Sum(Q, enow, Qnew, hh);
Q = Qnew;
if (hh != 0.0) {
h[hindex++] = hh;
}
}
if ((Q != 0.0) || (hindex == 0)) {
h[hindex++] = Q;
}
return hindex;
}
/*****************************************************************************/
/* */
/* expansion_sum() Sum two expansions. */
/* */
/* Sets h = e + f. See the long version of my paper for details. */
/* */
/* Maintains the nonoverlapping property. If round-to-even is used (as */
/* with IEEE 754), maintains the nonadjacent property as well. (That is, */
/* if e has one of these properties, so will h.) Does NOT maintain the */
/* strongly nonoverlapping property. */
/* */
/*****************************************************************************/
static int expansion_sum(elen, e, flen, f, h)
/* e and h can be the same, but f and h cannot. */
int elen;
REAL *e;
int flen;
REAL *f;
REAL *h;
{
REAL Q;
INEXACT REAL Qnew;
int findex, hindex, hlast;
REAL hnow;
INEXACT REAL bvirt;
REAL avirt, bround, around;
Q = f[0];
for (hindex = 0; hindex < elen; hindex++) {
hnow = e[hindex];
Two_Sum(Q, hnow, Qnew, h[hindex]);
Q = Qnew;
}
h[hindex] = Q;
hlast = hindex;
for (findex = 1; findex < flen; findex++) {
Q = f[findex];
for (hindex = findex; hindex <= hlast; hindex++) {
hnow = h[hindex];
Two_Sum(Q, hnow, Qnew, h[hindex]);
Q = Qnew;
}
h[++hlast] = Q;
}
return hlast + 1;
}
/*****************************************************************************/
/* */
/* expansion_sum_zeroelim1() Sum two expansions, eliminating zero */
/* components from the output expansion. */
/* */
/* Sets h = e + f. See the long version of my paper for details. */
/* */
/* Maintains the nonoverlapping property. If round-to-even is used (as */
/* with IEEE 754), maintains the nonadjacent property as well. (That is, */
/* if e has one of these properties, so will h.) Does NOT maintain the */
/* strongly nonoverlapping property. */
/* */
/*****************************************************************************/
static int expansion_sum_zeroelim1(elen, e, flen, f, h)
/* e and h can be the same, but f and h cannot. */
int elen;
REAL *e;
int flen;
REAL *f;
REAL *h;
{
REAL Q;
INEXACT REAL Qnew;
int index, findex, hindex, hlast;
REAL hnow;
INEXACT REAL bvirt;
REAL avirt, bround, around;
Q = f[0];
for (hindex = 0; hindex < elen; hindex++) {
hnow = e[hindex];
Two_Sum(Q, hnow, Qnew, h[hindex]);
Q = Qnew;
}
h[hindex] = Q;
hlast = hindex;
for (findex = 1; findex < flen; findex++) {
Q = f[findex];
for (hindex = findex; hindex <= hlast; hindex++) {
hnow = h[hindex];
Two_Sum(Q, hnow, Qnew, h[hindex]);
Q = Qnew;
}
h[++hlast] = Q;
}
hindex = -1;
for (index = 0; index <= hlast; index++) {
hnow = h[index];
if (hnow != 0.0) {
h[++hindex] = hnow;
}
}
if (hindex == -1) {
return 1;
} else {
return hindex + 1;
}
}
/*****************************************************************************/
/* */
/* expansion_sum_zeroelim2() Sum two expansions, eliminating zero */
/* components from the output expansion. */
/* */
/* Sets h = e + f. See the long version of my paper for details. */
/* */
/* Maintains the nonoverlapping property. If round-to-even is used (as */
/* with IEEE 754), maintains the nonadjacent property as well. (That is, */
/* if e has one of these properties, so will h.) Does NOT maintain the */
/* strongly nonoverlapping property. */
/* */
/*****************************************************************************/
static int expansion_sum_zeroelim2(elen, e, flen, f, h)
/* e and h can be the same, but f and h cannot. */
int elen;
REAL *e;
int flen;
REAL *f;
REAL *h;
{
REAL Q, hh;
INEXACT REAL Qnew;
int eindex, findex, hindex, hlast;
REAL enow;
INEXACT REAL bvirt;
REAL avirt, bround, around;
hindex = 0;
Q = f[0];
for (eindex = 0; eindex < elen; eindex++) {
enow = e[eindex];
Two_Sum(Q, enow, Qnew, hh);
Q = Qnew;
if (hh != 0.0) {
h[hindex++] = hh;
}
}
h[hindex] = Q;
hlast = hindex;
for (findex = 1; findex < flen; findex++) {
hindex = 0;
Q = f[findex];
for (eindex = 0; eindex <= hlast; eindex++) {
enow = h[eindex];
Two_Sum(Q, enow, Qnew, hh);
Q = Qnew;
if (hh != 0) {
h[hindex++] = hh;
}
}
h[hindex] = Q;
hlast = hindex;
}
return hlast + 1;
}
/*****************************************************************************/
/* */
/* fast_expansion_sum() Sum two expansions. */
/* */
/* Sets h = e + f. See the long version of my paper for details. */
/* */
/* If round-to-even is used (as with IEEE 754), maintains the strongly */
/* nonoverlapping property. (That is, if e is strongly nonoverlapping, h */
/* will be also.) Does NOT maintain the nonoverlapping or nonadjacent */
/* properties. */
/* */
/*****************************************************************************/
static int fast_expansion_sum(elen, e, flen, f, h) /* h cannot be e or f. */
int elen;
REAL *e;
int flen;
REAL *f;
REAL *h;
{
REAL Q;
INEXACT REAL Qnew;
INEXACT REAL bvirt;
REAL avirt, bround, around;
int eindex, findex, hindex;
REAL enow, fnow;
enow = e[0];
fnow = f[0];
eindex = findex = 0;
if ((fnow > enow) == (fnow > -enow)) {
Q = enow;
enow = e[++eindex];
} else {
Q = fnow;
fnow = f[++findex];
}
hindex = 0;
if ((eindex < elen) && (findex < flen)) {
if ((fnow > enow) == (fnow > -enow)) {
Fast_Two_Sum(enow, Q, Qnew, h[0]);
enow = e[++eindex];
} else {
Fast_Two_Sum(fnow, Q, Qnew, h[0]);
fnow = f[++findex];
}
Q = Qnew;
hindex = 1;
while ((eindex < elen) && (findex < flen)) {
if ((fnow > enow) == (fnow > -enow)) {
Two_Sum(Q, enow, Qnew, h[hindex]);
enow = e[++eindex];
} else {
Two_Sum(Q, fnow, Qnew, h[hindex]);
fnow = f[++findex];
}
Q = Qnew;
hindex++;
}
}