Extend support for float/double tolerances in testing
[alexxy/gromacs.git] / src / testutils / testasserts.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2013,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 /*! \libinternal \file
36  * \brief
37  * Extra assertions for unit tests.
38  *
39  * This file provides assertion macros that extend/replace Google Test
40  * assertions for:
41  *  - exceptions
42  *  - floating-point comparison
43  *  - comparison against NULL
44  *
45  * \if internal
46  * \todo
47  * The implementation is somewhat ugly, and accesses some Google Test
48  * internals.  Could be nice to clean it up a bit.
49  * \endif
50  *
51  * \author Teemu Murtola <teemu.murtola@gmail.com>
52  * \inlibraryapi
53  * \ingroup module_testutils
54  */
55 #ifndef GMX_TESTUTILS_TESTASSERTS_H
56 #define GMX_TESTUTILS_TESTASSERTS_H
57
58 #include <string>
59
60 #include <gtest/gtest.h>
61
62 #include "gromacs/utility/basedefinitions.h"
63 #include "gromacs/utility/exceptions.h"
64 #include "gromacs/utility/real.h"
65
66 namespace gmx
67 {
68 namespace test
69 {
70
71 //! \libinternal \addtogroup module_testutils
72 //! \{
73
74 /*! \name Assertions for exceptions
75  *
76  * These macros replace `(ASSERT|EXPECT)(_NO)?_THROW` from Google Test.
77  * They are used exactly like the Google Test ones, but also print details of
78  * any unexpected exceptions using \Gromacs-specific routines.
79  * This makes it much easier to see at one glance what went wrong.
80  * See Google Test documentation for details on how to use the macros.
81  */
82 //! \{
83
84 //! \cond internal
85 /*! \brief
86  * Internal implementation macro for exception assertations.
87  *
88  * \param statement          Statements to execute.
89  * \param expected_exception Exception type that \p statement should throw.
90  * \param fail               Function/macro to call on failure.
91  *
92  * The implementation is copied and adjusted from
93  * include/gtest/internal/gtest-internal.h in Google Test 1.6.0.
94  */
95 #define GMX_TEST_THROW_(statement, expected_exception, fail) \
96     GTEST_AMBIGUOUS_ELSE_BLOCKER_ \
97     if (::testing::AssertionResult gmx_ar = ::testing::AssertionSuccess()) { \
98         bool gmx_caught_expected = false; \
99         try { \
100             GTEST_SUPPRESS_UNREACHABLE_CODE_WARNING_BELOW_(statement); \
101         } \
102         catch (expected_exception const &) { \
103             gmx_caught_expected = true; \
104         } \
105         catch (std::exception const &ex) { \
106             gmx_ar << "Expected: " #statement " throws an exception of type " \
107             << #expected_exception ".\n  Actual: it throws a different type.\n" \
108             << "Exception details:\n" << ::gmx::formatExceptionMessageToString(ex); \
109             goto GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__); \
110         } \
111         catch (...) { \
112             gmx_ar << "Expected: " #statement " throws an exception of type " \
113             << #expected_exception ".\n  Actual: it throws a different type."; \
114             goto GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__); \
115         } \
116         if (!gmx_caught_expected) { \
117             gmx_ar << "Expected: " #statement " throws an exception of type " \
118             << #expected_exception ".\n  Actual: it throws nothing."; \
119             goto GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__); \
120         } \
121     } else \
122         GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__) : \
123             fail(gmx_ar.message())
124
125 /*! \brief
126  * Internal implementation macro for exception assertations.
127  *
128  * \param statement          Statements to execute.
129  * \param fail               Function/macro to call on failure.
130  *
131  * The implementation is copied and adjusted from
132  * include/gtest/internal/gtest-internal.h in Google Test 1.6.0.
133  */
134 #define GMX_TEST_NO_THROW_(statement, fail) \
135     GTEST_AMBIGUOUS_ELSE_BLOCKER_ \
136     if (::testing::AssertionResult gmx_ar = ::testing::AssertionSuccess()) { \
137         try { \
138             GTEST_SUPPRESS_UNREACHABLE_CODE_WARNING_BELOW_(statement); \
139         } \
140         catch (std::exception const &ex) { \
141             gmx_ar << "Expected: " #statement " doesn't throw an exception.\n" \
142             << "  Actual: it throws.\n" \
143             << "Exception details:\n" << ::gmx::formatExceptionMessageToString(ex); \
144             goto GTEST_CONCAT_TOKEN_(gmx_label_testnothrow_, __LINE__); \
145         } \
146         catch (...) { \
147             gmx_ar << "Expected: " #statement " doesn't throw an exception.\n" \
148             << "  Actual: it throws."; \
149             goto GTEST_CONCAT_TOKEN_(gmx_label_testnothrow_, __LINE__); \
150         } \
151     } else \
152         GTEST_CONCAT_TOKEN_(gmx_label_testnothrow_, __LINE__) : \
153             fail(gmx_ar.message())
154 //! \endcond
155
156 /*! \brief
157  * Asserts that a statement throws a given exception.
158  *
159  * \hideinitializer
160  */
161 #define EXPECT_THROW_GMX(statement, expected_exception) \
162     GMX_TEST_THROW_(statement, expected_exception, GTEST_NONFATAL_FAILURE_)
163 /*! \brief
164  * Asserts that a statement does not throw.
165  *
166  * \hideinitializer
167  */
168 #define EXPECT_NO_THROW_GMX(statement) \
169     GMX_TEST_NO_THROW_(statement, GTEST_NONFATAL_FAILURE_)
170 /*! \brief
171  * Asserts that a statement throws a given exception.
172  *
173  * \hideinitializer
174  */
175 #define ASSERT_THROW_GMX(statement, expected_exception) \
176     GMX_TEST_THROW_(statement, expected_exception, GTEST_FATAL_FAILURE_)
177 /*! \brief
178  * Asserts that a statement does not throw.
179  *
180  * \hideinitializer
181  */
182 #define ASSERT_NO_THROW_GMX(statement) \
183     GMX_TEST_NO_THROW_(statement, GTEST_FATAL_FAILURE_)
184
185 //! \}
186
187 /*! \libinternal \brief
188  * Computes and represents a floating-point difference value.
189  *
190  * Methods in this class do not throw, except for toString(), which may throw
191  * std::bad_alloc.
192  *
193  * \see FloatingPointTolerance
194  */
195 class FloatingPointDifference
196 {
197     public:
198         //! Initializes a single-precision difference.
199         FloatingPointDifference(float value1, float value2);
200         //! Initializes a double-precision difference.
201         FloatingPointDifference(double value1, double value2);
202
203         /*! \brief
204          * Whether one or both of the compared values were NaN.
205          *
206          * If this returns `true`, other accessors return meaningless values.
207          */
208         bool isNaN() const;
209         //! Returns the difference as an absolute number (always non-negative).
210         double asAbsolute() const { return absoluteDifference_; }
211         /*! \brief
212          * Returns the difference as ULPs (always non-negative).
213          *
214          * The ULPs are calculated for the type that corresponds to the
215          * constructor used to initialize the difference.
216          * The ULP difference between 0.0 and -0.0 is zero.
217          */
218         gmx_uint64_t asUlps() const { return ulpDifference_; }
219         /*! \brief
220          * Whether the compared values were of different sign.
221          *
222          * 0.0 and -0.0 are treated as positive and negative, respectively.
223          */
224         bool signsDiffer() const { return bSignDifference_; }
225         /*! \brief
226          * Whether the difference is between single- or double-precision
227          * numbers.
228          */
229         bool isDouble() const { return bDouble_; }
230         //! Formats the difference as a string for assertion failure messages.
231         std::string toString() const;
232
233     private:
234         //! Stores the absolute difference, or NaN if one or both values were NaN.
235         double       absoluteDifference_;
236         gmx_uint64_t ulpDifference_;
237         bool         bSignDifference_;
238         /*! \brief
239          * Whether the difference was computed for single or double precision.
240          *
241          * This sets the units for `ulpDifference_`.
242          */
243         bool         bDouble_;
244 };
245
246 /*! \libinternal \brief
247  * Specifies a floating-point comparison tolerance and checks whether a
248  * difference is within the tolerance.
249  *
250  * The related functions section lists methods that can be construct methods
251  * using less parameters than the full constructor, and with more obvious
252  * semantics.  These should be preferred over using the constructor directly.
253  *
254  * Several types of tolerances are possible:
255  *  - _absolute tolerance_: difference between the values must be smaller than
256  *    the given tolerance for the check to pass.
257  *    Setting the absolute tolerance to zero disables the absolute tolerance
258  *    check.
259  *  - _ULP tolerance_: ULP (units of least precision) difference between the
260  *    values must be smaller than the given tolerance for the check to pass.
261  *    Setting the ULP tolerance to zero requires exact match.
262  *    Setting the ULP tolerance to GMX_UINT64_MAX disables the ULP check.
263  *    `0.0` and `-0.0` are treated as equal for the ULP check.
264  *  - _sign check_: if set, any values that are of different signs fail the
265  *    check (note that this also applies to `0.0` and `-0.0`: a value with a
266  *    different sign than the zero will fail the check).
267  *
268  * Either an absolute or a ULP tolerance must always be specified.
269  * If both are specified, then the check passes if either of the tolerances is
270  * satisfied.
271  *
272  * Any combination of absolute and ULP tolerance can be combined with the sign
273  * check.  In this case, the sign check must succeed for the check to pass,
274  * even if other tolerances are satisfied.
275  *
276  * The tolerances can be specified separately for single and double precision
277  * comparison.  Different initialization functions have different semantics on
278  * how the provided tolerance values are interpreted; check their
279  * documentation.
280  *
281  * Methods in this class do not throw, except for toString(), which may throw
282  * std::bad_alloc.
283  *
284  * \todo
285  * The factory methods that take ULP difference could be better formulated as
286  * methods that take the acceptable number of incorrect bits and/or the number
287  * of accurate bits.
288  *
289  * \see FloatingPointDifference
290  */
291 class FloatingPointTolerance
292 {
293     public:
294         /*! \brief
295          * Creates a tolerance with the specified values.
296          *
297          * \param[in]  singleAbsoluteTolerance
298          *     Allowed absolute difference in a single-precision number.
299          * \param[in]  doubleAbsoluteTolerance
300          *     Allowed absolute difference in a double-precision number.
301          * \param[in]  singleUlpTolerance
302          *     Allowed ULP difference in a single-precision number.
303          * \param[in]  doubleUlpTolerance
304          *     Allowed ULP difference in a double-precision number.
305          * \param[in]  bSignMustMatch
306          *     Whether sign mismatch fails the comparison.
307          */
308         FloatingPointTolerance(float        singleAbsoluteTolerance,
309                                double       doubleAbsoluteTolerance,
310                                gmx_uint64_t singleUlpTolerance,
311                                gmx_uint64_t doubleUlpTolerance,
312                                bool         bSignMustMatch)
313             : singleAbsoluteTolerance_(singleAbsoluteTolerance),
314               doubleAbsoluteTolerance_(doubleAbsoluteTolerance),
315               singleUlpTolerance_(singleUlpTolerance),
316               doubleUlpTolerance_(doubleUlpTolerance),
317               bSignMustMatch_(bSignMustMatch)
318         {
319         }
320
321         /*! \brief
322          * Checks whether a difference is within the specified tolerance.
323          *
324          * NaNs are always treated outside the tolerance.
325          */
326         bool isWithin(const FloatingPointDifference &difference) const;
327
328         //! Formats the tolerance as a string for assertion failure messages.
329         std::string toString(const FloatingPointDifference &difference) const;
330
331     private:
332         float        singleAbsoluteTolerance_;
333         double       doubleAbsoluteTolerance_;
334         gmx_uint64_t singleUlpTolerance_;
335         gmx_uint64_t doubleUlpTolerance_;
336         bool         bSignMustMatch_;
337 };
338
339 /*! \brief
340  * Creates a tolerance that only allows a specified ULP difference.
341  *
342  * The tolerance uses the given ULP value for both precisions, i.e., double
343  * precision will have much stricter tolerance.
344  *
345  * \related FloatingPointTolerance
346  */
347 static inline FloatingPointTolerance
348 ulpTolerance(gmx_uint64_t ulpDiff)
349 {
350     return FloatingPointTolerance(0.0, 0.0, ulpDiff, ulpDiff, false);
351 }
352
353 /*! \brief
354  * Creates a tolerance that allows a difference in two compared values that is
355  * relative to the given magnitude.
356  *
357  * \param[in] magnitude  Magnitude of the numbers the computation operates in.
358  * \param[in] tolerance  Relative tolerance permitted (e.g. 1e-4).
359  *
360  * In addition to setting an ULP tolerance equivalent to \p tolerance for both
361  * precisions, this sets the absolute tolerance such that values close to zero
362  * (in general, smaller than \p magnitude) do not fail the check if they
363  * differ by less than \p tolerance evaluated at \p magnitude.  This accounts
364  * for potential loss of precision for small values, and should be used when
365  * accuracy of values much less than \p magnitude do not matter for
366  * correctness.
367  *
368  * The ULP tolerance for different precisions will be different to make them
369  * both match \p tolerance.
370  *
371  * \related FloatingPointTolerance
372  */
373 FloatingPointTolerance
374     relativeToleranceAsFloatingPoint(double magnitude, double tolerance);
375
376 /*! \brief
377  * Creates a tolerance that allows a precision-dependent relative difference in
378  * a complex computation.
379  *
380  * \param[in] magnitude      Magnitude of the numbers the computation operates in.
381  * \param[in] singleUlpDiff  Expected accuracy of single-precision
382  *     computation (in ULPs).
383  * \param[in] doubleUlpDiff  Expected accuracy of double-precision
384  *     computation (in ULPs).
385  *
386  * This works as relativeToleranceAsUlp(), but allows setting the ULP
387  * difference separately for the different precisions.  This supports
388  * cases where the double-precision calculation can acceptably has a higher ULP
389  * difference, but relaxing the single-precision tolerance would lead to an
390  * unnecessarily loose test.
391  *
392  * \related FloatingPointTolerance
393  */
394 static inline FloatingPointTolerance
395 relativeToleranceAsPrecisionDependentUlp(double       magnitude,
396                                          gmx_uint64_t singleUlpDiff,
397                                          gmx_uint64_t doubleUlpDiff)
398 {
399     return FloatingPointTolerance(magnitude*singleUlpDiff*GMX_FLOAT_EPS,
400                                   magnitude*doubleUlpDiff*GMX_DOUBLE_EPS,
401                                   singleUlpDiff, doubleUlpDiff, false);
402 }
403
404 /*! \brief
405  * Creates a tolerance that allows a relative difference in a complex
406  * computation.
407  *
408  * \param[in] magnitude  Magnitude of the numbers the computation operates in.
409  * \param[in] ulpDiff    Expected accuracy of the computation (in ULPs).
410  *
411  * In addition to setting the ULP tolerance as ulpTolerance(), this sets the
412  * absolute tolerance such that values close to zero (in general, smaller than
413  * \p magnitude) do not fail the check if they differ by less than \p ulpDiff
414  * evaluated at \p magnitude.  This accounts for potential loss of precision
415  * for small values, and should be used when accuracy of values much less than
416  * \p magnitude do not matter for correctness.
417  *
418  * \related FloatingPointTolerance
419  */
420 static inline FloatingPointTolerance
421 relativeToleranceAsUlp(double magnitude, gmx_uint64_t ulpDiff)
422 {
423     return relativeToleranceAsPrecisionDependentUlp(magnitude, ulpDiff, ulpDiff);
424 }
425
426 /*! \brief
427  * Returns the default tolerance for comparing `real` numbers.
428  *
429  * \related FloatingPointTolerance
430  */
431 static inline FloatingPointTolerance defaultRealTolerance()
432 {
433     return relativeToleranceAsUlp(1.0, 4);
434 }
435
436 /*! \name Assertions for floating-point comparison
437  *
438  * These routines extend `(EXPECT|ASSERT)_(FLOAT|DOUBLE)_EQ` and
439  * `(EXPECT|ASSERT)_NEAR` from Google Test to provide more flexible assertions
440  * for floating-point values.
441  *
442  * See gmx::test::FloatingPointTolerance for the possible ways to specify the
443  * tolerance, and gmx::test::FloatingPointDifference for some additional
444  * details of the difference calculation.
445  */
446 //! \{
447
448 //! \cond internal
449 /*! \internal \brief
450  * Assertion predicate formatter for comparing two floating-point values.
451  */
452 template <typename FloatType>
453 static inline ::testing::AssertionResult assertEqualWithinTolerance(
454         const char *expr1, const char *expr2, const char * /*exprTolerance*/,
455         FloatType value1, FloatType value2,
456         const FloatingPointTolerance &tolerance)
457 {
458     FloatingPointDifference diff(value1, value2);
459     if (tolerance.isWithin(diff))
460     {
461         return ::testing::AssertionSuccess();
462     }
463     return ::testing::AssertionFailure()
464            << "  Value of: " << expr2 << std::endl
465            << "    Actual: " << value2 << std::endl
466            << "  Expected: " << expr1 << std::endl
467            << "  Which is: " << value1 << std::endl
468            << "Difference: " << diff.toString() << std::endl
469            << " Tolerance: " << tolerance.toString(diff);
470 }
471 //! \endcond
472
473 /*! \brief
474  * Asserts that two single-precision values are within the given tolerance.
475  *
476  * \hideinitializer
477  */
478 #define EXPECT_FLOAT_EQ_TOL(value1, value2, tolerance) \
479     EXPECT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<float>, \
480                         value1, value2, tolerance)
481 /*! \brief
482  * Asserts that two double-precision values are within the given tolerance.
483  *
484  * \hideinitializer
485  */
486 #define EXPECT_DOUBLE_EQ_TOL(value1, value2, tolerance) \
487     EXPECT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<double>, \
488                         value1, value2, tolerance)
489 /*! \def EXPECT_REAL_EQ_TOL
490  * \brief
491  * Asserts that two `real` values are within the given tolerance.
492  *
493  * \hideinitializer
494  */
495 /*! \brief
496  * Asserts that two single-precision values are within the given tolerance.
497  *
498  * \hideinitializer
499  */
500 #define ASSERT_FLOAT_EQ_TOL(value1, value2, tolerance) \
501     ASSERT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<float>, \
502                         value1, value2, tolerance)
503 /*! \brief
504  * Asserts that two double-precision values are within the given tolerance.
505  *
506  * \hideinitializer
507  */
508 #define ASSERT_DOUBLE_EQ_TOL(value1, value2, tolerance) \
509     ASSERT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<double>, \
510                         value1, value2, tolerance)
511 /*! \def ASSERT_REAL_EQ_TOL
512  * \brief
513  * Asserts that two `real` values are within the given tolerance.
514  *
515  * \hideinitializer
516  */
517
518 #ifdef GMX_DOUBLE
519 #define EXPECT_REAL_EQ_TOL(value1, value2, tolerance) \
520     EXPECT_DOUBLE_EQ_TOL(value1, value2, tolerance)
521 #define ASSERT_REAL_EQ_TOL(value1, value2, tolerance) \
522     ASSERT_DOUBLE_EQ_TOL(value1, value2, tolerance)
523 #else
524 #define EXPECT_REAL_EQ_TOL(value1, value2, tolerance) \
525     EXPECT_FLOAT_EQ_TOL(value1, value2, tolerance)
526 #define ASSERT_REAL_EQ_TOL(value1, value2, tolerance) \
527     ASSERT_FLOAT_EQ_TOL(value1, value2, tolerance)
528 #endif
529
530 //! \}
531
532 /*! \name Assertions for NULL comparison
533  *
534  * These macros should be used instead of `(EXPECT|ASSERT)_EQ(NULL, ...)`,
535  * because Google Test doesn't support the NULL comparison with xlC++ 12.1 on
536  * BG/Q.
537  */
538 //! \{
539
540 /*! \brief
541  * Asserts that a pointer is null.
542  *
543  * Works exactly like EXPECT_EQ comparing with a null pointer. */
544 #define EXPECT_NULL(val) EXPECT_EQ((void *) NULL, val)
545 /*! \brief
546  * Asserts that a pointer is null.
547  *
548  * Works exactly like ASSERT_EQ comparing with a null pointer. */
549 #define ASSERT_NULL(val) ASSERT_EQ((void *) NULL, val)
550
551 //! \}
552
553 } // namespace test
554 } // namespace gmx
555
556 #endif