Extend support for float/double tolerances in testing
[alexxy/gromacs.git] / src / testutils / testasserts.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \internal \file
36  * \brief
37  * Implements floating-point comparison routines from testasserts.h.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_testutils
41  */
42 #include "gmxpre.h"
43
44 #include "testasserts.h"
45
46 #include <cmath>
47
48 #include <limits>
49
50 #include <gtest/gtest.h>
51
52 #include "gromacs/utility/basedefinitions.h"
53 #include "gromacs/utility/stringutil.h"
54
55 namespace gmx
56 {
57
58 namespace test
59 {
60
61 namespace
62 {
63
64 using ::testing::internal::FloatingPoint;
65
66 //! \internal \addtogroup module_testutils
67 //! \{
68
69 /*! \name Helper functions for computing floating-point differences
70  *
71  * These routines are used to initialize FloatingPointDifference.
72  * They peek into some internal types from Google Test (gtest-internal.h),
73  * and duplicate some other functionality from there, but that is likely
74  * a better alternative than just copying all that code here.
75  */
76 //! \{
77
78 /*! \brief
79  * Computes biased integer representation for a floating-point value.
80  *
81  * This moves the integer representation from a sign-and-magnitude
82  * representation to a biased representation where the 0x8000... represents
83  * zero, and the order of the integer values matches the order of the
84  * floating-point values.
85  */
86 template <typename FloatType>
87 typename FloatingPoint<FloatType>::Bits
88 floatingPointToBiasedInteger(const FloatingPoint<FloatType> &value)
89 {
90     if (value.sign_bit())
91     {
92         return ~value.bits() + 1;
93     }
94     else
95     {
96         return value.bits() | FloatingPoint<FloatType>::kSignBitMask;
97     }
98 }
99
100 /*! \brief
101  * Computes the magnitude of the difference in ULPs between two numbers,
102  * treating also values of different sign.
103  */
104 template <typename FloatType>
105 gmx_uint64_t calculateUlpDifference(const FloatingPoint<FloatType> &value1,
106                                     const FloatingPoint<FloatType> &value2)
107 {
108     typename FloatingPoint<FloatType>::Bits biased1
109         = floatingPointToBiasedInteger(value1);
110     typename FloatingPoint<FloatType>::Bits biased2
111         = floatingPointToBiasedInteger(value2);
112     return biased1 > biased2 ? biased1 - biased2 : biased2 - biased1;
113 }
114
115 /*! \brief
116  * Helper to implement the constructors for FloatingPointDifference.
117  */
118 template <typename FloatType>
119 void initDifference(FloatType raw1, FloatType raw2, double *absoluteDifference,
120                     gmx_uint64_t *ulpDifference, bool *bSignDifference)
121 {
122     FloatingPoint<FloatType> value1(raw1);
123     FloatingPoint<FloatType> value2(raw2);
124
125     if (value1.is_nan() || value2.is_nan())
126     {
127         *absoluteDifference = std::numeric_limits<double>::quiet_NaN();
128         *bSignDifference    = false;
129         *ulpDifference      = 0;
130         return;
131     }
132     *absoluteDifference = std::fabs(raw1 - raw2);
133     *bSignDifference    = (value1.sign_bit() != value2.sign_bit());
134     *ulpDifference      = calculateUlpDifference(value1, value2);
135 }
136
137 /*! \brief
138  * Converts a relative tolerance into an ULP difference.
139  */
140 template <typename FloatType>
141 gmx_uint64_t relativeToleranceToUlp(FloatType tolerance)
142 {
143     FloatingPoint<FloatType> m(1.0);
144     FloatingPoint<FloatType> t(1.0 + tolerance);
145     return calculateUlpDifference<FloatType>(m, t);
146 }
147
148 //! \}
149 //! \}
150
151 }       // namespace
152
153 /********************************************************************
154  * FloatingPointDifference
155  */
156
157 FloatingPointDifference::FloatingPointDifference(float value1, float value2)
158 {
159     initDifference(value1, value2,
160                    &absoluteDifference_, &ulpDifference_, &bSignDifference_);
161     bDouble_ = false;
162 }
163
164 FloatingPointDifference::FloatingPointDifference(double value1, double value2)
165 {
166     initDifference(value1, value2,
167                    &absoluteDifference_, &ulpDifference_, &bSignDifference_);
168     bDouble_ = true;
169 }
170
171 bool FloatingPointDifference::isNaN() const
172 {
173     return FloatingPoint<double>(absoluteDifference_).is_nan();
174 }
175
176 std::string FloatingPointDifference::toString() const
177 {
178     const double eps = isDouble() ? GMX_DOUBLE_EPS : GMX_FLOAT_EPS;
179     return formatString("%g (%" GMX_PRIu64 " %s-prec. ULPs, rel. %.3g)%s",
180                         absoluteDifference_, ulpDifference_,
181                         isDouble() ? "double" : "single",
182                         ulpDifference_ * eps,
183                         bSignDifference_ ? ", signs differ" : "");
184 }
185
186 /********************************************************************
187  * FloatingPointTolerance
188  */
189
190 bool FloatingPointTolerance::isWithin(
191         const FloatingPointDifference &difference) const
192 {
193     if (difference.isNaN())
194     {
195         return false;
196     }
197
198     if (bSignMustMatch_ && difference.signsDiffer())
199     {
200         return false;
201     }
202
203     const double absoluteTolerance
204         = difference.isDouble() ? doubleAbsoluteTolerance_ : singleAbsoluteTolerance_;
205     if (difference.asAbsolute() < absoluteTolerance)
206     {
207         return true;
208     }
209
210     const gmx_uint64_t ulpTolerance
211         = difference.isDouble() ? doubleUlpTolerance_ : singleUlpTolerance_;
212     if (ulpTolerance < GMX_UINT64_MAX && difference.asUlps() <= ulpTolerance)
213     {
214         return true;
215     }
216     return false;
217 }
218
219 std::string FloatingPointTolerance::toString(const FloatingPointDifference &difference) const
220 {
221     std::string        result;
222     const double       absoluteTolerance
223         = difference.isDouble() ? doubleAbsoluteTolerance_ : singleAbsoluteTolerance_;
224     const gmx_uint64_t ulpTolerance
225         = difference.isDouble() ? doubleUlpTolerance_ : singleUlpTolerance_;
226     const double       eps
227         = difference.isDouble() ? GMX_DOUBLE_EPS : GMX_FLOAT_EPS;
228
229     if (absoluteTolerance > 0.0)
230     {
231         result.append(formatString("abs. %g", absoluteTolerance));
232     }
233     if (ulpTolerance < GMX_UINT64_MAX)
234     {
235         if (!result.empty())
236         {
237             result.append(", ");
238         }
239         result.append(formatString("%" GMX_PRIu64 " ULPs (rel. %.3g)", ulpTolerance, ulpTolerance * eps));
240     }
241     if (bSignMustMatch_)
242     {
243         if (!result.empty())
244         {
245             result.append(", ");
246         }
247         result.append("sign must match");
248     }
249     return result;
250 }
251
252 // Doxygen does not recognize this as the same function as in the header...
253 //! \cond
254 FloatingPointTolerance
255 relativeToleranceAsFloatingPoint(double magnitude, double tolerance)
256 {
257     const double absoluteTolerance = magnitude * tolerance;
258     return FloatingPointTolerance(absoluteTolerance, absoluteTolerance,
259                                   relativeToleranceToUlp<float>(tolerance),
260                                   relativeToleranceToUlp<double>(tolerance),
261                                   false);
262 }
263 //! \endcond
264
265 } // namespace test
266 } // namespace gmx