Fix parallel 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,2015,2016,2017,2018,2019,2020, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*! \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  *  - death tests
45  *
46  * \if internal
47  * \todo
48  * The implementation is somewhat ugly, and accesses some Google Test
49  * internals.  Could be nice to clean it up a bit.
50  * \endif
51  *
52  * \author Teemu Murtola <teemu.murtola@gmail.com>
53  * \inlibraryapi
54  * \ingroup module_testutils
55  */
56 #ifndef GMX_TESTUTILS_TESTASSERTS_H
57 #define GMX_TESTUTILS_TESTASSERTS_H
58
59 #include <string>
60
61 #include <gtest/gtest.h>
62
63 #include "gromacs/utility/basedefinitions.h"
64 #include "gromacs/utility/exceptions.h"
65 #include "gromacs/utility/real.h"
66
67 namespace gmx
68 {
69 namespace test
70 {
71
72 namespace internal
73 {
74 //! \cond internal
75 /*! \internal
76  * \brief
77  * Called for an expected exception from EXPECT_THROW_GMX().
78  *
79  * \param[in] ex  Exception that was thrown.
80  */
81 void processExpectedException(const std::exception& ex);
82 //! \endcond
83 } // namespace internal
84
85 //! \libinternal \addtogroup module_testutils
86 //! \{
87
88 /*! \name Assertions for exceptions
89  *
90  * These macros replace `(ASSERT|EXPECT)(_NO)?_THROW` from Google Test.
91  * They are used exactly like the Google Test ones, but also print details of
92  * any unexpected exceptions using \Gromacs-specific routines.
93  * This makes it much easier to see at one glance what went wrong.
94  * See Google Test documentation for details on how to use the macros.
95  */
96 //! \{
97
98 //! \cond internal
99 /*! \brief
100  * Internal implementation macro for exception assertations.
101  *
102  * \param statement          Statements to execute.
103  * \param expected_exception Exception type that \p statement should throw.
104  * \param fail               Function/macro to call on failure.
105  *
106  * The implementation is copied and adjusted from
107  * include/gtest/internal/gtest-internal.h in Google Test 1.6.0.
108  */
109 #define GMX_TEST_THROW_(statement, expected_exception, fail)                           \
110     GTEST_AMBIGUOUS_ELSE_BLOCKER_                                                      \
111     if (::testing::AssertionResult gmx_ar = ::testing::AssertionSuccess())             \
112     {                                                                                  \
113         bool gmx_caught_expected = false;                                              \
114         try                                                                            \
115         {                                                                              \
116             GTEST_SUPPRESS_UNREACHABLE_CODE_WARNING_BELOW_(statement);                 \
117         }                                                                              \
118         catch (expected_exception const& ex)                                           \
119         {                                                                              \
120             gmx_caught_expected = true;                                                \
121             ::gmx::test::internal::processExpectedException(ex);                       \
122         }                                                                              \
123         catch (std::exception const& ex)                                               \
124         {                                                                              \
125             gmx_ar << "Expected: " #statement " throws an exception of type "          \
126                    << #expected_exception ".\n  Actual: it throws a different type.\n" \
127                    << "Exception details:\n"                                           \
128                    << ::gmx::formatExceptionMessageToString(ex);                       \
129             goto GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__);                  \
130         }                                                                              \
131         catch (...)                                                                    \
132         {                                                                              \
133             gmx_ar << "Expected: " #statement " throws an exception of type "          \
134                    << #expected_exception ".\n  Actual: it throws a different type.";  \
135             goto GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__);                  \
136         }                                                                              \
137         if (!gmx_caught_expected)                                                      \
138         {                                                                              \
139             gmx_ar << "Expected: " #statement " throws an exception of type "          \
140                    << #expected_exception ".\n  Actual: it throws nothing.";           \
141             goto GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__);                  \
142         }                                                                              \
143     }                                                                                  \
144     else                                                                               \
145         GTEST_CONCAT_TOKEN_(gmx_label_testthrow_, __LINE__) : fail(gmx_ar.message())
146
147 /*! \brief
148  * Internal implementation macro for exception assertations.
149  *
150  * \param statement          Statements to execute.
151  * \param fail               Function/macro to call on failure.
152  *
153  * The implementation is copied and adjusted from
154  * include/gtest/internal/gtest-internal.h in Google Test 1.6.0.
155  */
156 #define GMX_TEST_NO_THROW_(statement, fail)                                    \
157     GTEST_AMBIGUOUS_ELSE_BLOCKER_                                              \
158     if (::testing::AssertionResult gmx_ar = ::testing::AssertionSuccess())     \
159     {                                                                          \
160         try                                                                    \
161         {                                                                      \
162             GTEST_SUPPRESS_UNREACHABLE_CODE_WARNING_BELOW_(statement);         \
163         }                                                                      \
164         catch (std::exception const& ex)                                       \
165         {                                                                      \
166             gmx_ar << "Expected: " #statement " doesn't throw an exception.\n" \
167                    << "  Actual: it throws.\n"                                 \
168                    << "Exception details:\n"                                   \
169                    << ::gmx::formatExceptionMessageToString(ex);               \
170             goto GTEST_CONCAT_TOKEN_(gmx_label_testnothrow_, __LINE__);        \
171         }                                                                      \
172         catch (...)                                                            \
173         {                                                                      \
174             gmx_ar << "Expected: " #statement " doesn't throw an exception.\n" \
175                    << "  Actual: it throws.";                                  \
176             goto GTEST_CONCAT_TOKEN_(gmx_label_testnothrow_, __LINE__);        \
177         }                                                                      \
178     }                                                                          \
179     else                                                                       \
180         GTEST_CONCAT_TOKEN_(gmx_label_testnothrow_, __LINE__) : fail(gmx_ar.message())
181 //! \endcond
182
183 /*! \brief
184  * Asserts that a statement throws a given exception.
185  *
186  * \hideinitializer
187  */
188 #define EXPECT_THROW_GMX(statement, expected_exception) \
189     GMX_TEST_THROW_(statement, expected_exception, GTEST_NONFATAL_FAILURE_)
190 /*! \brief
191  * Asserts that a statement does not throw.
192  *
193  * \hideinitializer
194  */
195 #define EXPECT_NO_THROW_GMX(statement) GMX_TEST_NO_THROW_(statement, GTEST_NONFATAL_FAILURE_)
196 /*! \brief
197  * Asserts that a statement throws a given exception.
198  *
199  * \hideinitializer
200  */
201 #define ASSERT_THROW_GMX(statement, expected_exception) \
202     GMX_TEST_THROW_(statement, expected_exception, GTEST_FATAL_FAILURE_)
203 /*! \brief
204  * Asserts that a statement does not throw.
205  *
206  * \hideinitializer
207  */
208 #define ASSERT_NO_THROW_GMX(statement) GMX_TEST_NO_THROW_(statement, GTEST_FATAL_FAILURE_)
209
210 /*! \brief
211  * Wrapper around EXPECT_DEATH_IF_SUPPORTED gtest macro for thread safe execution.
212  *
213  * Makes sure that tests that are supposed to trigger an assertion are only executed
214  * in a threadsafe environment, and not when when running under e.g. MPI.
215  *
216  * \hideinitializer
217  */
218 #define GMX_EXPECT_DEATH_IF_SUPPORTED(expr, msg) \
219     if (!GMX_LIB_MPI)                            \
220     {                                            \
221         EXPECT_DEATH_IF_SUPPORTED(expr, msg);    \
222     }
223
224 /*! \brief
225  * Wrapper around ASSERT_DEATH_IF_SUPPORTED gtest macro for thread safe execution.
226  *
227  * Makes sure that tests that are supposed to trigger an assertion are only executed
228  * in a threadsafe environment when running under e.g. MPI.
229  *
230  * \hideinitializer
231  */
232 #define GMX_ASSERT_DEATH_IF_SUPPORTED(expr, msg) \
233     if (!GMX_LIB_MPI)                            \
234     {                                            \
235         ASSERT_DEATH_IF_SUPPORTED(expr, msg);    \
236     }
237
238 //! \}
239
240 /*! \libinternal \brief
241  * Computes and represents a floating-point difference value.
242  *
243  * Methods in this class do not throw, except for toString(), which may throw
244  * std::bad_alloc.
245  *
246  * \see FloatingPointTolerance
247  */
248 class FloatingPointDifference
249 {
250 public:
251     /*! \brief Initializes a single-precision difference.
252      *
253      *  \param ref    First term in difference
254      *  \param value  Second term in difference
255      *
256      *  For absolute and ULP differences the two parameters are equivalent,
257      *  since the difference is symmetric. For relative differences
258      *  the first term is interpreted as the reference value, from which
259      *  we extract the magnitude to compare with.
260      */
261     FloatingPointDifference(float ref, float value);
262
263     /*! \brief Initializes a double-precision difference.
264      *
265      *  \param ref    First term in difference
266      *  \param value  Second term in difference
267      *
268      *  For absolute and ULP differences the two parameters are equivalent,
269      *  since the difference is symmetric. For relative differences
270      *  the first term is interpreted as the reference value, from which
271      *  we extract the magnitude to compare with.
272      */
273     FloatingPointDifference(double ref, double value);
274
275     /*! \brief
276      * Whether one or both of the compared values were NaN.
277      *
278      * If this returns `true`, other accessors return meaningless values.
279      */
280     bool isNaN() const;
281     //! Returns the difference as an absolute number (always non-negative).
282     double asAbsolute() const { return absoluteDifference_; }
283     /*! \brief
284      * Returns the difference as ULPs (always non-negative).
285      *
286      * The ULPs are calculated for the type that corresponds to the
287      * constructor used to initialize the difference.
288      * The ULP difference between 0.0 and -0.0 is zero.
289      */
290     uint64_t asUlps() const { return ulpDifference_; }
291     /*! \brief
292      * Whether the compared values were of different sign.
293      *
294      * 0.0 and -0.0 are treated as positive and negative, respectively.
295      */
296     bool signsDiffer() const { return bSignDifference_; }
297     /*! \brief
298      * Whether the difference is between single- or double-precision
299      * numbers.
300      */
301     bool isDouble() const { return bDouble_; }
302     //! Formats the difference as a string for assertion failure messages.
303     std::string toString() const;
304
305     //! Returns the magnitude of the original second term of the difference.
306     double termMagnitude() const { return termMagnitude_; }
307
308 private:
309     //! Save the magnitude of the reference value for relative (i.e., not ULP) tolerance
310     double termMagnitude_;
311     //! Stores the absolute difference, or NaN if one or both values were NaN.
312     double   absoluteDifference_;
313     uint64_t ulpDifference_;
314     bool     bSignDifference_;
315     /*! \brief
316      * Whether the difference was computed for single or double precision.
317      *
318      * This sets the units for `ulpDifference_`.
319      */
320     bool bDouble_;
321 };
322
323 /*! \libinternal \brief
324  * Specifies a floating-point comparison tolerance and checks whether a
325  * difference is within the tolerance.
326  *
327  * The related functions section lists methods that can be construct methods
328  * using less parameters than the full constructor, and with more obvious
329  * semantics.  These should be preferred over using the constructor directly.
330  *
331  * Several types of tolerances are possible:
332  *  - _absolute tolerance_: difference between the values must be smaller than
333  *    the given tolerance for the check to pass.
334  *    Setting the absolute tolerance to zero disables the absolute tolerance
335  *    check.
336  *  - _relative tolerance_: the absolute difference between the numbers must
337  *    be smaller than the tolerance multiplied by the first number. Setting
338  *    the relative tolerance to zero disables this check.
339  *  - _ULP tolerance_: ULP (units of least precision) difference between the
340  *    values must be smaller than the given tolerance for the check to pass.
341  *    Setting the ULP tolerance to zero requires exact match.
342  *    Setting the ULP tolerance to UINT64_MAX disables the ULP check.
343  *    `0.0` and `-0.0` are treated as equal for the ULP check.
344  *  - _sign check_: if set, any values that are of different signs fail the
345  *    check (note that this also applies to `0.0` and `-0.0`: a value with a
346  *    different sign than the zero will fail the check).
347  *
348  * Either an absolute, relative, or ULP tolerance must always be specified.
349  * If several of them are specified, then the check passes if either of the
350  * tolerances is satisfied.
351  *
352  * Any combination of absolute, relative, and ULP tolerance can be combined with
353  * the sign check.  In this case, the sign check must succeed for the check to
354  * pass, even if other tolerances are satisfied.
355  *
356  * The tolerances can be specified separately for single and double precision
357  * comparison.  Different initialization functions have different semantics on
358  * how the provided tolerance values are interpreted; check their
359  * documentation.
360  *
361  * Methods in this class do not throw, except for toString(), which may throw
362  * std::bad_alloc.
363  *
364  * \todo
365  * The factory methods that take ULP difference could be better formulated as
366  * methods that take the acceptable number of incorrect bits and/or the number
367  * of accurate bits.
368  *
369  * \see FloatingPointDifference
370  */
371 class FloatingPointTolerance
372 {
373 public:
374     /*! \brief
375      * Creates a tolerance with the specified values.
376      *
377      * \param[in]  singleAbsoluteTolerance
378      *     Allowed absolute difference in a single-precision number.
379      * \param[in]  doubleAbsoluteTolerance
380      *     Allowed absolute difference in a double-precision number.
381      * \param[in]  singleRelativeTolerance
382      *     Allowed relative difference in a single-precision number.
383      * \param[in]  doubleRelativeTolerance
384      *     Allowed relative difference in a double-precision number.
385      * \param[in]  singleUlpTolerance
386      *     Allowed ULP difference in a single-precision number.
387      * \param[in]  doubleUlpTolerance
388      *     Allowed ULP difference in a double-precision number.
389      * \param[in]  bSignMustMatch
390      *     Whether sign mismatch fails the comparison.
391      */
392     FloatingPointTolerance(float    singleAbsoluteTolerance,
393                            double   doubleAbsoluteTolerance,
394                            float    singleRelativeTolerance,
395                            double   doubleRelativeTolerance,
396                            uint64_t singleUlpTolerance,
397                            uint64_t doubleUlpTolerance,
398                            bool     bSignMustMatch) :
399         singleAbsoluteTolerance_(singleAbsoluteTolerance),
400         doubleAbsoluteTolerance_(doubleAbsoluteTolerance),
401         singleRelativeTolerance_(singleRelativeTolerance),
402         doubleRelativeTolerance_(doubleRelativeTolerance),
403         singleUlpTolerance_(singleUlpTolerance),
404         doubleUlpTolerance_(doubleUlpTolerance),
405         bSignMustMatch_(bSignMustMatch)
406     {
407     }
408
409     /*! \brief
410      * Checks whether a difference is within the specified tolerance.
411      *
412      * NaNs are always treated outside the tolerance.
413      */
414     bool isWithin(const FloatingPointDifference& difference) const;
415
416     //! Formats the tolerance as a string for assertion failure messages.
417     std::string toString(const FloatingPointDifference& difference) const;
418
419 private:
420     float    singleAbsoluteTolerance_;
421     double   doubleAbsoluteTolerance_;
422     float    singleRelativeTolerance_;
423     double   doubleRelativeTolerance_;
424     uint64_t singleUlpTolerance_;
425     uint64_t doubleUlpTolerance_;
426     bool     bSignMustMatch_;
427 };
428
429 /*! \brief
430  * Creates a tolerance that only allows a specified ULP difference.
431  *
432  * The tolerance uses the given ULP value for both precisions, i.e., double
433  * precision will have much stricter tolerance.
434  *
435  * \related FloatingPointTolerance
436  */
437 static inline FloatingPointTolerance ulpTolerance(uint64_t ulpDiff)
438 {
439     return FloatingPointTolerance(0.0, 0.0, 0.0, 0.0, ulpDiff, ulpDiff, false);
440 }
441
442 /*! \brief
443  * Creates a tolerance that allows a difference in two compared values that is
444  * relative to the given magnitude.
445  *
446  * \param[in] magnitude  Magnitude of the numbers the computation operates in.
447  * \param[in] tolerance  Relative tolerance permitted (e.g. 1e-4).
448  *
449  * In addition to setting an relative tolerance for both
450  * precisions, this sets the absolute tolerance such that values close to zero
451  * (in general, smaller than \p magnitude) do not fail the check if they
452  * differ by less than \p tolerance evaluated at \p magnitude.  This accounts
453  * for potential loss of precision for small values, and should be used when
454  * accuracy of values much less than \p magnitude do not matter for
455  * correctness.
456  *
457  * \related FloatingPointTolerance
458  */
459 FloatingPointTolerance relativeToleranceAsFloatingPoint(double magnitude, double tolerance);
460
461 /*! \brief
462  * Creates a tolerance that allows a precision-dependent difference in two
463  * compared values that is relative to the given magnitude.
464  *
465  * \param[in] magnitude        Magnitude of the numbers the computation
466  *     operates in.
467  * \param[in] singleTolerance  Relative tolerance permitted (e.g. 1e-4)
468  *     in single precision.
469  * \param[in] doubleTolerance  Relative tolerance permitted (e.g. 1e-4)
470  *     in double precision.
471  *
472  * In addition to setting an relative tolerance for both
473  * precisions, this sets the absolute tolerance such that values close to zero
474  * (in general, smaller than \p magnitude) do not fail the check if they
475  * differ by less than \p tolerance evaluated at \p magnitude.  This accounts
476  * for potential loss of precision for small values, and should be used when
477  * accuracy of values much less than \p magnitude do not matter for
478  * correctness.
479  *
480  * \related FloatingPointTolerance
481  */
482 FloatingPointTolerance relativeToleranceAsPrecisionDependentFloatingPoint(double magnitude,
483                                                                           float  singleTolerance,
484                                                                           double doubleTolerance);
485
486 /*! \brief
487  * Creates a tolerance that allows a precision-dependent relative difference in
488  * a complex computation.
489  *
490  * \param[in] magnitude      Magnitude of the numbers the computation operates in.
491  * \param[in] singleUlpDiff  Expected accuracy of single-precision
492  *     computation (in ULPs).
493  * \param[in] doubleUlpDiff  Expected accuracy of double-precision
494  *     computation (in ULPs).
495  *
496  * This works as relativeToleranceAsUlp(), but allows setting the ULP
497  * difference separately for the different precisions.  This supports
498  * cases where the double-precision calculation can acceptably has a higher ULP
499  * difference, but relaxing the single-precision tolerance would lead to an
500  * unnecessarily loose test.
501  *
502  * \related FloatingPointTolerance
503  */
504 static inline FloatingPointTolerance relativeToleranceAsPrecisionDependentUlp(double   magnitude,
505                                                                               uint64_t singleUlpDiff,
506                                                                               uint64_t doubleUlpDiff)
507 {
508     return FloatingPointTolerance(float(magnitude) * singleUlpDiff * GMX_FLOAT_EPS,
509                                   magnitude * doubleUlpDiff * GMX_DOUBLE_EPS, 0.0, 0.0,
510                                   singleUlpDiff, doubleUlpDiff, false);
511 }
512
513 /*! \brief
514  * Creates a tolerance that allows a specified absolute difference.
515  *
516  * \related FloatingPointTolerance
517  */
518 static inline FloatingPointTolerance absoluteTolerance(double tolerance)
519 {
520     return FloatingPointTolerance(float(tolerance), tolerance, 0.0, 0.0, UINT64_MAX, UINT64_MAX, false);
521 }
522
523 /*! \brief
524  * Creates a tolerance that allows a relative difference in a complex
525  * computation.
526  *
527  * \param[in] magnitude  Magnitude of the numbers the computation operates in.
528  * \param[in] ulpDiff    Expected accuracy of the computation (in ULPs).
529  *
530  * In addition to setting the ULP tolerance as ulpTolerance(), this sets the
531  * absolute tolerance such that values close to zero (in general, smaller than
532  * \p magnitude) do not fail the check if they differ by less than \p ulpDiff
533  * evaluated at \p magnitude.  This accounts for potential loss of precision
534  * for small values, and should be used when accuracy of values much less than
535  * \p magnitude do not matter for correctness.
536  *
537  * \related FloatingPointTolerance
538  */
539 static inline FloatingPointTolerance relativeToleranceAsUlp(double magnitude, uint64_t ulpDiff)
540 {
541     return relativeToleranceAsPrecisionDependentUlp(magnitude, ulpDiff, ulpDiff);
542 }
543
544 namespace detail
545 {
546 //! Default tolerance in ULPs for two floating-point values to compare equal.
547 constexpr uint64_t g_defaultUlpTolerance = 4;
548 } // namespace detail
549
550 /*! \brief
551  * Returns the default tolerance for comparing `real` numbers.
552  *
553  * \related FloatingPointTolerance
554  */
555 static inline FloatingPointTolerance defaultRealTolerance()
556 {
557     return relativeToleranceAsUlp(1.0, detail::g_defaultUlpTolerance);
558 }
559
560
561 /*! \brief
562  * Returns the default tolerance for comparing single-precision numbers when
563  * compared by \Gromacs built in either precision mode.
564  *
565  * This permits a checker compiled with any \Gromacs precision to compare
566  * equal or not in the same way.
567  *
568  * \related FloatingPointTolerance
569  */
570 static inline FloatingPointTolerance defaultFloatTolerance()
571 {
572     return relativeToleranceAsPrecisionDependentUlp(
573             1.0, detail::g_defaultUlpTolerance,
574             static_cast<uint64_t>(detail::g_defaultUlpTolerance * (GMX_FLOAT_EPS / GMX_DOUBLE_EPS)));
575 }
576
577 /*! \name Assertions for floating-point comparison
578  *
579  * These routines extend `(EXPECT|ASSERT)_(FLOAT|DOUBLE)_EQ` and
580  * `(EXPECT|ASSERT)_NEAR` from Google Test to provide more flexible assertions
581  * for floating-point values.
582  *
583  * See gmx::test::FloatingPointTolerance for the possible ways to specify the
584  * tolerance, and gmx::test::FloatingPointDifference for some additional
585  * details of the difference calculation.
586  */
587 //! \{
588
589 //! \cond internal
590 /*! \internal \brief
591  * Assertion predicate formatter for comparing two floating-point values.
592  */
593 template<typename FloatType>
594 static inline ::testing::AssertionResult assertEqualWithinTolerance(const char* expr1,
595                                                                     const char* expr2,
596                                                                     const char* /*exprTolerance*/,
597                                                                     FloatType value1,
598                                                                     FloatType value2,
599                                                                     const FloatingPointTolerance& tolerance)
600 {
601     FloatingPointDifference diff(value1, value2);
602     if (tolerance.isWithin(diff))
603     {
604         return ::testing::AssertionSuccess();
605     }
606     return ::testing::AssertionFailure() << "  Value of: " << expr2 << std::endl
607                                          << "    Actual: " << value2 << std::endl
608                                          << "  Expected: " << expr1 << std::endl
609                                          << "  Which is: " << value1 << std::endl
610                                          << "Difference: " << diff.toString() << std::endl
611                                          << " Tolerance: " << tolerance.toString(diff);
612 }
613 //! \endcond
614
615 /*! \brief
616  * Asserts that two single-precision values are within the given tolerance.
617  *
618  * \hideinitializer
619  */
620 #define EXPECT_FLOAT_EQ_TOL(value1, value2, tolerance) \
621     EXPECT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<float>, value1, value2, tolerance)
622 /*! \brief
623  * Asserts that two double-precision values are within the given tolerance.
624  *
625  * \hideinitializer
626  */
627 #define EXPECT_DOUBLE_EQ_TOL(value1, value2, tolerance) \
628     EXPECT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<double>, value1, value2, tolerance)
629 /*! \def EXPECT_REAL_EQ_TOL
630  * \brief
631  * Asserts that two `real` values are within the given tolerance.
632  *
633  * \hideinitializer
634  */
635 /*! \brief
636  * Asserts that two single-precision values are within the given tolerance.
637  *
638  * \hideinitializer
639  */
640 #define ASSERT_FLOAT_EQ_TOL(value1, value2, tolerance) \
641     ASSERT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<float>, value1, value2, tolerance)
642 /*! \brief
643  * Asserts that two double-precision values are within the given tolerance.
644  *
645  * \hideinitializer
646  */
647 #define ASSERT_DOUBLE_EQ_TOL(value1, value2, tolerance) \
648     ASSERT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<double>, value1, value2, tolerance)
649 /*! \def ASSERT_REAL_EQ_TOL
650  * \brief
651  * Asserts that two `real` values are within the given tolerance.
652  *
653  * \hideinitializer
654  */
655
656 #if GMX_DOUBLE
657 #    define EXPECT_REAL_EQ_TOL(value1, value2, tolerance) \
658         EXPECT_DOUBLE_EQ_TOL(value1, value2, tolerance)
659 #    define ASSERT_REAL_EQ_TOL(value1, value2, tolerance) \
660         ASSERT_DOUBLE_EQ_TOL(value1, value2, tolerance)
661 #else
662 #    define EXPECT_REAL_EQ_TOL(value1, value2, tolerance) \
663         EXPECT_FLOAT_EQ_TOL(value1, value2, tolerance)
664 #    define ASSERT_REAL_EQ_TOL(value1, value2, tolerance) \
665         ASSERT_FLOAT_EQ_TOL(value1, value2, tolerance)
666 #endif
667
668 //! EXPECT_REAL_EQ_TOL with default tolerance
669 #define EXPECT_REAL_EQ(value1, value2) \
670     EXPECT_REAL_EQ_TOL(value1, value2, ::gmx::test::defaultRealTolerance())
671 //! ASSERT_REAL_EQ_TOL with default tolerance
672 #define ASSERT_REAL_EQ(value1, value2) \
673     ASSERT_REAL_EQ_TOL(value1, value2, ::gmx::test::defaultRealTolerance())
674 //! \}
675
676 //! \cond internal
677 /*! \internal \brief
678  * Helper method for `(EXPECT|ASSERT)_PLAIN`.
679  */
680 static inline ::testing::AssertionResult plainAssertHelper(const char* /*expr*/,
681                                                            const ::testing::AssertionResult& expr)
682 {
683     return expr;
684 }
685 //! \endcond
686
687 /*! \brief
688  * Assert for predicates that return AssertionResult and produce a full failure
689  * message.
690  *
691  * `expr` should evaluate to AssertionResult, and on failure the message from
692  * the result is used as-is, unlike in EXPECT_TRUE().
693  *
694  * \hideinitializer
695  */
696 #define EXPECT_PLAIN(expr) EXPECT_PRED_FORMAT1(plainAssertHelper, expr)
697 /*! \brief
698  * Assert for predicates that return AssertionResult and produce a full failure
699  * message.
700  *
701  * `expr` should evaluate to AssertionResult, and on failure the message from
702  * the result is used as-is, unlike in ASSERT_TRUE().
703  *
704  * \hideinitializer
705  */
706 #define ASSERT_PLAIN(expr) ASSERT_PRED_FORMAT1(plainAssertHelper, expr)
707
708 //! \}
709
710 } // namespace test
711 } // namespace gmx
712
713 #endif