===================================================
Errors Moderate Dataset
gtr eq current b1div b2div new
====== ======== ======== ======== ========
1 bit 0.24707% 0.92986% 0.24707% 0.24707%
2 bits 0.01762% 0.01770% 0.01762% 0.01762%
8 bits 0.00026% 0.00026% 0.00026% 0.00026%
16 bits 0.00000% 0.00000% 0.00000% 0.00000%
24 bits 0% 0% 0% 0%
52 bits 0% 0% 0% 0%
===================================================
Table 1: Errors with Moderate Dataset (Double Precision)
Note in Table 1 that both the old and new methods give identical error
rates for data with moderate exponents. Errors exceeding 16 bits are
exceedingly rare. There are substantial increases in the 1 bit error
rates for b1div (the 1 divide/2 multiplys method) as compared to b2div
(the 2 divides method). These differences are minimal for 2 bits and
larger error measurements.
===================================================
Errors Full Dataset
gtr eq current b1div b2div new
====== ======== ======== ======== ========
1 bit 2.05% 1.23842% 0.67130% 0.16664%
2 bits 1.88% 0.51615% 0.50354% 0.00900%
8 bits 1.77% 0.42856% 0.42168% 0.00011%
16 bits 1.63% 0.33840% 0.32879% 0.00001%
24 bits 1.51% 0.25583% 0.24405% 0.00000%
52 bits 1.13% 0.01886% 0.00350% 0.00000%
===================================================
Table 2: Errors with Full Dataset (Double Precision)
Table 2 shows significant differences in error rates. First, the
difference between b1div and b2div show a significantly higher error
rate for the b1div method both for single bit errros and well
beyond. Even for 52 bits, we see the b1div method gets completely
wrong answers more than 5 times as often as b2div. To retain
comparable accuracy with current complex divide results for small
exponents and due to the increase in errors for large exponents, I
choose to use the more accurate method of two divides.
The current method has more 1.6% of cases where it is getting results
where the low 24 bits of the mantissa differ from the correct
answer. More than 1.1% of cases where the answer is completely wrong.
The new method shows less than one case in 10,000 with greater than
two bits of error and only one case in 10 million with greater than
16 bits of errors. The new patch reduces 8 bit errors by
a factor of 16,000 and virtually eliminates completely wrong
answers.
As noted above, for architectures with double precision
hardware, the new method uses that hardware for the
intermediate calculations before returning the
result in float precision. Testing of the new patch
has shown zero errors found as seen in Tables 3 and 4.
Correctness for float
=============================
Errors Moderate Dataset
gtr eq current new
====== ======== ========
1 bit 28.68070% 0%
2 bits 0.64386% 0%
8 bits 0.00401% 0%
16 bits 0.00001% 0%
24 bits 0% 0%
=============================
Table 3: Errors with Moderate Dataset (float)
=============================
Errors Full Dataset
gtr eq current new
====== ======== ========
1 bit 19.98% 0%
2 bits 3.20% 0%
8 bits 1.97% 0%
16 bits 1.08% 0%
24 bits 0.55% 0%
=============================
Table 4: Errors with Full Dataset (float)
As before, the current method shows an troubling rate of extreme
errors.
There very mincr changes in accuracy for half-precision since the code
changes from Smith's method to the simple method. 5 out of 1 million
test cases show correct answers instead of 1 or 2 bit errors.
libgcc computes half-precision functions in float precision
allowing the existing methods to avoid overflow/underflow issues
for the allowed range of exponents for half-precision.
Extended precision (using x87 80-bit format on x86) and Long double
(using IEEE-754 128-bit on x86 and aarch64) both have 15-bit exponents
as compared to 11-bit exponents in double precision. We note that the
C standard also allows Long Double to be implemented in the equivalent
range of Double. The RMIN2 and RMINSCAL constants are selected to work
within the Double range as well as with extended and 128-bit ranges.
We will limit our performance and accurancy discussions to the 80-bit
and 128-bit formats as seen on x86 here.
The extended and long double precision investigations were more
limited. Aarch64 does not support extended precision but does support
the software implementation of 128-bit long double precision. For x86,
long double defaults to the 80-bit precision but using the
-mlong-double-128 flag switches to using the software implementation
of 128-bit precision. Both 80-bit and 128-bit precisions have the same
exponent range, with the 128-bit precision has extended mantissas.
Since this change is only aimed at avoiding underflow/overflow for
extreme exponents, I studied the extended precision results on x86 for
100,000 values. The limited exponent dataset showed no differences.
For the dataset with full exponent range, the current and new values
showed major differences (greater than 32 bits) in 567 cases out of
100,000 (0.56%). In every one of these cases, the ratio of c/d or d/c
(as appropriate) was zero or subnormal, indicating the advantage of
the new method and its continued correctness where needed.
PERFORMANCE Test results
In order for a library change to be practical, it is necessary to show
the slowdown is tolerable. The slowdowns observed are much less than
would be seen by (for example) switching from hardware double precison
to a software quad precision, which on the tested machines causes a
slowdown of around 100x).
The actual slowdown depends on the machine architecture. It also
depends on the nature of the input data. If underflow/overflow is
rare, then implementations that have strong branch prediction will
only slowdown by a few cycles. If underflow/overflow is common, then
the branch predictors will be less accurate and the cost will be
higher.
Results from two machines are presented as examples of the overhead
for the new method. The one labeled x86 is a 5 year old Intel x86
processor and the one labeled aarch64 is a 3 year old arm64 processor.
In the following chart, the times are averaged over a one million
value data set. All values are scaled to set the time of the current
method to be 1.0. Lower values are better. A value of less than 1.0
would be faster than the current method and a value greater than 1.0
would be slower than the current method.
================================================
Moderate set full set
x86 aarch64 x86 aarch64
======== =============== ===============
float 0.59 0.79 0.45 0.81
double 1.04 1.24 1.38 1.56
long double 1.13 1.24 1.29 1.25
================================================
Table 5: Performance Comparisons (ratio new/current)
The above tables omit the timing for the 1 divide and 2 multiply
comparison with the 2 divide approach.
The float results show clear performance improvement due to using the
simple method with double precision for intermediate calculations.
The double results with the newer method show less overhead for the
moderate dataset than for the full dataset. That's because the moderate
dataset does not ever take the new branches which protect from
under/overflow. The better the branch predictor, the lower the cost
for these untaken branches. Both platforms are somewhat dated, with
the x86 having a better branch predictor which reduces the cost of the
additional branches in the new code. Of course, the relative slowdown
may be greater for some architectures, especially those with limited
branch prediction combined with a high cost of misprediction.
The long double results are fairly consistent in showing the moderate
additional cost of the extra branches and calculations for all cases.
The observed cost for all precisions is claimed to be tolerable on the
grounds that:
(a) the cost is worthwhile considering the accuracy improvement shown.
(b) most applications will only spend a small fraction of their time
calculating complex divide.
(c) it is much less than the cost of extended precision
(d) users are not forced to use it (as described below)
Those users who find this degree of slowdown unsatisfactory may use
the gcc switch -fcx-fortran-rules which does not use the library
routine, instead inlining Smith's method without the C99 requirement
for dealing with NaN results. The proposed patch for libgcc complex
divide does not affect the code generated by -fcx-fortran-rules.
SUMMARY
When input data to complex divide has exponents whose absolute value
is less than half of *_MAX_EXP, this patch makes no changes in
accuracy and has only a modest effect on performance. When input data
contains values outside those ranges, the patch eliminates more than
99.9% of major errors with a tolerable cost in performance.
In comparison to Elen Kalda's method, this patch introduces more
performance overhead but reduces major errors by a factor of
greater than 4000.
REFERENCES
[1] Nelson H.F. Beebe, "The Mathematical-Function Computation Handbook.
Springer International Publishing AG, 2017.
[2] Robert L. Smith. Algorithm 116: Complex division. Commun. ACM,
5(8):435, 1962.
[3] Michael Baudin and Robert L. Smith. "A robust complex division in
Scilab," October 2012, available at http://arxiv.org/abs/1210.4539.
[4] Elen Kalda: Complex division improvements in libgcc
https://gcc.gnu.org/legacy-ml/gcc-patches/2019-08/msg01629.html
2020-12-08 Patrick McGehearty <patrick.mcgehearty@oracle.com>
gcc/c-family/
* c-cppbuiltin.c (c_cpp_builtins): Add supporting macros for new
complex divide
libgcc/
* libgcc2.c (XMTYPE, XCTYPE, RBIG, RMIN, RMIN2, RMINSCAL, RMAX2):
Define.
(__divsc3, __divdc3, __divxc3, __divtc3): Improve complex divide.
* config/rs6000/_divkc3.c (RBIG, RMIN, RMIN2, RMINSCAL, RMAX2):
Define.
(__divkc3): Improve complex divide.
gcc/testsuite/
* gcc.c-torture/execute/ieee/cdivchkd.c - New test.
* gcc.c-torture/execute/ieee/cdivchkf.c - New test.
* gcc.c-torture/execute/ieee/cdivchkld.c - New test.
---
gcc/c-family/c-cppbuiltin.c | 58 ++++++--
.../gcc.c-torture/execute/ieee/cdivchkd.c | 123 +++++++++++++++++
.../gcc.c-torture/execute/ieee/cdivchkf.c | 122 +++++++++++++++++
.../gcc.c-torture/execute/ieee/cdivchkld.c | 123 +++++++++++++++++
libgcc/config/rs6000/_divkc3.c | 105 +++++++++++++--
libgcc/libgcc2.c | 150 +++++++++++++++++++--
6 files changed, 650 insertions(+), 31 deletions(-)
create mode 100644 gcc/testsuite/gcc.c-torture/execute/ieee/cdivchkd.c
create mode 100644 gcc/testsuite/gcc.c-torture/execute/ieee/cdivchkf.c
create mode 100644 gcc/testsuite/gcc.c-torture/execute/ieee/cdivchkld.c
@@ -1026,10 +1026,10 @@ c_cpp_builtins (cpp_reader *pfile)
cpp_define (pfile, "__cpp_using_enum=201907L");
}
if (cxx_dialect > cxx20)
- {
- /* Set feature test macros for C++23. */
- cpp_define (pfile, "__cpp_size_t_suffix=202011L");
- }
+ {
+ /* Set feature test macros for C++23. */
+ cpp_define (pfile, "__cpp_size_t_suffix=202011L");
+ }
if (flag_concepts)
{
if (cxx_dialect >= cxx20)
@@ -1277,8 +1277,10 @@ c_cpp_builtins (cpp_reader *pfile)
{
scalar_float_mode mode = mode_iter.require ();
const char *name = GET_MODE_NAME (mode);
+ const int name_len = strlen (name);
+ char float_h_prefix[16] = "";
char *macro_name
- = (char *) alloca (strlen (name)
+ = (char *) alloca (name_len
+ sizeof ("__LIBGCC__MANT_DIG__"));
sprintf (macro_name, "__LIBGCC_%s_MANT_DIG__", name);
builtin_define_with_int_value (macro_name,
@@ -1286,20 +1288,29 @@ c_cpp_builtins (cpp_reader *pfile)
if (!targetm.scalar_mode_supported_p (mode)
|| !targetm.libgcc_floating_mode_supported_p (mode))
continue;
- macro_name = (char *) alloca (strlen (name)
+ macro_name = (char *) alloca (name_len
+ sizeof ("__LIBGCC_HAS__MODE__"));
sprintf (macro_name, "__LIBGCC_HAS_%s_MODE__", name);
cpp_define (pfile, macro_name);
- macro_name = (char *) alloca (strlen (name)
+ macro_name = (char *) alloca (name_len
+ sizeof ("__LIBGCC__FUNC_EXT__"));
sprintf (macro_name, "__LIBGCC_%s_FUNC_EXT__", name);
char suffix[20] = "";
if (mode == TYPE_MODE (double_type_node))
- ; /* Empty suffix correct. */
+ {
+ ; /* Empty suffix correct. */
+ memcpy (float_h_prefix, "DBL", 4);
+ }
else if (mode == TYPE_MODE (float_type_node))
- suffix[0] = 'f';
+ {
+ suffix[0] = 'f';
+ memcpy (float_h_prefix, "FLT", 4);
+ }
else if (mode == TYPE_MODE (long_double_type_node))
- suffix[0] = 'l';
+ {
+ suffix[0] = 'l';
+ memcpy (float_h_prefix, "LDBL", 4);
+ }
else
{
bool found_suffix = false;
@@ -1310,6 +1321,8 @@ c_cpp_builtins (cpp_reader *pfile)
sprintf (suffix, "f%d%s", floatn_nx_types[i].n,
floatn_nx_types[i].extended ? "x" : "");
found_suffix = true;
+ sprintf (float_h_prefix, "FLT%d%s", floatn_nx_types[i].n,
+ floatn_nx_types[i].extended ? "X" : "");
break;
}
gcc_assert (found_suffix);
@@ -1347,11 +1360,34 @@ c_cpp_builtins (cpp_reader *pfile)
default:
gcc_unreachable ();
}
- macro_name = (char *) alloca (strlen (name)
+ macro_name = (char *) alloca (name_len
+ sizeof ("__LIBGCC__EXCESS_"
"PRECISION__"));
sprintf (macro_name, "__LIBGCC_%s_EXCESS_PRECISION__", name);
builtin_define_with_int_value (macro_name, excess_precision);
+
+ char val_name[64];
+
+ macro_name = (char *) alloca (name_len
+ + sizeof ("__LIBGCC_EPSILON__"));
+ sprintf (macro_name, "__LIBGCC_%s_EPSILON__", name);
+ sprintf (val_name, "__%s_EPSILON__", float_h_prefix);
+ builtin_define_with_value (macro_name, val_name, 0);
+
+ macro_name = (char *) alloca (name_len + sizeof ("__LIBGCC_MAX__"));
+ sprintf (macro_name, "__LIBGCC_%s_MAX__", name);
+ sprintf (val_name, "__%s_MAX__", float_h_prefix);
+ builtin_define_with_value (macro_name, val_name, 0);
+
+ macro_name = (char *) alloca (name_len + sizeof ("__LIBGCC_MIN__"));
+ sprintf (macro_name, "__LIBGCC_%s_MIN__", name);
+ sprintf (val_name, "__%s_MIN__", float_h_prefix);
+ builtin_define_with_value (macro_name, val_name, 0);
+
+#ifdef HAVE_adddf3
+ builtin_define_with_int_value ("__LIBGCC_HAVE_HWDBL__",
+ HAVE_adddf3);
+#endif
}
/* For libgcc crtstuff.c and libgcc2.c. */
new file mode 100644
@@ -0,0 +1,123 @@
+/*
+ Program to test complex divide for correct results on selected values.
+ Checking known failure points.
+*/
+
+extern void abort();
+extern void exit();
+
+extern double log2();
+
+#define SMALL (0x1.0p-1024)
+#define MAXBIT 52
+#define ERRLIM 6
+
+/*
+ Compare c (computed value) with z (expected value).
+ Return 0 if within allowed range. Return 1 if not.
+*/
+int match(_Complex c, _Complex z)
+{
+ double rz, iz, rc, ic;
+ double rerr, ierr, rmax;
+ double biterr;
+ rz = __real__ z;
+ iz = __imag__ z;
+ rc = __real__ c;
+ ic = __imag__ c;
+
+ if (__builtin_fabs(rz) > SMALL)
+ {
+ rerr = __builtin_fabs(rz-rc)/__builtin_fabs(rz);
+ }
+ else if (__builtin_fabs(rz) == 0.0)
+ {
+ rerr = __builtin_fabs(rc);
+ }
+ else
+ {
+ rerr = __builtin_fabs(rz-rc)/SMALL;
+ }
+
+ if (__builtin_fabs(iz) > SMALL)
+ {
+ ierr = __builtin_fabs(iz-ic)/__builtin_fabs(iz);
+ }
+ else if (__builtin_fabs(iz) == 0.0)
+ {
+ ierr = __builtin_fabs(ic);
+ }
+ else
+ {
+ ierr = __builtin_fabs(iz-ic)/SMALL;
+ }
+ rmax = __builtin_fmax(rerr, ierr);
+ biterr = 0;
+ if ( rmax != 0.0)
+ {
+ biterr = log2(rmax) + MAXBIT + 1;
+ }
+
+ if (biterr >= ERRLIM)
+ return 0;
+ else
+ return 1;
+}
+
+
+int main(int argc, char** argv)
+{
+ double _Complex a,b,c,z;
+ double xr[4], xi[4], yr[4], yi[4], zr[4], zi[4];
+ double cr, ci;
+ int i;
+ int ok = 1;
+ xr[0] = -0x1.16e7fad79e45ep+651;
+ xi[0] = -0x1.f7f75b94c6c6ap-860;
+ yr[0] = -0x1.2f40d8ff7e55ep+245;
+ yi[0] = -0x0.0000000004ebcp-1022;
+ zr[0] = 0x1.d6e4b0e282869p+405;
+ zi[0] = -0x1.e9095e311e706p-900;
+
+ xr[1] = -0x1.21ff587f953d3p-310;
+ xi[1] = -0x1.5a526dcc59960p+837;
+ yr[1] = 0x1.b88b8b552eaadp+735;
+ yi[1] = -0x1.873e2d6544d92p-327;
+ zr[1] = 0x1.65734a88b2de0p-961;
+ zi[1] = -0x1.927e85b8b5770p+101;
+
+ xr[2] = 0x1.4612e41aa8080p-846;
+ xi[2] = -0x0.0000000613e07p-1022;
+ yr[2] = 0x1.df9cd0d58caafp-820;
+ yi[2] = -0x1.e47051a9036dbp-584;
+ zr[2] = 0x1.9b194f3fffa32p-469;
+ zi[2] = 0x1.58a00ab740a6bp-263;
+
+ xr[3] = 0x1.cb27eece7c585p-355;
+ xi[3] = 0x0.000000223b8a8p-1022;
+ yr[3] = -0x1.74e7ed2b9189fp-22;
+ yi[3] = 0x1.3d80439e9a119p-731;
+ zr[3] = -0x1.3b35ed806ae5ap-333;
+ zi[3] = -0x0.05e01bcbfd9f6p-1022;
+
+
+ for (i = 0; i < 4; i++)
+ {
+ __real__ a = xr[i];
+ __imag__ a = xi[i];
+ __real__ b = yr[i];
+ __imag__ b = yi[i];
+ __real__ z = zr[i];
+ __imag__ z = zi[i];
+ c = a/b;
+ cr = __real__ c;
+ ci = __imag__ c;
+
+ if (!match(c,z)){
+ ok = 0;
+ }
+ }
+ if (!ok)
+ abort();
+ exit(0);
+}
new file mode 100644
@@ -0,0 +1,122 @@
+/*
+ Program to test complex divide for correct results on selected values.
+ Checking known failure points.
+*/
+
+extern void abort();
+extern void exit();
+
+extern float log2f(float);
+
+#define SMALL (0x1.0p-126)
+#define MAXBIT 24
+#define ERRLIM 6
+
+/*
+ Compare c (computed value) with z (expected value).
+ Return 0 if within allowed range. Return 1 if not.
+*/
+int match(float _Complex c, float _Complex z)
+{
+ float rz, iz, rc, ic;
+ float rerr, ierr, rmax;
+ float biterr;
+ rz = __real__ z;
+ iz = __imag__ z;
+ rc = __real__ c;
+ ic = __imag__ c;
+
+ if (__builtin_fabsf(rz) > SMALL)
+ {
+ rerr = __builtin_fabsf(rz-rc)/__builtin_fabsf(rz);
+ }
+ else if (__builtin_fabsf(rz) == 0.0)
+ {
+ rerr = __builtin_fabsf(rc);
+ }
+ else
+ {
+ rerr = __builtin_fabsf(rz-rc)/SMALL;
+ }
+
+ if (__builtin_fabsf(iz) > SMALL)
+ {
+ ierr = __builtin_fabsf(iz-ic)/__builtin_fabsf(iz);
+ }
+ else if (__builtin_fabsf(iz) == 0.0)
+ {
+ ierr = __builtin_fabsf(ic);
+ }
+ else
+ {
+ ierr = __builtin_fabsf(iz-ic)/SMALL;
+ }
+ rmax = __builtin_fmaxf(rerr, ierr);
+ biterr = 0;
+ if ( rmax != 0.0)
+ {
+ biterr = log2f(rmax) + MAXBIT + 1;
+ }
+
+ if (biterr >= ERRLIM)
+ return 0;
+ else
+ return 1;
+}
+
+
+int main(int argc, char** argv)
+{
+ float _Complex a,b,c,z;
+ float xr[4], xi[4], yr[4], yi[4], zr[4], zi[4];
+ float cr, ci;
+ int i;
+ int ok = 1;
+ xr[0] = 0x1.0b1600p-133;
+ xi[0] = 0x1.5e1c28p+54;
+ yr[0] = -0x1.cdec8cp-119;
+ yi[0] = 0x1.1e72ccp+32;
+ zr[0] = 0x1.38e502p+22;
+ zi[0] = -0x1.f89220p-129;
+
+ xr[1] = -0x1.b1bee2p+121;
+ xi[1] = -0x1.cb403ep-59;
+ yr[1] = 0x1.480000p-144;
+ yi[1] = -0x1.c66fc4p+5;
+ zr[1] = -0x1.60b8cap-34;
+ zi[1] = -0x1.e8b02ap+115;
+
+ xr[2] = -0x1.3f6e00p-97;
+ xi[2] = -0x1.c00000p-146;
+ yr[2] = 0x1.000000p-148;
+ yi[2] = -0x1.0c4e70p-91;
+ zr[2] = 0x1.aa50d0p-55;
+ zi[2] = -0x1.30c746p-6;
+
+ xr[3] = 0x1.000000p-148;
+ xi[3] = 0x1.f4bc04p-84;
+ yr[3] = 0x1.00ad74p-20;
+ yi[3] = 0x1.2ad02ep-85;
+ zr[3] = 0x1.1102ccp-127;
+ zi[3] = 0x1.f369a4p-64;
+
+ for (i = 0; i < 4; i++)
+ {
+ __real__ a = xr[i];
+ __imag__ a = xi[i];
+ __real__ b = yr[i];
+ __imag__ b = yi[i];
+ __real__ z = zr[i];
+ __imag__ z = zi[i];
+ c = a/b;
+ cr = __real__ c;
+ ci = __imag__ c;
+
+ if (!match(c,z)){
+ ok = 0;
+ }
+ }
+ if (!ok)
+ abort();
+ exit(0);
+}
new file mode 100644
@@ -0,0 +1,123 @@
+/*
+ Program to test complex divide for correct results on selected values.
+ Checking known failure points.
+*/
+
+extern void abort();
+extern void exit();
+
+extern long double log2l(long double);
+
+#define SMALL (0x1.0p-16384L)
+#define MAXBIT 63
+#define ERRLIM 6
+
+/*
+ Compare c (computed value) with z (expected value).
+ Return 0 if within allowed range. Return 1 if not.
+*/
+int match(_Complex c, _Complex z)
+{
+ long double rz, iz, rc, ic;
+ long double rerr, ierr, rmax;
+ long double biterr;
+ rz = __real__ z;
+ iz = __imag__ z;
+ rc = __real__ c;
+ ic = __imag__ c;
+
+ if (__builtin_fabsl(rz) > SMALL)
+ {
+ rerr = __builtin_fabsl(rz-rc)/__builtin_fabsl(rz);
+ }
+ else if (__builtin_fabsl(rz) == 0.0)
+ {
+ rerr = __builtin_fabsl(rc);
+ }
+ else
+ {
+ rerr = __builtin_fabsl(rz-rc)/SMALL;
+ }
+
+ if (__builtin_fabsl(iz) > SMALL)
+ {
+ ierr = __builtin_fabsl(iz-ic)/__builtin_fabsl(iz);
+ }
+ else if (__builtin_fabsl(iz) == 0.0)
+ {
+ ierr = __builtin_fabsl(ic);
+ }
+ else
+ {
+ ierr = __builtin_fabsl(iz-ic)/SMALL;
+ }
+ rmax = __builtin_fmaxl(rerr, ierr);
+ biterr = 0;
+ if ( rmax != 0.0)
+ {
+ biterr = log2l(rmax) + MAXBIT + 1;
+ }
+
+ if (biterr >= ERRLIM)
+ return 0;
+ else
+ return 1;
+}
+
+
+int main(int argc, char** argv)
+{
+ long double _Complex a,b,c,z;
+ long double xr[4], xi[4], yr[4], yi[4], zr[4], zi[4];
+ long double cr, ci;
+ int i;
+ int ok = 1;
+ xr[0] = -0x9.c793985b7d029d90p-8480L;
+ xi[0] = 0x8.018745ffa61a8fe0p+16329L;
+ yr[0] = -0xe.d5bee9c523a35ad0p-15599L;
+ yi[0] = -0xa.8c93c5a4f94128f0p+869L;
+ zr[0] = -0xc.248bc228e01adcb0p+15456L;
+ zi[0] = -0x8.89baf6960fac5cf0p-1011L;
+
+ xr[1] = 0xb.68e44bc6d0b91a30p+16026L;
+ xi[1] = 0xb.ab10f5453e972f30p-14239L;
+ yr[1] = 0x8.8cbd470705428ff0p-16350L;
+ yi[1] = -0xa.0c1cbeae4e4b69f0p+347L;
+ zr[1] = 0xf.76204243c2f28070p-1022L;
+ zi[1] = 0x9.15b5abcc93d1f920p+15676L;
+
+ xr[2] = -0x9.e8c093a43b546a90p+15983L;
+ xi[2] = 0xc.95b18274208311e0p-2840L;
+ yr[2] = -0x8.dedb729b5c1b2ec0p+8L;
+ yi[2] = 0xa.a49fb81b24738370p-16385L;
+ zr[2] = 0x8.efccf744dd88c7a0p+15972L;
+ zi[2] = 0xa.b8f3f7c826b5b780p-421L;
+
+ xr[3] = 0xc.4687f251c0f48bd0p-3940L;
+ xi[3] = -0xe.a3f2138992d85fa0p+15598L;
+ yr[3] = 0xe.4b0c25c3d5ebb830p-16344L;
+ yi[3] = -0xa.6cbf1ba80f7b97a0p+78L;
+ zr[3] = 0xb.3c2dd11dfdba2670p+15517L;
+ zi[3] = -0xf.6773dc63dee9b110p-905L;
+
+
+ for (i = 0; i < 4; i++)
+ {
+ __real__ a = xr[i];
+ __imag__ a = xi[i];
+ __real__ b = yr[i];
+ __imag__ b = yi[i];
+ __real__ z = zr[i];
+ __imag__ z = zi[i];
+ c = a/b;
+ cr = __real__ c;
+ ci = __imag__ c;
+
+ if (!match(c,z)){
+ ok = 0;
+ }
+ }
+ if (!ok)
+ abort();
+ exit(0);
+}
@@ -37,31 +37,118 @@ see the files COPYING3 and COPYING.RUNTIME respectively. If not, see
#define __divkc3 __divkc3_sw
#endif
+# define RBIG ((__LIBGCC_TF_MAX__)/2)
+# define RMIN (__LIBGCC_TF_MIN__)
+# define RMIN2 (__LIBGCC_TF_EPSILON__)
+# define RMINSCAL (1/__LIBGCC_TF_EPSILON__)
+# define RMAX2 ((RBIG)*(RMIN2))
+
+
TCtype
__divkc3 (TFtype a, TFtype b, TFtype c, TFtype d)
{
TFtype denom, ratio, x, y;
TCtype res;
- /* ??? We can get better behavior from logarithmic scaling instead of
- the division. But that would mean starting to link libgcc against
- libm. We could implement something akin to ldexp/frexp as gcc builtins
- fairly easily... */
- if (FABS (c) < FABS (d))
+ /* double, extended, long double have significant potential
+ underflow/overflow errors than can be greatly reduced with
+ a limited number of tests and adjustments. float is handled
+ the same way when no HW double is available.
+ */
+
+ /* Scale by max(c,d) to reduce chances of denominator overflowing. */
+ if (FABS(c) < FABS(d))
{
+ /* Prevent underflow when denominator is near max representable. */
+ if (FABS(d) >= RBIG)
+ {
+ a = a / 2;
+ b = b / 2;
+ c = c / 2;
+ d = d / 2;
+ }
+ /* Avoid overflow/underflow issues when c and d are small.
+ Scaling up helps avoid some underflows.
+ No new overflow possible since c&d < RMIN2. */
+ if (FABS(d) < RMIN2)
+ {
+ a = a * RMINSCAL;
+ b = b * RMINSCAL;
+ c = c * RMINSCAL;
+ d = d * RMINSCAL;
+ }
+ else
+ {
+ if (((FABS(a) < RMIN) && (FABS(b) < RMAX2) && (FABS(d) < RMAX2))
+ || ((FABS(b) < RMIN) && (FABS(a) < RMAX2) && (FABS(d) < RMAX2)))
+ {
+ a = a * RMINSCAL;
+ b = b * RMINSCAL;
+ c = c * RMINSCAL;
+ d = d * RMINSCAL;
+ }
+ }
ratio = c / d;
denom = (c * ratio) + d;
- x = ((a * ratio) + b) / denom;
- y = ((b * ratio) - a) / denom;
+ /* Choose alternate order of computation if ratio is subnormal. */
+ if (FABS(ratio) > RMIN)
+ {
+ x = ((a * ratio) + b) / denom;
+ y = ((b * ratio) - a) / denom;
+ }
+ else
+ {
+ x = ((c * (a / d)) + b) / denom;
+ y = ((c * (b / d)) - a) / denom;
+ }
}
else
{
+ /* Prevent underflow when denominator is near max representable. */
+ if (FABS(c) >= RBIG)
+ {
+ a = a / 2;
+ b = b / 2;
+ c = c / 2;
+ d = d / 2;
+ }
+ /* Avoid overflow/underflow issues when both c and d are small.
+ Scaling up helps avoid some underflows.
+ No new overflow possible since both c&d are less than RMIN2. */
+ if (FABS(c) < RMIN2)
+ {
+ a = a * RMINSCAL;
+ b = b * RMINSCAL;
+ c = c * RMINSCAL;
+ d = d * RMINSCAL;
+ }
+ else
+ {
+ if (((FABS(a) < RMIN) && (FABS(b) < RMAX2) && (FABS(c) < RMAX2))
+ || ((FABS(b) < RMIN) && (FABS(a) < RMAX2) && (FABS(c) < RMAX2)))
+ {
+ a = a * RMINSCAL;
+ b = b * RMINSCAL;
+ c = c * RMINSCAL;
+ d = d * RMINSCAL;
+ }
+ }
ratio = d / c;
denom = (d * ratio) + c;
- x = ((b * ratio) + a) / denom;
- y = (b - (a * ratio)) / denom;
+ /* Choose alternate order of computation if ratio is subnormal. */
+ if (FABS(ratio) > RMIN)
+ {
+ x = ((b * ratio) + a) / denom;
+ y = (b - (a * ratio)) / denom;
+ }
+ else
+ {
+ x = (a + (d * (b / c))) / denom;
+ y = (b - (d * (a / c))) / denom;
+ }
}
+
/* Recover infinities and zeros that computed as NaN+iNaN; the only cases
are nonzero/zero, infinite/finite, and finite/infinite. */
if (isnan (x) && isnan (y))
@@ -1860,33 +1860,57 @@ NAME (TYPE x, int m)
#if defined(L_mulhc3) || defined(L_divhc3)
# define MTYPE HFtype
# define CTYPE HCtype
+# define XMTYPE SFtype
+# define XCTYPE SCtype
# define MODE hc
# define CEXT __LIBGCC_HF_FUNC_EXT__
# define NOTRUNC (!__LIBGCC_HF_EXCESS_PRECISION__)
#elif defined(L_mulsc3) || defined(L_divsc3)
# define MTYPE SFtype
# define CTYPE SCtype
+# define XMTYPE DFtype
+# define XCTYPE DCtype
# define MODE sc
# define CEXT __LIBGCC_SF_FUNC_EXT__
# define NOTRUNC (!__LIBGCC_SF_EXCESS_PRECISION__)
+# define RBIG ((__LIBGCC_SF_MAX__)/2)
+# define RMIN (__LIBGCC_SF_MIN__)
+# define RMIN2 (__LIBGCC_SF_EPSILON__)
+# define RMINSCAL (1/__LIBGCC_SF_EPSILON__)
+# define RMAX2 ((RBIG)*(RMIN2))
#elif defined(L_muldc3) || defined(L_divdc3)
# define MTYPE DFtype
# define CTYPE DCtype
# define MODE dc
# define CEXT __LIBGCC_DF_FUNC_EXT__
# define NOTRUNC (!__LIBGCC_DF_EXCESS_PRECISION__)
+# define RBIG ((__LIBGCC_DF_MAX__)/2)
+# define RMIN (__LIBGCC_DF_MIN__)
+# define RMIN2 (__LIBGCC_DF_EPSILON__)
+# define RMINSCAL (1/__LIBGCC_DF_EPSILON__)
+# define RMAX2 ((RBIG)*(RMIN2))
#elif defined(L_mulxc3) || defined(L_divxc3)
# define MTYPE XFtype
# define CTYPE XCtype
# define MODE xc
# define CEXT __LIBGCC_XF_FUNC_EXT__
# define NOTRUNC (!__LIBGCC_XF_EXCESS_PRECISION__)
+# define RBIG ((__LIBGCC_XF_MAX__)/2)
+# define RMIN (__LIBGCC_XF_MIN__)
+# define RMIN2 (__LIBGCC_XF_EPSILON__)
+# define RMINSCAL (1/__LIBGCC_XF_EPSILON__)
+# define RMAX2 ((RBIG)*(RMIN2))
#elif defined(L_multc3) || defined(L_divtc3)
# define MTYPE TFtype
# define CTYPE TCtype
# define MODE tc
# define CEXT __LIBGCC_TF_FUNC_EXT__
# define NOTRUNC (!__LIBGCC_TF_EXCESS_PRECISION__)
+# define RBIG ((__LIBGCC_TF_MAX__)/2)
+# define RMIN (__LIBGCC_TF_MIN__)
+# define RMIN2 (__LIBGCC_TF_EPSILON__)
+# define RMINSCAL (1/__LIBGCC_TF_EPSILON__)
+# define RMAX2 ((RBIG)*(RMIN2))
#else
# error
#endif
@@ -1994,30 +2018,134 @@ CONCAT3(__mul,MODE,3) (MTYPE a, MTYPE b, MTYPE c, MTYPE d)
CTYPE
CONCAT3(__div,MODE,3) (MTYPE a, MTYPE b, MTYPE c, MTYPE d)
{
+#if defined(L_divhc3) \
+ || (defined(L_divsc3) && defined(__LIBGCC_HAVE_HWDBL__) )
+
+ /* Half precision is handled with float precision.
+ float is handled with double precision when double precision
+ hardware is available.
+ Due to the additional precision, the simple complex divide
+ method (without Smith's method) is sufficient to get accurate
+ answers and runs slightly faster than Smith's method. */
+
+ XMTYPE aa, bb, cc, dd;
+ XMTYPE denom;
+ MTYPE x, y;
+ CTYPE res;
+ aa = a;
+ bb = b;
+ cc = c;
+ dd = d;
+
+ denom = (cc * cc) + (dd * dd);
+ x = ((aa * cc) + (bb * dd)) / denom;
+ y = ((bb * cc) - (aa * dd)) / denom;
+
+#else
MTYPE denom, ratio, x, y;
CTYPE res;
- /* ??? We can get better behavior from logarithmic scaling instead of
- the division. But that would mean starting to link libgcc against
- libm. We could implement something akin to ldexp/frexp as gcc builtins
- fairly easily... */
- if (FABS (c) < FABS (d))
+ /* double, extended, long double have significant potential
+ underflow/overflow errors than can be greatly reduced with
+ a limited number of tests and adjustments. float is handled
+ the same way when no HW double is available.
+ */
+
+ /* Scale by max(c,d) to reduce chances of denominator overflowing. */
+ if (FABS(c) < FABS(d))
{
+ /* Prevent underflow when denominator is near max representable. */
+ if (FABS(d) >= RBIG)
+ {
+ a = a / 2;
+ b = b / 2;
+ c = c / 2;
+ d = d / 2;
+ }
+ /* Avoid overflow/underflow issues when c and d are small.
+ Scaling up helps avoid some underflows.
+ No new overflow possible since c&d < RMIN2. */
+ if (FABS(d) < RMIN2)
+ {
+ a = a * RMINSCAL;
+ b = b * RMINSCAL;
+ c = c * RMINSCAL;
+ d = d * RMINSCAL;
+ }
+ else
+ {
+ if (((FABS(a) < RMIN) && (FABS(b) < RMAX2) && (FABS(d) < RMAX2))
+ || ((FABS(b) < RMIN) && (FABS(a) < RMAX2) && (FABS(d) < RMAX2)))
+ {
+ a = a * RMINSCAL;
+ b = b * RMINSCAL;
+ c = c * RMINSCAL;
+ d = d * RMINSCAL;
+ }
+ }
ratio = c / d;
denom = (c * ratio) + d;
- x = ((a * ratio) + b) / denom;
- y = ((b * ratio) - a) / denom;
+ /* Choose alternate order of computation if ratio is subnormal. */
+ if (FABS(ratio) > RMIN)
+ {
+ x = ((a * ratio) + b) / denom;
+ y = ((b * ratio) - a) / denom;
+ }
+ else
+ {
+ x = ((c * (a / d)) + b) / denom;
+ y = ((c * (b / d)) - a) / denom;
+ }
}
else
{
+ /* Prevent underflow when denominator is near max representable. */
+ if (FABS(c) >= RBIG)
+ {
+ a = a / 2;
+ b = b / 2;
+ c = c / 2;
+ d = d / 2;
+ }
+ /* Avoid overflow/underflow issues when both c and d are small.
+ Scaling up helps avoid some underflows.
+ No new overflow possible since both c&d are less than RMIN2. */
+ if (FABS(c) < RMIN2)
+ {
+ a = a * RMINSCAL;
+ b = b * RMINSCAL;
+ c = c * RMINSCAL;
+ d = d * RMINSCAL;
+ }
+ else
+ {
+ if (((FABS(a) < RMIN) && (FABS(b) < RMAX2) && (FABS(c) < RMAX2))
+ || ((FABS(b) < RMIN) && (FABS(a) < RMAX2) && (FABS(c) < RMAX2)))
+ {
+ a = a * RMINSCAL;
+ b = b * RMINSCAL;
+ c = c * RMINSCAL;
+ d = d * RMINSCAL;
+ }
+ }
ratio = d / c;
denom = (d * ratio) + c;
- x = ((b * ratio) + a) / denom;
- y = (b - (a * ratio)) / denom;
+ /* Choose alternate order of computation if ratio is subnormal. */
+ if (FABS(ratio) > RMIN)
+ {
+ x = ((b * ratio) + a) / denom;
+ y = (b - (a * ratio)) / denom;
+ }
+ else
+ {
+ x = (a + (d * (b / c))) / denom;
+ y = (b - (d * (a / c))) / denom;
+ }
}
+#endif
- /* Recover infinities and zeros that computed as NaN+iNaN; the only cases
- are nonzero/zero, infinite/finite, and finite/infinite. */
+ /* Recover infinities and zeros that computed as NaN+iNaN; the only
+ cases are nonzero/zero, infinite/finite, and finite/infinite. */
if (isnan (x) && isnan (y))
{
if (c == 0.0 && d == 0.0 && (!isnan (a) || !isnan (b)))