Apply clang-format to source tree
[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,2015,2016,2017,2018,2019, 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 namespace internal
72 {
73 //! \cond internal
74 /*! \internal
75  * \brief
76  * Called for an expected exception from EXPECT_THROW_GMX().
77  *
78  * \param[in] ex  Exception that was thrown.
79  */
80 void processExpectedException(const std::exception& ex);
81 //! \endcond
82 } // namespace internal
83
84 //! \libinternal \addtogroup module_testutils
85 //! \{
86
87 /*! \name Assertions for exceptions
88  *
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.
94  */
95 //! \{
96
97 //! \cond internal
98 /*! \brief
99  * Internal implementation macro for exception assertations.
100  *
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.
104  *
105  * The implementation is copied and adjusted from
106  * include/gtest/internal/gtest-internal.h in Google Test 1.6.0.
107  */
108 #define GMX_TEST_THROW_(statement, expected_exception, fail)                           \
109     GTEST_AMBIGUOUS_ELSE_BLOCKER_                                                      \
110     if (::testing::AssertionResult gmx_ar = ::testing::AssertionSuccess())             \
111     {                                                                                  \
112         bool gmx_caught_expected = false;                                              \
113         try                                                                            \
114         {                                                                              \
115             GTEST_SUPPRESS_UNREACHABLE_CODE_WARNING_BELOW_(statement);                 \
116         }                                                                              \
117         catch (expected_exception const& ex)                                           \
118         {                                                                              \
119             gmx_caught_expected = true;                                                \
120             ::gmx::test::internal::processExpectedException(ex);                       \
121         }                                                                              \
122         catch (std::exception const& ex)                                               \
123         {                                                                              \
124             gmx_ar << "Expected: " #statement " throws an exception of type "          \
125                    << #expected_exception ".\n  Actual: it throws a different type.\n" \
126                    << "Exception details:\n"                                           \
127                    << ::gmx::formatExceptionMessageToString(ex);                       \
128             goto GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__);                  \
129         }                                                                              \
130         catch (...)                                                                    \
131         {                                                                              \
132             gmx_ar << "Expected: " #statement " throws an exception of type "          \
133                    << #expected_exception ".\n  Actual: it throws a different type.";  \
134             goto GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__);                  \
135         }                                                                              \
136         if (!gmx_caught_expected)                                                      \
137         {                                                                              \
138             gmx_ar << "Expected: " #statement " throws an exception of type "          \
139                    << #expected_exception ".\n  Actual: it throws nothing.";           \
140             goto GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__);                  \
141         }                                                                              \
142     }                                                                                  \
143     else                                                                               \
144         GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__) : fail(gmx_ar.message())
145
146 /*! \brief
147  * Internal implementation macro for exception assertations.
148  *
149  * \param statement          Statements to execute.
150  * \param fail               Function/macro to call on failure.
151  *
152  * The implementation is copied and adjusted from
153  * include/gtest/internal/gtest-internal.h in Google Test 1.6.0.
154  */
155 #define GMX_TEST_NO_THROW_(statement, fail)                                    \
156     GTEST_AMBIGUOUS_ELSE_BLOCKER_                                              \
157     if (::testing::AssertionResult gmx_ar = ::testing::AssertionSuccess())     \
158     {                                                                          \
159         try                                                                    \
160         {                                                                      \
161             GTEST_SUPPRESS_UNREACHABLE_CODE_WARNING_BELOW_(statement);         \
162         }                                                                      \
163         catch (std::exception const& ex)                                       \
164         {                                                                      \
165             gmx_ar << "Expected: " #statement " doesn't throw an exception.\n" \
166                    << "  Actual: it throws.\n"                                 \
167                    << "Exception details:\n"                                   \
168                    << ::gmx::formatExceptionMessageToString(ex);               \
169             goto GTEST_CONCAT_TOKEN_(gmx_label_testnothrow_, __LINE__);        \
170         }                                                                      \
171         catch (...)                                                            \
172         {                                                                      \
173             gmx_ar << "Expected: " #statement " doesn't throw an exception.\n" \
174                    << "  Actual: it throws.";                                  \
175             goto GTEST_CONCAT_TOKEN_(gmx_label_testnothrow_, __LINE__);        \
176         }                                                                      \
177     }                                                                          \
178     else                                                                       \
179         GTEST_CONCAT_TOKEN_(gmx_label_testnothrow_, __LINE__) : fail(gmx_ar.message())
180 //! \endcond
181
182 /*! \brief
183  * Asserts that a statement throws a given exception.
184  *
185  * \hideinitializer
186  */
187 #define EXPECT_THROW_GMX(statement, expected_exception) \
188     GMX_TEST_THROW_(statement, expected_exception, GTEST_NONFATAL_FAILURE_)
189 /*! \brief
190  * Asserts that a statement does not throw.
191  *
192  * \hideinitializer
193  */
194 #define EXPECT_NO_THROW_GMX(statement) GMX_TEST_NO_THROW_(statement, GTEST_NONFATAL_FAILURE_)
195 /*! \brief
196  * Asserts that a statement throws a given exception.
197  *
198  * \hideinitializer
199  */
200 #define ASSERT_THROW_GMX(statement, expected_exception) \
201     GMX_TEST_THROW_(statement, expected_exception, GTEST_FATAL_FAILURE_)
202 /*! \brief
203  * Asserts that a statement does not throw.
204  *
205  * \hideinitializer
206  */
207 #define ASSERT_NO_THROW_GMX(statement) GMX_TEST_NO_THROW_(statement, GTEST_FATAL_FAILURE_)
208
209 //! \}
210
211 /*! \libinternal \brief
212  * Computes and represents a floating-point difference value.
213  *
214  * Methods in this class do not throw, except for toString(), which may throw
215  * std::bad_alloc.
216  *
217  * \see FloatingPointTolerance
218  */
219 class FloatingPointDifference
220 {
221 public:
222     /*! \brief Initializes a single-precision difference.
223      *
224      *  \param ref    First term in difference
225      *  \param value  Second term in difference
226      *
227      *  For absolute and ULP differences the two parameters are equivalent,
228      *  since the difference is symmetric. For relative differences
229      *  the first term is interpreted as the reference value, from which
230      *  we extract the magnitude to compare with.
231      */
232     FloatingPointDifference(float ref, float value);
233
234     /*! \brief Initializes a double-precision difference.
235      *
236      *  \param ref    First term in difference
237      *  \param value  Second term in difference
238      *
239      *  For absolute and ULP differences the two parameters are equivalent,
240      *  since the difference is symmetric. For relative differences
241      *  the first term is interpreted as the reference value, from which
242      *  we extract the magnitude to compare with.
243      */
244     FloatingPointDifference(double ref, double value);
245
246     /*! \brief
247      * Whether one or both of the compared values were NaN.
248      *
249      * If this returns `true`, other accessors return meaningless values.
250      */
251     bool isNaN() const;
252     //! Returns the difference as an absolute number (always non-negative).
253     double asAbsolute() const { return absoluteDifference_; }
254     /*! \brief
255      * Returns the difference as ULPs (always non-negative).
256      *
257      * The ULPs are calculated for the type that corresponds to the
258      * constructor used to initialize the difference.
259      * The ULP difference between 0.0 and -0.0 is zero.
260      */
261     uint64_t asUlps() const { return ulpDifference_; }
262     /*! \brief
263      * Whether the compared values were of different sign.
264      *
265      * 0.0 and -0.0 are treated as positive and negative, respectively.
266      */
267     bool signsDiffer() const { return bSignDifference_; }
268     /*! \brief
269      * Whether the difference is between single- or double-precision
270      * numbers.
271      */
272     bool isDouble() const { return bDouble_; }
273     //! Formats the difference as a string for assertion failure messages.
274     std::string toString() const;
275
276     //! Returns the magnitude of the original second term of the difference.
277     double termMagnitude() const { return termMagnitude_; }
278
279 private:
280     //! Save the magnitude of the reference value for relative (i.e., not ULP) tolerance
281     double termMagnitude_;
282     //! Stores the absolute difference, or NaN if one or both values were NaN.
283     double   absoluteDifference_;
284     uint64_t ulpDifference_;
285     bool     bSignDifference_;
286     /*! \brief
287      * Whether the difference was computed for single or double precision.
288      *
289      * This sets the units for `ulpDifference_`.
290      */
291     bool bDouble_;
292 };
293
294 /*! \libinternal \brief
295  * Specifies a floating-point comparison tolerance and checks whether a
296  * difference is within the tolerance.
297  *
298  * The related functions section lists methods that can be construct methods
299  * using less parameters than the full constructor, and with more obvious
300  * semantics.  These should be preferred over using the constructor directly.
301  *
302  * Several types of tolerances are possible:
303  *  - _absolute tolerance_: difference between the values must be smaller than
304  *    the given tolerance for the check to pass.
305  *    Setting the absolute tolerance to zero disables the absolute tolerance
306  *    check.
307  *  - _relative tolerance_: the absolute difference between the numbers must
308  *    be smaller than the tolerance multiplied by the first number. Setting
309  *    the relative tolerance to zero disables this check.
310  *  - _ULP tolerance_: ULP (units of least precision) difference between the
311  *    values must be smaller than the given tolerance for the check to pass.
312  *    Setting the ULP tolerance to zero requires exact match.
313  *    Setting the ULP tolerance to UINT64_MAX disables the ULP check.
314  *    `0.0` and `-0.0` are treated as equal for the ULP check.
315  *  - _sign check_: if set, any values that are of different signs fail the
316  *    check (note that this also applies to `0.0` and `-0.0`: a value with a
317  *    different sign than the zero will fail the check).
318  *
319  * Either an absolute, relative, or ULP tolerance must always be specified.
320  * If several of them are specified, then the check passes if either of the
321  * tolerances is satisfied.
322  *
323  * Any combination of absolute, relative, and ULP tolerance can be combined with
324  * the sign check.  In this case, the sign check must succeed for the check to
325  * pass, even if other tolerances are satisfied.
326  *
327  * The tolerances can be specified separately for single and double precision
328  * comparison.  Different initialization functions have different semantics on
329  * how the provided tolerance values are interpreted; check their
330  * documentation.
331  *
332  * Methods in this class do not throw, except for toString(), which may throw
333  * std::bad_alloc.
334  *
335  * \todo
336  * The factory methods that take ULP difference could be better formulated as
337  * methods that take the acceptable number of incorrect bits and/or the number
338  * of accurate bits.
339  *
340  * \see FloatingPointDifference
341  */
342 class FloatingPointTolerance
343 {
344 public:
345     /*! \brief
346      * Creates a tolerance with the specified values.
347      *
348      * \param[in]  singleAbsoluteTolerance
349      *     Allowed absolute difference in a single-precision number.
350      * \param[in]  doubleAbsoluteTolerance
351      *     Allowed absolute difference in a double-precision number.
352      * \param[in]  singleRelativeTolerance
353      *     Allowed relative difference in a single-precision number.
354      * \param[in]  doubleRelativeTolerance
355      *     Allowed relative difference in a double-precision number.
356      * \param[in]  singleUlpTolerance
357      *     Allowed ULP difference in a single-precision number.
358      * \param[in]  doubleUlpTolerance
359      *     Allowed ULP difference in a double-precision number.
360      * \param[in]  bSignMustMatch
361      *     Whether sign mismatch fails the comparison.
362      */
363     FloatingPointTolerance(float    singleAbsoluteTolerance,
364                            double   doubleAbsoluteTolerance,
365                            float    singleRelativeTolerance,
366                            double   doubleRelativeTolerance,
367                            uint64_t singleUlpTolerance,
368                            uint64_t doubleUlpTolerance,
369                            bool     bSignMustMatch) :
370         singleAbsoluteTolerance_(singleAbsoluteTolerance),
371         doubleAbsoluteTolerance_(doubleAbsoluteTolerance),
372         singleRelativeTolerance_(singleRelativeTolerance),
373         doubleRelativeTolerance_(doubleRelativeTolerance),
374         singleUlpTolerance_(singleUlpTolerance),
375         doubleUlpTolerance_(doubleUlpTolerance),
376         bSignMustMatch_(bSignMustMatch)
377     {
378     }
379
380     /*! \brief
381      * Checks whether a difference is within the specified tolerance.
382      *
383      * NaNs are always treated outside the tolerance.
384      */
385     bool isWithin(const FloatingPointDifference& difference) const;
386
387     //! Formats the tolerance as a string for assertion failure messages.
388     std::string toString(const FloatingPointDifference& difference) const;
389
390 private:
391     float    singleAbsoluteTolerance_;
392     double   doubleAbsoluteTolerance_;
393     float    singleRelativeTolerance_;
394     double   doubleRelativeTolerance_;
395     uint64_t singleUlpTolerance_;
396     uint64_t doubleUlpTolerance_;
397     bool     bSignMustMatch_;
398 };
399
400 /*! \brief
401  * Creates a tolerance that only allows a specified ULP difference.
402  *
403  * The tolerance uses the given ULP value for both precisions, i.e., double
404  * precision will have much stricter tolerance.
405  *
406  * \related FloatingPointTolerance
407  */
408 static inline FloatingPointTolerance ulpTolerance(uint64_t ulpDiff)
409 {
410     return FloatingPointTolerance(0.0, 0.0, 0.0, 0.0, ulpDiff, ulpDiff, false);
411 }
412
413 /*! \brief
414  * Creates a tolerance that allows a difference in two compared values that is
415  * relative to the given magnitude.
416  *
417  * \param[in] magnitude  Magnitude of the numbers the computation operates in.
418  * \param[in] tolerance  Relative tolerance permitted (e.g. 1e-4).
419  *
420  * In addition to setting an relative tolerance for both
421  * precisions, this sets the absolute tolerance such that values close to zero
422  * (in general, smaller than \p magnitude) do not fail the check if they
423  * differ by less than \p tolerance evaluated at \p magnitude.  This accounts
424  * for potential loss of precision for small values, and should be used when
425  * accuracy of values much less than \p magnitude do not matter for
426  * correctness.
427  *
428  * \related FloatingPointTolerance
429  */
430 FloatingPointTolerance relativeToleranceAsFloatingPoint(double magnitude, double tolerance);
431
432 /*! \brief
433  * Creates a tolerance that allows a precision-dependent difference in two
434  * compared values that is relative to the given magnitude.
435  *
436  * \param[in] magnitude        Magnitude of the numbers the computation
437  *     operates in.
438  * \param[in] singleTolerance  Relative tolerance permitted (e.g. 1e-4)
439  *     in single precision.
440  * \param[in] doubleTolerance  Relative tolerance permitted (e.g. 1e-4)
441  *     in double precision.
442  *
443  * In addition to setting an relative tolerance for both
444  * precisions, this sets the absolute tolerance such that values close to zero
445  * (in general, smaller than \p magnitude) do not fail the check if they
446  * differ by less than \p tolerance evaluated at \p magnitude.  This accounts
447  * for potential loss of precision for small values, and should be used when
448  * accuracy of values much less than \p magnitude do not matter for
449  * correctness.
450  *
451  * \related FloatingPointTolerance
452  */
453 FloatingPointTolerance relativeToleranceAsPrecisionDependentFloatingPoint(double magnitude,
454                                                                           float  singleTolerance,
455                                                                           double doubleTolerance);
456
457 /*! \brief
458  * Creates a tolerance that allows a precision-dependent relative difference in
459  * a complex computation.
460  *
461  * \param[in] magnitude      Magnitude of the numbers the computation operates in.
462  * \param[in] singleUlpDiff  Expected accuracy of single-precision
463  *     computation (in ULPs).
464  * \param[in] doubleUlpDiff  Expected accuracy of double-precision
465  *     computation (in ULPs).
466  *
467  * This works as relativeToleranceAsUlp(), but allows setting the ULP
468  * difference separately for the different precisions.  This supports
469  * cases where the double-precision calculation can acceptably has a higher ULP
470  * difference, but relaxing the single-precision tolerance would lead to an
471  * unnecessarily loose test.
472  *
473  * \related FloatingPointTolerance
474  */
475 static inline FloatingPointTolerance relativeToleranceAsPrecisionDependentUlp(double   magnitude,
476                                                                               uint64_t singleUlpDiff,
477                                                                               uint64_t doubleUlpDiff)
478 {
479     return FloatingPointTolerance(float(magnitude) * singleUlpDiff * GMX_FLOAT_EPS,
480                                   magnitude * doubleUlpDiff * GMX_DOUBLE_EPS, 0.0, 0.0,
481                                   singleUlpDiff, doubleUlpDiff, false);
482 }
483
484 /*! \brief
485  * Creates a tolerance that allows a specified absolute difference.
486  *
487  * \related FloatingPointTolerance
488  */
489 static inline FloatingPointTolerance absoluteTolerance(double tolerance)
490 {
491     return FloatingPointTolerance(float(tolerance), tolerance, 0.0, 0.0, UINT64_MAX, UINT64_MAX, false);
492 }
493
494 /*! \brief
495  * Creates a tolerance that allows a relative difference in a complex
496  * computation.
497  *
498  * \param[in] magnitude  Magnitude of the numbers the computation operates in.
499  * \param[in] ulpDiff    Expected accuracy of the computation (in ULPs).
500  *
501  * In addition to setting the ULP tolerance as ulpTolerance(), this sets the
502  * absolute tolerance such that values close to zero (in general, smaller than
503  * \p magnitude) do not fail the check if they differ by less than \p ulpDiff
504  * evaluated at \p magnitude.  This accounts for potential loss of precision
505  * for small values, and should be used when accuracy of values much less than
506  * \p magnitude do not matter for correctness.
507  *
508  * \related FloatingPointTolerance
509  */
510 static inline FloatingPointTolerance relativeToleranceAsUlp(double magnitude, uint64_t ulpDiff)
511 {
512     return relativeToleranceAsPrecisionDependentUlp(magnitude, ulpDiff, ulpDiff);
513 }
514
515 namespace detail
516 {
517 //! Default tolerance in ULPs for two floating-point values to compare equal.
518 constexpr uint64_t g_defaultUlpTolerance = 4;
519 } // namespace detail
520
521 /*! \brief
522  * Returns the default tolerance for comparing `real` numbers.
523  *
524  * \related FloatingPointTolerance
525  */
526 static inline FloatingPointTolerance defaultRealTolerance()
527 {
528     return relativeToleranceAsUlp(1.0, detail::g_defaultUlpTolerance);
529 }
530
531
532 /*! \brief
533  * Returns the default tolerance for comparing single-precision numbers when
534  * compared by \Gromacs built in either precision mode.
535  *
536  * This permits a checker compiled with any \Gromacs precision to compare
537  * equal or not in the same way.
538  *
539  * \related FloatingPointTolerance
540  */
541 static inline FloatingPointTolerance defaultFloatTolerance()
542 {
543     return relativeToleranceAsPrecisionDependentUlp(
544             1.0, detail::g_defaultUlpTolerance,
545             static_cast<uint64_t>(detail::g_defaultUlpTolerance * (GMX_FLOAT_EPS / GMX_DOUBLE_EPS)));
546 }
547
548 /*! \name Assertions for floating-point comparison
549  *
550  * These routines extend `(EXPECT|ASSERT)_(FLOAT|DOUBLE)_EQ` and
551  * `(EXPECT|ASSERT)_NEAR` from Google Test to provide more flexible assertions
552  * for floating-point values.
553  *
554  * See gmx::test::FloatingPointTolerance for the possible ways to specify the
555  * tolerance, and gmx::test::FloatingPointDifference for some additional
556  * details of the difference calculation.
557  */
558 //! \{
559
560 //! \cond internal
561 /*! \internal \brief
562  * Assertion predicate formatter for comparing two floating-point values.
563  */
564 template<typename FloatType>
565 static inline ::testing::AssertionResult assertEqualWithinTolerance(const char* expr1,
566                                                                     const char* expr2,
567                                                                     const char* /*exprTolerance*/,
568                                                                     FloatType value1,
569                                                                     FloatType value2,
570                                                                     const FloatingPointTolerance& tolerance)
571 {
572     FloatingPointDifference diff(value1, value2);
573     if (tolerance.isWithin(diff))
574     {
575         return ::testing::AssertionSuccess();
576     }
577     return ::testing::AssertionFailure() << "  Value of: " << expr2 << std::endl
578                                          << "    Actual: " << value2 << std::endl
579                                          << "  Expected: " << expr1 << std::endl
580                                          << "  Which is: " << value1 << std::endl
581                                          << "Difference: " << diff.toString() << std::endl
582                                          << " Tolerance: " << tolerance.toString(diff);
583 }
584 //! \endcond
585
586 /*! \brief
587  * Asserts that two single-precision values are within the given tolerance.
588  *
589  * \hideinitializer
590  */
591 #define EXPECT_FLOAT_EQ_TOL(value1, value2, tolerance) \
592     EXPECT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<float>, value1, value2, tolerance)
593 /*! \brief
594  * Asserts that two double-precision values are within the given tolerance.
595  *
596  * \hideinitializer
597  */
598 #define EXPECT_DOUBLE_EQ_TOL(value1, value2, tolerance) \
599     EXPECT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<double>, value1, value2, tolerance)
600 /*! \def EXPECT_REAL_EQ_TOL
601  * \brief
602  * Asserts that two `real` values are within the given tolerance.
603  *
604  * \hideinitializer
605  */
606 /*! \brief
607  * Asserts that two single-precision values are within the given tolerance.
608  *
609  * \hideinitializer
610  */
611 #define ASSERT_FLOAT_EQ_TOL(value1, value2, tolerance) \
612     ASSERT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<float>, value1, value2, tolerance)
613 /*! \brief
614  * Asserts that two double-precision values are within the given tolerance.
615  *
616  * \hideinitializer
617  */
618 #define ASSERT_DOUBLE_EQ_TOL(value1, value2, tolerance) \
619     ASSERT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<double>, value1, value2, tolerance)
620 /*! \def ASSERT_REAL_EQ_TOL
621  * \brief
622  * Asserts that two `real` values are within the given tolerance.
623  *
624  * \hideinitializer
625  */
626
627 #if GMX_DOUBLE
628 #    define EXPECT_REAL_EQ_TOL(value1, value2, tolerance) \
629         EXPECT_DOUBLE_EQ_TOL(value1, value2, tolerance)
630 #    define ASSERT_REAL_EQ_TOL(value1, value2, tolerance) \
631         ASSERT_DOUBLE_EQ_TOL(value1, value2, tolerance)
632 #else
633 #    define EXPECT_REAL_EQ_TOL(value1, value2, tolerance) \
634         EXPECT_FLOAT_EQ_TOL(value1, value2, tolerance)
635 #    define ASSERT_REAL_EQ_TOL(value1, value2, tolerance) \
636         ASSERT_FLOAT_EQ_TOL(value1, value2, tolerance)
637 #endif
638
639 //! EXPECT_REAL_EQ_TOL with default tolerance
640 #define EXPECT_REAL_EQ(value1, value2) \
641     EXPECT_REAL_EQ_TOL(value1, value2, ::gmx::test::defaultRealTolerance())
642 //! ASSERT_REAL_EQ_TOL with default tolerance
643 #define ASSERT_REAL_EQ(value1, value2) \
644     ASSERT_REAL_EQ_TOL(value1, value2, ::gmx::test::defaultRealTolerance())
645 //! \}
646
647 //! \cond internal
648 /*! \internal \brief
649  * Helper method for `(EXPECT|ASSERT)_PLAIN`.
650  */
651 static inline ::testing::AssertionResult plainAssertHelper(const char* /*expr*/,
652                                                            const ::testing::AssertionResult& expr)
653 {
654     return expr;
655 }
656 //! \endcond
657
658 /*! \brief
659  * Assert for predicates that return AssertionResult and produce a full failure
660  * message.
661  *
662  * `expr` should evaluate to AssertionResult, and on failure the message from
663  * the result is used as-is, unlike in EXPECT_TRUE().
664  *
665  * \hideinitializer
666  */
667 #define EXPECT_PLAIN(expr) EXPECT_PRED_FORMAT1(plainAssertHelper, expr)
668 /*! \brief
669  * Assert for predicates that return AssertionResult and produce a full failure
670  * message.
671  *
672  * `expr` should evaluate to AssertionResult, and on failure the message from
673  * the result is used as-is, unlike in ASSERT_TRUE().
674  *
675  * \hideinitializer
676  */
677 #define ASSERT_PLAIN(expr) ASSERT_PRED_FORMAT1(plainAssertHelper, expr)
678
679 //! \}
680
681 } // namespace test
682 } // namespace gmx
683
684 #endif