f1baab8be797c28ccd514d6807bcf0e532466bd0
[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 by the GROMACS development team.
5  * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \libinternal \file
37  * \brief
38  * Extra assertions for unit tests.
39  *
40  * This file provides assertion macros that extend/replace Google Test
41  * assertions for:
42  *  - exceptions
43  *  - floating-point comparison
44  *  - comparison against NULL
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 //! \}
211
212 /*! \libinternal \brief
213  * Computes and represents a floating-point difference value.
214  *
215  * Methods in this class do not throw, except for toString(), which may throw
216  * std::bad_alloc.
217  *
218  * \see FloatingPointTolerance
219  */
220 class FloatingPointDifference
221 {
222 public:
223     /*! \brief Initializes a single-precision difference.
224      *
225      *  \param ref    First term in difference
226      *  \param value  Second term in difference
227      *
228      *  For absolute and ULP differences the two parameters are equivalent,
229      *  since the difference is symmetric. For relative differences
230      *  the first term is interpreted as the reference value, from which
231      *  we extract the magnitude to compare with.
232      */
233     FloatingPointDifference(float ref, float value);
234
235     /*! \brief Initializes a double-precision difference.
236      *
237      *  \param ref    First term in difference
238      *  \param value  Second term in difference
239      *
240      *  For absolute and ULP differences the two parameters are equivalent,
241      *  since the difference is symmetric. For relative differences
242      *  the first term is interpreted as the reference value, from which
243      *  we extract the magnitude to compare with.
244      */
245     FloatingPointDifference(double ref, double value);
246
247     /*! \brief
248      * Whether one or both of the compared values were NaN.
249      *
250      * If this returns `true`, other accessors return meaningless values.
251      */
252     bool isNaN() const;
253     //! Returns the difference as an absolute number (always non-negative).
254     double asAbsolute() const { return absoluteDifference_; }
255     /*! \brief
256      * Returns the difference as ULPs (always non-negative).
257      *
258      * The ULPs are calculated for the type that corresponds to the
259      * constructor used to initialize the difference.
260      * The ULP difference between 0.0 and -0.0 is zero.
261      */
262     uint64_t asUlps() const { return ulpDifference_; }
263     /*! \brief
264      * Whether the compared values were of different sign.
265      *
266      * 0.0 and -0.0 are treated as positive and negative, respectively.
267      */
268     bool signsDiffer() const { return bSignDifference_; }
269     /*! \brief
270      * Whether the difference is between single- or double-precision
271      * numbers.
272      */
273     bool isDouble() const { return bDouble_; }
274     //! Formats the difference as a string for assertion failure messages.
275     std::string toString() const;
276
277     //! Returns the magnitude of the original second term of the difference.
278     double termMagnitude() const { return termMagnitude_; }
279
280 private:
281     //! Save the magnitude of the reference value for relative (i.e., not ULP) tolerance
282     double termMagnitude_;
283     //! Stores the absolute difference, or NaN if one or both values were NaN.
284     double   absoluteDifference_;
285     uint64_t ulpDifference_;
286     bool     bSignDifference_;
287     /*! \brief
288      * Whether the difference was computed for single or double precision.
289      *
290      * This sets the units for `ulpDifference_`.
291      */
292     bool bDouble_;
293 };
294
295 /*! \libinternal \brief
296  * Specifies a floating-point comparison tolerance and checks whether a
297  * difference is within the tolerance.
298  *
299  * The related functions section lists methods that can be construct methods
300  * using less parameters than the full constructor, and with more obvious
301  * semantics.  These should be preferred over using the constructor directly.
302  *
303  * Several types of tolerances are possible:
304  *  - _absolute tolerance_: difference between the values must be smaller than
305  *    the given tolerance for the check to pass.
306  *    Setting the absolute tolerance to zero disables the absolute tolerance
307  *    check.
308  *  - _relative tolerance_: the absolute difference between the numbers must
309  *    be smaller than the tolerance multiplied by the first number. Setting
310  *    the relative tolerance to zero disables this check.
311  *  - _ULP tolerance_: ULP (units of least precision) difference between the
312  *    values must be smaller than the given tolerance for the check to pass.
313  *    Setting the ULP tolerance to zero requires exact match.
314  *    Setting the ULP tolerance to UINT64_MAX disables the ULP check.
315  *    `0.0` and `-0.0` are treated as equal for the ULP check.
316  *  - _sign check_: if set, any values that are of different signs fail the
317  *    check (note that this also applies to `0.0` and `-0.0`: a value with a
318  *    different sign than the zero will fail the check).
319  *
320  * Either an absolute, relative, or ULP tolerance must always be specified.
321  * If several of them are specified, then the check passes if either of the
322  * tolerances is satisfied.
323  *
324  * Any combination of absolute, relative, and ULP tolerance can be combined with
325  * the sign check.  In this case, the sign check must succeed for the check to
326  * pass, even if other tolerances are satisfied.
327  *
328  * The tolerances can be specified separately for single and double precision
329  * comparison.  Different initialization functions have different semantics on
330  * how the provided tolerance values are interpreted; check their
331  * documentation.
332  *
333  * Methods in this class do not throw, except for toString(), which may throw
334  * std::bad_alloc.
335  *
336  * \todo
337  * The factory methods that take ULP difference could be better formulated as
338  * methods that take the acceptable number of incorrect bits and/or the number
339  * of accurate bits.
340  *
341  * \see FloatingPointDifference
342  */
343 class FloatingPointTolerance
344 {
345 public:
346     /*! \brief
347      * Creates a tolerance with the specified values.
348      *
349      * \param[in]  singleAbsoluteTolerance
350      *     Allowed absolute difference in a single-precision number.
351      * \param[in]  doubleAbsoluteTolerance
352      *     Allowed absolute difference in a double-precision number.
353      * \param[in]  singleRelativeTolerance
354      *     Allowed relative difference in a single-precision number.
355      * \param[in]  doubleRelativeTolerance
356      *     Allowed relative difference in a double-precision number.
357      * \param[in]  singleUlpTolerance
358      *     Allowed ULP difference in a single-precision number.
359      * \param[in]  doubleUlpTolerance
360      *     Allowed ULP difference in a double-precision number.
361      * \param[in]  bSignMustMatch
362      *     Whether sign mismatch fails the comparison.
363      */
364     FloatingPointTolerance(float    singleAbsoluteTolerance,
365                            double   doubleAbsoluteTolerance,
366                            float    singleRelativeTolerance,
367                            double   doubleRelativeTolerance,
368                            uint64_t singleUlpTolerance,
369                            uint64_t doubleUlpTolerance,
370                            bool     bSignMustMatch) :
371         singleAbsoluteTolerance_(singleAbsoluteTolerance),
372         doubleAbsoluteTolerance_(doubleAbsoluteTolerance),
373         singleRelativeTolerance_(singleRelativeTolerance),
374         doubleRelativeTolerance_(doubleRelativeTolerance),
375         singleUlpTolerance_(singleUlpTolerance),
376         doubleUlpTolerance_(doubleUlpTolerance),
377         bSignMustMatch_(bSignMustMatch)
378     {
379     }
380
381     /*! \brief
382      * Checks whether a difference is within the specified tolerance.
383      *
384      * NaNs are always treated outside the tolerance.
385      */
386     bool isWithin(const FloatingPointDifference& difference) const;
387
388     //! Formats the tolerance as a string for assertion failure messages.
389     std::string toString(const FloatingPointDifference& difference) const;
390
391 private:
392     float    singleAbsoluteTolerance_;
393     double   doubleAbsoluteTolerance_;
394     float    singleRelativeTolerance_;
395     double   doubleRelativeTolerance_;
396     uint64_t singleUlpTolerance_;
397     uint64_t doubleUlpTolerance_;
398     bool     bSignMustMatch_;
399 };
400
401 /*! \brief
402  * Creates a tolerance that only allows a specified ULP difference.
403  *
404  * The tolerance uses the given ULP value for both precisions, i.e., double
405  * precision will have much stricter tolerance.
406  *
407  * \related FloatingPointTolerance
408  */
409 static inline FloatingPointTolerance ulpTolerance(uint64_t ulpDiff)
410 {
411     return FloatingPointTolerance(0.0, 0.0, 0.0, 0.0, ulpDiff, ulpDiff, false);
412 }
413
414 /*! \brief
415  * Creates a tolerance that allows a difference in two compared values that is
416  * relative to the given magnitude.
417  *
418  * \param[in] magnitude  Magnitude of the numbers the computation operates in.
419  * \param[in] tolerance  Relative tolerance permitted (e.g. 1e-4).
420  *
421  * In addition to setting an relative tolerance for both
422  * precisions, this sets the absolute tolerance such that values close to zero
423  * (in general, smaller than \p magnitude) do not fail the check if they
424  * differ by less than \p tolerance evaluated at \p magnitude.  This accounts
425  * for potential loss of precision for small values, and should be used when
426  * accuracy of values much less than \p magnitude do not matter for
427  * correctness.
428  *
429  * \related FloatingPointTolerance
430  */
431 FloatingPointTolerance relativeToleranceAsFloatingPoint(double magnitude, double tolerance);
432
433 /*! \brief
434  * Creates a tolerance that allows a precision-dependent difference in two
435  * compared values that is relative to the given magnitude.
436  *
437  * \param[in] magnitude        Magnitude of the numbers the computation
438  *     operates in.
439  * \param[in] singleTolerance  Relative tolerance permitted (e.g. 1e-4)
440  *     in single precision.
441  * \param[in] doubleTolerance  Relative tolerance permitted (e.g. 1e-4)
442  *     in double precision.
443  *
444  * In addition to setting an relative tolerance for both
445  * precisions, this sets the absolute tolerance such that values close to zero
446  * (in general, smaller than \p magnitude) do not fail the check if they
447  * differ by less than \p tolerance evaluated at \p magnitude.  This accounts
448  * for potential loss of precision for small values, and should be used when
449  * accuracy of values much less than \p magnitude do not matter for
450  * correctness.
451  *
452  * \related FloatingPointTolerance
453  */
454 FloatingPointTolerance relativeToleranceAsPrecisionDependentFloatingPoint(double magnitude,
455                                                                           float  singleTolerance,
456                                                                           double doubleTolerance);
457
458 /*! \brief
459  * Creates a tolerance that allows a precision-dependent relative difference in
460  * a complex computation.
461  *
462  * \param[in] magnitude      Magnitude of the numbers the computation operates in.
463  * \param[in] singleUlpDiff  Expected accuracy of single-precision
464  *     computation (in ULPs).
465  * \param[in] doubleUlpDiff  Expected accuracy of double-precision
466  *     computation (in ULPs).
467  *
468  * This works as relativeToleranceAsUlp(), but allows setting the ULP
469  * difference separately for the different precisions.  This supports
470  * cases where the double-precision calculation can acceptably has a higher ULP
471  * difference, but relaxing the single-precision tolerance would lead to an
472  * unnecessarily loose test.
473  *
474  * \related FloatingPointTolerance
475  */
476 static inline FloatingPointTolerance relativeToleranceAsPrecisionDependentUlp(double   magnitude,
477                                                                               uint64_t singleUlpDiff,
478                                                                               uint64_t doubleUlpDiff)
479 {
480     return FloatingPointTolerance(float(magnitude) * singleUlpDiff * GMX_FLOAT_EPS,
481                                   magnitude * doubleUlpDiff * GMX_DOUBLE_EPS, 0.0, 0.0,
482                                   singleUlpDiff, doubleUlpDiff, false);
483 }
484
485 /*! \brief
486  * Creates a tolerance that allows a specified absolute difference.
487  *
488  * \related FloatingPointTolerance
489  */
490 static inline FloatingPointTolerance absoluteTolerance(double tolerance)
491 {
492     return FloatingPointTolerance(float(tolerance), tolerance, 0.0, 0.0, UINT64_MAX, UINT64_MAX, false);
493 }
494
495 /*! \brief
496  * Creates a tolerance that allows a relative difference in a complex
497  * computation.
498  *
499  * \param[in] magnitude  Magnitude of the numbers the computation operates in.
500  * \param[in] ulpDiff    Expected accuracy of the computation (in ULPs).
501  *
502  * In addition to setting the ULP tolerance as ulpTolerance(), this sets the
503  * absolute tolerance such that values close to zero (in general, smaller than
504  * \p magnitude) do not fail the check if they differ by less than \p ulpDiff
505  * evaluated at \p magnitude.  This accounts for potential loss of precision
506  * for small values, and should be used when accuracy of values much less than
507  * \p magnitude do not matter for correctness.
508  *
509  * \related FloatingPointTolerance
510  */
511 static inline FloatingPointTolerance relativeToleranceAsUlp(double magnitude, uint64_t ulpDiff)
512 {
513     return relativeToleranceAsPrecisionDependentUlp(magnitude, ulpDiff, ulpDiff);
514 }
515
516 namespace detail
517 {
518 //! Default tolerance in ULPs for two floating-point values to compare equal.
519 constexpr uint64_t g_defaultUlpTolerance = 4;
520 } // namespace detail
521
522 /*! \brief
523  * Returns the default tolerance for comparing `real` numbers.
524  *
525  * \related FloatingPointTolerance
526  */
527 static inline FloatingPointTolerance defaultRealTolerance()
528 {
529     return relativeToleranceAsUlp(1.0, detail::g_defaultUlpTolerance);
530 }
531
532
533 /*! \brief
534  * Returns the default tolerance for comparing single-precision numbers when
535  * compared by \Gromacs built in either precision mode.
536  *
537  * This permits a checker compiled with any \Gromacs precision to compare
538  * equal or not in the same way.
539  *
540  * \related FloatingPointTolerance
541  */
542 static inline FloatingPointTolerance defaultFloatTolerance()
543 {
544     return relativeToleranceAsPrecisionDependentUlp(
545             1.0, detail::g_defaultUlpTolerance,
546             static_cast<uint64_t>(detail::g_defaultUlpTolerance * (GMX_FLOAT_EPS / GMX_DOUBLE_EPS)));
547 }
548
549 /*! \name Assertions for floating-point comparison
550  *
551  * These routines extend `(EXPECT|ASSERT)_(FLOAT|DOUBLE)_EQ` and
552  * `(EXPECT|ASSERT)_NEAR` from Google Test to provide more flexible assertions
553  * for floating-point values.
554  *
555  * See gmx::test::FloatingPointTolerance for the possible ways to specify the
556  * tolerance, and gmx::test::FloatingPointDifference for some additional
557  * details of the difference calculation.
558  */
559 //! \{
560
561 //! \cond internal
562 /*! \internal \brief
563  * Assertion predicate formatter for comparing two floating-point values.
564  */
565 template<typename FloatType>
566 static inline ::testing::AssertionResult assertEqualWithinTolerance(const char* expr1,
567                                                                     const char* expr2,
568                                                                     const char* /*exprTolerance*/,
569                                                                     FloatType value1,
570                                                                     FloatType value2,
571                                                                     const FloatingPointTolerance& tolerance)
572 {
573     FloatingPointDifference diff(value1, value2);
574     if (tolerance.isWithin(diff))
575     {
576         return ::testing::AssertionSuccess();
577     }
578     return ::testing::AssertionFailure() << "  Value of: " << expr2 << std::endl
579                                          << "    Actual: " << value2 << std::endl
580                                          << "  Expected: " << expr1 << std::endl
581                                          << "  Which is: " << value1 << std::endl
582                                          << "Difference: " << diff.toString() << std::endl
583                                          << " Tolerance: " << tolerance.toString(diff);
584 }
585 //! \endcond
586
587 /*! \brief
588  * Asserts that two single-precision values are within the given tolerance.
589  *
590  * \hideinitializer
591  */
592 #define EXPECT_FLOAT_EQ_TOL(value1, value2, tolerance) \
593     EXPECT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<float>, value1, value2, tolerance)
594 /*! \brief
595  * Asserts that two double-precision values are within the given tolerance.
596  *
597  * \hideinitializer
598  */
599 #define EXPECT_DOUBLE_EQ_TOL(value1, value2, tolerance) \
600     EXPECT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<double>, value1, value2, tolerance)
601 /*! \def EXPECT_REAL_EQ_TOL
602  * \brief
603  * Asserts that two `real` values are within the given tolerance.
604  *
605  * \hideinitializer
606  */
607 /*! \brief
608  * Asserts that two single-precision values are within the given tolerance.
609  *
610  * \hideinitializer
611  */
612 #define ASSERT_FLOAT_EQ_TOL(value1, value2, tolerance) \
613     ASSERT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<float>, value1, value2, tolerance)
614 /*! \brief
615  * Asserts that two double-precision values are within the given tolerance.
616  *
617  * \hideinitializer
618  */
619 #define ASSERT_DOUBLE_EQ_TOL(value1, value2, tolerance) \
620     ASSERT_PRED_FORMAT3(::gmx::test::assertEqualWithinTolerance<double>, value1, value2, tolerance)
621 /*! \def ASSERT_REAL_EQ_TOL
622  * \brief
623  * Asserts that two `real` values are within the given tolerance.
624  *
625  * \hideinitializer
626  */
627
628 #if GMX_DOUBLE
629 #    define EXPECT_REAL_EQ_TOL(value1, value2, tolerance) \
630         EXPECT_DOUBLE_EQ_TOL(value1, value2, tolerance)
631 #    define ASSERT_REAL_EQ_TOL(value1, value2, tolerance) \
632         ASSERT_DOUBLE_EQ_TOL(value1, value2, tolerance)
633 #else
634 #    define EXPECT_REAL_EQ_TOL(value1, value2, tolerance) \
635         EXPECT_FLOAT_EQ_TOL(value1, value2, tolerance)
636 #    define ASSERT_REAL_EQ_TOL(value1, value2, tolerance) \
637         ASSERT_FLOAT_EQ_TOL(value1, value2, tolerance)
638 #endif
639
640 //! EXPECT_REAL_EQ_TOL with default tolerance
641 #define EXPECT_REAL_EQ(value1, value2) \
642     EXPECT_REAL_EQ_TOL(value1, value2, ::gmx::test::defaultRealTolerance())
643 //! ASSERT_REAL_EQ_TOL with default tolerance
644 #define ASSERT_REAL_EQ(value1, value2) \
645     ASSERT_REAL_EQ_TOL(value1, value2, ::gmx::test::defaultRealTolerance())
646 //! \}
647
648 //! \cond internal
649 /*! \internal \brief
650  * Helper method for `(EXPECT|ASSERT)_PLAIN`.
651  */
652 static inline ::testing::AssertionResult plainAssertHelper(const char* /*expr*/,
653                                                            const ::testing::AssertionResult& expr)
654 {
655     return expr;
656 }
657 //! \endcond
658
659 /*! \brief
660  * Assert for predicates that return AssertionResult and produce a full failure
661  * message.
662  *
663  * `expr` should evaluate to AssertionResult, and on failure the message from
664  * the result is used as-is, unlike in EXPECT_TRUE().
665  *
666  * \hideinitializer
667  */
668 #define EXPECT_PLAIN(expr) EXPECT_PRED_FORMAT1(plainAssertHelper, expr)
669 /*! \brief
670  * Assert for predicates that return AssertionResult and produce a full failure
671  * message.
672  *
673  * `expr` should evaluate to AssertionResult, and on failure the message from
674  * the result is used as-is, unlike in ASSERT_TRUE().
675  *
676  * \hideinitializer
677  */
678 #define ASSERT_PLAIN(expr) ASSERT_PRED_FORMAT1(plainAssertHelper, expr)
679
680 //! \}
681
682 } // namespace test
683 } // namespace gmx
684
685 #endif