2 * This file is part of the GROMACS molecular simulation package.
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.
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"
71 //! \libinternal \addtogroup module_testutils
74 /*! \name Assertions for exceptions
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.
86 * Internal implementation macro for exception assertations.
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.
92 * The implementation is copied and adjusted from
93 * include/gtest/internal/gtest-internal.h in Google Test 1.6.0.
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; \
100 GTEST_SUPPRESS_UNREACHABLE_CODE_WARNING_BELOW_(statement); \
102 catch (expected_exception const &) { \
103 gmx_caught_expected = true; \
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__); \
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__); \
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__); \
122 GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__) : \
123 fail(gmx_ar.message())
126 * Internal implementation macro for exception assertations.
128 * \param statement Statements to execute.
129 * \param fail Function/macro to call on failure.
131 * The implementation is copied and adjusted from
132 * include/gtest/internal/gtest-internal.h in Google Test 1.6.0.
134 #define GMX_TEST_NO_THROW_(statement, fail) \
135 GTEST_AMBIGUOUS_ELSE_BLOCKER_ \
136 if (::testing::AssertionResult gmx_ar = ::testing::AssertionSuccess()) { \
138 GTEST_SUPPRESS_UNREACHABLE_CODE_WARNING_BELOW_(statement); \
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__); \
147 gmx_ar << "Expected: " #statement " doesn't throw an exception.\n" \
148 << " Actual: it throws."; \
149 goto GTEST_CONCAT_TOKEN_(gmx_label_testnothrow_, __LINE__); \
152 GTEST_CONCAT_TOKEN_(gmx_label_testnothrow_, __LINE__) : \
153 fail(gmx_ar.message())
157 * Asserts that a statement throws a given exception.
161 #define EXPECT_THROW_GMX(statement, expected_exception) \
162 GMX_TEST_THROW_(statement, expected_exception, GTEST_NONFATAL_FAILURE_)
164 * Asserts that a statement does not throw.
168 #define EXPECT_NO_THROW_GMX(statement) \
169 GMX_TEST_NO_THROW_(statement, GTEST_NONFATAL_FAILURE_)
171 * Asserts that a statement throws a given exception.
175 #define ASSERT_THROW_GMX(statement, expected_exception) \
176 GMX_TEST_THROW_(statement, expected_exception, GTEST_FATAL_FAILURE_)
178 * Asserts that a statement does not throw.
182 #define ASSERT_NO_THROW_GMX(statement) \
183 GMX_TEST_NO_THROW_(statement, GTEST_FATAL_FAILURE_)
187 /*! \libinternal \brief
188 * Computes and represents a floating-point difference value.
190 * Methods in this class do not throw, except for toString(), which may throw
193 * \see FloatingPointTolerance
195 class FloatingPointDifference
198 //! Initializes a single-precision difference.
199 FloatingPointDifference(float value1
, float value2
);
200 //! Initializes a double-precision difference.
201 FloatingPointDifference(double value1
, double value2
);
204 * Whether one or both of the compared values were NaN.
206 * If this returns `true`, other accessors return meaningless values.
209 //! Returns the difference as an absolute number (always non-negative).
210 double asAbsolute() const { return absoluteDifference_
; }
212 * Returns the difference as ULPs (always non-negative).
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.
218 gmx_uint64_t
asUlps() const { return ulpDifference_
; }
220 * Whether the compared values were of different sign.
222 * 0.0 and -0.0 are treated as positive and negative, respectively.
224 bool signsDiffer() const { return bSignDifference_
; }
226 * Whether the difference is between single- or double-precision
229 bool isDouble() const { return bDouble_
; }
230 //! Formats the difference as a string for assertion failure messages.
231 std::string
toString() const;
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_
;
239 * Whether the difference was computed for single or double precision.
241 * This sets the units for `ulpDifference_`.
246 /*! \libinternal \brief
247 * Specifies a floating-point comparison tolerance and checks whether a
248 * difference is within the tolerance.
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.
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
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).
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
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.
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
281 * Methods in this class do not throw, except for toString(), which may throw
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
289 * \see FloatingPointDifference
291 class FloatingPointTolerance
295 * Creates a tolerance with the specified values.
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.
308 FloatingPointTolerance(float singleAbsoluteTolerance
,
309 double doubleAbsoluteTolerance
,
310 gmx_uint64_t singleUlpTolerance
,
311 gmx_uint64_t doubleUlpTolerance
,
313 : singleAbsoluteTolerance_(singleAbsoluteTolerance
),
314 doubleAbsoluteTolerance_(doubleAbsoluteTolerance
),
315 singleUlpTolerance_(singleUlpTolerance
),
316 doubleUlpTolerance_(doubleUlpTolerance
),
317 bSignMustMatch_(bSignMustMatch
)
322 * Checks whether a difference is within the specified tolerance.
324 * NaNs are always treated outside the tolerance.
326 bool isWithin(const FloatingPointDifference
&difference
) const;
328 //! Formats the tolerance as a string for assertion failure messages.
329 std::string
toString(const FloatingPointDifference
&difference
) const;
332 float singleAbsoluteTolerance_
;
333 double doubleAbsoluteTolerance_
;
334 gmx_uint64_t singleUlpTolerance_
;
335 gmx_uint64_t doubleUlpTolerance_
;
336 bool bSignMustMatch_
;
340 * Creates a tolerance that only allows a specified ULP difference.
342 * The tolerance uses the given ULP value for both precisions, i.e., double
343 * precision will have much stricter tolerance.
345 * \related FloatingPointTolerance
347 static inline FloatingPointTolerance
348 ulpTolerance(gmx_uint64_t ulpDiff
)
350 return FloatingPointTolerance(0.0, 0.0, ulpDiff
, ulpDiff
, false);
354 * Creates a tolerance that allows a difference in two compared values that is
355 * relative to the given magnitude.
357 * \param[in] magnitude Magnitude of the numbers the computation operates in.
358 * \param[in] tolerance Relative tolerance permitted (e.g. 1e-4).
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
368 * The ULP tolerance for different precisions will be different to make them
369 * both match \p tolerance.
371 * \related FloatingPointTolerance
373 FloatingPointTolerance
374 relativeToleranceAsFloatingPoint(double magnitude
, double tolerance
);
377 * Creates a tolerance that allows a precision-dependent relative difference in
378 * a complex computation.
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).
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.
392 * \related FloatingPointTolerance
394 static inline FloatingPointTolerance
395 relativeToleranceAsPrecisionDependentUlp(double magnitude
,
396 gmx_uint64_t singleUlpDiff
,
397 gmx_uint64_t doubleUlpDiff
)
399 return FloatingPointTolerance(magnitude
*singleUlpDiff
*GMX_FLOAT_EPS
,
400 magnitude
*doubleUlpDiff
*GMX_DOUBLE_EPS
,
401 singleUlpDiff
, doubleUlpDiff
, false);
405 * Creates a tolerance that allows a specified absolute difference.
407 * \related FloatingPointTolerance
409 static inline FloatingPointTolerance
410 absoluteTolerance(double tolerance
)
412 return FloatingPointTolerance(tolerance
, tolerance
,
413 GMX_UINT64_MAX
, GMX_UINT64_MAX
, false);
417 * Creates a tolerance that allows a relative difference in a complex
420 * \param[in] magnitude Magnitude of the numbers the computation operates in.
421 * \param[in] ulpDiff Expected accuracy of the computation (in ULPs).
423 * In addition to setting the ULP tolerance as ulpTolerance(), this sets the
424 * absolute tolerance such that values close to zero (in general, smaller than
425 * \p magnitude) do not fail the check if they differ by less than \p ulpDiff
426 * evaluated at \p magnitude. This accounts for potential loss of precision
427 * for small values, and should be used when accuracy of values much less than
428 * \p magnitude do not matter for correctness.
430 * \related FloatingPointTolerance
432 static inline FloatingPointTolerance
433 relativeToleranceAsUlp(double magnitude
, gmx_uint64_t ulpDiff
)
435 return relativeToleranceAsPrecisionDependentUlp(magnitude
, ulpDiff
, ulpDiff
);
439 * Returns the default tolerance for comparing `real` numbers.
441 * \related FloatingPointTolerance
443 static inline FloatingPointTolerance
defaultRealTolerance()
445 return relativeToleranceAsUlp(1.0, 4);
448 /*! \name Assertions for floating-point comparison
450 * These routines extend `(EXPECT|ASSERT)_(FLOAT|DOUBLE)_EQ` and
451 * `(EXPECT|ASSERT)_NEAR` from Google Test to provide more flexible assertions
452 * for floating-point values.
454 * See gmx::test::FloatingPointTolerance for the possible ways to specify the
455 * tolerance, and gmx::test::FloatingPointDifference for some additional
456 * details of the difference calculation.
462 * Assertion predicate formatter for comparing two floating-point values.
464 template <typename FloatType
>
465 static inline ::testing::AssertionResult
assertEqualWithinTolerance(
466 const char *expr1
, const char *expr2
, const char * /*exprTolerance*/,
467 FloatType value1
, FloatType value2
,
468 const FloatingPointTolerance
&tolerance
)
470 FloatingPointDifference
diff(value1
, value2
);
471 if (tolerance
.isWithin(diff
))
473 return ::testing::AssertionSuccess();
475 return ::testing::AssertionFailure()
476 << " Value of: " << expr2
<< std::endl
477 << " Actual: " << value2
<< std::endl
478 << " Expected: " << expr1
<< std::endl
479 << " Which is: " << value1
<< std::endl
480 << "Difference: " << diff
.toString() << std::endl
481 << " Tolerance: " << tolerance
.toString(diff
);
486 * Asserts that two single-precision values are within the given tolerance.
490 #define EXPECT_FLOAT_EQ_TOL(value1, value2, tolerance) \
491 EXPECT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<float>, \
492 value1, value2, tolerance)
494 * Asserts that two double-precision values are within the given tolerance.
498 #define EXPECT_DOUBLE_EQ_TOL(value1, value2, tolerance) \
499 EXPECT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<double>, \
500 value1, value2, tolerance)
501 /*! \def EXPECT_REAL_EQ_TOL
503 * Asserts that two `real` values are within the given tolerance.
508 * Asserts that two single-precision values are within the given tolerance.
512 #define ASSERT_FLOAT_EQ_TOL(value1, value2, tolerance) \
513 ASSERT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<float>, \
514 value1, value2, tolerance)
516 * Asserts that two double-precision values are within the given tolerance.
520 #define ASSERT_DOUBLE_EQ_TOL(value1, value2, tolerance) \
521 ASSERT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<double>, \
522 value1, value2, tolerance)
523 /*! \def ASSERT_REAL_EQ_TOL
525 * Asserts that two `real` values are within the given tolerance.
531 #define EXPECT_REAL_EQ_TOL(value1, value2, tolerance) \
532 EXPECT_DOUBLE_EQ_TOL(value1, value2, tolerance)
533 #define ASSERT_REAL_EQ_TOL(value1, value2, tolerance) \
534 ASSERT_DOUBLE_EQ_TOL(value1, value2, tolerance)
536 #define EXPECT_REAL_EQ_TOL(value1, value2, tolerance) \
537 EXPECT_FLOAT_EQ_TOL(value1, value2, tolerance)
538 #define ASSERT_REAL_EQ_TOL(value1, value2, tolerance) \
539 ASSERT_FLOAT_EQ_TOL(value1, value2, tolerance)
544 /*! \name Assertions for NULL comparison
546 * These macros should be used instead of `(EXPECT|ASSERT)_EQ(NULL, ...)`,
547 * because Google Test doesn't support the NULL comparison with xlC++ 12.1 on
553 * Asserts that a pointer is null.
555 * Works exactly like EXPECT_EQ comparing with a null pointer. */
556 #define EXPECT_NULL(val) EXPECT_EQ((void *) NULL, val)
558 * Asserts that a pointer is null.
560 * Works exactly like ASSERT_EQ comparing with a null pointer. */
561 #define ASSERT_NULL(val) ASSERT_EQ((void *) NULL, val)