2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2013,2014,2015,2016,2017,2018, 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.
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.
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.
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.
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.
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.
35 /*! \libinternal \file
37 * Extra assertions for unit tests.
39 * This file provides assertion macros that extend/replace Google Test
42 * - floating-point comparison
43 * - comparison against NULL
47 * The implementation is somewhat ugly, and accesses some Google Test
48 * internals. Could be nice to clean it up a bit.
51 * \author Teemu Murtola <teemu.murtola@gmail.com>
53 * \ingroup module_testutils
55 #ifndef GMX_TESTUTILS_TESTASSERTS_H
56 #define GMX_TESTUTILS_TESTASSERTS_H
60 #include <gtest/gtest.h>
62 #include "gromacs/utility/basedefinitions.h"
63 #include "gromacs/utility/exceptions.h"
64 #include "gromacs/utility/real.h"
76 * Called for an expected exception from EXPECT_THROW_GMX().
78 * \param[in] ex Exception that was thrown.
80 void processExpectedException(const std::exception
&ex
);
82 } // namespace internal
84 //! \libinternal \addtogroup module_testutils
87 /*! \name Assertions for exceptions
89 * These macros replace `(ASSERT|EXPECT)(_NO)?_THROW` from Google Test.
90 * They are used exactly like the Google Test ones, but also print details of
91 * any unexpected exceptions using \Gromacs-specific routines.
92 * This makes it much easier to see at one glance what went wrong.
93 * See Google Test documentation for details on how to use the macros.
99 * Internal implementation macro for exception assertations.
101 * \param statement Statements to execute.
102 * \param expected_exception Exception type that \p statement should throw.
103 * \param fail Function/macro to call on failure.
105 * The implementation is copied and adjusted from
106 * include/gtest/internal/gtest-internal.h in Google Test 1.6.0.
108 #define GMX_TEST_THROW_(statement, expected_exception, fail) \
109 GTEST_AMBIGUOUS_ELSE_BLOCKER_ \
110 if (::testing::AssertionResult gmx_ar = ::testing::AssertionSuccess()) { \
111 bool gmx_caught_expected = false; \
113 GTEST_SUPPRESS_UNREACHABLE_CODE_WARNING_BELOW_(statement); \
115 catch (expected_exception const &ex) { \
116 gmx_caught_expected = true; \
117 ::gmx::test::internal::processExpectedException(ex); \
119 catch (std::exception const &ex) { \
120 gmx_ar << "Expected: " #statement " throws an exception of type " \
121 << #expected_exception ".\n Actual: it throws a different type.\n" \
122 << "Exception details:\n" << ::gmx::formatExceptionMessageToString(ex); \
123 goto GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__); \
126 gmx_ar << "Expected: " #statement " throws an exception of type " \
127 << #expected_exception ".\n Actual: it throws a different type."; \
128 goto GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__); \
130 if (!gmx_caught_expected) { \
131 gmx_ar << "Expected: " #statement " throws an exception of type " \
132 << #expected_exception ".\n Actual: it throws nothing."; \
133 goto GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__); \
136 GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__) : \
137 fail(gmx_ar.message())
140 * Internal implementation macro for exception assertations.
142 * \param statement Statements to execute.
143 * \param fail Function/macro to call on failure.
145 * The implementation is copied and adjusted from
146 * include/gtest/internal/gtest-internal.h in Google Test 1.6.0.
148 #define GMX_TEST_NO_THROW_(statement, fail) \
149 GTEST_AMBIGUOUS_ELSE_BLOCKER_ \
150 if (::testing::AssertionResult gmx_ar = ::testing::AssertionSuccess()) { \
152 GTEST_SUPPRESS_UNREACHABLE_CODE_WARNING_BELOW_(statement); \
154 catch (std::exception const &ex) { \
155 gmx_ar << "Expected: " #statement " doesn't throw an exception.\n" \
156 << " Actual: it throws.\n" \
157 << "Exception details:\n" << ::gmx::formatExceptionMessageToString(ex); \
158 goto GTEST_CONCAT_TOKEN_(gmx_label_testnothrow_, __LINE__); \
161 gmx_ar << "Expected: " #statement " doesn't throw an exception.\n" \
162 << " Actual: it throws."; \
163 goto GTEST_CONCAT_TOKEN_(gmx_label_testnothrow_, __LINE__); \
166 GTEST_CONCAT_TOKEN_(gmx_label_testnothrow_, __LINE__) : \
167 fail(gmx_ar.message())
171 * Asserts that a statement throws a given exception.
175 #define EXPECT_THROW_GMX(statement, expected_exception) \
176 GMX_TEST_THROW_(statement, expected_exception, GTEST_NONFATAL_FAILURE_)
178 * Asserts that a statement does not throw.
182 #define EXPECT_NO_THROW_GMX(statement) \
183 GMX_TEST_NO_THROW_(statement, GTEST_NONFATAL_FAILURE_)
185 * Asserts that a statement throws a given exception.
189 #define ASSERT_THROW_GMX(statement, expected_exception) \
190 GMX_TEST_THROW_(statement, expected_exception, GTEST_FATAL_FAILURE_)
192 * Asserts that a statement does not throw.
196 #define ASSERT_NO_THROW_GMX(statement) \
197 GMX_TEST_NO_THROW_(statement, GTEST_FATAL_FAILURE_)
201 /*! \libinternal \brief
202 * Computes and represents a floating-point difference value.
204 * Methods in this class do not throw, except for toString(), which may throw
207 * \see FloatingPointTolerance
209 class FloatingPointDifference
213 /*! \brief Initializes a single-precision difference.
215 * \param ref First term in difference
216 * \param value Second term in difference
218 * For absolute and ULP differences the two parameters are equivalent,
219 * since the difference is symmetric. For relative differences
220 * the first term is interpreted as the reference value, from which
221 * we extract the magnitude to compare with.
223 FloatingPointDifference(float ref
, float value
);
225 /*! \brief Initializes a double-precision difference.
227 * \param ref First term in difference
228 * \param value Second term in difference
230 * For absolute and ULP differences the two parameters are equivalent,
231 * since the difference is symmetric. For relative differences
232 * the first term is interpreted as the reference value, from which
233 * we extract the magnitude to compare with.
235 FloatingPointDifference(double ref
, double value
);
238 * Whether one or both of the compared values were NaN.
240 * If this returns `true`, other accessors return meaningless values.
243 //! Returns the difference as an absolute number (always non-negative).
244 double asAbsolute() const { return absoluteDifference_
; }
246 * Returns the difference as ULPs (always non-negative).
248 * The ULPs are calculated for the type that corresponds to the
249 * constructor used to initialize the difference.
250 * The ULP difference between 0.0 and -0.0 is zero.
252 uint64_t asUlps() const { return ulpDifference_
; }
254 * Whether the compared values were of different sign.
256 * 0.0 and -0.0 are treated as positive and negative, respectively.
258 bool signsDiffer() const { return bSignDifference_
; }
260 * Whether the difference is between single- or double-precision
263 bool isDouble() const { return bDouble_
; }
264 //! Formats the difference as a string for assertion failure messages.
265 std::string
toString() const;
267 //! Returns the magnitude of the original second term of the difference.
268 double termMagnitude() const { return termMagnitude_
; }
272 //! Save the magnitude of the reference value for relative (i.e., not ULP) tolerance
273 double termMagnitude_
;
274 //! Stores the absolute difference, or NaN if one or both values were NaN.
275 double absoluteDifference_
;
276 uint64_t ulpDifference_
;
277 bool bSignDifference_
;
279 * Whether the difference was computed for single or double precision.
281 * This sets the units for `ulpDifference_`.
286 /*! \libinternal \brief
287 * Specifies a floating-point comparison tolerance and checks whether a
288 * difference is within the tolerance.
290 * The related functions section lists methods that can be construct methods
291 * using less parameters than the full constructor, and with more obvious
292 * semantics. These should be preferred over using the constructor directly.
294 * Several types of tolerances are possible:
295 * - _absolute tolerance_: difference between the values must be smaller than
296 * the given tolerance for the check to pass.
297 * Setting the absolute tolerance to zero disables the absolute tolerance
299 * - _relative tolerance_: the absolute difference between the numbers must
300 * be smaller than the tolerance multiplied by the first number. Setting
301 * the relative tolerance to zero disables this check.
302 * - _ULP tolerance_: ULP (units of least precision) difference between the
303 * values must be smaller than the given tolerance for the check to pass.
304 * Setting the ULP tolerance to zero requires exact match.
305 * Setting the ULP tolerance to UINT64_MAX disables the ULP check.
306 * `0.0` and `-0.0` are treated as equal for the ULP check.
307 * - _sign check_: if set, any values that are of different signs fail the
308 * check (note that this also applies to `0.0` and `-0.0`: a value with a
309 * different sign than the zero will fail the check).
311 * Either an absolute, relative, or ULP tolerance must always be specified.
312 * If several of them are specified, then the check passes if either of the
313 * tolerances is satisfied.
315 * Any combination of absolute, relative, and ULP tolerance can be combined with
316 * the sign check. In this case, the sign check must succeed for the check to
317 * pass, even if other tolerances are satisfied.
319 * The tolerances can be specified separately for single and double precision
320 * comparison. Different initialization functions have different semantics on
321 * how the provided tolerance values are interpreted; check their
324 * Methods in this class do not throw, except for toString(), which may throw
328 * The factory methods that take ULP difference could be better formulated as
329 * methods that take the acceptable number of incorrect bits and/or the number
332 * \see FloatingPointDifference
334 class FloatingPointTolerance
338 * Creates a tolerance with the specified values.
340 * \param[in] singleAbsoluteTolerance
341 * Allowed absolute difference in a single-precision number.
342 * \param[in] doubleAbsoluteTolerance
343 * Allowed absolute difference in a double-precision number.
344 * \param[in] singleRelativeTolerance
345 * Allowed relative difference in a single-precision number.
346 * \param[in] doubleRelativeTolerance
347 * Allowed relative difference in a double-precision number.
348 * \param[in] singleUlpTolerance
349 * Allowed ULP difference in a single-precision number.
350 * \param[in] doubleUlpTolerance
351 * Allowed ULP difference in a double-precision number.
352 * \param[in] bSignMustMatch
353 * Whether sign mismatch fails the comparison.
355 FloatingPointTolerance(float singleAbsoluteTolerance
,
356 double doubleAbsoluteTolerance
,
357 float singleRelativeTolerance
,
358 double doubleRelativeTolerance
,
359 uint64_t singleUlpTolerance
,
360 uint64_t doubleUlpTolerance
,
362 : singleAbsoluteTolerance_(singleAbsoluteTolerance
),
363 doubleAbsoluteTolerance_(doubleAbsoluteTolerance
),
364 singleRelativeTolerance_(singleRelativeTolerance
),
365 doubleRelativeTolerance_(doubleRelativeTolerance
),
366 singleUlpTolerance_(singleUlpTolerance
),
367 doubleUlpTolerance_(doubleUlpTolerance
),
368 bSignMustMatch_(bSignMustMatch
)
373 * Checks whether a difference is within the specified tolerance.
375 * NaNs are always treated outside the tolerance.
377 bool isWithin(const FloatingPointDifference
&difference
) const;
379 //! Formats the tolerance as a string for assertion failure messages.
380 std::string
toString(const FloatingPointDifference
&difference
) const;
383 float singleAbsoluteTolerance_
;
384 double doubleAbsoluteTolerance_
;
385 float singleRelativeTolerance_
;
386 double doubleRelativeTolerance_
;
387 uint64_t singleUlpTolerance_
;
388 uint64_t doubleUlpTolerance_
;
389 bool bSignMustMatch_
;
393 * Creates a tolerance that only allows a specified ULP difference.
395 * The tolerance uses the given ULP value for both precisions, i.e., double
396 * precision will have much stricter tolerance.
398 * \related FloatingPointTolerance
400 static inline FloatingPointTolerance
401 ulpTolerance(uint64_t ulpDiff
)
403 return FloatingPointTolerance(0.0, 0.0, 0.0, 0.0, ulpDiff
, ulpDiff
, false);
407 * Creates a tolerance that allows a difference in two compared values that is
408 * relative to the given magnitude.
410 * \param[in] magnitude Magnitude of the numbers the computation operates in.
411 * \param[in] tolerance Relative tolerance permitted (e.g. 1e-4).
413 * In addition to setting an relative tolerance for both
414 * precisions, this sets the absolute tolerance such that values close to zero
415 * (in general, smaller than \p magnitude) do not fail the check if they
416 * differ by less than \p tolerance evaluated at \p magnitude. This accounts
417 * for potential loss of precision for small values, and should be used when
418 * accuracy of values much less than \p magnitude do not matter for
421 * \related FloatingPointTolerance
423 FloatingPointTolerance
424 relativeToleranceAsFloatingPoint(double magnitude
, double tolerance
);
427 * Creates a tolerance that allows a precision-dependent relative difference in
428 * a complex computation.
430 * \param[in] magnitude Magnitude of the numbers the computation operates in.
431 * \param[in] singleUlpDiff Expected accuracy of single-precision
432 * computation (in ULPs).
433 * \param[in] doubleUlpDiff Expected accuracy of double-precision
434 * computation (in ULPs).
436 * This works as relativeToleranceAsUlp(), but allows setting the ULP
437 * difference separately for the different precisions. This supports
438 * cases where the double-precision calculation can acceptably has a higher ULP
439 * difference, but relaxing the single-precision tolerance would lead to an
440 * unnecessarily loose test.
442 * \related FloatingPointTolerance
444 static inline FloatingPointTolerance
445 relativeToleranceAsPrecisionDependentUlp(double magnitude
,
446 uint64_t singleUlpDiff
,
447 uint64_t doubleUlpDiff
)
449 return FloatingPointTolerance(magnitude
*singleUlpDiff
*GMX_FLOAT_EPS
,
450 magnitude
*doubleUlpDiff
*GMX_DOUBLE_EPS
,
452 singleUlpDiff
, doubleUlpDiff
, false);
456 * Creates a tolerance that allows a specified absolute difference.
458 * \related FloatingPointTolerance
460 static inline FloatingPointTolerance
461 absoluteTolerance(double tolerance
)
463 return FloatingPointTolerance(tolerance
, tolerance
, 0.0, 0.0,
464 UINT64_MAX
, UINT64_MAX
, false);
468 * Creates a tolerance that allows a relative difference in a complex
471 * \param[in] magnitude Magnitude of the numbers the computation operates in.
472 * \param[in] ulpDiff Expected accuracy of the computation (in ULPs).
474 * In addition to setting the ULP tolerance as ulpTolerance(), this sets the
475 * absolute tolerance such that values close to zero (in general, smaller than
476 * \p magnitude) do not fail the check if they differ by less than \p ulpDiff
477 * evaluated at \p magnitude. This accounts for potential loss of precision
478 * for small values, and should be used when accuracy of values much less than
479 * \p magnitude do not matter for correctness.
481 * \related FloatingPointTolerance
483 static inline FloatingPointTolerance
484 relativeToleranceAsUlp(double magnitude
, uint64_t ulpDiff
)
486 return relativeToleranceAsPrecisionDependentUlp(magnitude
, ulpDiff
, ulpDiff
);
491 //! Default tolerance in ULPs for two floating-point values to compare equal.
492 constexpr uint64_t g_defaultUlpTolerance
= 4;
496 * Returns the default tolerance for comparing `real` numbers.
498 * \related FloatingPointTolerance
500 static inline FloatingPointTolerance
defaultRealTolerance()
502 return relativeToleranceAsUlp(1.0, detail::g_defaultUlpTolerance
);
507 * Returns the default tolerance for comparing single-precision numbers when
508 * compared by \Gromacs built in either precision mode.
510 * This permits a checker compiled with any \Gromacs precision to compare
511 * equal or not in the same way.
513 * \related FloatingPointTolerance
515 static inline FloatingPointTolerance
defaultFloatTolerance()
517 return relativeToleranceAsPrecisionDependentUlp
518 (1.0, detail::g_defaultUlpTolerance
,
519 static_cast<uint64_t>(detail::g_defaultUlpTolerance
* (GMX_FLOAT_EPS
/ GMX_DOUBLE_EPS
)));
522 /*! \name Assertions for floating-point comparison
524 * These routines extend `(EXPECT|ASSERT)_(FLOAT|DOUBLE)_EQ` and
525 * `(EXPECT|ASSERT)_NEAR` from Google Test to provide more flexible assertions
526 * for floating-point values.
528 * See gmx::test::FloatingPointTolerance for the possible ways to specify the
529 * tolerance, and gmx::test::FloatingPointDifference for some additional
530 * details of the difference calculation.
536 * Assertion predicate formatter for comparing two floating-point values.
538 template <typename FloatType
>
539 static inline ::testing::AssertionResult
assertEqualWithinTolerance(
540 const char *expr1
, const char *expr2
, const char * /*exprTolerance*/,
541 FloatType value1
, FloatType value2
,
542 const FloatingPointTolerance
&tolerance
)
544 FloatingPointDifference
diff(value1
, value2
);
545 if (tolerance
.isWithin(diff
))
547 return ::testing::AssertionSuccess();
549 return ::testing::AssertionFailure()
550 << " Value of: " << expr2
<< std::endl
551 << " Actual: " << value2
<< std::endl
552 << " Expected: " << expr1
<< std::endl
553 << " Which is: " << value1
<< std::endl
554 << "Difference: " << diff
.toString() << std::endl
555 << " Tolerance: " << tolerance
.toString(diff
);
560 * Asserts that two single-precision values are within the given tolerance.
564 #define EXPECT_FLOAT_EQ_TOL(value1, value2, tolerance) \
565 EXPECT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<float>, \
566 value1, value2, tolerance)
568 * Asserts that two double-precision values are within the given tolerance.
572 #define EXPECT_DOUBLE_EQ_TOL(value1, value2, tolerance) \
573 EXPECT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<double>, \
574 value1, value2, tolerance)
575 /*! \def EXPECT_REAL_EQ_TOL
577 * Asserts that two `real` values are within the given tolerance.
582 * Asserts that two single-precision values are within the given tolerance.
586 #define ASSERT_FLOAT_EQ_TOL(value1, value2, tolerance) \
587 ASSERT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<float>, \
588 value1, value2, tolerance)
590 * Asserts that two double-precision values are within the given tolerance.
594 #define ASSERT_DOUBLE_EQ_TOL(value1, value2, tolerance) \
595 ASSERT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<double>, \
596 value1, value2, tolerance)
597 /*! \def ASSERT_REAL_EQ_TOL
599 * Asserts that two `real` values are within the given tolerance.
605 #define EXPECT_REAL_EQ_TOL(value1, value2, tolerance) \
606 EXPECT_DOUBLE_EQ_TOL(value1, value2, tolerance)
607 #define ASSERT_REAL_EQ_TOL(value1, value2, tolerance) \
608 ASSERT_DOUBLE_EQ_TOL(value1, value2, tolerance)
610 #define EXPECT_REAL_EQ_TOL(value1, value2, tolerance) \
611 EXPECT_FLOAT_EQ_TOL(value1, value2, tolerance)
612 #define ASSERT_REAL_EQ_TOL(value1, value2, tolerance) \
613 ASSERT_FLOAT_EQ_TOL(value1, value2, tolerance)
616 //! EXPECT_REAL_EQ_TOL with default tolerance
617 #define EXPECT_REAL_EQ(value1, value2) EXPECT_REAL_EQ_TOL(value1, value2, ::gmx::test::defaultRealTolerance())
618 //! ASSERT_REAL_EQ_TOL with default tolerance
619 #define ASSERT_REAL_EQ(value1, value2) ASSERT_REAL_EQ_TOL(value1, value2, ::gmx::test::defaultRealTolerance())
624 * Helper method for `(EXPECT|ASSERT)_PLAIN`.
626 static inline ::testing::AssertionResult
627 plainAssertHelper(const char * /*expr*/, const ::testing::AssertionResult
&expr
)
634 * Assert for predicates that return AssertionResult and produce a full failure
637 * `expr` should evaluate to AssertionResult, and on failure the message from
638 * the result is used as-is, unlike in EXPECT_TRUE().
642 #define EXPECT_PLAIN(expr) EXPECT_PRED_FORMAT1(plainAssertHelper, expr)
644 * Assert for predicates that return AssertionResult and produce a full failure
647 * `expr` should evaluate to AssertionResult, and on failure the message from
648 * the result is used as-is, unlike in ASSERT_TRUE().
652 #define ASSERT_PLAIN(expr) ASSERT_PRED_FORMAT1(plainAssertHelper, expr)