1#![cfg_attr(feature = "no_std", no_std)]
2#![doc(html_logo_url = "https://raw.githubusercontent.com/georust/meta/master/logo/logo.png")]
3#![allow(non_snake_case)]
12
13#[cfg(test)]
34mod tests;
35
36#[derive(Copy, Clone, Debug, PartialEq)]
38pub struct Coord<T: Into<f64>> {
39 pub x: T,
40 pub y: T,
41}
42
43#[derive(Copy, Clone, Debug, PartialEq)]
45pub struct Coord3D<T: Into<f64>> {
46 pub x: T,
47 pub y: T,
48 pub z: T,
49}
50
51const SPLITTER: f64 = 134_217_729f64;
54const EPSILON: f64 = 0.000_000_000_000_000_111_022_302_462_515_65;
55const RESULTERRBOUND: f64 = (3.0 + 8.0 * EPSILON) * EPSILON;
56const CCWERRBOUND_B: f64 = (2.0 + 12.0 * EPSILON) * EPSILON;
57const CCWERRBOUND_C: f64 = (9.0 + 64.0 * EPSILON) * EPSILON * EPSILON;
58const O3DERRBOUND_A: f64 = (7.0 + 56.0 * EPSILON) * EPSILON;
59const O3DERRBOUND_B: f64 = (3.0 + 28.0 * EPSILON) * EPSILON;
60const O3DERRBOUND_C: f64 = (26.0 + 288.0 * EPSILON) * EPSILON * EPSILON;
61const ICCERRBOUND_A: f64 = (10.0 + 96.0 * EPSILON) * EPSILON;
62const ICCERRBOUND_B: f64 = (4.0 + 48.0 * EPSILON) * EPSILON;
63const ICCERRBOUND_C: f64 = (44.0 + 576.0 * EPSILON) * EPSILON * EPSILON;
64const ISPERRBOUND_A: f64 = (16.0 + 224.0 * EPSILON) * EPSILON;
65const ISPERRBOUND_B: f64 = (5.0 + 72.0 * EPSILON) * EPSILON;
66const ISPERRBOUND_C: f64 = (71.0 + 1408.0 * EPSILON) * EPSILON * EPSILON;
67
68#[inline(always)]
69fn abs(x: f64) -> f64 {
70 x.abs()
71}
72
73pub fn orient2d<T: Into<f64>>(pa: Coord<T>, pb: Coord<T>, pc: Coord<T>) -> f64 {
78 let pa = Coord {
79 x: pa.x.into(),
80 y: pa.y.into(),
81 };
82 let pb = Coord {
83 x: pb.x.into(),
84 y: pb.y.into(),
85 };
86 let pc = Coord {
87 x: pc.x.into(),
88 y: pc.y.into(),
89 };
90
91 let detleft = (pa.x - pc.x) * (pb.y - pc.y);
92 let detright = (pa.y - pc.y) * (pb.x - pc.x);
93 let det = detleft - detright;
94
95 let detsum = abs(detleft + detright);
101 const THETA: f64 = 3.3306690621773722e-16;
102 let errbound = THETA * detsum;
103 if det >= errbound || -det >= errbound {
104 det
105 } else {
106 orient2dadapt(pa, pb, pc, detsum)
107 }
108}
109
110fn orient2dadapt(pa: Coord<f64>, pb: Coord<f64>, pc: Coord<f64>, detsum: f64) -> f64 {
111 let acx = pa.x - pc.x;
112 let bcx = pb.x - pc.x;
113 let acy = pa.y - pc.y;
114 let bcy = pb.y - pc.y;
115
116 let (detleft, detlefttail) = two_product(acx, bcy);
117 let (detright, detrighttail) = two_product(acy, bcx);
118
119 let (B3, B2, B1, B0) = two_two_diff(detleft, detlefttail, detright, detrighttail);
120 let B = [B0, B1, B2, B3];
121
122 let mut det = estimate(&B);
123 let errbound = CCWERRBOUND_B * detsum;
124 if det >= errbound || (-det >= errbound) {
125 return det;
126 }
127
128 let acxtail = two_diff_tail(pa.x, pc.x, acx);
129 let bcxtail = two_diff_tail(pb.x, pc.x, bcx);
130 let acytail = two_diff_tail(pa.y, pc.y, acy);
131 let bcytail = two_diff_tail(pb.y, pc.y, bcy);
132
133 if acxtail == 0.0 && acytail == 0.0 && bcxtail == 0.0 && bcytail == 0.0 {
134 return det;
135 }
136
137 let errbound = CCWERRBOUND_C * detsum + RESULTERRBOUND * abs(det);
138 det += (acx * bcytail + bcy * acxtail) - (acy * bcxtail + bcx * acytail);
139
140 if det >= errbound || -det >= errbound {
141 return det;
142 }
143
144 let (s1, s0) = two_product(acxtail, bcy);
145 let (t1, t0) = two_product(acytail, bcx);
146 let (u3, u2, u1, u0) = two_two_diff(s1, s0, t1, t0);
147 let U = [u0, u1, u2, u3];
148
149 let mut C1 = [0.0f64; 8];
150 let c1length = fast_expansion_sum_zeroelim(&B, &U, &mut C1);
151
152 let (s1, s0) = two_product(acx, bcytail);
153 let (t1, t0) = two_product(acy, bcxtail);
154 let (u3, u2, u1, u0) = two_two_diff(s1, s0, t1, t0);
155 let U = [u0, u1, u2, u3];
156
157 let mut C2 = [0.0f64; 12];
158 let c2length = fast_expansion_sum_zeroelim(&C1[..c1length], &U, &mut C2);
159
160 let (s1, s0) = two_product(acxtail, bcytail);
161 let (t1, t0) = two_product(acytail, bcxtail);
162 let (u3, u2, u1, u0) = two_two_diff(s1, s0, t1, t0);
163 let U = [u0, u1, u2, u3];
164 let mut D = [0.0f64; 16];
165 let dlength = fast_expansion_sum_zeroelim(&C2[..c2length], &U, &mut D);
166 D[dlength - 1]
167}
168
169pub fn orient3d<T: Into<f64>>(
174 pa: Coord3D<T>,
175 pb: Coord3D<T>,
176 pc: Coord3D<T>,
177 pd: Coord3D<T>,
178) -> f64 {
179 let pa = Coord3D {
180 x: pa.x.into(),
181 y: pa.y.into(),
182 z: pa.z.into(),
183 };
184 let pb = Coord3D {
185 x: pb.x.into(),
186 y: pb.y.into(),
187 z: pb.z.into(),
188 };
189 let pc = Coord3D {
190 x: pc.x.into(),
191 y: pc.y.into(),
192 z: pc.z.into(),
193 };
194 let pd = Coord3D {
195 x: pd.x.into(),
196 y: pd.y.into(),
197 z: pd.z.into(),
198 };
199
200 let adx = pa.x - pd.x;
201 let bdx = pb.x - pd.x;
202 let cdx = pc.x - pd.x;
203 let ady = pa.y - pd.y;
204 let bdy = pb.y - pd.y;
205 let cdy = pc.y - pd.y;
206 let adz = pa.z - pd.z;
207 let bdz = pb.z - pd.z;
208 let cdz = pc.z - pd.z;
209
210 let bdxcdy = bdx * cdy;
211 let cdxbdy = cdx * bdy;
212
213 let cdxady = cdx * ady;
214 let adxcdy = adx * cdy;
215
216 let adxbdy = adx * bdy;
217 let bdxady = bdx * ady;
218
219 let det = adz * (bdxcdy - cdxbdy) + bdz * (cdxady - adxcdy) + cdz * (adxbdy - bdxady);
220
221 let permanent = (abs(bdxcdy) + abs(cdxbdy)) * abs(adz)
222 + (abs(cdxady) + abs(adxcdy)) * abs(bdz)
223 + (abs(adxbdy) + abs(bdxady)) * abs(cdz);
224
225 let errbound = O3DERRBOUND_A * permanent;
226 if det > errbound || -det > errbound {
227 return det;
228 }
229
230 orient3dadapt(pa, pb, pc, pd, permanent)
231}
232
233fn orient3dadapt(
234 pa: Coord3D<f64>,
235 pb: Coord3D<f64>,
236 pc: Coord3D<f64>,
237 pd: Coord3D<f64>,
238 permanent: f64,
239) -> f64 {
240 let adx = pa.x - pd.x;
241 let bdx = pb.x - pd.x;
242 let cdx = pc.x - pd.x;
243 let ady = pa.y - pd.y;
244 let bdy = pb.y - pd.y;
245 let cdy = pc.y - pd.y;
246 let adz = pa.z - pd.z;
247 let bdz = pb.z - pd.z;
248 let cdz = pc.z - pd.z;
249
250 let (bdxcdy1, bdxcdy0) = two_product(bdx, cdy);
251 let (cdxbdy1, cdxbdy0) = two_product(cdx, bdy);
252 let (bc3, bc2, bc1, bc0) = two_two_diff(bdxcdy1, bdxcdy0, cdxbdy1, cdxbdy0);
253 let bc = [bc0, bc1, bc2, bc3];
254 let mut adet = [0f64; 8];
255 let alen = scale_expansion_zeroelim(&bc[..4], adz, &mut adet);
256
257 let (cdxady1, cdxady0) = two_product(cdx, ady);
258 let (adxcdy1, adxcdy0) = two_product(adx, cdy);
259 let (ca3, ca2, ca1, ca0) = two_two_diff(cdxady1, cdxady0, adxcdy1, adxcdy0);
260 let ca = [ca0, ca1, ca2, ca3];
261 let mut bdet = [0f64; 8];
262 let blen = scale_expansion_zeroelim(&ca[..4], bdz, &mut bdet);
263
264 let (adxbdy1, adxbdy0) = two_product(adx, bdy);
265 let (bdxady1, bdxady0) = two_product(bdx, ady);
266 let (ab3, ab2, ab1, ab0) = two_two_diff(adxbdy1, adxbdy0, bdxady1, bdxady0);
267 let ab = [ab0, ab1, ab2, ab3];
268 let mut cdet = [0f64; 8];
269 let clen = scale_expansion_zeroelim(&ab[..4], cdz, &mut cdet);
270
271 let mut abdet = [0f64; 16];
272 let ablen = fast_expansion_sum_zeroelim(&adet[..alen], &bdet[..blen], &mut abdet);
273 let mut fin1 = [0f64; 192];
274 let finlength = fast_expansion_sum_zeroelim(&abdet[..ablen], &cdet[..clen], &mut fin1);
275
276 let mut det = estimate(&fin1[..finlength]);
277 let mut errbound = O3DERRBOUND_B * permanent;
278 if det >= errbound || -det >= errbound {
279 return det;
280 }
281
282 let adxtail = two_diff_tail(pa.x, pd.x, adx);
283 let bdxtail = two_diff_tail(pb.x, pd.x, bdx);
284 let cdxtail = two_diff_tail(pc.x, pd.x, cdx);
285 let adytail = two_diff_tail(pa.y, pd.y, ady);
286 let bdytail = two_diff_tail(pb.y, pd.y, bdy);
287 let cdytail = two_diff_tail(pc.y, pd.y, cdy);
288 let adztail = two_diff_tail(pa.z, pd.z, adz);
289 let bdztail = two_diff_tail(pb.z, pd.z, bdz);
290 let cdztail = two_diff_tail(pc.z, pd.z, cdz);
291
292 if (adxtail == 0.0)
293 && (bdxtail == 0.0)
294 && (cdxtail == 0.0)
295 && (adytail == 0.0)
296 && (bdytail == 0.0)
297 && (cdytail == 0.0)
298 && (adztail == 0.0)
299 && (bdztail == 0.0)
300 && (cdztail == 0.0)
301 {
302 return det;
303 }
304
305 errbound = O3DERRBOUND_C * permanent + RESULTERRBOUND * abs(det);
306 det += (adz * ((bdx * cdytail + cdy * bdxtail) - (bdy * cdxtail + cdx * bdytail))
307 + adztail * (bdx * cdy - bdy * cdx))
308 + (bdz * ((cdx * adytail + ady * cdxtail) - (cdy * adxtail + adx * cdytail))
309 + bdztail * (cdx * ady - cdy * adx))
310 + (cdz * ((adx * bdytail + bdy * adxtail) - (ady * bdxtail + bdx * adytail))
311 + cdztail * (adx * bdy - ady * bdx));
312 if det >= errbound || -det >= errbound {
313 return det;
314 }
315
316 let mut finnow = fin1;
317 let mut finother = [0f64; 192];
318
319 let mut at_b = [0f64; 4];
320 let mut at_c = [0f64; 4];
321 let mut bt_c = [0f64; 4];
322 let mut bt_a = [0f64; 4];
323 let mut ct_a = [0f64; 4];
324 let mut ct_b = [0f64; 4];
325 let at_blen: usize;
326 let at_clen: usize;
327 let bt_clen: usize;
328 let bt_alen: usize;
329 let ct_alen: usize;
330 let ct_blen: usize;
331 if adxtail == 0.0 {
332 if adytail == 0.0 {
333 at_b[0] = 0.0;
334 at_blen = 1;
335 at_c[0] = 0.0;
336 at_clen = 1;
337 } else {
338 let negate = -adytail;
339 (at_b[1], at_b[0]) = two_product(negate, bdx);
340 at_blen = 2;
341 (at_c[1], at_c[0]) = two_product(adytail, cdx);
342 at_clen = 2;
343 }
344 } else if adytail == 0.0 {
345 (at_b[1], at_b[0]) = two_product(adxtail, bdy);
346 at_blen = 2;
347 let negate = -adxtail;
348 (at_c[1], at_c[0]) = two_product(negate, cdy);
349 at_clen = 2;
350 } else {
351 let (adxt_bdy1, adxt_bdy0) = two_product(adxtail, bdy);
352 let (adyt_bdx1, adyt_bdx0) = two_product(adytail, bdx);
353 (at_b[3], at_b[2], at_b[1], at_b[0]) =
354 two_two_diff(adxt_bdy1, adxt_bdy0, adyt_bdx1, adyt_bdx0);
355 at_blen = 4;
356 let (adyt_cdx1, adyt_cdx0) = two_product(adytail, cdx);
357 let (adxt_cdy1, adxt_cdy0) = two_product(adxtail, cdy);
358 (at_c[3], at_c[2], at_c[1], at_c[0]) =
359 two_two_diff(adyt_cdx1, adyt_cdx0, adxt_cdy1, adxt_cdy0);
360 at_clen = 4;
361 }
362 if bdxtail == 0.0 {
363 if bdytail == 0.0 {
364 bt_c[0] = 0.0;
365 bt_clen = 1;
366 bt_a[0] = 0.0;
367 bt_alen = 1;
368 } else {
369 let negate = -bdytail;
370 (bt_c[1], bt_c[0]) = two_product(negate, cdx);
371 bt_clen = 2;
372 (bt_a[1], bt_a[0]) = two_product(bdytail, adx);
373 bt_alen = 2;
374 }
375 } else if bdytail == 0.0 {
376 (bt_c[1], bt_c[0]) = two_product(bdxtail, cdy);
377 bt_clen = 2;
378 let negate = -bdxtail;
379 (bt_a[1], bt_a[0]) = two_product(negate, ady);
380 bt_alen = 2;
381 } else {
382 let (bdxt_cdy1, bdxt_cdy0) = two_product(bdxtail, cdy);
383 let (bdyt_cdx1, bdyt_cdx0) = two_product(bdytail, cdx);
384 (bt_c[3], bt_c[2], bt_c[1], bt_c[0]) =
385 two_two_diff(bdxt_cdy1, bdxt_cdy0, bdyt_cdx1, bdyt_cdx0);
386 bt_clen = 4;
387 let (bdyt_adx1, bdyt_adx0) = two_product(bdytail, adx);
388 let (bdxt_ady1, bdxt_ady0) = two_product(bdxtail, ady);
389 (bt_a[3], bt_a[2], bt_a[1], bt_a[0]) =
390 two_two_diff(bdyt_adx1, bdyt_adx0, bdxt_ady1, bdxt_ady0);
391 bt_alen = 4;
392 }
393 if cdxtail == 0.0 {
394 if cdytail == 0.0 {
395 ct_a[0] = 0.0;
396 ct_alen = 1;
397 ct_b[0] = 0.0;
398 ct_blen = 1;
399 } else {
400 let negate = -cdytail;
401 (ct_a[1], ct_a[0]) = two_product(negate, adx);
402 ct_alen = 2;
403 (ct_b[1], ct_b[0]) = two_product(cdytail, bdx);
404 ct_blen = 2;
405 }
406 } else if cdytail == 0.0 {
407 (ct_a[1], ct_a[0]) = two_product(cdxtail, ady);
408 ct_alen = 2;
409 let negate = -cdxtail;
410 (ct_b[1], ct_b[0]) = two_product(negate, bdy);
411 ct_blen = 2;
412 } else {
413 let (cdxt_ady1, cdxt_ady0) = two_product(cdxtail, ady);
414 let (cdyt_adx1, cdyt_adx0) = two_product(cdytail, adx);
415 (ct_a[3], ct_a[2], ct_a[1], ct_a[0]) =
416 two_two_diff(cdxt_ady1, cdxt_ady0, cdyt_adx1, cdyt_adx0);
417 ct_alen = 4;
418 let (cdyt_bdx1, cdyt_bdx0) = two_product(cdytail, bdx);
419 let (cdxt_bdy1, cdxt_bdy0) = two_product(cdxtail, bdy);
420 (ct_b[3], ct_b[2], ct_b[1], ct_b[0]) =
421 two_two_diff(cdyt_bdx1, cdyt_bdx0, cdxt_bdy1, cdxt_bdy0);
422 ct_blen = 4;
423 }
424
425 let mut bct = [0f64; 8];
426 let mut cat = [0f64; 8];
427 let mut abt = [0f64; 8];
428 let mut u = [0f64; 4];
429 let mut v = [0f64; 12];
430 let mut w = [0f64; 16];
431 let mut vlength: usize;
432 let mut wlength: usize;
433
434 let bctlen = fast_expansion_sum_zeroelim(&bt_c[..bt_clen], &ct_b[..ct_blen], &mut bct);
435 wlength = scale_expansion_zeroelim(&bct[..bctlen], adz, &mut w);
436 let mut finlength =
437 fast_expansion_sum_zeroelim(&finnow[..finlength], &w[..wlength], &mut finother);
438 ::core::mem::swap(&mut finnow, &mut finother);
439
440 let catlen = fast_expansion_sum_zeroelim(&ct_a[..ct_alen], &at_c[..at_clen], &mut cat);
441 wlength = scale_expansion_zeroelim(&cat[..catlen], bdz, &mut w);
442 finlength = fast_expansion_sum_zeroelim(&finnow[..finlength], &w[..wlength], &mut finother);
443 ::core::mem::swap(&mut finnow, &mut finother);
444
445 let abtlen = fast_expansion_sum_zeroelim(&at_b[..at_blen], &bt_a[..bt_alen], &mut abt);
446 wlength = scale_expansion_zeroelim(&abt[..abtlen], cdz, &mut w);
447 finlength = fast_expansion_sum_zeroelim(&finnow[..finlength], &w[..wlength], &mut finother);
448 ::core::mem::swap(&mut finnow, &mut finother);
449
450 if adztail != 0.0 {
451 vlength = scale_expansion_zeroelim(&bc[..4], adztail, &mut v);
452 finlength = fast_expansion_sum_zeroelim(&finnow[..finlength], &v[..vlength], &mut finother);
453 ::core::mem::swap(&mut finnow, &mut finother);
454 }
455 if bdztail != 0.0 {
456 vlength = scale_expansion_zeroelim(&ca[..4], bdztail, &mut v);
457 finlength = fast_expansion_sum_zeroelim(&finnow[..finlength], &v[..vlength], &mut finother);
458 ::core::mem::swap(&mut finnow, &mut finother);
459 }
460 if cdztail != 0.0 {
461 vlength = scale_expansion_zeroelim(&ab[..4], cdztail, &mut v);
462 finlength = fast_expansion_sum_zeroelim(&finnow[..finlength], &v[..vlength], &mut finother);
463 ::core::mem::swap(&mut finnow, &mut finother);
464 }
465
466 if adxtail != 0.0 {
467 if bdytail != 0.0 {
468 let (adxt_bdyt1, adxt_bdyt0) = two_product(adxtail, bdytail);
469 (u[3], u[2], u[1], u[0]) = two_one_product(adxt_bdyt1, adxt_bdyt0, cdz);
470 finlength = fast_expansion_sum_zeroelim(&finnow[..finlength], &u[..4], &mut finother);
471 ::core::mem::swap(&mut finnow, &mut finother);
472 if cdztail != 0.0 {
473 (u[3], u[2], u[1], u[0]) = two_one_product(adxt_bdyt1, adxt_bdyt0, cdztail);
474 finlength =
475 fast_expansion_sum_zeroelim(&finnow[..finlength], &u[..4], &mut finother);
476 ::core::mem::swap(&mut finnow, &mut finother);
477 }
478 }
479 if cdytail != 0.0 {
480 let negate = -adxtail;
481 let (adxt_cdyt1, adxt_cdyt0) = two_product(negate, cdytail);
482 (u[3], u[2], u[1], u[0]) = two_one_product(adxt_cdyt1, adxt_cdyt0, bdz);
483 finlength = fast_expansion_sum_zeroelim(&finnow[..finlength], &u[..4], &mut finother);
484 ::core::mem::swap(&mut finnow, &mut finother);
485 if bdztail != 0.0 {
486 (u[3], u[2], u[1], u[0]) = two_one_product(adxt_cdyt1, adxt_cdyt0, bdztail);
487 finlength =
488 fast_expansion_sum_zeroelim(&finnow[..finlength], &u[..4], &mut finother);
489 ::core::mem::swap(&mut finnow, &mut finother);
490 }
491 }
492 }
493 if bdxtail != 0.0 {
494 if cdytail != 0.0 {
495 let (bdxt_cdyt1, bdxt_cdyt0) = two_product(bdxtail, cdytail);
496 (u[3], u[2], u[1], u[0]) = two_one_product(bdxt_cdyt1, bdxt_cdyt0, adz);
497 finlength = fast_expansion_sum_zeroelim(&finnow[..finlength], &u[..4], &mut finother);
498 ::core::mem::swap(&mut finnow, &mut finother);
499 if adztail != 0.0 {
500 (u[3], u[2], u[1], u[0]) = two_one_product(bdxt_cdyt1, bdxt_cdyt0, adztail);
501 finlength =
502 fast_expansion_sum_zeroelim(&finnow[..finlength], &u[..4], &mut finother);
503 ::core::mem::swap(&mut finnow, &mut finother);
504 }
505 }
506 if adytail != 0.0 {
507 let negate = -bdxtail;
508 let (bdxt_adyt1, bdxt_adyt0) = two_product(negate, adytail);
509 (u[3], u[2], u[1], u[0]) = two_one_product(bdxt_adyt1, bdxt_adyt0, cdz);
510 finlength = fast_expansion_sum_zeroelim(&finnow[..finlength], &u[..4], &mut finother);
511 ::core::mem::swap(&mut finnow, &mut finother);
512 if cdztail != 0.0 {
513 (u[3], u[2], u[1], u[0]) = two_one_product(bdxt_adyt1, bdxt_adyt0, cdztail);
514 finlength =
515 fast_expansion_sum_zeroelim(&finnow[..finlength], &u[..4], &mut finother);
516 ::core::mem::swap(&mut finnow, &mut finother);
517 }
518 }
519 }
520 if cdxtail != 0.0 {
521 if adytail != 0.0 {
522 let (cdxt_adyt1, cdxt_adyt0) = two_product(cdxtail, adytail);
523 (u[3], u[2], u[1], u[0]) = two_one_product(cdxt_adyt1, cdxt_adyt0, bdz);
524 finlength = fast_expansion_sum_zeroelim(&finnow[..finlength], &u[..4], &mut finother);
525 ::core::mem::swap(&mut finnow, &mut finother);
526 if bdztail != 0.0 {
527 (u[3], u[2], u[1], u[0]) = two_one_product(cdxt_adyt1, cdxt_adyt0, bdztail);
528 finlength =
529 fast_expansion_sum_zeroelim(&finnow[..finlength], &u[..4], &mut finother);
530 ::core::mem::swap(&mut finnow, &mut finother);
531 }
532 }
533 if bdytail != 0.0 {
534 let negate = -cdxtail;
535 let (cdxt_bdyt1, cdxt_bdyt0) = two_product(negate, bdytail);
536 (u[3], u[2], u[1], u[0]) = two_one_product(cdxt_bdyt1, cdxt_bdyt0, adz);
537 finlength = fast_expansion_sum_zeroelim(&finnow[..finlength], &u[..4], &mut finother);
538 ::core::mem::swap(&mut finnow, &mut finother);
539 if adztail != 0.0 {
540 (u[3], u[2], u[1], u[0]) = two_one_product(cdxt_bdyt1, cdxt_bdyt0, adztail);
541 finlength =
542 fast_expansion_sum_zeroelim(&finnow[..finlength], &u[..4], &mut finother);
543 ::core::mem::swap(&mut finnow, &mut finother);
544 }
545 }
546 }
547
548 if adztail != 0.0 {
549 wlength = scale_expansion_zeroelim(&bct[..bctlen], adztail, &mut w);
550 finlength = fast_expansion_sum_zeroelim(&finnow[..finlength], &w[..wlength], &mut finother);
551 ::core::mem::swap(&mut finnow, &mut finother);
552 }
553 if bdztail != 0.0 {
554 wlength = scale_expansion_zeroelim(&cat[..catlen], bdztail, &mut w);
555 finlength = fast_expansion_sum_zeroelim(&finnow[..finlength], &w[..wlength], &mut finother);
556 ::core::mem::swap(&mut finnow, &mut finother);
557 }
558 if cdztail != 0.0 {
559 wlength = scale_expansion_zeroelim(&abt[..abtlen], cdztail, &mut w);
560 finlength = fast_expansion_sum_zeroelim(&finnow[..finlength], &w[..wlength], &mut finother);
561 ::core::mem::swap(&mut finnow, &mut finother);
562 }
563
564 finnow[finlength - 1]
565}
566
567pub fn incircle<T: Into<f64>>(pa: Coord<T>, pb: Coord<T>, pc: Coord<T>, pd: Coord<T>) -> f64 {
572 let pa = Coord {
573 x: pa.x.into(),
574 y: pa.y.into(),
575 };
576 let pb = Coord {
577 x: pb.x.into(),
578 y: pb.y.into(),
579 };
580 let pc = Coord {
581 x: pc.x.into(),
582 y: pc.y.into(),
583 };
584 let pd = Coord {
585 x: pd.x.into(),
586 y: pd.y.into(),
587 };
588
589 let adx = pa.x - pd.x;
590 let bdx = pb.x - pd.x;
591 let cdx = pc.x - pd.x;
592 let ady = pa.y - pd.y;
593 let bdy = pb.y - pd.y;
594 let cdy = pc.y - pd.y;
595
596 let bdxcdy = bdx * cdy;
597 let cdxbdy = cdx * bdy;
598 let alift = adx * adx + ady * ady;
599
600 let cdxady = cdx * ady;
601 let adxcdy = adx * cdy;
602 let blift = bdx * bdx + bdy * bdy;
603
604 let adxbdy = adx * bdy;
605 let bdxady = bdx * ady;
606 let clift = cdx * cdx + cdy * cdy;
607
608 let det = alift * (bdxcdy - cdxbdy) + blift * (cdxady - adxcdy) + clift * (adxbdy - bdxady);
609
610 let permanent = (abs(bdxcdy) + abs(cdxbdy)) * alift
611 + (abs(cdxady) + abs(adxcdy)) * blift
612 + (abs(adxbdy) + abs(bdxady)) * clift;
613 let errbound = ICCERRBOUND_A * permanent;
614 if det > errbound || -det > errbound {
615 return det;
616 }
617 incircleadapt(pa, pb, pc, pd, permanent)
618}
619
620fn incircleadapt(
621 pa: Coord<f64>,
622 pb: Coord<f64>,
623 pc: Coord<f64>,
624 pd: Coord<f64>,
625 permanent: f64,
626) -> f64 {
627 let mut temp8 = [0f64; 8];
628 let mut temp16a = [0f64; 16];
629 let mut temp16b = [0f64; 16];
630 let mut temp16c = [0f64; 16];
631 let mut temp32a = [0f64; 32];
632 let mut temp32b = [0f64; 32];
633 let mut temp48 = [0f64; 48];
634 let mut temp64 = [0f64; 64];
635
636 let adx = pa.x - pd.x;
637 let bdx = pb.x - pd.x;
638 let cdx = pc.x - pd.x;
639 let ady = pa.y - pd.y;
640 let bdy = pb.y - pd.y;
641 let cdy = pc.y - pd.y;
642
643 let (bdxcdy1, bdxcdy0) = two_product(bdx, cdy);
644 let (cdxbdy1, cdxbdy0) = two_product(cdx, bdy);
645 let (bc3, bc2, bc1, bc0) = two_two_diff(bdxcdy1, bdxcdy0, cdxbdy1, cdxbdy0);
646 let bc = [bc0, bc1, bc2, bc3];
647
648 let mut axbc = [0f64; 8];
649 let axbclen = scale_expansion_zeroelim(&bc, adx, &mut axbc);
650 let mut axxbc = [0f64; 16];
651 let axxbclen = scale_expansion_zeroelim(&axbc[..axbclen], adx, &mut axxbc);
652 let mut aybc = [0f64; 8];
653 let aybclen = scale_expansion_zeroelim(&bc, ady, &mut aybc);
654 let mut ayybc = [0f64; 16];
655 let ayybclen = scale_expansion_zeroelim(&aybc[..aybclen], ady, &mut ayybc);
656 let mut adet = [0f64; 32];
657 let alen = fast_expansion_sum_zeroelim(&axxbc[0..axxbclen], &ayybc[0..ayybclen], &mut adet);
658
659 let (cdxady1, cdxady0) = two_product(cdx, ady);
660 let (adxcdy1, adxcdy0) = two_product(adx, cdy);
661 let (c3, c2, c1, c0) = two_two_diff(cdxady1, cdxady0, adxcdy1, adxcdy0);
662 let ca = [c0, c1, c2, c3];
663
664 let mut bxca = [0f64; 8];
665 let bxcalen = scale_expansion_zeroelim(&ca, bdx, &mut bxca);
666 let mut bxxca = [0f64; 16];
667 let bxxcalen = scale_expansion_zeroelim(&bxca[..bxcalen], bdx, &mut bxxca);
668 let mut byca = [0f64; 8];
669 let bycalen = scale_expansion_zeroelim(&ca, bdy, &mut byca);
670 let mut byyca = [0f64; 16];
671 let byycalen = scale_expansion_zeroelim(&byca[..bycalen], bdy, &mut byyca);
672 let mut bdet = [0f64; 32];
673 let blen = fast_expansion_sum_zeroelim(&bxxca[..bxxcalen], &byyca[0..byycalen], &mut bdet);
674
675 let (adxbdy1, adxbdy0) = two_product(adx, bdy);
676 let (bdxady1, bdxady0) = two_product(bdx, ady);
677 let (ab3, ab2, ab1, ab0) = two_two_diff(adxbdy1, adxbdy0, bdxady1, bdxady0);
678 let ab = [ab0, ab1, ab2, ab3];
679
680 let mut cxab = [0f64; 8];
681 let cxablen = scale_expansion_zeroelim(&ab, cdx, &mut cxab);
682 let mut cxxab = [0f64; 16];
683 let cxxablen = scale_expansion_zeroelim(&cxab[..cxablen], cdx, &mut cxxab);
684 let mut cyab = [0f64; 8];
685 let cyablen = scale_expansion_zeroelim(&ab, cdy, &mut cyab);
686 let mut cyyab = [0f64; 16];
687 let cyyablen = scale_expansion_zeroelim(&cyab[..cyablen], cdy, &mut cyyab);
688 let mut cdet = [0f64; 32];
689 let clen = fast_expansion_sum_zeroelim(&cxxab[..cxxablen], &cyyab[..cyyablen], &mut cdet);
690
691 let mut abdet = [0f64; 64];
692 let ablen = fast_expansion_sum_zeroelim(&adet[..alen], &bdet[..blen], &mut abdet);
693 let mut fin1 = [0f64; 1152];
694 let mut finlength = fast_expansion_sum_zeroelim(&abdet[..ablen], &cdet[..clen], &mut fin1);
695
696 let mut det = estimate(&fin1[..finlength]);
697 let errbound = ICCERRBOUND_B * permanent;
698 if det >= errbound || -det >= errbound {
699 return det;
700 }
701
702 let adxtail = two_diff_tail(pa.x, pd.x, adx);
703 let adytail = two_diff_tail(pa.y, pd.y, ady);
704 let bdxtail = two_diff_tail(pb.x, pd.x, bdx);
705 let bdytail = two_diff_tail(pb.y, pd.y, bdy);
706 let cdxtail = two_diff_tail(pc.x, pd.x, cdx);
707 let cdytail = two_diff_tail(pc.y, pd.y, cdy);
708 if adxtail == 0.0
709 && bdxtail == 0.0
710 && cdxtail == 0.0
711 && adytail == 0.0
712 && bdytail == 0.0
713 && cdytail == 0.0
714 {
715 return det;
716 }
717
718 let errbound = ICCERRBOUND_C * permanent + RESULTERRBOUND * abs(det);
719 det += ((adx * adx + ady * ady)
720 * ((bdx * cdytail + cdy * bdxtail) - (bdy * cdxtail + cdx * bdytail))
721 + 2.0 * (adx * adxtail + ady * adytail) * (bdx * cdy - bdy * cdx))
722 + ((bdx * bdx + bdy * bdy)
723 * ((cdx * adytail + ady * cdxtail) - (cdy * adxtail + adx * cdytail))
724 + 2.0 * (bdx * bdxtail + bdy * bdytail) * (cdx * ady - cdy * adx))
725 + ((cdx * cdx + cdy * cdy)
726 * ((adx * bdytail + bdy * adxtail) - (ady * bdxtail + bdx * adytail))
727 + 2.0 * (cdx * cdxtail + cdy * cdytail) * (adx * bdy - ady * bdx));
728
729 if det >= errbound || -det >= errbound {
730 return det;
731 }
732
733 let mut fin2 = [0f64; 1152];
734
735 let aa = if bdxtail != 0.0 || bdytail != 0.0 || cdxtail != 0.0 || cdytail != 0.0 {
736 let (adxadx1, adxadx0) = square(adx);
737 let (adyady1, adyady0) = square(ady);
738 let (aa3, aa2, aa1, aa0) = two_two_sum(adxadx1, adxadx0, adyady1, adyady0);
739 [aa0, aa1, aa2, aa3]
740 } else {
741 [0f64; 4]
742 };
743
744 let bb = if cdxtail != 0.0 || cdytail != 0.0 || adxtail != 0.0 || adytail != 0.0 {
745 let (bdxbdx1, bdxbdx0) = square(bdx);
746 let (bdybdy1, bdybdy0) = square(bdy);
747 let (bb3, bb2, bb1, bb0) = two_two_sum(bdxbdx1, bdxbdx0, bdybdy1, bdybdy0);
748 [bb0, bb1, bb2, bb3]
749 } else {
750 [0f64; 4]
751 };
752
753 let cc = if adxtail != 0.0 || adytail != 0.0 || bdxtail != 0.0 || bdytail != 0.0 {
754 let (cdxcdx1, cdxcdx0) = square(cdx);
755 let (cdycdy1, cdycdy0) = square(cdy);
756 let (cc3, cc2, cc1, cc0) = two_two_sum(cdxcdx1, cdxcdx0, cdycdy1, cdycdy0);
757 [cc0, cc1, cc2, cc3]
758 } else {
759 [0f64; 4]
760 };
761
762 let mut axtbclen = 9;
763 let mut axtbc = [0f64; 8];
764 if adxtail != 0.0 {
765 axtbclen = scale_expansion_zeroelim(&bc, adxtail, &mut axtbc);
766 let mut temp16a = [0f64; 16];
767 let temp16alen = scale_expansion_zeroelim(&axtbc[..axtbclen], 2.0 * adx, &mut temp16a);
768
769 let mut axtcc = [0f64; 8];
770 let axtcclen = scale_expansion_zeroelim(&cc, adxtail, &mut axtcc);
771 let mut temp16b = [0f64; 16];
772 let temp16blen = scale_expansion_zeroelim(&axtcc[..axtcclen], bdy, &mut temp16b);
773
774 let mut axtbb = [0f64; 8];
775 let axtbblen = scale_expansion_zeroelim(&bb, adxtail, &mut axtbb);
776 let mut temp16c = [0f64; 16];
777 let temp16clen = scale_expansion_zeroelim(&axtbb[..axtbblen], -cdy, &mut temp16c);
778
779 let mut temp32a = [0f64; 32];
780 let temp32alen = fast_expansion_sum_zeroelim(
781 &temp16a[..temp16alen],
782 &temp16b[..temp16blen],
783 &mut temp32a,
784 );
785 let mut temp48 = [0f64; 48];
786 let temp48len = fast_expansion_sum_zeroelim(
787 &temp16c[..temp16clen],
788 &temp32a[..temp32alen],
789 &mut temp48,
790 );
791 finlength =
792 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp48[..temp48len], &mut fin2);
793 ::core::mem::swap(&mut fin1, &mut fin2)
794 }
795
796 let mut aytbclen = 9;
797 let mut aytbc = [0f64; 8];
798 if adytail != 0.0 {
799 aytbclen = scale_expansion_zeroelim(&bc, adytail, &mut aytbc);
800 let temp16alen = scale_expansion_zeroelim(&aytbc[..aytbclen], 2.0 * ady, &mut temp16a);
801
802 let mut aytbb = [0f64; 8];
803 let aytbblen = scale_expansion_zeroelim(&bb, adytail, &mut aytbb);
804 let temp16blen = scale_expansion_zeroelim(&aytbb[..aytbblen], cdx, &mut temp16b);
805
806 let mut aytcc = [0f64; 8];
807 let aytcclen = scale_expansion_zeroelim(&cc, adytail, &mut aytcc);
808 let temp16clen = scale_expansion_zeroelim(&aytcc[..aytcclen], -bdx, &mut temp16c);
809
810 let temp32alen = fast_expansion_sum_zeroelim(
811 &temp16a[..temp16alen],
812 &temp16b[..temp16blen],
813 &mut temp32a,
814 );
815
816 let temp48len = fast_expansion_sum_zeroelim(
817 &temp16c[..temp16clen],
818 &temp32a[..temp32alen],
819 &mut temp48,
820 );
821 finlength =
822 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp48[..temp48len], &mut fin2);
823 ::core::mem::swap(&mut fin1, &mut fin2)
824 }
825
826 let mut bxtcalen = 9;
827 let mut bxtca = [0f64; 8];
828 if bdxtail != 0.0 {
829 bxtcalen = scale_expansion_zeroelim(&ca, bdxtail, &mut bxtca);
830 let temp16alen = scale_expansion_zeroelim(&bxtca[..bxtcalen], 2.0 * bdx, &mut temp16a);
831
832 let mut bxtaa = [0f64; 8];
833 let bxtaalen = scale_expansion_zeroelim(&aa, bdxtail, &mut bxtaa);
834 let temp16blen = scale_expansion_zeroelim(&bxtaa[..bxtaalen], cdy, &mut temp16b);
835
836 let mut bxtcc = [0f64; 8];
837 let bxtcclen = scale_expansion_zeroelim(&cc, bdxtail, &mut bxtcc);
838 let temp16clen = scale_expansion_zeroelim(&bxtcc[..bxtcclen], -ady, &mut temp16c);
839
840 let temp32alen = fast_expansion_sum_zeroelim(
841 &temp16a[..temp16alen],
842 &temp16b[..temp16blen],
843 &mut temp32a,
844 );
845 let temp48len = fast_expansion_sum_zeroelim(
846 &temp16c[..temp16clen],
847 &temp32a[..temp32alen],
848 &mut temp48,
849 );
850 finlength =
851 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp48[..temp48len], &mut fin2);
852 ::core::mem::swap(&mut fin1, &mut fin2)
853 }
854
855 let mut bytcalen = 9;
856 let mut bytca = [0f64; 8];
857 if bdytail != 0.0 {
858 bytcalen = scale_expansion_zeroelim(&ca, bdytail, &mut bytca);
859 let temp16alen = scale_expansion_zeroelim(&bytca[..bytcalen], 2.0 * bdy, &mut temp16a);
860
861 let mut bytcc = [0f64; 8];
862 let bytcclen = scale_expansion_zeroelim(&cc, bdytail, &mut bytcc);
863 let temp16blen = scale_expansion_zeroelim(&bytcc[..bytcclen], adx, &mut temp16b);
864
865 let mut bytaa = [0f64; 8];
866 let bytaalen = scale_expansion_zeroelim(&aa, bdytail, &mut bytaa);
867 let temp16clen = scale_expansion_zeroelim(&bytaa[..bytaalen], -cdx, &mut temp16c);
868
869 let temp32alen = fast_expansion_sum_zeroelim(
870 &temp16a[..temp16alen],
871 &temp16b[..temp16blen],
872 &mut temp32a,
873 );
874 let temp48len = fast_expansion_sum_zeroelim(
875 &temp16c[..temp16clen],
876 &temp32a[..temp32alen],
877 &mut temp48,
878 );
879
880 finlength =
881 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp48[..temp48len], &mut fin2);
882 ::core::mem::swap(&mut fin1, &mut fin2)
883 }
884
885 let mut cxtab = [0f64; 8];
886 let mut cxtablen = 9;
887 if cdxtail != 0.0 {
888 cxtablen = scale_expansion_zeroelim(&ab, cdxtail, &mut cxtab);
889 let temp16alen = scale_expansion_zeroelim(&cxtab[..cxtablen], 2.0 * cdx, &mut temp16a);
890
891 let mut cxtbb = [0f64; 8];
892 let cxtbblen = scale_expansion_zeroelim(&bb, cdxtail, &mut cxtbb);
893 let temp16blen = scale_expansion_zeroelim(&cxtbb[..cxtbblen], ady, &mut temp16b);
894
895 let mut cxtaa = [0f64; 8];
896 let cxtaalen = scale_expansion_zeroelim(&aa, cdxtail, &mut cxtaa);
897 let temp16clen = scale_expansion_zeroelim(&cxtaa[..cxtaalen], -bdy, &mut temp16c);
898
899 let temp32alen = fast_expansion_sum_zeroelim(
900 &temp16a[..temp16alen],
901 &temp16b[..temp16blen],
902 &mut temp32a,
903 );
904 let temp48len = fast_expansion_sum_zeroelim(
905 &temp16c[..temp16clen],
906 &temp32a[..temp32alen],
907 &mut temp48,
908 );
909 finlength =
910 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp48[..temp48len], &mut fin2);
911 ::core::mem::swap(&mut fin1, &mut fin2);
912 }
913
914 let mut cytab = [0f64; 8];
915 let mut cytablen = 9;
916 if cdytail != 0.0 {
917 cytablen = scale_expansion_zeroelim(&ab, cdytail, &mut cytab);
918 let temp16alen = scale_expansion_zeroelim(&cytab[..cytablen], 2.0 * cdy, &mut temp16a);
919
920 let mut cytaa = [0f64; 8];
921 let cytaalen = scale_expansion_zeroelim(&aa, cdytail, &mut cytaa);
922 let temp16blen = scale_expansion_zeroelim(&cytaa[..cytaalen], bdx, &mut temp16b);
923
924 let mut cytbb = [0f64; 8];
925 let cytbblen = scale_expansion_zeroelim(&bb, cdytail, &mut cytbb);
926 let temp16clen = scale_expansion_zeroelim(&cytbb[..cytbblen], -adx, &mut temp16c);
927
928 let temp32alen = fast_expansion_sum_zeroelim(
929 &temp16a[..temp16alen],
930 &temp16b[..temp16blen],
931 &mut temp32a,
932 );
933 let temp48len = fast_expansion_sum_zeroelim(
934 &temp16c[..temp16clen],
935 &temp32a[..temp32alen],
936 &mut temp48,
937 );
938 finlength =
939 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp48[..temp48len], &mut fin2);
940 ::core::mem::swap(&mut fin1, &mut fin2);
941 }
942
943 if adxtail != 0.0 || adytail != 0.0 {
944 let mut bctt = [0f64; 4];
945 let mut bct = [0f64; 8];
946 let bcttlen;
947 let bctlen;
948 if bdxtail != 0.0 || bdytail != 0.0 || cdxtail != 0.0 || cdytail != 0.0 {
949 let (ti1, ti0) = two_product(bdxtail, cdy);
950 let (tj1, tj0) = two_product(bdx, cdytail);
951 let (u3, u2, u1, u0) = two_two_sum(ti1, ti0, tj1, tj0);
952 let u = [u0, u1, u2, u3];
953 let negate = -bdy;
954 let (ti1, ti0) = two_product(cdxtail, negate);
955 let negate = -bdytail;
956 let (tj1, tj0) = two_product(cdx, negate);
957 let (v3, v2, v1, v0) = two_two_sum(ti1, ti0, tj1, tj0);
958 let v = [v0, v1, v2, v3];
959 bctlen = fast_expansion_sum_zeroelim(&u, &v, &mut bct);
960 let (ti1, ti0) = two_product(bdxtail, cdytail);
961 let (tj1, tj0) = two_product(cdxtail, bdytail);
962 let (bctt3, bctt2, bctt1, bctt0) = two_two_diff(ti1, ti0, tj1, tj0);
963 bctt = [bctt0, bctt1, bctt2, bctt3];
964 bcttlen = 4;
965 } else {
966 bct[0] = 0.0;
967 bctlen = 1;
968 bctt[0] = 0.0;
969 bcttlen = 1;
970 }
971
972 if adxtail != 0.0 {
973 let temp16alen = scale_expansion_zeroelim(&axtbc[..axtbclen], adxtail, &mut temp16a);
974 let mut axtbct = [0f64; 16];
975 let axtbctlen = scale_expansion_zeroelim(&bct[..bctlen], adxtail, &mut axtbct);
976 let temp32alen =
977 scale_expansion_zeroelim(&axtbct[..axtbctlen], 2.0 * adx, &mut temp32a);
978 let temp48len = fast_expansion_sum_zeroelim(
979 &temp16a[..temp16alen],
980 &temp32a[..temp32alen],
981 &mut temp48,
982 );
983 finlength =
984 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp48[..temp48len], &mut fin2);
985 ::core::mem::swap(&mut fin1, &mut fin2);
986
987 if bdytail != 0.0 {
988 let temp8len = scale_expansion_zeroelim(&cc, adxtail, &mut temp8);
989 let temp16alen =
990 scale_expansion_zeroelim(&temp8[..temp8len], bdytail, &mut temp16a);
991 finlength = fast_expansion_sum_zeroelim(
992 &fin1[..finlength],
993 &temp16a[..temp16alen],
994 &mut fin2,
995 );
996 ::core::mem::swap(&mut fin1, &mut fin2);
997 }
998 if cdytail != 0.0 {
999 let temp8len = scale_expansion_zeroelim(&bb, -adxtail, &mut temp8);
1000 let temp16alen =
1001 scale_expansion_zeroelim(&temp8[..temp8len], cdytail, &mut temp16a);
1002 finlength = fast_expansion_sum_zeroelim(
1003 &fin1[..finlength],
1004 &temp16a[..temp16alen],
1005 &mut fin2,
1006 );
1007 ::core::mem::swap(&mut fin1, &mut fin2);
1008 }
1009
1010 let temp32alen = scale_expansion_zeroelim(&axtbct[..axtbctlen], adxtail, &mut temp32a);
1011 let mut axtbctt = [0f64; 8];
1012 let axtbcttlen = scale_expansion_zeroelim(&bctt[..bcttlen], adxtail, &mut axtbctt);
1013 let temp16alen =
1014 scale_expansion_zeroelim(&axtbctt[..axtbcttlen], 2.0 * adx, &mut temp16a);
1015 let temp16blen =
1016 scale_expansion_zeroelim(&axtbctt[..axtbcttlen], adxtail, &mut temp16b);
1017 let temp32blen = fast_expansion_sum_zeroelim(
1018 &temp16a[..temp16alen],
1019 &temp16b[..temp16blen],
1020 &mut temp32b,
1021 );
1022 let temp64len = fast_expansion_sum_zeroelim(
1023 &temp32a[..temp32alen],
1024 &temp32b[..temp32blen],
1025 &mut temp64,
1026 );
1027 finlength =
1028 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp64[..temp64len], &mut fin2);
1029 ::core::mem::swap(&mut fin1, &mut fin2);
1030 }
1031
1032 if adytail != 0.0 {
1033 let temp16alen = scale_expansion_zeroelim(&aytbc[..aytbclen], adytail, &mut temp16a);
1034 let mut aytbct = [0f64; 16];
1035 let aytbctlen = scale_expansion_zeroelim(&bct[..bctlen], adytail, &mut aytbct);
1036 let temp32alen =
1037 scale_expansion_zeroelim(&aytbct[..aytbctlen], 2.0 * ady, &mut temp32a);
1038 let temp48len = fast_expansion_sum_zeroelim(
1039 &temp16a[..temp16alen],
1040 &temp32a[..temp32alen],
1041 &mut temp48,
1042 );
1043 finlength =
1044 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp48[..temp48len], &mut fin2);
1045 ::core::mem::swap(&mut fin1, &mut fin2);
1046
1047 let temp32alen = scale_expansion_zeroelim(&aytbct[..aytbctlen], adytail, &mut temp32a);
1048 let mut aytbctt = [0f64; 8];
1049 let aytbcttlen = scale_expansion_zeroelim(&bctt[..bcttlen], adytail, &mut aytbctt);
1050 let temp16alen =
1051 scale_expansion_zeroelim(&aytbctt[..aytbcttlen], 2.0 * ady, &mut temp16a);
1052 let temp16blen =
1053 scale_expansion_zeroelim(&aytbctt[..aytbcttlen], adytail, &mut temp16b);
1054 let temp32blen = fast_expansion_sum_zeroelim(
1055 &temp16a[..temp16alen],
1056 &temp16b[..temp16blen],
1057 &mut temp32b,
1058 );
1059 let temp64len = fast_expansion_sum_zeroelim(
1060 &temp32a[..temp32alen],
1061 &temp32b[..temp32blen],
1062 &mut temp64,
1063 );
1064 finlength =
1065 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp64[..temp64len], &mut fin2);
1066 ::core::mem::swap(&mut fin1, &mut fin2);
1067 }
1068 }
1069
1070 if bdxtail != 0.0 || bdytail != 0.0 {
1071 let mut catt = [0f64; 4];
1072 let mut cat = [0f64; 8];
1073 let cattlen;
1074 let catlen;
1075
1076 if cdxtail != 0.0 || cdytail != 0.0 || adxtail != 0.0 || adytail != 0.0 {
1077 let (ti1, ti0) = two_product(cdxtail, ady);
1078 let (tj1, tj0) = two_product(cdx, adytail);
1079 let (u3, u2, u1, u0) = two_two_sum(ti1, ti0, tj1, tj0);
1080 let u = [u0, u1, u2, u3];
1081 let negate = -cdy;
1082 let (ti1, ti0) = two_product(adxtail, negate);
1083 let negate = -cdytail;
1084 let (tj1, tj0) = two_product(adx, negate);
1085 let (v3, v2, v1, v0) = two_two_sum(ti1, ti0, tj1, tj0);
1086 let v = [v0, v1, v2, v3];
1087 catlen = fast_expansion_sum_zeroelim(&u, &v, &mut cat);
1088
1089 let (ti1, ti0) = two_product(cdxtail, adytail);
1090 let (tj1, tj0) = two_product(adxtail, cdytail);
1091 let (catt3, catt2, catt1, catt0) = two_two_diff(ti1, ti0, tj1, tj0);
1092 catt = [catt0, catt1, catt2, catt3];
1093 cattlen = 4;
1094 } else {
1095 cat[0] = 0.0;
1096 catlen = 1;
1097 catt[0] = 0.0;
1098 cattlen = 1;
1099 }
1100
1101 if bdxtail != 0.0 {
1102 let temp16alen = scale_expansion_zeroelim(&bxtca[..bxtcalen], bdxtail, &mut temp16a);
1103 let mut bxtcat = [0f64; 16];
1104 let bxtcatlen = scale_expansion_zeroelim(&cat[..catlen], bdxtail, &mut bxtcat);
1105 let temp32alen =
1106 scale_expansion_zeroelim(&bxtcat[..bxtcatlen], 2.0 * bdx, &mut temp32a);
1107 let temp48len = fast_expansion_sum_zeroelim(
1108 &temp16a[..temp16alen],
1109 &temp32a[..temp32alen],
1110 &mut temp48,
1111 );
1112 finlength =
1113 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp48[..temp48len], &mut fin2);
1114 ::core::mem::swap(&mut fin1, &mut fin2);
1115
1116 if cdytail != 0.0 {
1117 let temp8len = scale_expansion_zeroelim(&aa, bdxtail, &mut temp8);
1118 let temp16alen =
1119 scale_expansion_zeroelim(&temp8[..temp8len], cdytail, &mut temp16a);
1120 finlength = fast_expansion_sum_zeroelim(
1121 &fin1[..finlength],
1122 &temp16a[..temp16alen],
1123 &mut fin2,
1124 );
1125 ::core::mem::swap(&mut fin1, &mut fin2);
1126 }
1127 if adytail != 0.0 {
1128 let temp8len = scale_expansion_zeroelim(&cc, -bdxtail, &mut temp8);
1129 let temp16alen =
1130 scale_expansion_zeroelim(&temp8[..temp8len], adytail, &mut temp16a);
1131 finlength = fast_expansion_sum_zeroelim(
1132 &fin1[..finlength],
1133 &temp16a[..temp16alen],
1134 &mut fin2,
1135 );
1136 ::core::mem::swap(&mut fin1, &mut fin2);
1137 }
1138
1139 let temp32alen = scale_expansion_zeroelim(&bxtcat[..bxtcatlen], bdxtail, &mut temp32a);
1140 let mut bxtcatt = [0f64; 8];
1141 let bxtcattlen = scale_expansion_zeroelim(&catt[..cattlen], bdxtail, &mut bxtcatt);
1142 let temp16alen =
1143 scale_expansion_zeroelim(&bxtcatt[..bxtcattlen], 2.0 * bdx, &mut temp16a);
1144 let temp16blen =
1145 scale_expansion_zeroelim(&bxtcatt[..bxtcattlen], bdxtail, &mut temp16b);
1146 let temp32blen = fast_expansion_sum_zeroelim(
1147 &temp16a[..temp16alen],
1148 &temp16b[..temp16blen],
1149 &mut temp32b,
1150 );
1151 let temp64len = fast_expansion_sum_zeroelim(
1152 &temp32a[..temp32alen],
1153 &temp32b[..temp32blen],
1154 &mut temp64,
1155 );
1156 finlength =
1157 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp64[..temp64len], &mut fin2);
1158 ::core::mem::swap(&mut fin1, &mut fin2);
1159 }
1160 if bdytail != 0.0 {
1161 let temp16alen = scale_expansion_zeroelim(&bytca[..bytcalen], bdytail, &mut temp16a);
1162 let mut bytcat = [0f64; 16];
1163 let bytcatlen = scale_expansion_zeroelim(&cat[..catlen], bdytail, &mut bytcat);
1164 let temp32alen =
1165 scale_expansion_zeroelim(&bytcat[..bytcatlen], 2.0 * bdy, &mut temp32a);
1166 let temp48len = fast_expansion_sum_zeroelim(
1167 &temp16a[..temp16alen],
1168 &temp32a[..temp32alen],
1169 &mut temp48,
1170 );
1171 finlength =
1172 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp48[..temp48len], &mut fin2);
1173 ::core::mem::swap(&mut fin1, &mut fin2);
1174
1175 let temp32alen = scale_expansion_zeroelim(&bytcat[..bytcatlen], bdytail, &mut temp32a);
1176 let mut bytcatt = [0f64; 8];
1177 let bytcattlen = scale_expansion_zeroelim(&catt[..cattlen], bdytail, &mut bytcatt);
1178 let temp16alen =
1179 scale_expansion_zeroelim(&bytcatt[..bytcattlen], 2.0 * bdy, &mut temp16a);
1180 let temp16blen =
1181 scale_expansion_zeroelim(&bytcatt[..bytcattlen], bdytail, &mut temp16b);
1182 let temp32blen = fast_expansion_sum_zeroelim(
1183 &temp16a[..temp16alen],
1184 &temp16b[..temp16blen],
1185 &mut temp32b,
1186 );
1187 let temp64len = fast_expansion_sum_zeroelim(
1188 &temp32a[..temp32alen],
1189 &temp32b[..temp32blen],
1190 &mut temp64,
1191 );
1192 finlength =
1193 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp64[..temp64len], &mut fin2);
1194 ::core::mem::swap(&mut fin1, &mut fin2);
1195 }
1196 }
1197
1198 if cdxtail != 0.0 || cdytail != 0.0 {
1199 let mut abtt = [0f64; 4];
1200 let mut abt = [0f64; 8];
1201 let abttlen;
1202 let abtlen;
1203
1204 if adxtail != 0.0 || adytail != 0.0 || bdxtail != 0.0 || bdytail != 0.0 {
1205 let (ti1, ti0) = two_product(adxtail, bdy);
1206 let (tj1, tj0) = two_product(adx, bdytail);
1207 let (u3, u2, u1, u0) = two_two_sum(ti1, ti0, tj1, tj0);
1208 let u = [u0, u1, u2, u3];
1209 let negate = -ady;
1210 let (ti1, ti0) = two_product(bdxtail, negate);
1211 let negate = -adytail;
1212 let (tj1, tj0) = two_product(bdx, negate);
1213 let (v3, v2, v1, v0) = two_two_sum(ti1, ti0, tj1, tj0);
1214 let v = [v0, v1, v2, v3];
1215 abtlen = fast_expansion_sum_zeroelim(&u, &v, &mut abt);
1216
1217 let (ti1, ti0) = two_product(adxtail, bdytail);
1218 let (tj1, tj0) = two_product(bdxtail, adytail);
1219 let (abtt3, abtt2, abtt1, abtt0) = two_two_diff(ti1, ti0, tj1, tj0);
1220 abtt = [abtt0, abtt1, abtt2, abtt3];
1221 abttlen = 4;
1222 } else {
1223 abt[0] = 0.0;
1224 abtlen = 1;
1225 abtt[0] = 0.0;
1226 abttlen = 1;
1227 }
1228
1229 if cdxtail != 0.0 {
1230 let temp16alen = scale_expansion_zeroelim(&cxtab[..cxtablen], cdxtail, &mut temp16a);
1231 let mut cxtabt = [0f64; 16];
1232 let cxtabtlen = scale_expansion_zeroelim(&abt[..abtlen], cdxtail, &mut cxtabt);
1233 let temp32alen =
1234 scale_expansion_zeroelim(&cxtabt[..cxtabtlen], 2.0 * cdx, &mut temp32a);
1235 let temp48len = fast_expansion_sum_zeroelim(
1236 &temp16a[..temp16alen],
1237 &temp32a[..temp32alen],
1238 &mut temp48,
1239 );
1240 finlength =
1241 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp48[..temp48len], &mut fin2);
1242 ::core::mem::swap(&mut fin1, &mut fin2);
1243
1244 if adytail != 0.0 {
1245 let temp8len = scale_expansion_zeroelim(&bb, cdxtail, &mut temp8);
1246 let temp16alen =
1247 scale_expansion_zeroelim(&temp8[..temp8len], adytail, &mut temp16a);
1248 finlength = fast_expansion_sum_zeroelim(
1249 &fin1[..finlength],
1250 &temp16a[..temp16alen],
1251 &mut fin2,
1252 );
1253 ::core::mem::swap(&mut fin1, &mut fin2);
1254 }
1255 if bdytail != 0.0 {
1256 let temp8len = scale_expansion_zeroelim(&aa, -cdxtail, &mut temp8);
1257 let temp16alen =
1258 scale_expansion_zeroelim(&temp8[..temp8len], bdytail, &mut temp16a);
1259 finlength = fast_expansion_sum_zeroelim(
1260 &fin1[..finlength],
1261 &temp16a[..temp16alen],
1262 &mut fin2,
1263 );
1264 ::core::mem::swap(&mut fin1, &mut fin2);
1265 }
1266
1267 let temp32alen = scale_expansion_zeroelim(&cxtabt[..cxtabtlen], cdxtail, &mut temp32a);
1268 let mut cxtabtt = [0f64; 8];
1269 let cxtabttlen = scale_expansion_zeroelim(&abtt[..abttlen], cdxtail, &mut cxtabtt);
1270 let temp16alen =
1271 scale_expansion_zeroelim(&cxtabtt[..cxtabttlen], 2.0 * cdx, &mut temp16a);
1272 let temp16blen =
1273 scale_expansion_zeroelim(&cxtabtt[..cxtabttlen], cdxtail, &mut temp16b);
1274 let temp32blen = fast_expansion_sum_zeroelim(
1275 &temp16a[..temp16alen],
1276 &temp16b[..temp16blen],
1277 &mut temp32b,
1278 );
1279 let temp64len = fast_expansion_sum_zeroelim(
1280 &temp32a[..temp32alen],
1281 &temp32b[..temp32blen],
1282 &mut temp64,
1283 );
1284 finlength =
1285 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp64[..temp64len], &mut fin2);
1286 ::core::mem::swap(&mut fin1, &mut fin2);
1287 }
1288 if cdytail != 0.0 {
1289 let temp16alen = scale_expansion_zeroelim(&cytab[..cytablen], cdytail, &mut temp16a);
1290 let mut cytabt = [0f64; 16];
1291 let cytabtlen = scale_expansion_zeroelim(&abt[..abtlen], cdytail, &mut cytabt);
1292 let temp32alen =
1293 scale_expansion_zeroelim(&cytabt[..cytabtlen], 2.0 * cdy, &mut temp32a);
1294 let temp48len = fast_expansion_sum_zeroelim(
1295 &temp16a[..temp16alen],
1296 &temp32a[..temp32alen],
1297 &mut temp48,
1298 );
1299 finlength =
1300 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp48[..temp48len], &mut fin2);
1301 ::core::mem::swap(&mut fin1, &mut fin2);
1302
1303 let temp32alen = scale_expansion_zeroelim(&cytabt[..cytabtlen], cdytail, &mut temp32a);
1304 let mut cytabtt = [0f64; 8];
1305 let cytabttlen = scale_expansion_zeroelim(&abtt[..abttlen], cdytail, &mut cytabtt);
1306 let temp16alen =
1307 scale_expansion_zeroelim(&cytabtt[..cytabttlen], 2.0 * cdy, &mut temp16a);
1308 let temp16blen =
1309 scale_expansion_zeroelim(&cytabtt[..cytabttlen], cdytail, &mut temp16b);
1310 let temp32blen = fast_expansion_sum_zeroelim(
1311 &temp16a[..temp16alen],
1312 &temp16b[..temp16blen],
1313 &mut temp32b,
1314 );
1315 let temp64len = fast_expansion_sum_zeroelim(
1316 &temp32a[..temp32alen],
1317 &temp32b[..temp32blen],
1318 &mut temp64,
1319 );
1320 finlength =
1321 fast_expansion_sum_zeroelim(&fin1[..finlength], &temp64[..temp64len], &mut fin2);
1322 ::core::mem::swap(&mut fin1, &mut fin2);
1323 }
1324 }
1325 fin1[finlength - 1]
1326}
1327
1328pub fn insphere<T: Into<f64>>(
1333 pa: Coord3D<T>,
1334 pb: Coord3D<T>,
1335 pc: Coord3D<T>,
1336 pd: Coord3D<T>,
1337 pe: Coord3D<T>,
1338) -> f64 {
1339 let pa = Coord3D {
1340 x: pa.x.into(),
1341 y: pa.y.into(),
1342 z: pa.z.into(),
1343 };
1344 let pb = Coord3D {
1345 x: pb.x.into(),
1346 y: pb.y.into(),
1347 z: pb.z.into(),
1348 };
1349 let pc = Coord3D {
1350 x: pc.x.into(),
1351 y: pc.y.into(),
1352 z: pc.z.into(),
1353 };
1354 let pd = Coord3D {
1355 x: pd.x.into(),
1356 y: pd.y.into(),
1357 z: pd.z.into(),
1358 };
1359 let pe = Coord3D {
1360 x: pe.x.into(),
1361 y: pe.y.into(),
1362 z: pe.z.into(),
1363 };
1364
1365 let aex = pa.x - pe.x;
1366 let bex = pb.x - pe.x;
1367 let cex = pc.x - pe.x;
1368 let dex = pd.x - pe.x;
1369 let aey = pa.y - pe.y;
1370 let bey = pb.y - pe.y;
1371 let cey = pc.y - pe.y;
1372 let dey = pd.y - pe.y;
1373 let aez = pa.z - pe.z;
1374 let bez = pb.z - pe.z;
1375 let cez = pc.z - pe.z;
1376 let dez = pd.z - pe.z;
1377
1378 let aexbey = aex * bey;
1379 let bexaey = bex * aey;
1380 let ab = aexbey - bexaey;
1381 let bexcey = bex * cey;
1382 let cexbey = cex * bey;
1383 let bc = bexcey - cexbey;
1384 let cexdey = cex * dey;
1385 let dexcey = dex * cey;
1386 let cd = cexdey - dexcey;
1387 let dexaey = dex * aey;
1388 let aexdey = aex * dey;
1389 let da = dexaey - aexdey;
1390
1391 let aexcey = aex * cey;
1392 let cexaey = cex * aey;
1393 let ac = aexcey - cexaey;
1394 let bexdey = bex * dey;
1395 let dexbey = dex * bey;
1396 let bd = bexdey - dexbey;
1397
1398 let abc = aez * bc - bez * ac + cez * ab;
1399 let bcd = bez * cd - cez * bd + dez * bc;
1400 let cda = cez * da + dez * ac + aez * cd;
1401 let dab = dez * ab + aez * bd + bez * da;
1402
1403 let alift = aex * aex + aey * aey + aez * aez;
1404 let blift = bex * bex + bey * bey + bez * bez;
1405 let clift = cex * cex + cey * cey + cez * cez;
1406 let dlift = dex * dex + dey * dey + dez * dez;
1407
1408 let det = (dlift * abc - clift * dab) + (blift * cda - alift * bcd);
1409
1410 let aezplus = abs(aez);
1411 let bezplus = abs(bez);
1412 let cezplus = abs(cez);
1413 let dezplus = abs(dez);
1414 let aexbeyplus = abs(aexbey);
1415 let bexaeyplus = abs(bexaey);
1416 let bexceyplus = abs(bexcey);
1417 let cexbeyplus = abs(cexbey);
1418 let cexdeyplus = abs(cexdey);
1419 let dexceyplus = abs(dexcey);
1420 let dexaeyplus = abs(dexaey);
1421 let aexdeyplus = abs(aexdey);
1422 let aexceyplus = abs(aexcey);
1423 let cexaeyplus = abs(cexaey);
1424 let bexdeyplus = abs(bexdey);
1425 let dexbeyplus = abs(dexbey);
1426 let permanent = ((cexdeyplus + dexceyplus) * bezplus
1427 + (dexbeyplus + bexdeyplus) * cezplus
1428 + (bexceyplus + cexbeyplus) * dezplus)
1429 * alift
1430 + ((dexaeyplus + aexdeyplus) * cezplus
1431 + (aexceyplus + cexaeyplus) * dezplus
1432 + (cexdeyplus + dexceyplus) * aezplus)
1433 * blift
1434 + ((aexbeyplus + bexaeyplus) * dezplus
1435 + (bexdeyplus + dexbeyplus) * aezplus
1436 + (dexaeyplus + aexdeyplus) * bezplus)
1437 * clift
1438 + ((bexceyplus + cexbeyplus) * aezplus
1439 + (cexaeyplus + aexceyplus) * bezplus
1440 + (aexbeyplus + bexaeyplus) * cezplus)
1441 * dlift;
1442 let errbound = ISPERRBOUND_A * permanent;
1443 if det > errbound || -det > errbound {
1444 return det;
1445 }
1446
1447 insphereadapt(pa, pb, pc, pd, pe, permanent)
1448}
1449
1450fn insphereadapt<T: Into<f64>>(
1451 pa: Coord3D<T>,
1452 pb: Coord3D<T>,
1453 pc: Coord3D<T>,
1454 pd: Coord3D<T>,
1455 pe: Coord3D<T>,
1456 permanent: f64,
1457) -> f64 {
1458 let pa = Coord3D {
1459 x: pa.x.into(),
1460 y: pa.y.into(),
1461 z: pa.z.into(),
1462 };
1463 let pb = Coord3D {
1464 x: pb.x.into(),
1465 y: pb.y.into(),
1466 z: pb.z.into(),
1467 };
1468 let pc = Coord3D {
1469 x: pc.x.into(),
1470 y: pc.y.into(),
1471 z: pc.z.into(),
1472 };
1473 let pd = Coord3D {
1474 x: pd.x.into(),
1475 y: pd.y.into(),
1476 z: pd.z.into(),
1477 };
1478 let pe = Coord3D {
1479 x: pe.x.into(),
1480 y: pe.y.into(),
1481 z: pe.z.into(),
1482 };
1483
1484 let aex = pa.x - pe.x;
1485 let bex = pb.x - pe.x;
1486 let cex = pc.x - pe.x;
1487 let dex = pd.x - pe.x;
1488 let aey = pa.y - pe.y;
1489 let bey = pb.y - pe.y;
1490 let cey = pc.y - pe.y;
1491 let dey = pd.y - pe.y;
1492 let aez = pa.z - pe.z;
1493 let bez = pb.z - pe.z;
1494 let cez = pc.z - pe.z;
1495 let dez = pd.z - pe.z;
1496
1497 let (aexbey1, aexbey0) = two_product(aex, bey);
1498 let (bexaey1, bexaey0) = two_product(bex, aey);
1499 let (ab3, ab2, ab1, ab0) = two_two_diff(aexbey1, aexbey0, bexaey1, bexaey0);
1500 let ab = [ab0, ab1, ab2, ab3];
1501
1502 let (bexcey1, bexcey0) = two_product(bex, cey);
1503 let (cexbey1, cexbey0) = two_product(cex, bey);
1504 let (bc3, bc2, bc1, bc0) = two_two_diff(bexcey1, bexcey0, cexbey1, cexbey0);
1505 let bc = [bc0, bc1, bc2, bc3];
1506
1507 let (cexdey1, cexdey0) = two_product(cex, dey);
1508 let (dexcey1, dexcey0) = two_product(dex, cey);
1509 let (cd3, cd2, cd1, cd0) = two_two_diff(cexdey1, cexdey0, dexcey1, dexcey0);
1510 let cd = [cd0, cd1, cd2, cd3];
1511
1512 let (dexaey1, dexaey0) = two_product(dex, aey);
1513 let (aexdey1, aexdey0) = two_product(aex, dey);
1514 let (da3, da2, da1, da0) = two_two_diff(dexaey1, dexaey0, aexdey1, aexdey0);
1515 let da = [da0, da1, da2, da3];
1516
1517 let (aexcey1, aexcey0) = two_product(aex, cey);
1518 let (cexaey1, cexaey0) = two_product(cex, aey);
1519 let (ac3, ac2, ac1, ac0) = two_two_diff(aexcey1, aexcey0, cexaey1, cexaey0);
1520 let ac = [ac0, ac1, ac2, ac3];
1521
1522 let (bexdey1, bexdey0) = two_product(bex, dey);
1523 let (dexbey1, dexbey0) = two_product(dex, bey);
1524 let (bd3, bd2, bd1, bd0) = two_two_diff(bexdey1, bexdey0, dexbey1, dexbey0);
1525 let bd = [bd0, bd1, bd2, bd3];
1526
1527 let mut temp8a = [0f64; 8];
1528 let mut temp8b = [0f64; 8];
1529 let mut temp8c = [0f64; 8];
1530 let mut temp16 = [0f64; 16];
1531 let mut temp24 = [0f64; 24];
1532 let mut temp48 = [0f64; 48];
1533 let mut xdet = [0f64; 96];
1534 let mut ydet = [0f64; 96];
1535 let mut zdet = [0f64; 96];
1536 let mut xydet = [0f64; 192];
1537 let mut adet = [0f64; 288];
1538 let mut bdet = [0f64; 288];
1539 let mut cdet = [0f64; 288];
1540 let mut ddet = [0f64; 288];
1541 let mut abdet = [0f64; 576];
1542 let mut cddet = [0f64; 576];
1543 let mut fin1 = [0f64; 1152];
1544
1545 let mut temp8alen = scale_expansion_zeroelim(&cd[..4], bez, &mut temp8a);
1546 let mut temp8blen = scale_expansion_zeroelim(&bd[..4], -cez, &mut temp8b);
1547 let mut temp8clen = scale_expansion_zeroelim(&bc[..4], dez, &mut temp8c);
1548 let mut temp16len =
1549 fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp8b[..temp8blen], &mut temp16);
1550 let mut temp24len =
1551 fast_expansion_sum_zeroelim(&temp8c[..temp8clen], &temp16[..temp16len], &mut temp24);
1552 let temp48len = scale_expansion_zeroelim(&temp24[..temp24len], aex, &mut temp48);
1553 let mut xlen = scale_expansion_zeroelim(&temp48[..temp48len], -aex, &mut xdet);
1554 let temp48len = scale_expansion_zeroelim(&temp24[..temp24len], aey, &mut temp48);
1555 let mut ylen = scale_expansion_zeroelim(&temp48[..temp48len], -aey, &mut ydet);
1556 let mut temp48len = scale_expansion_zeroelim(&temp24[..temp24len], aez, &mut temp48);
1557 let mut zlen = scale_expansion_zeroelim(&temp48[..temp48len], -aez, &mut zdet);
1558 let mut xylen = fast_expansion_sum_zeroelim(&xdet[..xlen], &ydet[..ylen], &mut xydet);
1559 let alen = fast_expansion_sum_zeroelim(&xydet[..xylen], &zdet[..zlen], &mut adet);
1560
1561 temp8alen = scale_expansion_zeroelim(&da[..4], cez, &mut temp8a);
1562 temp8blen = scale_expansion_zeroelim(&ac[..4], dez, &mut temp8b);
1563 temp8clen = scale_expansion_zeroelim(&cd[..4], aez, &mut temp8c);
1564 temp16len =
1565 fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp8b[..temp8blen], &mut temp16);
1566 temp24len =
1567 fast_expansion_sum_zeroelim(&temp8c[..temp8clen], &temp16[..temp16len], &mut temp24);
1568 temp48len = scale_expansion_zeroelim(&temp24[..temp24len], bex, &mut temp48);
1569 xlen = scale_expansion_zeroelim(&temp48[..temp48len], bex, &mut xdet);
1570 temp48len = scale_expansion_zeroelim(&temp24[..temp24len], bey, &mut temp48);
1571 ylen = scale_expansion_zeroelim(&temp48[..temp48len], bey, &mut ydet);
1572 temp48len = scale_expansion_zeroelim(&temp24[..temp24len], bez, &mut temp48);
1573 zlen = scale_expansion_zeroelim(&temp48[..temp48len], bez, &mut zdet);
1574 xylen = fast_expansion_sum_zeroelim(&xdet[..xlen], &ydet[..ylen], &mut xydet);
1575 let blen = fast_expansion_sum_zeroelim(&xydet[..xylen], &zdet[..zlen], &mut bdet);
1576
1577 temp8alen = scale_expansion_zeroelim(&ab[..4], dez, &mut temp8a);
1578 temp8blen = scale_expansion_zeroelim(&bd[..4], aez, &mut temp8b);
1579 temp8clen = scale_expansion_zeroelim(&da[..4], bez, &mut temp8c);
1580 temp16len =
1581 fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp8b[..temp8blen], &mut temp16);
1582 temp24len =
1583 fast_expansion_sum_zeroelim(&temp8c[..temp8clen], &temp16[..temp16len], &mut temp24);
1584 temp48len = scale_expansion_zeroelim(&temp24[..temp24len], cex, &mut temp48);
1585 xlen = scale_expansion_zeroelim(&temp48[..temp48len], -cex, &mut xdet);
1586 temp48len = scale_expansion_zeroelim(&temp24[..temp24len], cey, &mut temp48);
1587 ylen = scale_expansion_zeroelim(&temp48[..temp48len], -cey, &mut ydet);
1588 temp48len = scale_expansion_zeroelim(&temp24[..temp24len], cez, &mut temp48);
1589 zlen = scale_expansion_zeroelim(&temp48[..temp48len], -cez, &mut zdet);
1590 xylen = fast_expansion_sum_zeroelim(&xdet[..xlen], &ydet[..ylen], &mut xydet);
1591 let clen = fast_expansion_sum_zeroelim(&xydet[..xylen], &zdet[..zlen], &mut cdet);
1592
1593 temp8alen = scale_expansion_zeroelim(&bc[..4], aez, &mut temp8a);
1594 temp8blen = scale_expansion_zeroelim(&ac[..4], -bez, &mut temp8b);
1595 temp8clen = scale_expansion_zeroelim(&ab[..4], cez, &mut temp8c);
1596 temp16len =
1597 fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp8b[..temp8blen], &mut temp16);
1598 temp24len =
1599 fast_expansion_sum_zeroelim(&temp8c[..temp8clen], &temp16[..temp16len], &mut temp24);
1600 temp48len = scale_expansion_zeroelim(&temp24[..temp24len], dex, &mut temp48);
1601 xlen = scale_expansion_zeroelim(&temp48[..temp48len], dex, &mut xdet);
1602 temp48len = scale_expansion_zeroelim(&temp24[..temp24len], dey, &mut temp48);
1603 ylen = scale_expansion_zeroelim(&temp48[..temp48len], dey, &mut ydet);
1604 temp48len = scale_expansion_zeroelim(&temp24[..temp24len], dez, &mut temp48);
1605 zlen = scale_expansion_zeroelim(&temp48[..temp48len], dez, &mut zdet);
1606 xylen = fast_expansion_sum_zeroelim(&xdet[..xlen], &ydet[..ylen], &mut xydet);
1607 let dlen = fast_expansion_sum_zeroelim(&xydet[..xylen], &zdet[..zlen], &mut ddet);
1608
1609 let ablen = fast_expansion_sum_zeroelim(&adet[..alen], &bdet[..blen], &mut abdet);
1610 let cdlen = fast_expansion_sum_zeroelim(&cdet[..clen], &ddet[..dlen], &mut cddet);
1611 let finlength = fast_expansion_sum_zeroelim(&abdet[..ablen], &cddet[..cdlen], &mut fin1);
1612
1613 let mut det = estimate(&fin1[..finlength]);
1614 let mut errbound = ISPERRBOUND_B * permanent;
1615 if det >= errbound || -det >= errbound {
1616 return det;
1617 }
1618
1619 let aextail = two_diff_tail(pa.x, pe.x, aex);
1620 let aeytail = two_diff_tail(pa.y, pe.y, aey);
1621 let aeztail = two_diff_tail(pa.z, pe.z, aez);
1622 let bextail = two_diff_tail(pb.x, pe.x, bex);
1623 let beytail = two_diff_tail(pb.y, pe.y, bey);
1624 let beztail = two_diff_tail(pb.z, pe.z, bez);
1625 let cextail = two_diff_tail(pc.x, pe.x, cex);
1626 let ceytail = two_diff_tail(pc.y, pe.y, cey);
1627 let ceztail = two_diff_tail(pc.z, pe.z, cez);
1628 let dextail = two_diff_tail(pd.x, pe.x, dex);
1629 let deytail = two_diff_tail(pd.y, pe.y, dey);
1630 let deztail = two_diff_tail(pd.z, pe.z, dez);
1631 if (aextail == 0.0)
1632 && (aeytail == 0.0)
1633 && (aeztail == 0.0)
1634 && (bextail == 0.0)
1635 && (beytail == 0.0)
1636 && (beztail == 0.0)
1637 && (cextail == 0.0)
1638 && (ceytail == 0.0)
1639 && (ceztail == 0.0)
1640 && (dextail == 0.0)
1641 && (deytail == 0.0)
1642 && (deztail == 0.0)
1643 {
1644 return det;
1645 }
1646
1647 errbound = ISPERRBOUND_C * permanent + RESULTERRBOUND * abs(det);
1648 let abeps = (aex * beytail + bey * aextail) - (aey * bextail + bex * aeytail);
1649 let bceps = (bex * ceytail + cey * bextail) - (bey * cextail + cex * beytail);
1650 let cdeps = (cex * deytail + dey * cextail) - (cey * dextail + dex * ceytail);
1651 let daeps = (dex * aeytail + aey * dextail) - (dey * aextail + aex * deytail);
1652 let aceps = (aex * ceytail + cey * aextail) - (aey * cextail + cex * aeytail);
1653 let bdeps = (bex * deytail + dey * bextail) - (bey * dextail + dex * beytail);
1654 det += (((bex * bex + bey * bey + bez * bez)
1655 * ((cez * daeps + dez * aceps + aez * cdeps)
1656 + (ceztail * da3 + deztail * ac3 + aeztail * cd3))
1657 + (dex * dex + dey * dey + dez * dez)
1658 * ((aez * bceps - bez * aceps + cez * abeps)
1659 + (aeztail * bc3 - beztail * ac3 + ceztail * ab3)))
1660 - ((aex * aex + aey * aey + aez * aez)
1661 * ((bez * cdeps - cez * bdeps + dez * bceps)
1662 + (beztail * cd3 - ceztail * bd3 + deztail * bc3))
1663 + (cex * cex + cey * cey + cez * cez)
1664 * ((dez * abeps + aez * bdeps + bez * daeps)
1665 + (deztail * ab3 + aeztail * bd3 + beztail * da3))))
1666 + 2.0
1667 * (((bex * bextail + bey * beytail + bez * beztail)
1668 * (cez * da3 + dez * ac3 + aez * cd3)
1669 + (dex * dextail + dey * deytail + dez * deztail)
1670 * (aez * bc3 - bez * ac3 + cez * ab3))
1671 - ((aex * aextail + aey * aeytail + aez * aeztail)
1672 * (bez * cd3 - cez * bd3 + dez * bc3)
1673 + (cex * cextail + cey * ceytail + cez * ceztail)
1674 * (dez * ab3 + aez * bd3 + bez * da3)));
1675 if det >= errbound || -det >= errbound {
1676 return det;
1677 }
1678
1679 insphereexact(pa, pb, pc, pd, pe)
1680}
1681
1682fn insphereexact<T: Into<f64>>(
1683 pa: Coord3D<T>,
1684 pb: Coord3D<T>,
1685 pc: Coord3D<T>,
1686 pd: Coord3D<T>,
1687 pe: Coord3D<T>,
1688) -> f64 {
1689 let pa = Coord3D {
1690 x: pa.x.into(),
1691 y: pa.y.into(),
1692 z: pa.z.into(),
1693 };
1694 let pb = Coord3D {
1695 x: pb.x.into(),
1696 y: pb.y.into(),
1697 z: pb.z.into(),
1698 };
1699 let pc = Coord3D {
1700 x: pc.x.into(),
1701 y: pc.y.into(),
1702 z: pc.z.into(),
1703 };
1704 let pd = Coord3D {
1705 x: pd.x.into(),
1706 y: pd.y.into(),
1707 z: pd.z.into(),
1708 };
1709 let pe = Coord3D {
1710 x: pe.x.into(),
1711 y: pe.y.into(),
1712 z: pe.z.into(),
1713 };
1714
1715 let (axby1, axby0) = two_product(pa.x, pb.y);
1716 let (bxay1, bxay0) = two_product(pb.x, pa.y);
1717 let (ab3, ab2, ab1, ab0) = two_two_diff(axby1, axby0, bxay1, bxay0);
1718 let ab = [ab0, ab1, ab2, ab3];
1719
1720 let (bxcy1, bxcy0) = two_product(pb.x, pc.y);
1721 let (cxby1, cxby0) = two_product(pc.x, pb.y);
1722 let (bc3, bc2, bc1, bc0) = two_two_diff(bxcy1, bxcy0, cxby1, cxby0);
1723 let bc = [bc0, bc1, bc2, bc3];
1724
1725 let (cxdy1, cxdy0) = two_product(pc.x, pd.y);
1726 let (dxcy1, dxcy0) = two_product(pd.x, pc.y);
1727 let (cd3, cd2, cd1, cd0) = two_two_diff(cxdy1, cxdy0, dxcy1, dxcy0);
1728 let cd = [cd0, cd1, cd2, cd3];
1729
1730 let (dxey1, dxey0) = two_product(pd.x, pe.y);
1731 let (exdy1, exdy0) = two_product(pe.x, pd.y);
1732 let (de3, de2, de1, de0) = two_two_diff(dxey1, dxey0, exdy1, exdy0);
1733 let de = [de0, de1, de2, de3];
1734
1735 let (exay1, exay0) = two_product(pe.x, pa.y);
1736 let (axey1, axey0) = two_product(pa.x, pe.y);
1737 let (ea3, ea2, ea1, ea0) = two_two_diff(exay1, exay0, axey1, axey0);
1738 let ea = [ea0, ea1, ea2, ea3];
1739
1740 let (axcy1, axcy0) = two_product(pa.x, pc.y);
1741 let (cxay1, cxay0) = two_product(pc.x, pa.y);
1742 let (ac3, ac2, ac1, ac0) = two_two_diff(axcy1, axcy0, cxay1, cxay0);
1743 let ac = [ac0, ac1, ac2, ac3];
1744
1745 let (bxdy1, bxdy0) = two_product(pb.x, pd.y);
1746 let (dxby1, dxby0) = two_product(pd.x, pb.y);
1747 let (bd3, bd2, bd1, bd0) = two_two_diff(bxdy1, bxdy0, dxby1, dxby0);
1748 let bd = [bd0, bd1, bd2, bd3];
1749
1750 let (cxey1, cxey0) = two_product(pc.x, pe.y);
1751 let (excy1, excy0) = two_product(pe.x, pc.y);
1752 let (ce3, ce2, ce1, ce0) = two_two_diff(cxey1, cxey0, excy1, excy0);
1753 let ce = [ce0, ce1, ce2, ce3];
1754
1755 let (dxay1, dxay0) = two_product(pd.x, pa.y);
1756 let (axdy1, axdy0) = two_product(pa.x, pd.y);
1757 let (da3, da2, da1, da0) = two_two_diff(dxay1, dxay0, axdy1, axdy0);
1758 let da = [da0, da1, da2, da3];
1759
1760 let (exby1, exby0) = two_product(pe.x, pb.y);
1761 let (bxey1, bxey0) = two_product(pb.x, pe.y);
1762 let (eb3, eb2, eb1, eb0) = two_two_diff(exby1, exby0, bxey1, bxey0);
1763 let eb = [eb0, eb1, eb2, eb3];
1764
1765 let mut temp8a = [0f64; 8];
1766 let mut temp8b = [0f64; 8];
1767 let mut temp16 = [0f64; 16];
1768 let mut temp48a = [0f64; 48];
1769 let mut temp48b = [0f64; 48];
1770
1771 let mut abc = [0f64; 24];
1772 let mut bcd = [0f64; 24];
1773 let mut cde = [0f64; 24];
1774 let mut dea = [0f64; 24];
1775 let mut eab = [0f64; 24];
1776 let mut abd = [0f64; 24];
1777 let mut bce = [0f64; 24];
1778 let mut cda = [0f64; 24];
1779 let mut deb = [0f64; 24];
1780 let mut eac = [0f64; 24];
1781
1782 let mut abcd = [0f64; 96];
1783 let mut bcde = [0f64; 96];
1784 let mut cdea = [0f64; 96];
1785 let mut deab = [0f64; 96];
1786 let mut eabc = [0f64; 96];
1787
1788 let mut temp192 = [0f64; 192];
1789 let mut det384x = [0f64; 384];
1790 let mut det384y = [0f64; 384];
1791 let mut det384z = [0f64; 384];
1792
1793 let mut detxy = [0f64; 768];
1794
1795 let mut adet = [0f64; 1152];
1796 let mut bdet = [0f64; 1152];
1797 let mut cdet = [0f64; 1152];
1798 let mut ddet = [0f64; 1152];
1799 let mut edet = [0f64; 1152];
1800
1801 let mut abdet = [0f64; 2304];
1802 let mut cddet = [0f64; 2304];
1803 let mut cdedet = [0f64; 3456];
1804
1805 let mut deter = [0f64; 5760];
1806
1807 let mut temp8alen = scale_expansion_zeroelim(&bc[..4], pa.z, &mut temp8a);
1808 let mut temp8blen = scale_expansion_zeroelim(&ac[..4], -pb.z, &mut temp8b);
1809 let mut temp16len =
1810 fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp8b[..temp8blen], &mut temp16);
1811 temp8alen = scale_expansion_zeroelim(&ab[..4], pc.z, &mut temp8a);
1812 let abclen = fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp16[..temp16len], &mut abc);
1813
1814 temp8alen = scale_expansion_zeroelim(&cd[..4], pb.z, &mut temp8a);
1815 temp8blen = scale_expansion_zeroelim(&bd[..4], -pc.z, &mut temp8b);
1816 temp16len =
1817 fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp8b[..temp8blen], &mut temp16);
1818 temp8alen = scale_expansion_zeroelim(&bc[..4], pd.z, &mut temp8a);
1819 let bcdlen = fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp16[..temp16len], &mut bcd);
1820
1821 temp8alen = scale_expansion_zeroelim(&de[..4], pc.z, &mut temp8a);
1822 temp8blen = scale_expansion_zeroelim(&ce[..4], -pd.z, &mut temp8b);
1823 temp16len =
1824 fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp8b[..temp8blen], &mut temp16);
1825 temp8alen = scale_expansion_zeroelim(&cd[..4], pe.z, &mut temp8a);
1826 let cdelen = fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp16[..temp16len], &mut cde);
1827
1828 temp8alen = scale_expansion_zeroelim(&ea[..4], pd.z, &mut temp8a);
1829 temp8blen = scale_expansion_zeroelim(&da[..4], -pe.z, &mut temp8b);
1830 temp16len =
1831 fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp8b[..temp8blen], &mut temp16);
1832 temp8alen = scale_expansion_zeroelim(&de[..4], pa.z, &mut temp8a);
1833 let dealen = fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp16[..temp16len], &mut dea);
1834
1835 temp8alen = scale_expansion_zeroelim(&ab[..4], pe.z, &mut temp8a);
1836 temp8blen = scale_expansion_zeroelim(&eb[..4], -pa.z, &mut temp8b);
1837 temp16len =
1838 fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp8b[..temp8blen], &mut temp16);
1839 temp8alen = scale_expansion_zeroelim(&ea[..4], pb.z, &mut temp8a);
1840 let eablen = fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp16[..temp16len], &mut eab);
1841
1842 temp8alen = scale_expansion_zeroelim(&bd[..4], pa.z, &mut temp8a);
1843 temp8blen = scale_expansion_zeroelim(&da[..4], pb.z, &mut temp8b);
1844 temp16len =
1845 fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp8b[..temp8blen], &mut temp16);
1846 temp8alen = scale_expansion_zeroelim(&ab[..4], pd.z, &mut temp8a);
1847 let abdlen = fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp16[..temp16len], &mut abd);
1848
1849 temp8alen = scale_expansion_zeroelim(&ce[..4], pb.z, &mut temp8a);
1850 temp8blen = scale_expansion_zeroelim(&eb[..4], pc.z, &mut temp8b);
1851 temp16len =
1852 fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp8b[..temp8blen], &mut temp16);
1853 temp8alen = scale_expansion_zeroelim(&bc[..4], pe.z, &mut temp8a);
1854 let bcelen = fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp16[..temp16len], &mut bce);
1855
1856 temp8alen = scale_expansion_zeroelim(&da[..4], pc.z, &mut temp8a);
1857 temp8blen = scale_expansion_zeroelim(&ac[..4], pd.z, &mut temp8b);
1858 temp16len =
1859 fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp8b[..temp8blen], &mut temp16);
1860 temp8alen = scale_expansion_zeroelim(&cd[..4], pa.z, &mut temp8a);
1861 let cdalen = fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp16[..temp16len], &mut cda);
1862
1863 temp8alen = scale_expansion_zeroelim(&eb[..4], pd.z, &mut temp8a);
1864 temp8blen = scale_expansion_zeroelim(&bd[..4], pe.z, &mut temp8b);
1865 temp16len =
1866 fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp8b[..temp8blen], &mut temp16);
1867 temp8alen = scale_expansion_zeroelim(&de[..4], pb.z, &mut temp8a);
1868 let deblen = fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp16[..temp16len], &mut deb);
1869
1870 temp8alen = scale_expansion_zeroelim(&ac[..4], pe.z, &mut temp8a);
1871 temp8blen = scale_expansion_zeroelim(&ce[..4], pa.z, &mut temp8b);
1872 temp16len =
1873 fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp8b[..temp8blen], &mut temp16);
1874 let temp8alen = scale_expansion_zeroelim(&ea[..4], pc.z, &mut temp8a);
1875 let eaclen = fast_expansion_sum_zeroelim(&temp8a[..temp8alen], &temp16[..temp16len], &mut eac);
1876
1877 let mut temp48alen = fast_expansion_sum_zeroelim(&cde[..cdelen], &bce[..bcelen], &mut temp48a);
1878 let mut temp48blen = fast_expansion_sum_zeroelim(&deb[..deblen], &bcd[..bcdlen], &mut temp48b);
1879 for i in 0..temp48blen {
1880 temp48b[i] = -temp48b[i];
1881 }
1882
1883 let bcdelen =
1884 fast_expansion_sum_zeroelim(&temp48a[..temp48alen], &temp48b[..temp48blen], &mut bcde);
1885 let mut xlen = scale_expansion_zeroelim(&bcde[..bcdelen], pa.x, &mut temp192);
1886 xlen = scale_expansion_zeroelim(&temp192[..xlen], pa.x, &mut det384x);
1887 let mut ylen = scale_expansion_zeroelim(&bcde[..bcdelen], pa.y, &mut temp192);
1888 ylen = scale_expansion_zeroelim(&temp192[..ylen], pa.y, &mut det384y);
1889 let mut zlen = scale_expansion_zeroelim(&bcde[..bcdelen], pa.z, &mut temp192);
1890 zlen = scale_expansion_zeroelim(&temp192[..zlen], pa.z, &mut det384z);
1891 let mut xylen = fast_expansion_sum_zeroelim(&det384x[..xlen], &det384y[..ylen], &mut detxy);
1892 let alen = fast_expansion_sum_zeroelim(&detxy[..xylen], &det384z[..zlen], &mut adet);
1893
1894 temp48alen = fast_expansion_sum_zeroelim(&dea[..dealen], &cda[..cdalen], &mut temp48a);
1895 temp48blen = fast_expansion_sum_zeroelim(&eac[..eaclen], &cde[..cdelen], &mut temp48b);
1896 for i in 0..temp48blen {
1897 temp48b[i] = -temp48b[i];
1898 }
1899 let cdealen =
1900 fast_expansion_sum_zeroelim(&temp48a[..temp48alen], &temp48b[..temp48blen], &mut cdea);
1901 xlen = scale_expansion_zeroelim(&cdea[..cdealen], pb.x, &mut temp192);
1902 xlen = scale_expansion_zeroelim(&temp192[..xlen], pb.x, &mut det384x);
1903 ylen = scale_expansion_zeroelim(&cdea[..cdealen], pb.y, &mut temp192);
1904 ylen = scale_expansion_zeroelim(&temp192[..ylen], pb.y, &mut det384y);
1905 zlen = scale_expansion_zeroelim(&cdea[..cdealen], pb.z, &mut temp192);
1906 zlen = scale_expansion_zeroelim(&temp192[..zlen], pb.z, &mut det384z);
1907 xylen = fast_expansion_sum_zeroelim(&det384x[..xlen], &det384y[..ylen], &mut detxy);
1908 let blen = fast_expansion_sum_zeroelim(&detxy[..xylen], &det384z[..zlen], &mut bdet);
1909
1910 temp48alen = fast_expansion_sum_zeroelim(&eab[..eablen], &deb[..deblen], &mut temp48a);
1911 temp48blen = fast_expansion_sum_zeroelim(&abd[..abdlen], &dea[..dealen], &mut temp48b);
1912 for i in 0..temp48blen {
1913 temp48b[i] = -temp48b[i];
1914 }
1915 let deablen =
1916 fast_expansion_sum_zeroelim(&temp48a[..temp48alen], &temp48b[..temp48blen], &mut deab);
1917 xlen = scale_expansion_zeroelim(&deab[..deablen], pc.x, &mut temp192);
1918 xlen = scale_expansion_zeroelim(&temp192[..xlen], pc.x, &mut det384x);
1919 ylen = scale_expansion_zeroelim(&deab[..deablen], pc.y, &mut temp192);
1920 ylen = scale_expansion_zeroelim(&temp192[..ylen], pc.y, &mut det384y);
1921 zlen = scale_expansion_zeroelim(&deab[..deablen], pc.z, &mut temp192);
1922 zlen = scale_expansion_zeroelim(&temp192[..zlen], pc.z, &mut det384z);
1923 xylen = fast_expansion_sum_zeroelim(&det384x[..xlen], &det384y[..ylen], &mut detxy);
1924 let clen = fast_expansion_sum_zeroelim(&detxy[..xylen], &det384z[..zlen], &mut cdet);
1925
1926 temp48alen = fast_expansion_sum_zeroelim(&abc[..abclen], &eac[..eaclen], &mut temp48a);
1927 temp48blen = fast_expansion_sum_zeroelim(&bce[..bcelen], &eab[..eablen], &mut temp48b);
1928 for i in 0..temp48blen {
1929 temp48b[i] = -temp48b[i];
1930 }
1931 let eabclen =
1932 fast_expansion_sum_zeroelim(&temp48a[..temp48alen], &temp48b[..temp48blen], &mut eabc);
1933 xlen = scale_expansion_zeroelim(&eabc[..eabclen], pd.x, &mut temp192);
1934 xlen = scale_expansion_zeroelim(&temp192[..xlen], pd.x, &mut det384x);
1935 ylen = scale_expansion_zeroelim(&eabc[..eabclen], pd.y, &mut temp192);
1936 ylen = scale_expansion_zeroelim(&temp192[..ylen], pd.y, &mut det384y);
1937 zlen = scale_expansion_zeroelim(&eabc[..eabclen], pd.z, &mut temp192);
1938 zlen = scale_expansion_zeroelim(&temp192[..zlen], pd.z, &mut det384z);
1939 xylen = fast_expansion_sum_zeroelim(&det384x[..xlen], &det384y[..ylen], &mut detxy);
1940 let dlen = fast_expansion_sum_zeroelim(&detxy[..xylen], &det384z[..zlen], &mut ddet);
1941
1942 temp48alen = fast_expansion_sum_zeroelim(&bcd[..bcdlen], &abd[..abdlen], &mut temp48a);
1943 temp48blen = fast_expansion_sum_zeroelim(&cda[..cdalen], &abc[..abclen], &mut temp48b);
1944 for i in 0..temp48blen {
1945 temp48b[i] = -temp48b[i];
1946 }
1947 let abcdlen =
1948 fast_expansion_sum_zeroelim(&temp48a[..temp48alen], &temp48b[..temp48blen], &mut abcd);
1949 xlen = scale_expansion_zeroelim(&abcd[..abcdlen], pe.x, &mut temp192);
1950 xlen = scale_expansion_zeroelim(&temp192[..xlen], pe.x, &mut det384x);
1951 ylen = scale_expansion_zeroelim(&abcd[..abcdlen], pe.y, &mut temp192);
1952 ylen = scale_expansion_zeroelim(&temp192[..ylen], pe.y, &mut det384y);
1953 zlen = scale_expansion_zeroelim(&abcd[..abcdlen], pe.z, &mut temp192);
1954 zlen = scale_expansion_zeroelim(&temp192[..zlen], pe.z, &mut det384z);
1955 xylen = fast_expansion_sum_zeroelim(&det384x[..xlen], &det384y[..ylen], &mut detxy);
1956 let elen = fast_expansion_sum_zeroelim(&detxy[..xylen], &det384z[..zlen], &mut edet);
1957
1958 let ablen = fast_expansion_sum_zeroelim(&adet[..alen], &bdet[..blen], &mut abdet);
1959 let cdlen = fast_expansion_sum_zeroelim(&cdet[..clen], &ddet[..dlen], &mut cddet);
1960 let cdelen = fast_expansion_sum_zeroelim(&cddet[..cdlen], &edet[..elen], &mut cdedet);
1961 let deterlen = fast_expansion_sum_zeroelim(&abdet[..ablen], &cdedet[..cdelen], &mut deter);
1962
1963 deter[deterlen - 1]
1964}
1965
1966fn scale_expansion_zeroelim(e: &[f64], b: f64, h: &mut [f64]) -> usize {
1967 let (bhi, blo) = split(b);
1968 let (mut Q, hh) = two_product_presplit(e[0], b, bhi, blo);
1969 let mut hindex = 0;
1970 if hh != 0.0 {
1971 h[hindex] = hh;
1972 hindex += 1;
1973 }
1974 for eindex in 1..e.len() {
1975 let enow = e[eindex];
1976 let (product1, product0) = two_product_presplit(enow, b, bhi, blo);
1977 let (sum, hh) = two_sum(Q, product0);
1978 if hh != 0.0 {
1979 h[hindex] = hh;
1980 hindex += 1;
1981 }
1982 let (new_q, hh) = fast_two_sum(product1, sum);
1983 Q = new_q;
1984 if hh != 0.0 {
1985 h[hindex] = hh;
1986 hindex += 1;
1987 }
1988 }
1989 if Q != 0.0 || hindex == 0 {
1990 h[hindex] = Q;
1991 hindex += 1;
1992 }
1993 hindex
1994}
1995
1996#[inline]
1997fn two_product(a: f64, b: f64) -> (f64, f64) {
1998 let x = a * b;
1999 (x, two_product_tail(a, b, x))
2000}
2001
2002#[inline]
2003fn two_product_tail(a: f64, b: f64, x: f64) -> f64 {
2004 let (ahi, alo) = split(a);
2005 let (bhi, blo) = split(b);
2006 let err1 = x - (ahi * bhi);
2007 let err2 = err1 - (alo * bhi);
2008 let err3 = err2 - (ahi * blo);
2009 (alo * blo) - err3
2010}
2011
2012#[inline]
2013fn split(a: f64) -> (f64, f64) {
2014 let c = SPLITTER * a;
2015 let abig = c - a;
2016 let ahi = c - abig;
2017 let alo = a - ahi;
2018 (ahi, alo)
2019}
2020
2021#[inline]
2022fn two_product_presplit(a: f64, b: f64, bhi: f64, blo: f64) -> (f64, f64) {
2023 let x = a * b;
2024 let (ahi, alo) = split(a);
2025 let err1 = x - ahi * bhi;
2026 let err2 = err1 - alo * bhi;
2027 let err3 = err2 - ahi * blo;
2028 let y = alo * blo - err3;
2029 (x, y)
2030}
2031
2032#[inline]
2033fn two_two_diff(a1: f64, a0: f64, b1: f64, b0: f64) -> (f64, f64, f64, f64) {
2034 let (j, _r0, x0) = two_one_diff(a1, a0, b0);
2035 let (x3, x2, x1) = two_one_diff(j, _r0, b1);
2036 (x3, x2, x1, x0)
2037}
2038
2039#[inline]
2040fn two_one_diff(a1: f64, a0: f64, b: f64) -> (f64, f64, f64) {
2041 let (i, x0) = two_diff(a0, b);
2042 let (x2, x1) = two_sum(a1, i);
2043 (x2, x1, x0)
2044}
2045
2046#[inline]
2047fn two_diff(a: f64, b: f64) -> (f64, f64) {
2048 let x = a - b;
2049 (x, two_diff_tail(a, b, x))
2050}
2051
2052#[inline]
2053fn two_diff_tail(a: f64, b: f64, x: f64) -> f64 {
2054 let bvirt = a - x;
2055 let avirt = x + bvirt;
2056 let bround = bvirt - b;
2057 let around = a - avirt;
2058 around + bround
2059}
2060
2061#[inline]
2062fn two_sum(a: f64, b: f64) -> (f64, f64) {
2063 let x = a + b;
2064 (x, two_sum_tail(a, b, x))
2065}
2066
2067#[inline]
2068fn two_sum_tail(a: f64, b: f64, x: f64) -> f64 {
2069 let bvirt = x - a;
2070 let avirt = x - bvirt;
2071 let bround = b - bvirt;
2072 let around = a - avirt;
2073 around + bround
2074}
2075
2076fn estimate(e: &[f64]) -> f64 {
2077 let mut q = e[0];
2078 for cur in &e[1..] {
2079 q += *cur;
2080 }
2081 q
2082}
2083
2084fn fast_expansion_sum_zeroelim(e: &[f64], f: &[f64], h: &mut [f64]) -> usize {
2085 let mut enow = e[0];
2086 let mut fnow = f[0];
2087 let mut eindex = 0;
2088 let mut findex = 0;
2089 let mut Qnew;
2090 let mut hh;
2091 let mut Q;
2092 if (fnow > enow) == (fnow > -enow) {
2093 Q = enow;
2094 eindex += 1;
2095 } else {
2096 Q = fnow;
2097 findex += 1;
2098 }
2099
2100 let mut hindex = 0;
2101 if eindex < e.len() && findex < f.len() {
2102 enow = e[eindex];
2103 fnow = f[findex];
2104 if (fnow > enow) == (fnow > -enow) {
2105 let r = fast_two_sum(enow, Q);
2106 Qnew = r.0;
2107 hh = r.1;
2108 eindex += 1;
2109 } else {
2110 let r = fast_two_sum(fnow, Q);
2111 Qnew = r.0;
2112 hh = r.1;
2113 findex += 1;
2114 }
2115 Q = Qnew;
2116 if hh != 0.0 {
2117 h[hindex] = hh;
2118 hindex += 1;
2119 }
2120
2121 while eindex < e.len() && findex < f.len() {
2122 enow = e[eindex];
2123 fnow = f[findex];
2124 if (fnow > enow) == (fnow > -enow) {
2125 let r = two_sum(Q, enow);
2126 Qnew = r.0;
2127 hh = r.1;
2128 eindex += 1;
2129 } else {
2130 let r = two_sum(Q, fnow);
2131 Qnew = r.0;
2132 hh = r.1;
2133 findex += 1;
2134 };
2135 Q = Qnew;
2136 if hh != 0.0 {
2137 h[hindex] = hh;
2138 hindex += 1;
2139 }
2140 }
2141 }
2142
2143 while eindex < e.len() {
2144 enow = e[eindex];
2145 let r = two_sum(Q, enow);
2146 Qnew = r.0;
2147 hh = r.1;
2148 Q = Qnew;
2149 eindex += 1;
2150 if hh != 0.0 {
2151 h[hindex] = hh;
2152 hindex += 1
2153 }
2154 }
2155
2156 while findex < f.len() {
2157 fnow = f[findex];
2158 let r = two_sum(Q, fnow);
2159 Qnew = r.0;
2160 hh = r.1;
2161 Q = Qnew;
2162 findex += 1;
2163 if hh != 0.0 {
2164 h[hindex] = hh;
2165 hindex += 1
2166 }
2167 }
2168
2169 if Q != 0.0 || hindex == 0 {
2170 h[hindex] = Q;
2171 hindex += 1;
2172 }
2173 hindex
2174}
2175
2176#[inline]
2177fn fast_two_sum_tail(a: f64, b: f64, x: f64) -> f64 {
2178 let bvirt = x - a;
2179 b - bvirt
2180}
2181
2182#[inline]
2183fn fast_two_sum(a: f64, b: f64) -> (f64, f64) {
2184 let x = a + b;
2185 (x, fast_two_sum_tail(a, b, x))
2186}
2187
2188#[inline]
2189fn square_tail(a: f64, x: f64) -> f64 {
2190 let (ahi, alo) = split(a);
2191 let err1 = x - ahi * ahi;
2192 let err3 = err1 - (ahi + ahi) * alo;
2193 alo * alo - err3
2194}
2195
2196#[inline]
2197fn square(a: f64) -> (f64, f64) {
2198 let x = a * a;
2199 (x, square_tail(a, x))
2200}
2201
2202#[inline]
2203fn two_one_sum(a1: f64, a0: f64, b: f64) -> (f64, f64, f64) {
2204 let (_i, x0) = two_sum(a0, b);
2205 let (x2, x1) = two_sum(a1, _i);
2206 (x2, x1, x0)
2207}
2208
2209#[inline]
2210fn two_two_sum(a1: f64, a0: f64, b1: f64, b0: f64) -> (f64, f64, f64, f64) {
2211 let (_j, _r0, x0) = two_one_sum(a1, a0, b0);
2212 let (x3, x2, x1) = two_one_sum(_j, _r0, b1);
2213 (x3, x2, x1, x0)
2214}
2215
2216#[inline]
2217fn two_one_product(a1: f64, a0: f64, b: f64) -> (f64, f64, f64, f64) {
2218 let (bhi, blo) = split(b);
2219 let (mut _i, x0) = two_product_presplit(a0, b, bhi, blo);
2220 let (mut _j, _0) = two_product_presplit(a1, b, bhi, blo);
2221 let (_k, x1) = two_sum(_i, _0);
2222 let (x3, x2) = fast_two_sum(_j, _k);
2223 (x3, x2, x1, x0)
2224}