Bring PME GPU/CUDA internal structure names to CamelCase
[gromacs/AngularHB.git] / src / testutils / tests / testasserts_tests.cpp
blob01e9d6146da17abb74e72d9d68dc94605c583e89
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
35 /*! \internal \file
36 * \brief
37 * Tests \Gromacs-specific test assertions.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_testutils
42 #include "gmxpre.h"
44 #include "testutils/testasserts.h"
46 #include <gtest/gtest.h>
48 namespace
51 using ::testing::internal::FloatingPoint;
53 /*! \brief
54 * Helper to produce floating-point numbers with specified ULP difference.
56 * This doesn't work if the value would change sign.
58 * \ingroup module_testutils
60 template <typename FloatType>
61 FloatType addUlps(FloatType value, int ulps)
63 return FloatingPoint<FloatType>::ReinterpretBits(
64 FloatingPoint<FloatType>(value).bits() + ulps);
67 using ::gmx::test::FloatingPointDifference;
69 TEST(FloatingPointDifferenceTest, HandlesEqualValues)
71 FloatingPointDifference diff(1.2, 1.2);
72 EXPECT_TRUE(diff.isDouble());
73 EXPECT_FALSE(diff.isNaN());
74 EXPECT_EQ(0.0, diff.asAbsolute());
75 EXPECT_EQ(0U, diff.asUlps());
76 EXPECT_FALSE(diff.signsDiffer());
79 // TODO: Use typed tests to run all the tests for single and double.
80 TEST(FloatingPointDifferenceTest, HandlesFloatValues)
82 FloatingPointDifference diff(1.2f, 1.2f);
83 EXPECT_FALSE(diff.isDouble());
84 EXPECT_FALSE(diff.isNaN());
85 EXPECT_EQ(0.0, diff.asAbsolute());
86 EXPECT_EQ(0U, diff.asUlps());
87 EXPECT_FALSE(diff.signsDiffer());
90 TEST(FloatingPointDifferenceTest, HandlesZerosOfDifferentSign)
92 FloatingPointDifference diff(0.0, GMX_DOUBLE_NEGZERO);
93 EXPECT_FALSE(diff.isNaN());
94 EXPECT_EQ(0.0, diff.asAbsolute());
95 EXPECT_EQ(0U, diff.asUlps());
96 EXPECT_TRUE(diff.signsDiffer());
99 TEST(FloatingPointDifferenceTest, HandlesSignComparisonWithZero)
102 FloatingPointDifference diff(0.0, -1.2);
103 EXPECT_FALSE(diff.isNaN());
104 EXPECT_DOUBLE_EQ(1.2, diff.asAbsolute());
105 EXPECT_TRUE(diff.signsDiffer());
108 FloatingPointDifference diff(GMX_DOUBLE_NEGZERO, -1.2);
109 EXPECT_FALSE(diff.isNaN());
110 EXPECT_DOUBLE_EQ(1.2, diff.asAbsolute());
111 EXPECT_FALSE(diff.signsDiffer());
115 TEST(FloatingPointDifferenceTest, HandlesUlpDifferences)
117 const double first = 1.0;
118 const double second = addUlps(first, 2);
120 FloatingPointDifference diff(first, second);
121 EXPECT_FALSE(diff.isNaN());
122 EXPECT_DOUBLE_EQ(second - first, diff.asAbsolute());
123 EXPECT_EQ(2U, diff.asUlps());
124 EXPECT_FALSE(diff.signsDiffer());
127 FloatingPointDifference diff(second, first);
128 EXPECT_FALSE(diff.isNaN());
129 EXPECT_DOUBLE_EQ(second - first, diff.asAbsolute());
130 EXPECT_EQ(2U, diff.asUlps());
131 EXPECT_FALSE(diff.signsDiffer());
135 TEST(FloatingPointDifferenceTest, HandlesUlpDifferenceAcrossZero)
137 const double first = addUlps(GMX_DOUBLE_NEGZERO, 2);
138 const double second = addUlps( 0.0, 2);
139 FloatingPointDifference diff(first, second);
140 EXPECT_FALSE(diff.isNaN());
141 EXPECT_DOUBLE_EQ(second - first, diff.asAbsolute());
142 EXPECT_EQ(4U, diff.asUlps());
143 EXPECT_TRUE(diff.signsDiffer());
146 TEST(FloatingPointDifferenceTest, HandlesNaN)
148 const double first = std::numeric_limits<double>::quiet_NaN();
149 const double second = 2.0;
150 FloatingPointDifference diff(first, second);
151 EXPECT_TRUE(diff.isNaN());
154 TEST(FloatingPointToleranceTest, UlpTolerance)
156 using gmx::test::ulpTolerance;
158 FloatingPointDifference fequal(1.0, 1.0);
159 FloatingPointDifference fulp2(1.0f, addUlps(1.0f, 2));
160 EXPECT_TRUE(ulpTolerance(0).isWithin(fequal));
161 EXPECT_FALSE(ulpTolerance(1).isWithin(fulp2));
162 EXPECT_TRUE(ulpTolerance(2).isWithin(fulp2));
164 FloatingPointDifference dequal(1.0, 1.0);
165 FloatingPointDifference dulp2(1.0, addUlps(1.0, 2));
166 FloatingPointDifference dulp2f(1.0, static_cast<double>(addUlps(1.0f, 2)));
167 EXPECT_TRUE(ulpTolerance(0).isWithin(dequal));
168 EXPECT_TRUE(ulpTolerance(2).isWithin(dulp2));
169 EXPECT_FALSE(ulpTolerance(2).isWithin(dulp2f));
172 TEST(FloatingPointToleranceTest, RelativeToleranceAsFloatingPoint)
174 using gmx::test::relativeToleranceAsFloatingPoint;
176 FloatingPointDifference fequal(1.0f, 1.0f);
177 FloatingPointDifference fulp2(1.0f, addUlps(1.0f, 2));
178 FloatingPointDifference fdiff(1.0f, 1.011f);
179 FloatingPointDifference fsmall(0.1f, 0.111f);
180 FloatingPointDifference fsmall2(0.1f, 0.121f);
181 EXPECT_TRUE(relativeToleranceAsFloatingPoint(1.0, 1e-2).isWithin(fequal));
182 EXPECT_TRUE(relativeToleranceAsFloatingPoint(1.0, 1e-9).isWithin(fequal));
183 EXPECT_TRUE(relativeToleranceAsFloatingPoint(1.0, 1e-2).isWithin(fulp2));
184 EXPECT_FALSE(relativeToleranceAsFloatingPoint(1.0, 1e-9).isWithin(fulp2));
185 EXPECT_FALSE(relativeToleranceAsFloatingPoint(1.0, 1e-2).isWithin(fdiff));
186 EXPECT_TRUE(relativeToleranceAsFloatingPoint(1.0, 2e-2).isWithin(fdiff));
187 EXPECT_TRUE(relativeToleranceAsFloatingPoint(1.0, 2e-2).isWithin(fsmall));
188 EXPECT_FALSE(relativeToleranceAsFloatingPoint(0.1, 2e-2).isWithin(fsmall));
189 EXPECT_FALSE(relativeToleranceAsFloatingPoint(1.0, 2e-2).isWithin(fsmall2));
191 FloatingPointDifference dequal(1.0, 1.0);
192 FloatingPointDifference dulp2f(1.0, static_cast<double>(addUlps(1.0f, 2)));
193 FloatingPointDifference ddiff(1.0, 1.011);
194 FloatingPointDifference dsmall(0.1, 0.111);
195 FloatingPointDifference dsmall2(0.1, 0.121);
196 EXPECT_TRUE(relativeToleranceAsFloatingPoint(1.0, 1e-2).isWithin(dequal));
197 EXPECT_TRUE(relativeToleranceAsFloatingPoint(1.0, 1e-9).isWithin(dequal));
198 EXPECT_TRUE(relativeToleranceAsFloatingPoint(1.0, 1e-2).isWithin(dulp2f));
199 EXPECT_FALSE(relativeToleranceAsFloatingPoint(1.0, 1e-9).isWithin(dulp2f));
200 EXPECT_FALSE(relativeToleranceAsFloatingPoint(1.0, 1e-2).isWithin(ddiff));
201 EXPECT_TRUE(relativeToleranceAsFloatingPoint(1.0, 2e-2).isWithin(ddiff));
202 EXPECT_TRUE(relativeToleranceAsFloatingPoint(1.0, 2e-2).isWithin(dsmall));
203 EXPECT_FALSE(relativeToleranceAsFloatingPoint(0.1, 2e-2).isWithin(dsmall));
204 EXPECT_FALSE(relativeToleranceAsFloatingPoint(1.0, 2e-2).isWithin(dsmall2));
207 TEST(FloatingPointToleranceTest, RelativeToleranceAsUlp)
209 using gmx::test::relativeToleranceAsUlp;
211 FloatingPointDifference fequal(1.0f, 1.0f);
212 FloatingPointDifference fulp4(1.0f, addUlps(1.0f, 4));
213 FloatingPointDifference fsmall(0.1f, addUlps(1.0f, 2) - 0.9f);
214 FloatingPointDifference fsmall2(0.1f, addUlps(1.0f, 6) - 0.9f);
215 EXPECT_TRUE(relativeToleranceAsUlp(1.0, 2).isWithin(fequal));
216 EXPECT_FALSE(relativeToleranceAsUlp(1.0, 2).isWithin(fulp4));
217 EXPECT_TRUE(relativeToleranceAsUlp(1.0, 4).isWithin(fulp4));
218 EXPECT_TRUE(relativeToleranceAsUlp(1.0, 4).isWithin(fsmall));
219 EXPECT_FALSE(relativeToleranceAsUlp(0.1, 4).isWithin(fsmall));
220 EXPECT_FALSE(relativeToleranceAsUlp(1.0, 4).isWithin(fsmall2));
222 FloatingPointDifference dequal(1.0, 1.0);
223 FloatingPointDifference dulp4(1.0, addUlps(1.0, 4));
224 FloatingPointDifference dulp4f(1.0, static_cast<double>(addUlps(1.0f, 4)));
225 FloatingPointDifference dsmall(0.1, addUlps(1.0, 2) - 0.9);
226 FloatingPointDifference dsmall2(0.1, addUlps(1.0, 6) - 0.9);
227 EXPECT_TRUE(relativeToleranceAsUlp(1.0, 2).isWithin(dequal));
228 EXPECT_FALSE(relativeToleranceAsUlp(1.0, 2).isWithin(dulp4));
229 EXPECT_TRUE(relativeToleranceAsUlp(1.0, 4).isWithin(dulp4));
230 EXPECT_FALSE(relativeToleranceAsUlp(1.0, 4).isWithin(dulp4f));
231 EXPECT_TRUE(relativeToleranceAsUlp(1.0, 4).isWithin(dsmall));
232 EXPECT_FALSE(relativeToleranceAsUlp(0.1, 4).isWithin(dsmall));
233 EXPECT_FALSE(relativeToleranceAsUlp(1.0, 4).isWithin(dsmall2));
236 } // namespace