123 #include <fpu_control.h>
128 #ifdef USE_CGAL_PREDICATES
129 #include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
130 typedef CGAL::Exact_predicates_inexact_constructions_kernel cgalEpick;
131 typedef cgalEpick::Point_3 Point;
132 cgalEpick cgal_pred_obj;
133 #endif // #ifdef USE_CGAL_PREDICATES
148 #define REALPRINT doubleprint
149 #define REALRAND doublerand
150 #define NARROWRAND narrowdoublerand
151 #define UNIFORMRAND uniformdoublerand
160 #define Absolute(a) fabs(a)
175 #define Fast_Two_Sum_Tail(a, b, x, y) \
179 #define Fast_Two_Sum(a, b, x, y) \
180 x = (REAL) (a + b); \
181 Fast_Two_Sum_Tail(a, b, x, y)
183 #define Fast_Two_Diff_Tail(a, b, x, y) \
187 #define Fast_Two_Diff(a, b, x, y) \
188 x = (REAL) (a - b); \
189 Fast_Two_Diff_Tail(a, b, x, y)
191 #define Two_Sum_Tail(a, b, x, y) \
192 bvirt = (REAL) (x - a); \
194 bround = b - bvirt; \
195 around = a - avirt; \
198 #define Two_Sum(a, b, x, y) \
199 x = (REAL) (a + b); \
200 Two_Sum_Tail(a, b, x, y)
202 #define Two_Diff_Tail(a, b, x, y) \
203 bvirt = (REAL) (a - x); \
205 bround = bvirt - b; \
206 around = a - avirt; \
209 #define Two_Diff(a, b, x, y) \
210 x = (REAL) (a - b); \
211 Two_Diff_Tail(a, b, x, y)
213 #define Split(a, ahi, alo) \
214 c = (REAL) (splitter * a); \
215 abig = (REAL) (c - a); \
219 #define Two_Product_Tail(a, b, x, y) \
220 Split(a, ahi, alo); \
221 Split(b, bhi, blo); \
222 err1 = x - (ahi * bhi); \
223 err2 = err1 - (alo * bhi); \
224 err3 = err2 - (ahi * blo); \
225 y = (alo * blo) - err3
227 #define Two_Product(a, b, x, y) \
228 x = (REAL) (a * b); \
229 Two_Product_Tail(a, b, x, y)
234 #define Two_Product_Presplit(a, b, bhi, blo, x, y) \
235 x = (REAL) (a * b); \
236 Split(a, ahi, alo); \
237 err1 = x - (ahi * bhi); \
238 err2 = err1 - (alo * bhi); \
239 err3 = err2 - (ahi * blo); \
240 y = (alo * blo) - err3
245 #define Two_Product_2Presplit(a, ahi, alo, b, bhi, blo, x, y) \
246 x = (REAL) (a * b); \
247 err1 = x - (ahi * bhi); \
248 err2 = err1 - (alo * bhi); \
249 err3 = err2 - (ahi * blo); \
250 y = (alo * blo) - err3
254 #define Square_Tail(a, x, y) \
255 Split(a, ahi, alo); \
256 err1 = x - (ahi * ahi); \
257 err3 = err1 - ((ahi + ahi) * alo); \
258 y = (alo * alo) - err3
260 #define Square(a, x, y) \
261 x = (REAL) (a * a); \
267 #define Two_One_Sum(a1, a0, b, x2, x1, x0) \
268 Two_Sum(a0, b , _i, x0); \
269 Two_Sum(a1, _i, x2, x1)
271 #define Two_One_Diff(a1, a0, b, x2, x1, x0) \
272 Two_Diff(a0, b , _i, x0); \
273 Two_Sum( a1, _i, x2, x1)
275 #define Two_Two_Sum(a1, a0, b1, b0, x3, x2, x1, x0) \
276 Two_One_Sum(a1, a0, b0, _j, _0, x0); \
277 Two_One_Sum(_j, _0, b1, x3, x2, x1)
279 #define Two_Two_Diff(a1, a0, b1, b0, x3, x2, x1, x0) \
280 Two_One_Diff(a1, a0, b0, _j, _0, x0); \
281 Two_One_Diff(_j, _0, b1, x3, x2, x1)
283 #define Four_One_Sum(a3, a2, a1, a0, b, x4, x3, x2, x1, x0) \
284 Two_One_Sum(a1, a0, b , _j, x1, x0); \
285 Two_One_Sum(a3, a2, _j, x4, x3, x2)
287 #define Four_Two_Sum(a3, a2, a1, a0, b1, b0, x5, x4, x3, x2, x1, x0) \
288 Four_One_Sum(a3, a2, a1, a0, b0, _k, _2, _1, _0, x0); \
289 Four_One_Sum(_k, _2, _1, _0, b1, x5, x4, x3, x2, x1)
291 #define Four_Four_Sum(a3, a2, a1, a0, b4, b3, b1, b0, x7, x6, x5, x4, x3, x2, \
293 Four_Two_Sum(a3, a2, a1, a0, b1, b0, _l, _2, _1, _0, x1, x0); \
294 Four_Two_Sum(_l, _2, _1, _0, b4, b3, x7, x6, x5, x4, x3, x2)
296 #define Eight_One_Sum(a7, a6, a5, a4, a3, a2, a1, a0, b, x8, x7, x6, x5, x4, \
298 Four_One_Sum(a3, a2, a1, a0, b , _j, x3, x2, x1, x0); \
299 Four_One_Sum(a7, a6, a5, a4, _j, x8, x7, x6, x5, x4)
301 #define Eight_Two_Sum(a7, a6, a5, a4, a3, a2, a1, a0, b1, b0, x9, x8, x7, \
302 x6, x5, x4, x3, x2, x1, x0) \
303 Eight_One_Sum(a7, a6, a5, a4, a3, a2, a1, a0, b0, _k, _6, _5, _4, _3, _2, \
305 Eight_One_Sum(_k, _6, _5, _4, _3, _2, _1, _0, b1, x9, x8, x7, x6, x5, x4, \
308 #define Eight_Four_Sum(a7, a6, a5, a4, a3, a2, a1, a0, b4, b3, b1, b0, x11, \
309 x10, x9, x8, x7, x6, x5, x4, x3, x2, x1, x0) \
310 Eight_Two_Sum(a7, a6, a5, a4, a3, a2, a1, a0, b1, b0, _l, _6, _5, _4, _3, \
311 _2, _1, _0, x1, x0); \
312 Eight_Two_Sum(_l, _6, _5, _4, _3, _2, _1, _0, b4, b3, x11, x10, x9, x8, \
313 x7, x6, x5, x4, x3, x2)
317 #define Two_One_Product(a1, a0, b, x3, x2, x1, x0) \
318 Split(b, bhi, blo); \
319 Two_Product_Presplit(a0, b, bhi, blo, _i, x0); \
320 Two_Product_Presplit(a1, b, bhi, blo, _j, _0); \
321 Two_Sum(_i, _0, _k, x1); \
322 Fast_Two_Sum(_j, _k, x3, x2)
324 #define Four_One_Product(a3, a2, a1, a0, b, x7, x6, x5, x4, x3, x2, x1, x0) \
325 Split(b, bhi, blo); \
326 Two_Product_Presplit(a0, b, bhi, blo, _i, x0); \
327 Two_Product_Presplit(a1, b, bhi, blo, _j, _0); \
328 Two_Sum(_i, _0, _k, x1); \
329 Fast_Two_Sum(_j, _k, _i, x2); \
330 Two_Product_Presplit(a2, b, bhi, blo, _j, _0); \
331 Two_Sum(_i, _0, _k, x3); \
332 Fast_Two_Sum(_j, _k, _i, x4); \
333 Two_Product_Presplit(a3, b, bhi, blo, _j, _0); \
334 Two_Sum(_i, _0, _k, x5); \
335 Fast_Two_Sum(_j, _k, x7, x6)
337 #define Two_Two_Product(a1, a0, b1, b0, x7, x6, x5, x4, x3, x2, x1, x0) \
338 Split(a0, a0hi, a0lo); \
339 Split(b0, bhi, blo); \
340 Two_Product_2Presplit(a0, a0hi, a0lo, b0, bhi, blo, _i, x0); \
341 Split(a1, a1hi, a1lo); \
342 Two_Product_2Presplit(a1, a1hi, a1lo, b0, bhi, blo, _j, _0); \
343 Two_Sum(_i, _0, _k, _1); \
344 Fast_Two_Sum(_j, _k, _l, _2); \
345 Split(b1, bhi, blo); \
346 Two_Product_2Presplit(a0, a0hi, a0lo, b1, bhi, blo, _i, _0); \
347 Two_Sum(_1, _0, _k, x1); \
348 Two_Sum(_2, _k, _j, _1); \
349 Two_Sum(_l, _j, _m, _2); \
350 Two_Product_2Presplit(a1, a1hi, a1lo, b1, bhi, blo, _j, _0); \
351 Two_Sum(_i, _0, _n, _0); \
352 Two_Sum(_1, _0, _i, x2); \
353 Two_Sum(_2, _i, _k, _1); \
354 Two_Sum(_m, _k, _l, _2); \
355 Two_Sum(_j, _n, _k, _0); \
356 Two_Sum(_1, _0, _j, x3); \
357 Two_Sum(_2, _j, _i, _1); \
358 Two_Sum(_l, _i, _m, _2); \
359 Two_Sum(_1, _k, _i, x4); \
360 Two_Sum(_2, _i, _k, x5); \
361 Two_Sum(_m, _k, x7, x6)
367 #define Two_Square(a1, a0, x5, x4, x3, x2, x1, x0) \
368 Square(a0, _j, x0); \
370 Two_Product(a1, _0, _k, _1); \
371 Two_One_Sum(_k, _1, _j, _l, _2, x1); \
372 Square(a1, _j, _1); \
373 Two_Two_Sum(_j, _1, _l, _2, x5, x4, x3, x2)
405 x = (n < 0) ? ((
double)1.0/(double)2.0) : (double)2.0;
406 n = (n < 0) ? -n : n;
415 float fstore(
float x)
428 (
void)
printf(
" sizeof(float) = %2u\n", (
unsigned int)
sizeof(float));
429 #ifdef CPU86 // <float.h>
430 (
void)
printf(
" FLT_MANT_DIG = %2d\n", FLT_MANT_DIG);
435 while (fstore((
float)1.0 + x/(float)2.0) != (float)1.0)
438 (
void)
printf(
" machine epsilon = %13.5e ", x);
440 if (x == (
float)
fppow2(-23)) {
444 (
void)
printf(
"[not IEEE 754 conformant] !!\n");
449 while (fstore(x / (
float)2.0) != (float)0.0)
452 (
void)
printf(
" smallest positive number = %13.5e ", x);
454 if (x == (
float)
fppow2(-149)) {
456 (
void)
printf(
"[smallest 32-bit subnormal]\n");
457 }
else if (x == (
float)
fppow2(-126)) {
461 (
void)
printf(
"[not IEEE 754 conformant] !!\n");
482 (
void)
printf(
" sizeof(double) = %2u\n", (
unsigned int)
sizeof(double));
483 #ifdef CPU86 // <float.h>
484 (
void)
printf(
" DBL_MANT_DIG = %2d\n", DBL_MANT_DIG);
489 while (
dstore(1.0 + x/2.0) != 1.0)
492 (
void)
printf(
" machine epsilon = %13.5le ", x);
494 if (x == (
double)
fppow2(-52)) {
498 (
void)
printf(
"[not IEEE 754 conformant] !!\n");
503 while (
dstore(x / 2.0) != 0.0)
508 if (x == (
double)
fppow2(-1074)) {
511 }
else if (x == (
double)
fppow2(-1022)) {
515 (
void)
printf(
"[not IEEE 754 conformant] !!\n");
547 REAL check, lastcheck;
555 _control87(_PC_24, _MCW_PC);
557 _control87(_PC_53, _MCW_PC);
572 printf(
" Initializing robust predicates.\n");
575 #ifdef USE_CGAL_PREDICATES
576 if (cgal_pred_obj.Has_static_filters) {
577 printf(
" Use static filter.\n");
579 printf(
" No static filter.\n");
581 #endif // USE_CGAL_PREDICATES
604 every_other = !every_other;
606 }
while ((check != 1.0) && (check != lastcheck));
637 half = maxx; maxx =
maxz; maxz =
half;
640 half = maxy; maxy =
maxz; maxz =
half;
642 else if (maxy < maxx) {
643 half = maxy; maxy = maxx; maxx =
half;
672 REAL avirt, bround, around;
675 for (eindex = 0; eindex < elen; eindex++) {
677 Two_Sum(Q, enow, Qnew, h[eindex]);
706 REAL avirt, bround, around;
710 for (eindex = 0; eindex < elen; eindex++) {
718 if ((Q != 0.0) || (hindex == 0)) {
742 int findex, hindex, hlast;
745 REAL avirt, bround, around;
748 for (hindex = 0; hindex < elen; hindex++) {
750 Two_Sum(Q, hnow, Qnew, h[hindex]);
755 for (findex = 1; findex < flen; findex++) {
757 for (hindex = findex; hindex <= hlast; hindex++) {
759 Two_Sum(Q, hnow, Qnew, h[hindex]);
786 int index, findex, hindex, hlast;
789 REAL avirt, bround, around;
792 for (hindex = 0; hindex < elen; hindex++) {
794 Two_Sum(Q, hnow, Qnew, h[hindex]);
799 for (findex = 1; findex < flen; findex++) {
801 for (hindex = findex; hindex <= hlast; hindex++) {
803 Two_Sum(Q, hnow, Qnew, h[hindex]);
809 for (index = 0; index <= hlast; index++) {
841 int eindex, findex, hindex, hlast;
844 REAL avirt, bround, around;
848 for (eindex = 0; eindex < elen; eindex++) {
858 for (findex = 1; findex < flen; findex++) {
861 for (eindex = 0; eindex <= hlast; eindex++) {
894 REAL avirt, bround, around;
895 int eindex, findex, hindex;
901 if ((fnow > enow) == (fnow > -enow)) {
909 if ((eindex < elen) && (findex < flen)) {
910 if ((fnow > enow) == (fnow > -enow)) {
919 while ((eindex < elen) && (findex < flen)) {
920 if ((fnow > enow) == (fnow > -enow)) {
921 Two_Sum(Q, enow, Qnew, h[hindex]);
924 Two_Sum(Q, fnow, Qnew, h[hindex]);
931 while (eindex < elen) {
932 Two_Sum(Q, enow, Qnew, h[hindex]);
937 while (findex < flen) {
938 Two_Sum(Q, fnow, Qnew, h[hindex]);
968 REAL avirt, bround, around;
969 int eindex, findex, hindex;
975 if ((fnow > enow) == (fnow > -enow)) {
983 if ((eindex < elen) && (findex < flen)) {
984 if ((fnow > enow) == (fnow > -enow)) {
995 while ((eindex < elen) && (findex < flen)) {
996 if ((fnow > enow) == (fnow > -enow)) {
1009 while (eindex < elen) {
1017 while (findex < flen) {
1025 if ((Q != 0.0) || (hindex == 0)) {
1049 REAL avirt, bround, around;
1050 int eindex, findex, hindex;
1056 eindex = findex = 0;
1057 if ((fnow > enow) == (fnow > -enow)) {
1064 if ((eindex < elen) && ((findex >= flen)
1065 || ((fnow > enow) == (fnow > -enow)))) {
1073 for (hindex = 0; hindex < elen + flen - 2; hindex++) {
1074 if ((eindex < elen) && ((findex >= flen)
1075 || ((fnow > enow) == (fnow > -enow)))) {
1110 REAL avirt, bround, around;
1111 int eindex, findex, hindex;
1118 eindex = findex = 0;
1120 if ((fnow > enow) == (fnow > -enow)) {
1127 if ((eindex < elen) && ((findex >= flen)
1128 || ((fnow > enow) == (fnow > -enow)))) {
1136 for (count = 2; count < elen + flen; count++) {
1137 if ((eindex < elen) && ((findex >= flen)
1138 || ((fnow > enow) == (fnow > -enow)))) {
1154 if ((Q != 0.0) || (hindex == 0)) {
1183 REAL avirt, bround, around;
1186 REAL ahi, alo, bhi, blo;
1187 REAL err1, err2, err3;
1192 for (eindex = 1; eindex < elen; eindex++) {
1195 Two_Sum(Q, product0, sum, h[hindex]);
1197 Two_Sum(product1, sum, Q, h[hindex]);
1229 REAL avirt, bround, around;
1232 REAL ahi, alo, bhi, blo;
1233 REAL err1, err2, err3;
1241 for (eindex = 1; eindex < elen; eindex++) {
1244 Two_Sum(Q, product0, sum, hh);
1253 if ((Q != 0.0) || (hindex == 0)) {
1283 for (eindex = elen - 2; eindex >= 0; eindex--) {
1294 for (hindex = bottom + 1; hindex < elen; hindex++) {
1320 for (eindex = 1; eindex < elen; eindex++) {
1354 REAL acx, bcx, acy, bcy;
1356 acx = pa[0] - pc[0];
1357 bcx = pb[0] - pc[0];
1358 acy = pa[1] - pc[1];
1359 bcy = pb[1] - pc[1];
1360 return acx * bcy - acy * bcx;
1365 INEXACT REAL axby1, axcy1, bxcy1, bxay1, cxay1, cxby1;
1366 REAL axby0, axcy0, bxcy0, bxay0, cxay0, cxby0;
1367 REAL aterms[4], bterms[4], cterms[4];
1370 int vlength, wlength;
1373 REAL avirt, bround, around;
1376 REAL ahi, alo, bhi, blo;
1377 REAL err1, err2, err3;
1384 aterms3, aterms[2], aterms[1], aterms[0]);
1385 aterms[3] = aterms3;
1390 bterms3, bterms[2], bterms[1], bterms[0]);
1391 bterms[3] = bterms3;
1396 cterms3, cterms[2], cterms[1], cterms[0]);
1397 cterms[3] = cterms3;
1402 return w[wlength - 1];
1408 REAL acxtail, acytail;
1409 REAL bcxtail, bcytail;
1410 REAL negate, negatetail;
1411 REAL axby[8], bxay[8];
1417 REAL avirt, bround, around;
1420 REAL a0hi, a0lo, a1hi, a1lo, bhi, blo;
1421 REAL err1, err2, err3;
1425 Two_Diff(pa[0], pc[0], acx, acxtail);
1426 Two_Diff(pa[1], pc[1], acy, acytail);
1427 Two_Diff(pb[0], pc[0], bcx, bcxtail);
1428 Two_Diff(pb[1], pc[1], bcy, bcytail);
1431 axby7, axby[6], axby[5], axby[4],
1432 axby[3], axby[2], axby[1], axby[0]);
1435 negatetail = -acytail;
1437 bxay7, bxay[6], bxay[5], bxay[4],
1438 bxay[3], bxay[2], bxay[1], bxay[0]);
1443 return deter[deterlen - 1];
1449 REAL acxtail, acytail, bcxtail, bcytail;
1451 REAL detlefttail, detrighttail;
1453 REAL B[4], C1[8], C2[12], D[16];
1455 int C1length, C2length, Dlength;
1462 REAL avirt, bround, around;
1465 REAL ahi, alo, bhi, blo;
1466 REAL err1, err2, err3;
1470 acx = (
REAL) (pa[0] - pc[0]);
1471 bcx = (
REAL) (pb[0] - pc[0]);
1472 acy = (
REAL) (pa[1] - pc[1]);
1473 bcy = (
REAL) (pb[1] - pc[1]);
1478 Two_Two_Diff(detleft, detlefttail, detright, detrighttail,
1479 B3, B[2], B[1], B[0]);
1484 if ((det >= errbound) || (-det >= errbound)) {
1493 if ((acxtail == 0.0) && (acytail == 0.0)
1494 && (bcxtail == 0.0) && (bcytail == 0.0)) {
1499 det += (acx * bcytail + bcy * acxtail)
1500 - (acy * bcxtail + bcx * acytail);
1501 if ((det >= errbound) || (-det >= errbound)) {
1523 return(D[Dlength - 1]);
1528 REAL detleft, detright, det;
1529 REAL detsum, errbound;
1531 detleft = (pa[0] - pc[0]) * (pb[1] - pc[1]);
1532 detright = (pa[1] - pc[1]) * (pb[0] - pc[0]);
1533 det = detleft - detright;
1535 if (detleft > 0.0) {
1536 if (detright <= 0.0) {
1539 detsum = detleft + detright;
1541 }
else if (detleft < 0.0) {
1542 if (detright >= 0.0) {
1545 detsum = -detleft - detright;
1552 if ((det >= errbound) || (-det >= errbound)) {
1594 adx = pa[0] - pd[0];
1595 bdx = pb[0] - pd[0];
1596 cdx = pc[0] - pd[0];
1597 ady = pa[1] - pd[1];
1598 bdy = pb[1] - pd[1];
1599 cdy = pc[1] - pd[1];
1600 adz = pa[2] - pd[2];
1601 bdz = pb[2] - pd[2];
1602 cdz = pc[2] - pd[2];
1604 return adx * (bdy * cdz - bdz * cdy)
1605 + bdx * (cdy * adz - cdz * ady)
1606 + cdx * (ady * bdz - adz * bdy);
1611 INEXACT REAL axby1, bxcy1, cxdy1, dxay1, axcy1, bxdy1;
1612 INEXACT REAL bxay1, cxby1, dxcy1, axdy1, cxay1, dxby1;
1613 REAL axby0, bxcy0, cxdy0, dxay0, axcy0, bxdy0;
1614 REAL bxay0, cxby0, dxcy0, axdy0, cxay0, dxby0;
1615 REAL ab[4], bc[4], cd[4], da[4], ac[4], bd[4];
1618 REAL abc[12], bcd[12], cda[12], dab[12];
1619 int abclen, bcdlen, cdalen, dablen;
1620 REAL adet[24], bdet[24], cdet[24], ddet[24];
1621 int alen, blen, clen, dlen;
1622 REAL abdet[48], cddet[48];
1629 REAL avirt, bround, around;
1632 REAL ahi, alo, bhi, blo;
1633 REAL err1, err2, err3;
1639 Two_Two_Diff(axby1, axby0, bxay1, bxay0, ab[3], ab[2], ab[1], ab[0]);
1643 Two_Two_Diff(bxcy1, bxcy0, cxby1, cxby0, bc[3], bc[2], bc[1], bc[0]);
1647 Two_Two_Diff(cxdy1, cxdy0, dxcy1, dxcy0, cd[3], cd[2], cd[1], cd[0]);
1651 Two_Two_Diff(dxay1, dxay0, axdy1, axdy0, da[3], da[2], da[1], da[0]);
1655 Two_Two_Diff(axcy1, axcy0, cxay1, cxay0, ac[3], ac[2], ac[1], ac[0]);
1659 Two_Two_Diff(bxdy1, bxdy0, dxby1, dxby0, bd[3], bd[2], bd[1], bd[0]);
1665 for (i = 0; i < 4; i++) {
1683 return deter[deterlen - 1];
1688 INEXACT REAL adx, ady, adz, bdx, bdy, bdz, cdx, cdy, cdz;
1689 REAL adxtail, adytail, adztail;
1690 REAL bdxtail, bdytail, bdztail;
1691 REAL cdxtail, cdytail, cdztail;
1692 REAL negate, negatetail;
1693 INEXACT REAL axby7, bxcy7, axcy7, bxay7, cxby7, cxay7;
1694 REAL axby[8], bxcy[8], axcy[8], bxay[8], cxby[8], cxay[8];
1695 REAL temp16[16], temp32[32], temp32t[32];
1696 int temp16len, temp32len, temp32tlen;
1697 REAL adet[64], bdet[64], cdet[64];
1698 int alen, blen, clen;
1705 REAL avirt, bround, around;
1708 REAL a0hi, a0lo, a1hi, a1lo, bhi, blo;
1709 REAL err1, err2, err3;
1713 Two_Diff(pa[0], pd[0], adx, adxtail);
1714 Two_Diff(pa[1], pd[1], ady, adytail);
1715 Two_Diff(pa[2], pd[2], adz, adztail);
1716 Two_Diff(pb[0], pd[0], bdx, bdxtail);
1717 Two_Diff(pb[1], pd[1], bdy, bdytail);
1718 Two_Diff(pb[2], pd[2], bdz, bdztail);
1719 Two_Diff(pc[0], pd[0], cdx, cdxtail);
1720 Two_Diff(pc[1], pd[1], cdy, cdytail);
1721 Two_Diff(pc[2], pd[2], cdz, cdztail);
1724 axby7, axby[6], axby[5], axby[4],
1725 axby[3], axby[2], axby[1], axby[0]);
1728 negatetail = -adytail;
1730 bxay7, bxay[6], bxay[5], bxay[4],
1731 bxay[3], bxay[2], bxay[1], bxay[0]);
1734 bxcy7, bxcy[6], bxcy[5], bxcy[4],
1735 bxcy[3], bxcy[2], bxcy[1], bxcy[0]);
1738 negatetail = -bdytail;
1740 cxby7, cxby[6], cxby[5], cxby[4],
1741 cxby[3], cxby[2], cxby[1], cxby[0]);
1744 cxay7, cxay[6], cxay[5], cxay[4],
1745 cxay[3], cxay[2], cxay[1], cxay[0]);
1748 negatetail = -cdytail;
1750 axcy7, axcy[6], axcy[5], axcy[4],
1751 axcy[3], axcy[2], axcy[1], axcy[0]);
1775 return deter[deterlen - 1];
1780 INEXACT REAL adx, bdx, cdx, ady, bdy, cdy, adz, bdz, cdz;
1783 INEXACT REAL bdxcdy1, cdxbdy1, cdxady1, adxcdy1, adxbdy1, bdxady1;
1784 REAL bdxcdy0, cdxbdy0, cdxady0, adxcdy0, adxbdy0, bdxady0;
1785 REAL bc[4], ca[4], ab[4];
1787 REAL adet[8], bdet[8], cdet[8];
1788 int alen, blen, clen;
1791 REAL *finnow, *finother, *finswap;
1792 REAL fin1[192], fin2[192];
1796 REAL adxtail, bdxtail, cdxtail;
1797 REAL adytail, bdytail, cdytail;
1798 REAL adztail, bdztail, cdztail;
1802 REAL at_b[4], at_c[4], bt_c[4], bt_a[4], ct_a[4], ct_b[4];
1803 int at_blen, at_clen, bt_clen, bt_alen, ct_alen, ct_blen;
1806 REAL bdxt_cdy0, cdxt_bdy0, cdxt_ady0;
1807 REAL adxt_cdy0, adxt_bdy0, bdxt_ady0;
1810 REAL bdyt_cdx0, cdyt_bdx0, cdyt_adx0;
1811 REAL adyt_cdx0, adyt_bdx0, bdyt_adx0;
1812 REAL bct[8], cat[8], abt[8];
1813 int bctlen, catlen, abtlen;
1816 REAL bdxt_cdyt0, cdxt_bdyt0, cdxt_adyt0;
1817 REAL adxt_cdyt0, adxt_bdyt0, bdxt_adyt0;
1820 int vlength, wlength;
1824 REAL avirt, bround, around;
1827 REAL ahi, alo, bhi, blo;
1828 REAL err1, err2, err3;
1833 adx = (
REAL) (pa[0] - pd[0]);
1834 bdx = (
REAL) (pb[0] - pd[0]);
1835 cdx = (
REAL) (pc[0] - pd[0]);
1836 ady = (
REAL) (pa[1] - pd[1]);
1837 bdy = (
REAL) (pb[1] - pd[1]);
1838 cdy = (
REAL) (pc[1] - pd[1]);
1839 adz = (
REAL) (pa[2] - pd[2]);
1840 bdz = (
REAL) (pb[2] - pd[2]);
1841 cdz = (
REAL) (pc[2] - pd[2]);
1845 Two_Two_Diff(bdxcdy1, bdxcdy0, cdxbdy1, cdxbdy0, bc3, bc[2], bc[1], bc[0]);
1851 Two_Two_Diff(cdxady1, cdxady0, adxcdy1, adxcdy0, ca3, ca[2], ca[1], ca[0]);
1857 Two_Two_Diff(adxbdy1, adxbdy0, bdxady1, bdxady0, ab3, ab[2], ab[1], ab[0]);
1866 if ((det >= errbound) || (-det >= errbound)) {
1880 if ((adxtail == 0.0) && (bdxtail == 0.0) && (cdxtail == 0.0)
1881 && (adytail == 0.0) && (bdytail == 0.0) && (cdytail == 0.0)
1882 && (adztail == 0.0) && (bdztail == 0.0) && (cdztail == 0.0)) {
1887 det += (adz * ((bdx * cdytail + cdy * bdxtail)
1888 - (bdy * cdxtail + cdx * bdytail))
1889 + adztail * (bdx * cdy - bdy * cdx))
1890 + (bdz * ((cdx * adytail + ady * cdxtail)
1891 - (cdy * adxtail + adx * cdytail))
1892 + bdztail * (cdx * ady - cdy * adx))
1893 + (cdz * ((adx * bdytail + bdy * adxtail)
1894 - (ady * bdxtail + bdx * adytail))
1895 + cdztail * (adx * bdy - ady * bdx));
1896 if ((det >= errbound) || (-det >= errbound)) {
1903 if (adxtail == 0.0) {
1904 if (adytail == 0.0) {
1912 at_b[1] = at_blarge;
1915 at_c[1] = at_clarge;
1919 if (adytail == 0.0) {
1921 at_b[1] = at_blarge;
1925 at_c[1] = at_clarge;
1930 Two_Two_Diff(adxt_bdy1, adxt_bdy0, adyt_bdx1, adyt_bdx0,
1931 at_blarge, at_b[2], at_b[1], at_b[0]);
1932 at_b[3] = at_blarge;
1936 Two_Two_Diff(adyt_cdx1, adyt_cdx0, adxt_cdy1, adxt_cdy0,
1937 at_clarge, at_c[2], at_c[1], at_c[0]);
1938 at_c[3] = at_clarge;
1942 if (bdxtail == 0.0) {
1943 if (bdytail == 0.0) {
1951 bt_c[1] = bt_clarge;
1954 bt_a[1] = bt_alarge;
1958 if (bdytail == 0.0) {
1960 bt_c[1] = bt_clarge;
1964 bt_a[1] = bt_alarge;
1969 Two_Two_Diff(bdxt_cdy1, bdxt_cdy0, bdyt_cdx1, bdyt_cdx0,
1970 bt_clarge, bt_c[2], bt_c[1], bt_c[0]);
1971 bt_c[3] = bt_clarge;
1975 Two_Two_Diff(bdyt_adx1, bdyt_adx0, bdxt_ady1, bdxt_ady0,
1976 bt_alarge, bt_a[2], bt_a[1], bt_a[0]);
1977 bt_a[3] = bt_alarge;
1981 if (cdxtail == 0.0) {
1982 if (cdytail == 0.0) {
1990 ct_a[1] = ct_alarge;
1993 ct_b[1] = ct_blarge;
1997 if (cdytail == 0.0) {
1999 ct_a[1] = ct_alarge;
2003 ct_b[1] = ct_blarge;
2008 Two_Two_Diff(cdxt_ady1, cdxt_ady0, cdyt_adx1, cdyt_adx0,
2009 ct_alarge, ct_a[2], ct_a[1], ct_a[0]);
2010 ct_a[3] = ct_alarge;
2014 Two_Two_Diff(cdyt_bdx1, cdyt_bdx0, cdxt_bdy1, cdxt_bdy0,
2015 ct_blarge, ct_b[2], ct_b[1], ct_b[0]);
2016 ct_b[3] = ct_blarge;
2025 finswap = finnow; finnow = finother; finother = finswap;
2031 finswap = finnow; finnow = finother; finother = finswap;
2037 finswap = finnow; finnow = finother; finother = finswap;
2039 if (adztail != 0.0) {
2043 finswap = finnow; finnow = finother; finother = finswap;
2045 if (bdztail != 0.0) {
2049 finswap = finnow; finnow = finother; finother = finswap;
2051 if (cdztail != 0.0) {
2055 finswap = finnow; finnow = finother; finother = finswap;
2058 if (adxtail != 0.0) {
2059 if (bdytail != 0.0) {
2060 Two_Product(adxtail, bdytail, adxt_bdyt1, adxt_bdyt0);
2065 finswap = finnow; finnow = finother; finother = finswap;
2066 if (cdztail != 0.0) {
2067 Two_One_Product(adxt_bdyt1, adxt_bdyt0, cdztail, u3, u[2], u[1], u[0]);
2071 finswap = finnow; finnow = finother; finother = finswap;
2074 if (cdytail != 0.0) {
2076 Two_Product(negate, cdytail, adxt_cdyt1, adxt_cdyt0);
2081 finswap = finnow; finnow = finother; finother = finswap;
2082 if (bdztail != 0.0) {
2083 Two_One_Product(adxt_cdyt1, adxt_cdyt0, bdztail, u3, u[2], u[1], u[0]);
2087 finswap = finnow; finnow = finother; finother = finswap;
2091 if (bdxtail != 0.0) {
2092 if (cdytail != 0.0) {
2093 Two_Product(bdxtail, cdytail, bdxt_cdyt1, bdxt_cdyt0);
2098 finswap = finnow; finnow = finother; finother = finswap;
2099 if (adztail != 0.0) {
2100 Two_One_Product(bdxt_cdyt1, bdxt_cdyt0, adztail, u3, u[2], u[1], u[0]);
2104 finswap = finnow; finnow = finother; finother = finswap;
2107 if (adytail != 0.0) {
2109 Two_Product(negate, adytail, bdxt_adyt1, bdxt_adyt0);
2114 finswap = finnow; finnow = finother; finother = finswap;
2115 if (cdztail != 0.0) {
2116 Two_One_Product(bdxt_adyt1, bdxt_adyt0, cdztail, u3, u[2], u[1], u[0]);
2120 finswap = finnow; finnow = finother; finother = finswap;
2124 if (cdxtail != 0.0) {
2125 if (adytail != 0.0) {
2126 Two_Product(cdxtail, adytail, cdxt_adyt1, cdxt_adyt0);
2131 finswap = finnow; finnow = finother; finother = finswap;
2132 if (bdztail != 0.0) {
2133 Two_One_Product(cdxt_adyt1, cdxt_adyt0, bdztail, u3, u[2], u[1], u[0]);
2137 finswap = finnow; finnow = finother; finother = finswap;
2140 if (bdytail != 0.0) {
2142 Two_Product(negate, bdytail, cdxt_bdyt1, cdxt_bdyt0);
2147 finswap = finnow; finnow = finother; finother = finswap;
2148 if (adztail != 0.0) {
2149 Two_One_Product(cdxt_bdyt1, cdxt_bdyt0, adztail, u3, u[2], u[1], u[0]);
2153 finswap = finnow; finnow = finother; finother = finswap;
2158 if (adztail != 0.0) {
2162 finswap = finnow; finnow = finother; finother = finswap;
2164 if (bdztail != 0.0) {
2168 finswap = finnow; finnow = finother; finother = finswap;
2170 if (cdztail != 0.0) {
2174 finswap = finnow; finnow = finother; finother = finswap;
2177 return finnow[finlength - 1];
2180 #ifdef USE_CGAL_PREDICATES
2185 - cgal_pred_obj.orientation_3_object()
2186 (Point(pa[0], pa[1], pa[2]),
2187 Point(pb[0], pb[1], pb[2]),
2188 Point(pc[0], pc[1], pc[2]),
2189 Point(pd[0], pd[1], pd[2]));
2196 REAL adx, bdx, cdx, ady, bdy, cdy, adz, bdz, cdz;
2197 REAL bdxcdy, cdxbdy, cdxady, adxcdy, adxbdy, bdxady;
2201 adx = pa[0] - pd[0];
2202 ady = pa[1] - pd[1];
2203 adz = pa[2] - pd[2];
2204 bdx = pb[0] - pd[0];
2205 bdy = pb[1] - pd[1];
2206 bdz = pb[2] - pd[2];
2207 cdx = pc[0] - pd[0];
2208 cdy = pc[1] - pd[1];
2209 cdz = pc[2] - pd[2];
2220 det = adz * (bdxcdy - cdxbdy)
2221 + bdz * (cdxady - adxcdy)
2222 + cdz * (adxbdy - bdxady);
2235 REAL permanent, errbound;
2241 if ((det > errbound) || (-det > errbound)) {
2248 #endif // #ifdef USE_CGAL_PREDICATES
2278 REAL adx, ady, bdx, bdy, cdx, cdy;
2279 REAL abdet, bcdet, cadet;
2280 REAL alift, blift, clift;
2282 adx = pa[0] - pd[0];
2283 ady = pa[1] - pd[1];
2284 bdx = pb[0] - pd[0];
2285 bdy = pb[1] - pd[1];
2286 cdx = pc[0] - pd[0];
2287 cdy = pc[1] - pd[1];
2289 abdet = adx * bdy - bdx * ady;
2290 bcdet = bdx * cdy - cdx * bdy;
2291 cadet = cdx * ady - adx * cdy;
2292 alift = adx * adx + ady * ady;
2293 blift = bdx * bdx + bdy * bdy;
2294 clift = cdx * cdx + cdy * cdy;
2296 return alift * bcdet + blift * cadet + clift * abdet;
2301 INEXACT REAL axby1, bxcy1, cxdy1, dxay1, axcy1, bxdy1;
2302 INEXACT REAL bxay1, cxby1, dxcy1, axdy1, cxay1, dxby1;
2303 REAL axby0, bxcy0, cxdy0, dxay0, axcy0, bxdy0;
2304 REAL bxay0, cxby0, dxcy0, axdy0, cxay0, dxby0;
2305 REAL ab[4], bc[4], cd[4], da[4], ac[4], bd[4];
2308 REAL abc[12], bcd[12], cda[12], dab[12];
2309 int abclen, bcdlen, cdalen, dablen;
2310 REAL det24x[24], det24y[24], det48x[48], det48y[48];
2312 REAL adet[96], bdet[96], cdet[96], ddet[96];
2313 int alen, blen, clen, dlen;
2314 REAL abdet[192], cddet[192];
2321 REAL avirt, bround, around;
2324 REAL ahi, alo, bhi, blo;
2325 REAL err1, err2, err3;
2331 Two_Two_Diff(axby1, axby0, bxay1, bxay0, ab[3], ab[2], ab[1], ab[0]);
2335 Two_Two_Diff(bxcy1, bxcy0, cxby1, cxby0, bc[3], bc[2], bc[1], bc[0]);
2339 Two_Two_Diff(cxdy1, cxdy0, dxcy1, dxcy0, cd[3], cd[2], cd[1], cd[0]);
2343 Two_Two_Diff(dxay1, dxay0, axdy1, axdy0, da[3], da[2], da[1], da[0]);
2347 Two_Two_Diff(axcy1, axcy0, cxay1, cxay0, ac[3], ac[2], ac[1], ac[0]);
2351 Two_Two_Diff(bxdy1, bxdy0, dxby1, dxby0, bd[3], bd[2], bd[1], bd[0]);
2357 for (i = 0; i < 4; i++) {
2394 return deter[deterlen - 1];
2400 REAL adxtail, bdxtail, cdxtail;
2401 REAL adytail, bdytail, cdytail;
2402 REAL negate, negatetail;
2403 INEXACT REAL axby7, bxcy7, axcy7, bxay7, cxby7, cxay7;
2404 REAL axby[8], bxcy[8], axcy[8], bxay[8], cxby[8], cxay[8];
2407 REAL detx[32], detxx[64], detxt[32], detxxt[64], detxtxt[64];
2408 int xlen, xxlen, xtlen, xxtlen, xtxtlen;
2411 REAL dety[32], detyy[64], detyt[32], detyyt[64], detytyt[64];
2412 int ylen, yylen, ytlen, yytlen, ytytlen;
2413 REAL y1[128], y2[192];
2415 REAL adet[384], bdet[384], cdet[384], abdet[768], deter[1152];
2416 int alen, blen, clen, ablen, deterlen;
2420 REAL avirt, bround, around;
2423 REAL a0hi, a0lo, a1hi, a1lo, bhi, blo;
2424 REAL err1, err2, err3;
2428 Two_Diff(pa[0], pd[0], adx, adxtail);
2429 Two_Diff(pa[1], pd[1], ady, adytail);
2430 Two_Diff(pb[0], pd[0], bdx, bdxtail);
2431 Two_Diff(pb[1], pd[1], bdy, bdytail);
2432 Two_Diff(pc[0], pd[0], cdx, cdxtail);
2433 Two_Diff(pc[1], pd[1], cdy, cdytail);
2436 axby7, axby[6], axby[5], axby[4],
2437 axby[3], axby[2], axby[1], axby[0]);
2440 negatetail = -adytail;
2442 bxay7, bxay[6], bxay[5], bxay[4],
2443 bxay[3], bxay[2], bxay[1], bxay[0]);
2446 bxcy7, bxcy[6], bxcy[5], bxcy[4],
2447 bxcy[3], bxcy[2], bxcy[1], bxcy[0]);
2450 negatetail = -bdytail;
2452 cxby7, cxby[6], cxby[5], cxby[4],
2453 cxby[3], cxby[2], cxby[1], cxby[0]);
2456 cxay7, cxay[6], cxay[5], cxay[4],
2457 cxay[3], cxay[2], cxay[1], cxay[0]);
2460 negatetail = -cdytail;
2462 axcy7, axcy[6], axcy[5], axcy[4],
2463 axcy[3], axcy[2], axcy[1], axcy[0]);
2473 for (i = 0; i < xxtlen; i++) {
2484 for (i = 0; i < yytlen; i++) {
2500 for (i = 0; i < xxtlen; i++) {
2511 for (i = 0; i < yytlen; i++) {
2527 for (i = 0; i < xxtlen; i++) {
2538 for (i = 0; i < yytlen; i++) {
2550 return deter[deterlen - 1];
2558 INEXACT REAL bdxcdy1, cdxbdy1, cdxady1, adxcdy1, adxbdy1, bdxady1;
2559 REAL bdxcdy0, cdxbdy0, cdxady0, adxcdy0, adxbdy0, bdxady0;
2560 REAL bc[4], ca[4], ab[4];
2562 REAL axbc[8], axxbc[16], aybc[8], ayybc[16], adet[32];
2563 int axbclen, axxbclen, aybclen, ayybclen, alen;
2564 REAL bxca[8], bxxca[16], byca[8], byyca[16], bdet[32];
2565 int bxcalen, bxxcalen, bycalen, byycalen, blen;
2566 REAL cxab[8], cxxab[16], cyab[8], cyyab[16], cdet[32];
2567 int cxablen, cxxablen, cyablen, cyyablen, clen;
2570 REAL fin1[1152], fin2[1152];
2571 REAL *finnow, *finother, *finswap;
2574 REAL adxtail, bdxtail, cdxtail, adytail, bdytail, cdytail;
2575 INEXACT REAL adxadx1, adyady1, bdxbdx1, bdybdy1, cdxcdx1, cdycdy1;
2576 REAL adxadx0, adyady0, bdxbdx0, bdybdy0, cdxcdx0, cdycdy0;
2577 REAL aa[4], bb[4], cc[4];
2583 REAL temp8[8], temp16a[16], temp16b[16], temp16c[16];
2584 REAL temp32a[32], temp32b[32], temp48[48], temp64[64];
2585 int temp8len, temp16alen, temp16blen, temp16clen;
2586 int temp32alen, temp32blen, temp48len, temp64len;
2587 REAL axtbb[8], axtcc[8], aytbb[8], aytcc[8];
2588 int axtbblen, axtcclen, aytbblen, aytcclen;
2589 REAL bxtaa[8], bxtcc[8], bytaa[8], bytcc[8];
2590 int bxtaalen, bxtcclen, bytaalen, bytcclen;
2591 REAL cxtaa[8], cxtbb[8], cytaa[8], cytbb[8];
2592 int cxtaalen, cxtbblen, cytaalen, cytbblen;
2593 REAL axtbc[8], aytbc[8], bxtca[8], bytca[8], cxtab[8], cytab[8];
2594 int axtbclen, aytbclen, bxtcalen, bytcalen, cxtablen, cytablen;
2595 REAL axtbct[16], aytbct[16], bxtcat[16], bytcat[16], cxtabt[16], cytabt[16];
2596 int axtbctlen, aytbctlen, bxtcatlen, bytcatlen, cxtabtlen, cytabtlen;
2597 REAL axtbctt[8], aytbctt[8], bxtcatt[8];
2598 REAL bytcatt[8], cxtabtt[8], cytabtt[8];
2599 int axtbcttlen, aytbcttlen, bxtcattlen, bytcattlen, cxtabttlen, cytabttlen;
2600 REAL abt[8], bct[8], cat[8];
2601 int abtlen, bctlen, catlen;
2602 REAL abtt[4], bctt[4], catt[4];
2603 int abttlen, bcttlen, cattlen;
2608 REAL avirt, bround, around;
2611 REAL ahi, alo, bhi, blo;
2612 REAL err1, err2, err3;
2617 axtbclen = aytbclen = bxtcalen = bytcalen = cxtablen = cytablen = 0;
2619 adx = (
REAL) (pa[0] - pd[0]);
2620 bdx = (
REAL) (pb[0] - pd[0]);
2621 cdx = (
REAL) (pc[0] - pd[0]);
2622 ady = (
REAL) (pa[1] - pd[1]);
2623 bdy = (
REAL) (pb[1] - pd[1]);
2624 cdy = (
REAL) (pc[1] - pd[1]);
2628 Two_Two_Diff(bdxcdy1, bdxcdy0, cdxbdy1, cdxbdy0, bc3, bc[2], bc[1], bc[0]);
2638 Two_Two_Diff(cdxady1, cdxady0, adxcdy1, adxcdy0, ca3, ca[2], ca[1], ca[0]);
2648 Two_Two_Diff(adxbdy1, adxbdy0, bdxady1, bdxady0, ab3, ab[2], ab[1], ab[0]);
2661 if ((det >= errbound) || (-det >= errbound)) {
2671 if ((adxtail == 0.0) && (bdxtail == 0.0) && (cdxtail == 0.0)
2672 && (adytail == 0.0) && (bdytail == 0.0) && (cdytail == 0.0)) {
2677 det += ((adx * adx + ady * ady) * ((bdx * cdytail + cdy * bdxtail)
2678 - (bdy * cdxtail + cdx * bdytail))
2679 + 2.0 * (adx * adxtail + ady * adytail) * (bdx * cdy - bdy * cdx))
2680 + ((bdx * bdx + bdy * bdy) * ((cdx * adytail + ady * cdxtail)
2681 - (cdy * adxtail + adx * cdytail))
2682 + 2.0 * (bdx * bdxtail + bdy * bdytail) * (cdx * ady - cdy * adx))
2683 + ((cdx * cdx + cdy * cdy) * ((adx * bdytail + bdy * adxtail)
2684 - (ady * bdxtail + bdx * adytail))
2685 + 2.0 * (cdx * cdxtail + cdy * cdytail) * (adx * bdy - ady * bdx));
2686 if ((det >= errbound) || (-det >= errbound)) {
2693 if ((bdxtail != 0.0) || (bdytail != 0.0)
2694 || (cdxtail != 0.0) || (cdytail != 0.0)) {
2695 Square(adx, adxadx1, adxadx0);
2696 Square(ady, adyady1, adyady0);
2697 Two_Two_Sum(adxadx1, adxadx0, adyady1, adyady0, aa3, aa[2], aa[1], aa[0]);
2700 if ((cdxtail != 0.0) || (cdytail != 0.0)
2701 || (adxtail != 0.0) || (adytail != 0.0)) {
2702 Square(bdx, bdxbdx1, bdxbdx0);
2703 Square(bdy, bdybdy1, bdybdy0);
2704 Two_Two_Sum(bdxbdx1, bdxbdx0, bdybdy1, bdybdy0, bb3, bb[2], bb[1], bb[0]);
2707 if ((adxtail != 0.0) || (adytail != 0.0)
2708 || (bdxtail != 0.0) || (bdytail != 0.0)) {
2709 Square(cdx, cdxcdx1, cdxcdx0);
2710 Square(cdy, cdycdy1, cdycdy0);
2711 Two_Two_Sum(cdxcdx1, cdxcdx0, cdycdy1, cdycdy0, cc3, cc[2], cc[1], cc[0]);
2715 if (adxtail != 0.0) {
2727 temp16blen, temp16b, temp32a);
2729 temp32alen, temp32a, temp48);
2732 finswap = finnow; finnow = finother; finother = finswap;
2734 if (adytail != 0.0) {
2746 temp16blen, temp16b, temp32a);
2748 temp32alen, temp32a, temp48);
2751 finswap = finnow; finnow = finother; finother = finswap;
2753 if (bdxtail != 0.0) {
2765 temp16blen, temp16b, temp32a);
2767 temp32alen, temp32a, temp48);
2770 finswap = finnow; finnow = finother; finother = finswap;
2772 if (bdytail != 0.0) {
2784 temp16blen, temp16b, temp32a);
2786 temp32alen, temp32a, temp48);
2789 finswap = finnow; finnow = finother; finother = finswap;
2791 if (cdxtail != 0.0) {
2803 temp16blen, temp16b, temp32a);
2805 temp32alen, temp32a, temp48);
2808 finswap = finnow; finnow = finother; finother = finswap;
2810 if (cdytail != 0.0) {
2822 temp16blen, temp16b, temp32a);
2824 temp32alen, temp32a, temp48);
2827 finswap = finnow; finnow = finother; finother = finswap;
2830 if ((adxtail != 0.0) || (adytail != 0.0)) {
2831 if ((bdxtail != 0.0) || (bdytail != 0.0)
2832 || (cdxtail != 0.0) || (cdytail != 0.0)) {
2835 Two_Two_Sum(ti1, ti0, tj1, tj0, u3, u[2], u[1], u[0]);
2841 Two_Two_Sum(ti1, ti0, tj1, tj0, v3, v[2], v[1], v[0]);
2847 Two_Two_Diff(ti1, ti0, tj1, tj0, bctt3, bctt[2], bctt[1], bctt[0]);
2857 if (adxtail != 0.0) {
2863 temp32alen, temp32a, temp48);
2866 finswap = finnow; finnow = finother; finother = finswap;
2867 if (bdytail != 0.0) {
2873 finswap = finnow; finnow = finother; finother = finswap;
2875 if (cdytail != 0.0) {
2881 finswap = finnow; finnow = finother; finother = finswap;
2892 temp16blen, temp16b, temp32b);
2894 temp32blen, temp32b, temp64);
2897 finswap = finnow; finnow = finother; finother = finswap;
2899 if (adytail != 0.0) {
2905 temp32alen, temp32a, temp48);
2908 finswap = finnow; finnow = finother; finother = finswap;
2919 temp16blen, temp16b, temp32b);
2921 temp32blen, temp32b, temp64);
2924 finswap = finnow; finnow = finother; finother = finswap;
2927 if ((bdxtail != 0.0) || (bdytail != 0.0)) {
2928 if ((cdxtail != 0.0) || (cdytail != 0.0)
2929 || (adxtail != 0.0) || (adytail != 0.0)) {
2932 Two_Two_Sum(ti1, ti0, tj1, tj0, u3, u[2], u[1], u[0]);
2938 Two_Two_Sum(ti1, ti0, tj1, tj0, v3, v[2], v[1], v[0]);
2944 Two_Two_Diff(ti1, ti0, tj1, tj0, catt3, catt[2], catt[1], catt[0]);
2954 if (bdxtail != 0.0) {
2960 temp32alen, temp32a, temp48);
2963 finswap = finnow; finnow = finother; finother = finswap;
2964 if (cdytail != 0.0) {
2970 finswap = finnow; finnow = finother; finother = finswap;
2972 if (adytail != 0.0) {
2978 finswap = finnow; finnow = finother; finother = finswap;
2989 temp16blen, temp16b, temp32b);
2991 temp32blen, temp32b, temp64);
2994 finswap = finnow; finnow = finother; finother = finswap;
2996 if (bdytail != 0.0) {
3002 temp32alen, temp32a, temp48);
3005 finswap = finnow; finnow = finother; finother = finswap;
3016 temp16blen, temp16b, temp32b);
3018 temp32blen, temp32b, temp64);
3021 finswap = finnow; finnow = finother; finother = finswap;
3024 if ((cdxtail != 0.0) || (cdytail != 0.0)) {
3025 if ((adxtail != 0.0) || (adytail != 0.0)
3026 || (bdxtail != 0.0) || (bdytail != 0.0)) {
3029 Two_Two_Sum(ti1, ti0, tj1, tj0, u3, u[2], u[1], u[0]);
3035 Two_Two_Sum(ti1, ti0, tj1, tj0, v3, v[2], v[1], v[0]);
3041 Two_Two_Diff(ti1, ti0, tj1, tj0, abtt3, abtt[2], abtt[1], abtt[0]);
3051 if (cdxtail != 0.0) {
3057 temp32alen, temp32a, temp48);
3060 finswap = finnow; finnow = finother; finother = finswap;
3061 if (adytail != 0.0) {
3067 finswap = finnow; finnow = finother; finother = finswap;
3069 if (bdytail != 0.0) {
3075 finswap = finnow; finnow = finother; finother = finswap;
3086 temp16blen, temp16b, temp32b);
3088 temp32blen, temp32b, temp64);
3091 finswap = finnow; finnow = finother; finother = finswap;
3093 if (cdytail != 0.0) {
3099 temp32alen, temp32a, temp48);
3102 finswap = finnow; finnow = finother; finother = finswap;
3113 temp16blen, temp16b, temp32b);
3115 temp32blen, temp32b, temp64);
3118 finswap = finnow; finnow = finother; finother = finswap;
3122 return finnow[finlength - 1];
3127 REAL adx, bdx, cdx, ady, bdy, cdy;
3128 REAL bdxcdy, cdxbdy, cdxady, adxcdy, adxbdy, bdxady;
3129 REAL alift, blift, clift;
3131 REAL permanent, errbound;
3133 adx = pa[0] - pd[0];
3134 bdx = pb[0] - pd[0];
3135 cdx = pc[0] - pd[0];
3136 ady = pa[1] - pd[1];
3137 bdy = pb[1] - pd[1];
3138 cdy = pc[1] - pd[1];
3142 alift = adx * adx + ady * ady;
3146 blift = bdx * bdx + bdy * bdy;
3150 clift = cdx * cdx + cdy * cdy;
3152 det = alift * (bdxcdy - cdxbdy)
3153 + blift * (cdxady - adxcdy)
3154 + clift * (adxbdy - bdxady);
3160 if ((det > errbound) || (-det > errbound)) {
3196 REAL aex, bex, cex, dex;
3197 REAL aey, bey, cey, dey;
3198 REAL aez, bez, cez, dez;
3199 REAL alift, blift, clift, dlift;
3200 REAL ab, bc, cd, da, ac, bd;
3201 REAL abc, bcd, cda, dab;
3203 aex = pa[0] - pe[0];
3204 bex = pb[0] - pe[0];
3205 cex = pc[0] - pe[0];
3206 dex = pd[0] - pe[0];
3207 aey = pa[1] - pe[1];
3208 bey = pb[1] - pe[1];
3209 cey = pc[1] - pe[1];
3210 dey = pd[1] - pe[1];
3211 aez = pa[2] - pe[2];
3212 bez = pb[2] - pe[2];
3213 cez = pc[2] - pe[2];
3214 dez = pd[2] - pe[2];
3216 ab = aex * bey - bex * aey;
3217 bc = bex * cey - cex * bey;
3218 cd = cex * dey - dex * cey;
3219 da = dex * aey - aex * dey;
3221 ac = aex * cey - cex * aey;
3222 bd = bex * dey - dex * bey;
3224 abc = aez * bc - bez * ac + cez * ab;
3225 bcd = bez * cd - cez * bd + dez * bc;
3226 cda = cez * da + dez * ac + aez * cd;
3227 dab = dez * ab + aez * bd + bez * da;
3229 alift = aex * aex + aey * aey + aez * aez;
3230 blift = bex * bex + bey * bey + bez * bez;
3231 clift = cex * cex + cey * cey + cez * cez;
3232 dlift = dex * dex + dey * dey + dez * dez;
3234 return (dlift * abc - clift * dab) + (blift * cda - alift * bcd);
3243 REAL axby0, bxcy0, cxdy0, dxey0, exay0;
3244 REAL bxay0, cxby0, dxcy0, exdy0, axey0;
3245 REAL axcy0, bxdy0, cxey0, dxay0, exby0;
3246 REAL cxay0, dxby0, excy0, axdy0, bxey0;
3247 REAL ab[4], bc[4], cd[4], de[4], ea[4];
3248 REAL ac[4], bd[4], ce[4], da[4], eb[4];
3249 REAL temp8a[8], temp8b[8], temp16[16];
3250 int temp8alen, temp8blen, temp16len;
3251 REAL abc[24], bcd[24], cde[24], dea[24], eab[24];
3252 REAL abd[24], bce[24], cda[24], deb[24], eac[24];
3253 int abclen, bcdlen, cdelen, dealen, eablen;
3254 int abdlen, bcelen, cdalen, deblen, eaclen;
3255 REAL temp48a[48], temp48b[48];
3256 int temp48alen, temp48blen;
3257 REAL abcd[96], bcde[96], cdea[96], deab[96], eabc[96];
3258 int abcdlen, bcdelen, cdealen, deablen, eabclen;
3260 REAL det384x[384], det384y[384], det384z[384];
3261 int xlen, ylen, zlen;
3264 REAL adet[1152], bdet[1152], cdet[1152], ddet[1152], edet[1152];
3265 int alen, blen, clen, dlen, elen;
3266 REAL abdet[2304], cddet[2304], cdedet[3456];
3273 REAL avirt, bround, around;
3276 REAL ahi, alo, bhi, blo;
3277 REAL err1, err2, err3;
3284 Two_Two_Diff(axby1, axby0, bxay1, bxay0, ab[3], ab[2], ab[1], ab[0]);
3288 Two_Two_Diff(bxcy1, bxcy0, cxby1, cxby0, bc[3], bc[2], bc[1], bc[0]);
3292 Two_Two_Diff(cxdy1, cxdy0, dxcy1, dxcy0, cd[3], cd[2], cd[1], cd[0]);
3296 Two_Two_Diff(dxey1, dxey0, exdy1, exdy0, de[3], de[2], de[1], de[0]);
3300 Two_Two_Diff(exay1, exay0, axey1, axey0, ea[3], ea[2], ea[1], ea[0]);
3304 Two_Two_Diff(axcy1, axcy0, cxay1, cxay0, ac[3], ac[2], ac[1], ac[0]);
3308 Two_Two_Diff(bxdy1, bxdy0, dxby1, dxby0, bd[3], bd[2], bd[1], bd[0]);
3312 Two_Two_Diff(cxey1, cxey0, excy1, excy0, ce[3], ce[2], ce[1], ce[0]);
3316 Two_Two_Diff(dxay1, dxay0, axdy1, axdy0, da[3], da[2], da[1], da[0]);
3320 Two_Two_Diff(exby1, exby0, bxey1, bxey0, eb[3], eb[2], eb[1], eb[0]);
3404 for (i = 0; i < temp48blen; i++) {
3405 temp48b[i] = -temp48b[i];
3408 temp48blen, temp48b, bcde);
3420 for (i = 0; i < temp48blen; i++) {
3421 temp48b[i] = -temp48b[i];
3424 temp48blen, temp48b, cdea);
3436 for (i = 0; i < temp48blen; i++) {
3437 temp48b[i] = -temp48b[i];
3440 temp48blen, temp48b, deab);
3452 for (i = 0; i < temp48blen; i++) {
3453 temp48b[i] = -temp48b[i];
3456 temp48blen, temp48b, eabc);
3468 for (i = 0; i < temp48blen; i++) {
3469 temp48b[i] = -temp48b[i];
3472 temp48blen, temp48b, abcd);
3487 return deter[deterlen - 1];
3492 INEXACT REAL aex, bex, cex, dex, aey, bey, cey, dey, aez, bez, cez, dez;
3493 REAL aextail, bextail, cextail, dextail;
3494 REAL aeytail, beytail, ceytail, deytail;
3495 REAL aeztail, beztail, ceztail, deztail;
3496 REAL negate, negatetail;
3497 INEXACT REAL axby7, bxcy7, cxdy7, dxay7, axcy7, bxdy7;
3498 INEXACT REAL bxay7, cxby7, dxcy7, axdy7, cxay7, dxby7;
3499 REAL axby[8], bxcy[8], cxdy[8], dxay[8], axcy[8], bxdy[8];
3500 REAL bxay[8], cxby[8], dxcy[8], axdy[8], cxay[8], dxby[8];
3501 REAL ab[16], bc[16], cd[16], da[16], ac[16], bd[16];
3502 int ablen, bclen, cdlen, dalen, aclen, bdlen;
3503 REAL temp32a[32], temp32b[32], temp64a[64], temp64b[64], temp64c[64];
3504 int temp32alen, temp32blen, temp64alen, temp64blen, temp64clen;
3505 REAL temp128[128], temp192[192];
3506 int temp128len, temp192len;
3507 REAL detx[384], detxx[768], detxt[384], detxxt[768], detxtxt[768];
3508 int xlen, xxlen, xtlen, xxtlen, xtxtlen;
3511 REAL dety[384], detyy[768], detyt[384], detyyt[768], detytyt[768];
3512 int ylen, yylen, ytlen, yytlen, ytytlen;
3513 REAL y1[1536], y2[2304];
3515 REAL detz[384], detzz[768], detzt[384], detzzt[768], detztzt[768];
3516 int zlen, zzlen, ztlen, zztlen, ztztlen;
3521 REAL adet[6912], bdet[6912], cdet[6912], ddet[6912];
3522 int alen, blen, clen, dlen;
3523 REAL abdet[13824], cddet[13824], deter[27648];
3528 REAL avirt, bround, around;
3531 REAL a0hi, a0lo, a1hi, a1lo, bhi, blo;
3532 REAL err1, err2, err3;
3536 Two_Diff(pa[0], pe[0], aex, aextail);
3537 Two_Diff(pa[1], pe[1], aey, aeytail);
3538 Two_Diff(pa[2], pe[2], aez, aeztail);
3539 Two_Diff(pb[0], pe[0], bex, bextail);
3540 Two_Diff(pb[1], pe[1], bey, beytail);
3541 Two_Diff(pb[2], pe[2], bez, beztail);
3542 Two_Diff(pc[0], pe[0], cex, cextail);
3543 Two_Diff(pc[1], pe[1], cey, ceytail);
3544 Two_Diff(pc[2], pe[2], cez, ceztail);
3545 Two_Diff(pd[0], pe[0], dex, dextail);
3546 Two_Diff(pd[1], pe[1], dey, deytail);
3547 Two_Diff(pd[2], pe[2], dez, deztail);
3550 axby7, axby[6], axby[5], axby[4],
3551 axby[3], axby[2], axby[1], axby[0]);
3554 negatetail = -aeytail;
3556 bxay7, bxay[6], bxay[5], bxay[4],
3557 bxay[3], bxay[2], bxay[1], bxay[0]);
3561 bxcy7, bxcy[6], bxcy[5], bxcy[4],
3562 bxcy[3], bxcy[2], bxcy[1], bxcy[0]);
3565 negatetail = -beytail;
3567 cxby7, cxby[6], cxby[5], cxby[4],
3568 cxby[3], cxby[2], cxby[1], cxby[0]);
3572 cxdy7, cxdy[6], cxdy[5], cxdy[4],
3573 cxdy[3], cxdy[2], cxdy[1], cxdy[0]);
3576 negatetail = -ceytail;
3578 dxcy7, dxcy[6], dxcy[5], dxcy[4],
3579 dxcy[3], dxcy[2], dxcy[1], dxcy[0]);
3583 dxay7, dxay[6], dxay[5], dxay[4],
3584 dxay[3], dxay[2], dxay[1], dxay[0]);
3587 negatetail = -deytail;
3589 axdy7, axdy[6], axdy[5], axdy[4],
3590 axdy[3], axdy[2], axdy[1], axdy[0]);
3594 axcy7, axcy[6], axcy[5], axcy[4],
3595 axcy[3], axcy[2], axcy[1], axcy[0]);
3598 negatetail = -aeytail;
3600 cxay7, cxay[6], cxay[5], cxay[4],
3601 cxay[3], cxay[2], cxay[1], cxay[0]);
3605 bxdy7, bxdy[6], bxdy[5], bxdy[4],
3606 bxdy[3], bxdy[2], bxdy[1], bxdy[0]);
3609 negatetail = -beytail;
3611 dxby7, dxby[6], dxby[5], dxby[4],
3612 dxby[3], dxby[2], dxby[1], dxby[0]);
3619 temp32blen, temp32b, temp64a);
3623 temp32blen, temp32b, temp64b);
3627 temp32blen, temp32b, temp64c);
3629 temp64blen, temp64b, temp128);
3631 temp128len, temp128, temp192);
3636 for (i = 0; i < xxtlen; i++) {
3646 for (i = 0; i < yytlen; i++) {
3656 for (i = 0; i < zztlen; i++) {
3668 temp32blen, temp32b, temp64a);
3672 temp32blen, temp32b, temp64b);
3676 temp32blen, temp32b, temp64c);
3678 temp64blen, temp64b, temp128);
3680 temp128len, temp128, temp192);
3685 for (i = 0; i < xxtlen; i++) {
3695 for (i = 0; i < yytlen; i++) {
3705 for (i = 0; i < zztlen; i++) {
3717 temp32blen, temp32b, temp64a);
3721 temp32blen, temp32b, temp64b);
3725 temp32blen, temp32b, temp64c);
3727 temp64blen, temp64b, temp128);
3729 temp128len, temp128, temp192);
3734 for (i = 0; i < xxtlen; i++) {
3744 for (i = 0; i < yytlen; i++) {
3754 for (i = 0; i < zztlen; i++) {
3766 temp32blen, temp32b, temp64a);
3770 temp32blen, temp32b, temp64b);
3774 temp32blen, temp32b, temp64c);
3776 temp64blen, temp64b, temp128);
3778 temp128len, temp128, temp192);
3783 for (i = 0; i < xxtlen; i++) {
3793 for (i = 0; i < yytlen; i++) {
3803 for (i = 0; i < zztlen; i++) {
3816 return deter[deterlen - 1];
3822 INEXACT REAL aex, bex, cex, dex, aey, bey, cey, dey, aez, bez, cez, dez;
3828 REAL aexbey0, bexaey0, bexcey0, cexbey0;
3829 REAL cexdey0, dexcey0, dexaey0, aexdey0;
3830 REAL aexcey0, cexaey0, bexdey0, dexbey0;
3831 REAL ab[4], bc[4], cd[4], da[4], ac[4], bd[4];
3833 REAL abeps, bceps, cdeps, daeps, aceps, bdeps;
3834 REAL temp8a[8], temp8b[8], temp8c[8], temp16[16], temp24[24], temp48[48];
3835 int temp8alen, temp8blen, temp8clen, temp16len, temp24len, temp48len;
3836 REAL xdet[96], ydet[96], zdet[96], xydet[192];
3837 int xlen, ylen, zlen, xylen;
3838 REAL adet[288], bdet[288], cdet[288], ddet[288];
3839 int alen, blen, clen, dlen;
3840 REAL abdet[576], cddet[576];
3845 REAL aextail, bextail, cextail, dextail;
3846 REAL aeytail, beytail, ceytail, deytail;
3847 REAL aeztail, beztail, ceztail, deztail;
3850 REAL avirt, bround, around;
3853 REAL ahi, alo, bhi, blo;
3854 REAL err1, err2, err3;
3859 aex = (
REAL) (pa[0] - pe[0]);
3860 bex = (
REAL) (pb[0] - pe[0]);
3861 cex = (
REAL) (pc[0] - pe[0]);
3862 dex = (
REAL) (pd[0] - pe[0]);
3863 aey = (
REAL) (pa[1] - pe[1]);
3864 bey = (
REAL) (pb[1] - pe[1]);
3865 cey = (
REAL) (pc[1] - pe[1]);
3866 dey = (
REAL) (pd[1] - pe[1]);
3867 aez = (
REAL) (pa[2] - pe[2]);
3868 bez = (
REAL) (pb[2] - pe[2]);
3869 cez = (
REAL) (pc[2] - pe[2]);
3870 dez = (
REAL) (pd[2] - pe[2]);
3874 Two_Two_Diff(aexbey1, aexbey0, bexaey1, bexaey0, ab3, ab[2], ab[1], ab[0]);
3879 Two_Two_Diff(bexcey1, bexcey0, cexbey1, cexbey0, bc3, bc[2], bc[1], bc[0]);
3884 Two_Two_Diff(cexdey1, cexdey0, dexcey1, dexcey0, cd3, cd[2], cd[1], cd[0]);
3889 Two_Two_Diff(dexaey1, dexaey0, aexdey1, aexdey0, da3, da[2], da[1], da[0]);
3894 Two_Two_Diff(aexcey1, aexcey0, cexaey1, cexaey0, ac3, ac[2], ac[1], ac[0]);
3899 Two_Two_Diff(bexdey1, bexdey0, dexbey1, dexbey0, bd3, bd[2], bd[1], bd[0]);
3906 temp8blen, temp8b, temp16);
3908 temp16len, temp16, temp24);
3922 temp8blen, temp8b, temp16);
3924 temp16len, temp16, temp24);
3938 temp8blen, temp8b, temp16);
3940 temp16len, temp16, temp24);
3954 temp8blen, temp8b, temp16);
3956 temp16len, temp16, temp24);
3972 if ((det >= errbound) || (-det >= errbound)) {
3988 if ((aextail == 0.0) && (aeytail == 0.0) && (aeztail == 0.0)
3989 && (bextail == 0.0) && (beytail == 0.0) && (beztail == 0.0)
3990 && (cextail == 0.0) && (ceytail == 0.0) && (ceztail == 0.0)
3991 && (dextail == 0.0) && (deytail == 0.0) && (deztail == 0.0)) {
3996 abeps = (aex * beytail + bey * aextail)
3997 - (aey * bextail + bex * aeytail);
3998 bceps = (bex * ceytail + cey * bextail)
3999 - (bey * cextail + cex * beytail);
4000 cdeps = (cex * deytail + dey * cextail)
4001 - (cey * dextail + dex * ceytail);
4002 daeps = (dex * aeytail + aey * dextail)
4003 - (dey * aextail + aex * deytail);
4004 aceps = (aex * ceytail + cey * aextail)
4005 - (aey * cextail + cex * aeytail);
4006 bdeps = (bex * deytail + dey * bextail)
4007 - (bey * dextail + dex * beytail);
4008 det += (((bex * bex + bey * bey + bez * bez)
4009 * ((cez * daeps + dez * aceps + aez * cdeps)
4010 + (ceztail * da3 + deztail * ac3 + aeztail * cd3))
4011 + (dex * dex + dey * dey + dez * dez)
4012 * ((aez * bceps - bez * aceps + cez * abeps)
4013 + (aeztail * bc3 - beztail * ac3 + ceztail * ab3)))
4014 - ((aex * aex + aey * aey + aez * aez)
4015 * ((bez * cdeps - cez * bdeps + dez * bceps)
4016 + (beztail * cd3 - ceztail * bd3 + deztail * bc3))
4017 + (cex * cex + cey * cey + cez * cez)
4018 * ((dez * abeps + aez * bdeps + bez * daeps)
4019 + (deztail * ab3 + aeztail * bd3 + beztail * da3))))
4020 + 2.0 * (((bex * bextail + bey * beytail + bez * beztail)
4021 * (cez * da3 + dez * ac3 + aez * cd3)
4022 + (dex * dextail + dey * deytail + dez * deztail)
4023 * (aez * bc3 - bez * ac3 + cez * ab3))
4024 - ((aex * aextail + aey * aeytail + aez * aeztail)
4025 * (bez * cd3 - cez * bd3 + dez * bc3)
4026 + (cex * cextail + cey * ceytail + cez * ceztail)
4027 * (dez * ab3 + aez * bd3 + bez * da3)));
4028 if ((det >= errbound) || (-det >= errbound)) {
4035 #ifdef USE_CGAL_PREDICATES
4040 - cgal_pred_obj.side_of_oriented_sphere_3_object()
4041 (Point(pa[0], pa[1], pa[2]),
4042 Point(pb[0], pb[1], pb[2]),
4043 Point(pc[0], pc[1], pc[2]),
4044 Point(pd[0], pd[1], pd[2]),
4045 Point(pe[0], pe[1], pe[2]));
4052 REAL aex, bex, cex, dex;
4053 REAL aey, bey, cey, dey;
4054 REAL aez, bez, cez, dez;
4055 REAL aexbey, bexaey, bexcey, cexbey, cexdey, dexcey, dexaey, aexdey;
4056 REAL aexcey, cexaey, bexdey, dexbey;
4057 REAL alift, blift, clift, dlift;
4058 REAL ab, bc, cd, da, ac, bd;
4059 REAL abc, bcd, cda, dab;
4063 aex = pa[0] - pe[0];
4064 bex = pb[0] - pe[0];
4065 cex = pc[0] - pe[0];
4066 dex = pd[0] - pe[0];
4067 aey = pa[1] - pe[1];
4068 bey = pb[1] - pe[1];
4069 cey = pc[1] - pe[1];
4070 dey = pd[1] - pe[1];
4071 aez = pa[2] - pe[2];
4072 bez = pb[2] - pe[2];
4073 cez = pc[2] - pe[2];
4074 dez = pd[2] - pe[2];
4078 ab = aexbey - bexaey;
4081 bc = bexcey - cexbey;
4084 cd = cexdey - dexcey;
4087 da = dexaey - aexdey;
4091 ac = aexcey - cexaey;
4094 bd = bexdey - dexbey;
4096 abc = aez * bc - bez * ac + cez * ab;
4097 bcd = bez * cd - cez * bd + dez * bc;
4098 cda = cez * da + dez * ac + aez * cd;
4099 dab = dez * ab + aez * bd + bez * da;
4101 alift = aex * aex + aey * aey + aez * aez;
4102 blift = bex * bex + bey * bey + bez * bez;
4103 clift = cex * cex + cey * cey + cez * cez;
4104 dlift = dex * dex + dey * dey + dez * dez;
4106 det = (dlift * abc - clift * dab) + (blift * cda - alift * bcd);
4119 REAL aezplus, bezplus, cezplus, dezplus;
4120 REAL aexbeyplus, bexaeyplus, bexceyplus, cexbeyplus;
4121 REAL cexdeyplus, dexceyplus, dexaeyplus, aexdeyplus;
4122 REAL aexceyplus, cexaeyplus, bexdeyplus, dexbeyplus;
4123 REAL permanent, errbound;
4141 permanent = ((cexdeyplus + dexceyplus) * bezplus
4142 + (dexbeyplus + bexdeyplus) * cezplus
4143 + (bexceyplus + cexbeyplus) * dezplus)
4145 + ((dexaeyplus + aexdeyplus) * cezplus
4146 + (aexceyplus + cexaeyplus) * dezplus
4147 + (cexdeyplus + dexceyplus) * aezplus)
4149 + ((aexbeyplus + bexaeyplus) * dezplus
4150 + (bexdeyplus + dexbeyplus) * aezplus
4151 + (dexaeyplus + aexdeyplus) * bezplus)
4153 + ((bexceyplus + cexbeyplus) * aezplus
4154 + (cexaeyplus + aexceyplus) * bezplus
4155 + (aexbeyplus + bexaeyplus) * cezplus)
4158 if ((det > errbound) || (-det > errbound)) {
4165 #endif // #ifdef USE_CGAL_PREDICATES
4197 REAL axby0, bxcy0, cxdy0, dxey0, exay0;
4198 REAL bxay0, cxby0, dxcy0, exdy0, axey0;
4199 REAL axcy0, bxdy0, cxey0, dxay0, exby0;
4200 REAL cxay0, dxby0, excy0, axdy0, bxey0;
4201 REAL ab[4], bc[4], cd[4], de[4], ea[4];
4202 REAL ac[4], bd[4], ce[4], da[4], eb[4];
4203 REAL temp8a[8], temp8b[8], temp16[16];
4204 int temp8alen, temp8blen, temp16len;
4205 REAL abc[24], bcd[24], cde[24], dea[24], eab[24];
4206 REAL abd[24], bce[24], cda[24], deb[24], eac[24];
4207 int abclen, bcdlen, cdelen, dealen, eablen;
4208 int abdlen, bcelen, cdalen, deblen, eaclen;
4209 REAL temp48a[48], temp48b[48];
4210 int temp48alen, temp48blen;
4211 REAL abcd[96], bcde[96], cdea[96], deab[96], eabc[96];
4212 int abcdlen, bcdelen, cdealen, deablen, eabclen;
4213 REAL adet[192], bdet[192], cdet[192], ddet[192], edet[192];
4214 int alen, blen, clen, dlen, elen;
4215 REAL abdet[384], cddet[384], cdedet[576];
4222 REAL avirt, bround, around;
4225 REAL ahi, alo, bhi, blo;
4226 REAL err1, err2, err3;
4233 Two_Two_Diff(axby1, axby0, bxay1, bxay0, ab[3], ab[2], ab[1], ab[0]);
4237 Two_Two_Diff(bxcy1, bxcy0, cxby1, cxby0, bc[3], bc[2], bc[1], bc[0]);
4241 Two_Two_Diff(cxdy1, cxdy0, dxcy1, dxcy0, cd[3], cd[2], cd[1], cd[0]);
4245 Two_Two_Diff(dxey1, dxey0, exdy1, exdy0, de[3], de[2], de[1], de[0]);
4249 Two_Two_Diff(exay1, exay0, axey1, axey0, ea[3], ea[2], ea[1], ea[0]);
4253 Two_Two_Diff(axcy1, axcy0, cxay1, cxay0, ac[3], ac[2], ac[1], ac[0]);
4257 Two_Two_Diff(bxdy1, bxdy0, dxby1, dxby0, bd[3], bd[2], bd[1], bd[0]);
4261 Two_Two_Diff(cxey1, cxey0, excy1, excy0, ce[3], ce[2], ce[1], ce[0]);
4265 Two_Two_Diff(dxay1, dxay0, axdy1, axdy0, da[3], da[2], da[1], da[0]);
4269 Two_Two_Diff(exby1, exby0, bxey1, bxey0, eb[3], eb[2], eb[1], eb[0]);
4353 for (i = 0; i < temp48blen; i++) {
4354 temp48b[i] = -temp48b[i];
4357 temp48blen, temp48b, bcde);
4362 for (i = 0; i < temp48blen; i++) {
4363 temp48b[i] = -temp48b[i];
4366 temp48blen, temp48b, cdea);
4371 for (i = 0; i < temp48blen; i++) {
4372 temp48b[i] = -temp48b[i];
4375 temp48blen, temp48b, deab);
4380 for (i = 0; i < temp48blen; i++) {
4381 temp48b[i] = -temp48b[i];
4384 temp48blen, temp48b, eabc);
4389 for (i = 0; i < temp48blen; i++) {
4390 temp48b[i] = -temp48b[i];
4393 temp48blen, temp48b, abcd);
4401 return deter[deterlen - 1];
4408 INEXACT REAL aex, bex, cex, dex, aey, bey, cey, dey, aez, bez, cez, dez;
4409 INEXACT REAL aeheight, beheight, ceheight, deheight;
4415 REAL aexbey0, bexaey0, bexcey0, cexbey0;
4416 REAL cexdey0, dexcey0, dexaey0, aexdey0;
4417 REAL aexcey0, cexaey0, bexdey0, dexbey0;
4418 REAL ab[4], bc[4], cd[4], da[4], ac[4], bd[4];
4420 REAL abeps, bceps, cdeps, daeps, aceps, bdeps;
4421 REAL temp8a[8], temp8b[8], temp8c[8], temp16[16], temp24[24];
4422 int temp8alen, temp8blen, temp8clen, temp16len, temp24len;
4423 REAL adet[48], bdet[48], cdet[48], ddet[48];
4424 int alen, blen, clen, dlen;
4425 REAL abdet[96], cddet[96];
4430 REAL aextail, bextail, cextail, dextail;
4431 REAL aeytail, beytail, ceytail, deytail;
4432 REAL aeztail, beztail, ceztail, deztail;
4433 REAL aeheighttail, beheighttail, ceheighttail, deheighttail;
4436 REAL avirt, bround, around;
4439 REAL ahi, alo, bhi, blo;
4440 REAL err1, err2, err3;
4445 aex = (
REAL) (pa[0] - pe[0]);
4446 bex = (
REAL) (pb[0] - pe[0]);
4447 cex = (
REAL) (pc[0] - pe[0]);
4448 dex = (
REAL) (pd[0] - pe[0]);
4449 aey = (
REAL) (pa[1] - pe[1]);
4450 bey = (
REAL) (pb[1] - pe[1]);
4451 cey = (
REAL) (pc[1] - pe[1]);
4452 dey = (
REAL) (pd[1] - pe[1]);
4453 aez = (
REAL) (pa[2] - pe[2]);
4454 bez = (
REAL) (pb[2] - pe[2]);
4455 cez = (
REAL) (pc[2] - pe[2]);
4456 dez = (
REAL) (pd[2] - pe[2]);
4457 aeheight = (
REAL) (aheight - eheight);
4458 beheight = (
REAL) (bheight - eheight);
4459 ceheight = (
REAL) (cheight - eheight);
4460 deheight = (
REAL) (dheight - eheight);
4464 Two_Two_Diff(aexbey1, aexbey0, bexaey1, bexaey0, ab3, ab[2], ab[1], ab[0]);
4469 Two_Two_Diff(bexcey1, bexcey0, cexbey1, cexbey0, bc3, bc[2], bc[1], bc[0]);
4474 Two_Two_Diff(cexdey1, cexdey0, dexcey1, dexcey0, cd3, cd[2], cd[1], cd[0]);
4479 Two_Two_Diff(dexaey1, dexaey0, aexdey1, aexdey0, da3, da[2], da[1], da[0]);
4484 Two_Two_Diff(aexcey1, aexcey0, cexaey1, cexaey0, ac3, ac[2], ac[1], ac[0]);
4489 Two_Two_Diff(bexdey1, bexdey0, dexbey1, dexbey0, bd3, bd[2], bd[1], bd[0]);
4496 temp8blen, temp8b, temp16);
4498 temp16len, temp16, temp24);
4505 temp8blen, temp8b, temp16);
4507 temp16len, temp16, temp24);
4514 temp8blen, temp8b, temp16);
4516 temp16len, temp16, temp24);
4523 temp8blen, temp8b, temp16);
4525 temp16len, temp16, temp24);
4534 if ((det >= errbound) || (-det >= errbound)) {
4554 if ((aextail == 0.0) && (aeytail == 0.0) && (aeztail == 0.0)
4555 && (bextail == 0.0) && (beytail == 0.0) && (beztail == 0.0)
4556 && (cextail == 0.0) && (ceytail == 0.0) && (ceztail == 0.0)
4557 && (dextail == 0.0) && (deytail == 0.0) && (deztail == 0.0)
4558 && (aeheighttail == 0.0) && (beheighttail == 0.0)
4559 && (ceheighttail == 0.0) && (deheighttail == 0.0)) {
4564 abeps = (aex * beytail + bey * aextail)
4565 - (aey * bextail + bex * aeytail);
4566 bceps = (bex * ceytail + cey * bextail)
4567 - (bey * cextail + cex * beytail);
4568 cdeps = (cex * deytail + dey * cextail)
4569 - (cey * dextail + dex * ceytail);
4570 daeps = (dex * aeytail + aey * dextail)
4571 - (dey * aextail + aex * deytail);
4572 aceps = (aex * ceytail + cey * aextail)
4573 - (aey * cextail + cex * aeytail);
4574 bdeps = (bex * deytail + dey * bextail)
4575 - (bey * dextail + dex * beytail);
4577 * ((cez * daeps + dez * aceps + aez * cdeps)
4578 + (ceztail * da3 + deztail * ac3 + aeztail * cd3))
4580 * ((aez * bceps - bez * aceps + cez * abeps)
4581 + (aeztail * bc3 - beztail * ac3 + ceztail * ab3)))
4583 * ((bez * cdeps - cez * bdeps + dez * bceps)
4584 + (beztail * cd3 - ceztail * bd3 + deztail * bc3))
4586 * ((dez * abeps + aez * bdeps + bez * daeps)
4587 + (deztail * ab3 + aeztail * bd3 + beztail * da3))))
4588 + ((beheighttail * (cez * da3 + dez * ac3 + aez * cd3)
4589 + deheighttail * (aez * bc3 - bez * ac3 + cez * ab3))
4590 - (aeheighttail * (bez * cd3 - cez * bd3 + dez * bc3)
4591 + ceheighttail * (dez * ab3 + aez * bd3 + bez * da3)));
4592 if ((det >= errbound) || (-det >= errbound)) {
4597 aheight, bheight, cheight, dheight, eheight);
4604 REAL aex, bex, cex, dex;
4605 REAL aey, bey, cey, dey;
4606 REAL aez, bez, cez, dez;
4607 REAL aexbey, bexaey, bexcey, cexbey, cexdey, dexcey, dexaey, aexdey;
4608 REAL aexcey, cexaey, bexdey, dexbey;
4609 REAL aeheight, beheight, ceheight, deheight;
4610 REAL ab, bc, cd, da, ac, bd;
4611 REAL abc, bcd, cda, dab;
4612 REAL aezplus, bezplus, cezplus, dezplus;
4613 REAL aexbeyplus, bexaeyplus, bexceyplus, cexbeyplus;
4614 REAL cexdeyplus, dexceyplus, dexaeyplus, aexdeyplus;
4615 REAL aexceyplus, cexaeyplus, bexdeyplus, dexbeyplus;
4617 REAL permanent, errbound;
4620 aex = pa[0] - pe[0];
4621 bex = pb[0] - pe[0];
4622 cex = pc[0] - pe[0];
4623 dex = pd[0] - pe[0];
4624 aey = pa[1] - pe[1];
4625 bey = pb[1] - pe[1];
4626 cey = pc[1] - pe[1];
4627 dey = pd[1] - pe[1];
4628 aez = pa[2] - pe[2];
4629 bez = pb[2] - pe[2];
4630 cez = pc[2] - pe[2];
4631 dez = pd[2] - pe[2];
4632 aeheight = aheight - eheight;
4633 beheight = bheight - eheight;
4634 ceheight = cheight - eheight;
4635 deheight = dheight - eheight;
4639 ab = aexbey - bexaey;
4642 bc = bexcey - cexbey;
4645 cd = cexdey - dexcey;
4648 da = dexaey - aexdey;
4652 ac = aexcey - cexaey;
4655 bd = bexdey - dexbey;
4657 abc = aez * bc - bez * ac + cez * ab;
4658 bcd = bez * cd - cez * bd + dez * bc;
4659 cda = cez * da + dez * ac + aez * cd;
4660 dab = dez * ab + aez * bd + bez * da;
4662 det = (deheight * abc - ceheight * dab) + (beheight * cda - aeheight * bcd);
4680 permanent = ((cexdeyplus + dexceyplus) * bezplus
4681 + (dexbeyplus + bexdeyplus) * cezplus
4682 + (bexceyplus + cexbeyplus) * dezplus)
4684 + ((dexaeyplus + aexdeyplus) * cezplus
4685 + (aexceyplus + cexaeyplus) * dezplus
4686 + (cexdeyplus + dexceyplus) * aezplus)
4688 + ((aexbeyplus + bexaeyplus) * dezplus
4689 + (bexdeyplus + dexbeyplus) * aezplus
4690 + (dexaeyplus + aexdeyplus) * bezplus)
4692 + ((bexceyplus + cexbeyplus) * aezplus
4693 + (cexaeyplus + aexceyplus) * bezplus
4694 + (aexbeyplus + bexaeyplus) * cezplus)
4697 if ((det > errbound) || (-det > errbound)) {
4702 aheight, bheight, cheight, dheight, eheight, permanent);