Modernize syntax of enable_if and traits to use _t helpers
[gromacs.git] / src / gromacs / utility / tests / arrayref.cpp
blob6cdb84b0585403cafffd4127fe43ee4f98beecfc
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
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 Tests for gmx::ArrayRef.
38 * \author Mark Abraham <mark.j.abraham@gmail.com>
39 * \ingroup module_utility
41 #include "gmxpre.h"
43 #include "gromacs/utility/arrayref.h"
45 #include <vector>
47 #include <gtest/gtest.h>
49 #include "gromacs/utility/basedefinitions.h"
50 #include "gromacs/utility/real.h"
52 namespace gmx
55 namespace
58 TEST(EmptyArrayRefTest, IsEmpty)
60 ArrayRef<real> empty;
62 EXPECT_EQ(0U, empty.size());
63 EXPECT_TRUE(empty.empty());
66 TEST(EmptyConstArrayRefTest, IsEmpty)
68 ArrayRef<const real> empty;
70 EXPECT_EQ(0U, empty.size());
71 EXPECT_TRUE(empty.empty());
74 #ifdef GTEST_HAS_TYPED_TEST
76 //! Define the types that end up being available as TypeParam in the test cases for both kinds of ArrayRef
77 typedef ::testing::Types<
78 ArrayRef<char>,
79 ArrayRef<unsigned char>,
80 ArrayRef<int>,
81 ArrayRef<unsigned int>,
82 ArrayRef<long>,
83 ArrayRef<unsigned long>,
84 ArrayRef<int64_t>,
85 ArrayRef<uint64_t>,
86 ArrayRef<float>,
87 ArrayRef<double>,
88 ArrayRef<const char>,
89 ArrayRef<const unsigned char>,
90 ArrayRef<const int>,
91 ArrayRef<const unsigned int>,
92 ArrayRef<const long>,
93 ArrayRef<const unsigned long>,
94 ArrayRef<const int64_t>,
95 ArrayRef<const uint64_t>,
96 ArrayRef<const float>,
97 ArrayRef<const double>
98 > ArrayRefTypes;
100 constexpr index aSize = 3;
102 /*! \brief Permit all the tests to run on all kinds of ArrayRefs
104 * The main objective is to verify that all the different kinds of
105 * construction lead to the expected result. */
106 template <typename TypeParam>
107 class ArrayRefTest : public ::testing::Test
109 public:
110 typedef TypeParam ArrayRefType;
111 typedef typename ArrayRefType::value_type ValueType;
112 typedef std::remove_const_t<ValueType> NonConstValueType;
114 /*! \brief Run the same tests all the time
116 * Note that test cases must call this->runTests(), because
117 * that's how the derived-class templates that implement
118 * type-parameterized tests actually work. */
119 void runTests(ValueType *aData,
120 ArrayRefType &arrayRef)
122 ASSERT_EQ(aSize, arrayRef.size());
123 ASSERT_FALSE(arrayRef.empty());
124 EXPECT_EQ(aData, arrayRef.data());
125 EXPECT_EQ(a[0], arrayRef.front());
126 EXPECT_EQ(a[aSize-1], arrayRef.back());
127 for (index i = 0; i != aSize; ++i)
129 EXPECT_EQ(a[i], arrayRef[i]);
133 ValueType a[aSize] = { ValueType(1.2), ValueType(2.4), ValueType(3.1) };
134 NonConstValueType ma[aSize] = { ValueType(1.2), ValueType(2.4), ValueType(3.1) };
137 TYPED_TEST_CASE(ArrayRefTest, ArrayRefTypes);
140 TYPED_TEST(ArrayRefTest, MakeWithAssignmentWorks)
142 typename TestFixture::ArrayRefType arrayRef = this->a;
143 this->runTests(this->a, arrayRef);
146 TYPED_TEST(ArrayRefTest, MakeWithNonConstAssignmentWorks)
148 typename TestFixture::ArrayRefType arrayRef = this->ma;
149 this->runTests(this->ma, arrayRef);
152 TYPED_TEST(ArrayRefTest, ConstructWithTemplateConstructorWorks)
154 typename TestFixture::ArrayRefType arrayRef(this->a);
155 this->runTests(this->a, arrayRef);
158 TYPED_TEST(ArrayRefTest, ConstructWithNonConstTemplateConstructorWorks)
160 typename TestFixture::ArrayRefType arrayRef(this->ma);
161 this->runTests(this->ma, arrayRef);
164 TYPED_TEST(ArrayRefTest, ConstructFromPointersWorks)
166 typename TestFixture::ArrayRefType arrayRef(this->a, this->a + aSize);
167 this->runTests(this->a, arrayRef);
170 TYPED_TEST(ArrayRefTest, ConstructFromNonConstPointersWorks)
172 typename TestFixture::ArrayRefType arrayRef(this->ma, this->ma + aSize);
173 this->runTests(this->ma, arrayRef);
176 template<bool c, typename T>
177 using makeConstIf_t = std::conditional_t<c, const T, T>;
179 TYPED_TEST(ArrayRefTest, ConstructFromVectorWorks)
181 makeConstIf_t<std::is_const<typename TestFixture::ValueType>::value,
182 std::vector<typename TestFixture::NonConstValueType> > v(this->a, this->a + aSize);
183 typename TestFixture::ArrayRefType arrayRef(v);
184 this->runTests(v.data(), arrayRef);
187 TYPED_TEST(ArrayRefTest, ConstructFromNonConstVectorWorks)
189 std::vector<typename TestFixture::NonConstValueType> v(this->a, this->a + aSize);
190 typename TestFixture::ArrayRefType arrayRef(v);
191 this->runTests(v.data(), arrayRef);
194 //! Helper struct for the case actually used in mdrun signalling
195 template <typename T>
196 struct Helper
198 public:
199 T a[3];
200 int size;
203 /*! \brief Test of the case actually used in mdrun signalling
205 * There, we take a non-const struct-field array of static length and
206 * make an ArrayRef to it using the template constructor that is
207 * supposed to infer the length from the static size. This has
208 * been a problem (for a compiler that we no longer support),
209 * so we test it.
212 TYPED_TEST(ArrayRefTest, ConstructFromStructFieldWithTemplateConstructorWorks)
214 Helper<typename TestFixture::NonConstValueType> h;
215 h.size = aSize;
216 for (int i = 0; i != h.size; ++i)
218 h.a[i] = this->a[i];
220 typename TestFixture::ArrayRefType arrayRef(h.a);
221 this->runTests(h.a, arrayRef);
224 #else // GTEST_HAS_TYPED_TEST
226 /* A dummy test that at least signals that something is missing if one runs the
227 * unit test executable itself.
229 TEST(DISABLED_ArrayRefTest, GenericTests)
231 ADD_FAILURE()
232 << "Tests for generic ArrayRef functionality require support for "
233 << "Google Test typed tests, which was not available when the tests "
234 << "were compiled.";
237 #endif // GTEST_HAS_TYPED_TEST
239 } // namespace
241 } // namespace gmx