83fbf5568cba4478bb5ca501aab9e6ec8b8b39fd
[alexxy/gromacs.git] / src / testutils / testasserts.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2014,2015,2016,2018,2019,2020, 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 "testutils/testasserts.h"
45
46 #include <cmath>
47 #include <cstdio>
48
49 #include <limits>
50
51 #include <gtest/gtest.h>
52
53 #include "gromacs/options/basicoptions.h"
54 #include "gromacs/options/ioptionscontainer.h"
55 #include "gromacs/utility/basedefinitions.h"
56 #include "gromacs/utility/exceptions.h"
57 #include "gromacs/utility/stringutil.h"
58
59 #include "testutils/testoptions.h"
60
61 namespace gmx
62 {
63
64 namespace test
65 {
66
67 namespace
68 {
69
70 //! Whether to print the message from expected exceptions.
71 bool g_showExpectedExceptions = false;
72
73 //! \cond
74 GMX_TEST_OPTIONS(ExceptionOptions, options)
75 {
76     options->addOption(BooleanOption("show-error-messages")
77                                .store(&g_showExpectedExceptions)
78                                .description("Show error messages from expected "
79                                             "exceptions"));
80 }
81 //! \endcond
82 } // namespace
83
84 namespace internal
85 {
86
87 //! \cond internal
88 void processExpectedException(const std::exception& ex)
89 {
90     if (g_showExpectedExceptions)
91     {
92         std::printf("Exception message (from expected exception):\n");
93         formatExceptionMessageToFile(stdout, ex);
94     }
95 }
96 //! \endcond
97
98 } // namespace internal
99
100 namespace
101 {
102
103 using ::testing::internal::FloatingPoint;
104
105 //! \internal \addtogroup module_testutils
106 //! \{
107
108 /*! \name Helper functions for computing floating-point differences
109  *
110  * These routines are used to initialize FloatingPointDifference.
111  * They peek into some internal types from Google Test (gtest-internal.h),
112  * and duplicate some other functionality from there, but that is likely
113  * a better alternative than just copying all that code here.
114  */
115 //! \{
116
117 /*! \brief
118  * Computes biased integer representation for a floating-point value.
119  *
120  * This moves the integer representation from a sign-and-magnitude
121  * representation to a biased representation where the 0x8000... represents
122  * zero, and the order of the integer values matches the order of the
123  * floating-point values.
124  */
125 template<typename FloatType>
126 typename FloatingPoint<FloatType>::Bits floatingPointToBiasedInteger(const FloatingPoint<FloatType>& value)
127 {
128     if (value.sign_bit())
129     {
130         return ~value.bits() + 1;
131     }
132     else
133     {
134         return value.bits() | FloatingPoint<FloatType>::kSignBitMask;
135     }
136 }
137
138 /*! \brief
139  * Computes the magnitude of the difference in ULPs between two numbers,
140  * treating also values of different sign.
141  */
142 template<typename FloatType>
143 uint64_t calculateUlpDifference(const FloatingPoint<FloatType>& value1, const FloatingPoint<FloatType>& value2)
144 {
145     typename FloatingPoint<FloatType>::Bits biased1 = floatingPointToBiasedInteger(value1);
146     typename FloatingPoint<FloatType>::Bits biased2 = floatingPointToBiasedInteger(value2);
147     return biased1 > biased2 ? biased1 - biased2 : biased2 - biased1;
148 }
149
150 /*! \brief
151  * Helper to implement the constructors for FloatingPointDifference.
152  */
153 template<typename FloatType>
154 void initDifference(FloatType raw1, FloatType raw2, double* absoluteDifference, uint64_t* ulpDifference, bool* bSignDifference)
155 {
156     FloatingPoint<FloatType> value1(raw1);
157     FloatingPoint<FloatType> value2(raw2);
158
159     if (value1.is_nan() || value2.is_nan())
160     {
161         *absoluteDifference = std::numeric_limits<double>::quiet_NaN();
162         *bSignDifference    = false;
163         *ulpDifference      = 0;
164         return;
165     }
166     *absoluteDifference = std::fabs(raw1 - raw2);
167     *bSignDifference    = (value1.sign_bit() != value2.sign_bit());
168     *ulpDifference      = calculateUlpDifference(value1, value2);
169 }
170
171 /*! \brief
172  * Converts a relative tolerance into an ULP difference.
173  */
174 template<typename FloatType>
175 uint64_t relativeToleranceToUlp(FloatType tolerance)
176 {
177     FloatingPoint<FloatType> m(1.0);
178     FloatingPoint<FloatType> t(1.0 + tolerance);
179     return calculateUlpDifference<FloatType>(m, t);
180 }
181
182 //! \}
183 //! \}
184
185 } // namespace
186
187 /********************************************************************
188  * FloatingPointDifference
189  */
190
191 FloatingPointDifference::FloatingPointDifference(float ref, float value) :
192     termMagnitude_(std::abs(ref))
193 {
194     initDifference(ref, value, &absoluteDifference_, &ulpDifference_, &bSignDifference_);
195     bDouble_ = false;
196 }
197
198 FloatingPointDifference::FloatingPointDifference(double ref, double value) :
199     termMagnitude_(std::abs(ref))
200 {
201     initDifference(ref, value, &absoluteDifference_, &ulpDifference_, &bSignDifference_);
202     bDouble_ = true;
203 }
204
205 bool FloatingPointDifference::isNaN() const
206 {
207     return FloatingPoint<double>(absoluteDifference_).is_nan();
208 }
209
210 std::string FloatingPointDifference::toString() const
211 {
212     std::string relDiffStr;
213
214     if (termMagnitude_ > 0)
215     {
216         // If the reference value is finite we calculate the proper quotient
217         relDiffStr = formatString("%.3g", std::abs(absoluteDifference_ / termMagnitude_));
218     }
219     else if (absoluteDifference_ == 0.0)
220     {
221         // If the numbers are identical the quotient is strictly NaN here, but
222         // there no reason to worry when we have a perfect match.
223         relDiffStr = formatString("%.3g", 0.0);
224     }
225     else
226     {
227         // If the reference value is zero and numbers are non-identical, relative difference is infinite.
228         relDiffStr = formatString("Inf");
229     }
230
231     return formatString("%g (%" PRIu64 " %s-prec. ULPs, rel. %s)%s", absoluteDifference_,
232                         ulpDifference_, isDouble() ? "double" : "single", relDiffStr.c_str(),
233                         bSignDifference_ ? ", signs differ" : "");
234 }
235
236 /********************************************************************
237  * FloatingPointTolerance
238  */
239
240 bool FloatingPointTolerance::isWithin(const FloatingPointDifference& difference) const
241 {
242     if (difference.isNaN())
243     {
244         return false;
245     }
246
247     if (bSignMustMatch_ && difference.signsDiffer())
248     {
249         return false;
250     }
251
252     const double absoluteTolerance =
253             difference.isDouble() ? doubleAbsoluteTolerance_ : singleAbsoluteTolerance_;
254     if (difference.asAbsolute() < absoluteTolerance)
255     {
256         return true;
257     }
258
259     // By using smaller-than-or-equal below, we allow the test to pass if
260     // the numbers are identical, even if the term magnitude is 0, which seems
261     // a reasonable thing to do...
262     const double relativeTolerance =
263             difference.isDouble() ? doubleRelativeTolerance_ : singleRelativeTolerance_;
264
265     if (difference.asAbsolute() <= relativeTolerance * difference.termMagnitude())
266     {
267         return true;
268     }
269
270     const uint64_t ulpTolerance = difference.isDouble() ? doubleUlpTolerance_ : singleUlpTolerance_;
271     return ulpTolerance < UINT64_MAX && difference.asUlps() <= ulpTolerance;
272 }
273
274 std::string FloatingPointTolerance::toString(const FloatingPointDifference& difference) const
275 {
276     std::string  result;
277     const double absoluteTolerance =
278             difference.isDouble() ? doubleAbsoluteTolerance_ : singleAbsoluteTolerance_;
279     const double relativeTolerance =
280             difference.isDouble() ? doubleRelativeTolerance_ : singleRelativeTolerance_;
281     const uint64_t ulpTolerance = difference.isDouble() ? doubleUlpTolerance_ : singleUlpTolerance_;
282
283     if (absoluteTolerance > 0.0)
284     {
285         result.append(formatString("abs. %g", absoluteTolerance));
286     }
287     if (relativeTolerance > 0.0)
288     {
289         if (!result.empty())
290         {
291             result.append(", ");
292         }
293         result.append(formatString("rel. %.3g", relativeTolerance));
294     }
295     if (ulpTolerance < UINT64_MAX)
296     {
297         if (!result.empty())
298         {
299             result.append(", ");
300         }
301         result.append(formatString("%" PRIu64 " ULPs", ulpTolerance));
302     }
303     if (bSignMustMatch_)
304     {
305         if (!result.empty())
306         {
307             result.append(", ");
308         }
309         result.append("sign must match");
310     }
311     return result;
312 }
313
314 // Doxygen does not recognize this as the same function as in the header...
315 //! \cond
316 FloatingPointTolerance relativeToleranceAsFloatingPoint(double magnitude, double tolerance)
317 {
318     return relativeToleranceAsPrecisionDependentFloatingPoint(magnitude, float(tolerance), tolerance);
319 }
320
321 FloatingPointTolerance relativeToleranceAsPrecisionDependentFloatingPoint(double magnitude,
322                                                                           float  singleTolerance,
323                                                                           double doubleTolerance)
324 {
325     const float  absoluteSingleTolerance = std::abs(float(magnitude)) * singleTolerance;
326     const double absoluteDoubleTolerance = std::abs(magnitude) * doubleTolerance;
327     return {
328         absoluteSingleTolerance, absoluteDoubleTolerance, singleTolerance, doubleTolerance, UINT64_MAX, UINT64_MAX, false
329     };
330 }
331 //! \endcond
332
333 } // namespace test
334 } // namespace gmx