Use common device context header for OpenCL
[gromacs.git] / src / gromacs / gpu_utils / tests / hostallocator.cpp
blob173ca15ee300ccba1c15093351769c1aa55c289e
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
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 for GPU host allocator.
39 * \author Mark Abraham <mark.j.abraham@gmail.com>
41 #include "gmxpre.h"
43 #include "gromacs/gpu_utils/hostallocator.h"
45 #include "config.h"
47 #include <type_traits>
48 #include <vector>
50 #include <gtest/gtest.h>
52 #include "gromacs/gpu_utils/gpu_utils.h"
53 #include "gromacs/hardware/device_management.h"
54 #include "gromacs/math/vectypes.h"
55 #include "gromacs/utility/arrayref.h"
56 #include "gromacs/utility/real.h"
58 #include "gromacs/math/tests/testarrayrefs.h"
59 #include "testutils/test_hardware_environment.h"
61 #include "devicetransfers.h"
63 namespace gmx
66 namespace test
69 /*! \internal \brief Typed test fixture for infrastructure for
70 * host-side memory used for GPU transfers. */
71 template<typename T>
72 class HostMemoryTest : public ::testing::Test
74 public:
75 //! Convenience type
76 using ValueType = T;
79 /*! \brief Convenience function to transform a view into one with base
80 * type of (non-const) char.
82 * This transformation is useful for using containers with C APIs
83 * where the function signature is not declared const even where the
84 * semantics of the usage actually are const.
86 * \param[in] data The data pointer.
87 * \param[in] size The size of the data pointer (in T).
88 * \tparam T The base type of the container
89 * */
90 template<typename T>
91 ArrayRef<char> charArrayRefFromArray(T* data, size_t size)
93 // Make a type like T, but without its possible const qualifier.
94 using NonConstT = std::remove_const_t<T>;
95 return arrayRefFromArray<char>(reinterpret_cast<char*>(const_cast<NonConstT*>(data)),
96 size * sizeof(T));
99 //! Does a device transfer of \c input to the device in \c gpuInfo, and back to \c output.
100 template<typename T>
101 void runTest(const DeviceInformation& deviceInfo, ArrayRef<T> input, ArrayRef<T> output)
103 // Convert the views of input and output to flat non-const chars,
104 // so that there's no templating when we call doDeviceTransfers.
105 auto inputRef = charArrayRefFromArray(input.data(), input.size());
106 auto outputRef = charArrayRefFromArray(output.data(), output.size());
108 ASSERT_EQ(inputRef.size(), outputRef.size());
110 doDeviceTransfers(deviceInfo, inputRef, outputRef);
111 compareViews(input, output);
114 struct MoveOnly
116 MoveOnly(real x = 0) : x(x) {}
117 MoveOnly(const MoveOnly&) = delete;
118 MoveOnly(MoveOnly&&) = default;
119 MoveOnly& operator=(const MoveOnly&) = delete;
120 MoveOnly& operator=(MoveOnly&&) = default;
121 bool operator==(const MoveOnly& o) const { return x == o.x; }
122 real operator*=(int scaleFactor) { return x *= scaleFactor; }
123 real x;
126 } // namespace test
128 namespace detail
131 template<>
132 struct PaddingTraits<test::MoveOnly>
134 using SimdBaseType = real;
135 static constexpr int maxSimdWidthOfBaseType = GMX_REAL_MAX_SIMD_WIDTH;
138 } // namespace detail
140 namespace test
143 //! The types used in testing of all operations.
144 typedef ::testing::Types<int32_t, real, RVec, test::MoveOnly> TestTypes;
146 //! Typed test fixture
147 template<typename T>
148 struct HostAllocatorTest : HostMemoryTest<T>
150 using VectorType = PaddedHostVector<T>; //!< PaddedHostVector of type tested
152 TYPED_TEST_CASE(HostAllocatorTest, TestTypes);
154 //! Typed test fixture (no mem/gpu initializtion - much faster)
155 template<typename T>
156 struct HostAllocatorTestNoMem : ::testing::Test
158 using VectorType = PaddedHostVector<T>; //!< PaddedHostVector of type tested
160 TYPED_TEST_CASE(HostAllocatorTestNoMem, TestTypes);
162 //! Typed test fixture for tests requiring a copyable type
163 template<typename T>
164 struct HostAllocatorTestNoMemCopyable : HostAllocatorTestNoMem<T>
167 //! The types used in testing minus move only types
168 using TestTypesCopyable = ::testing::Types<int32_t, real, RVec>;
170 TYPED_TEST_CASE(HostAllocatorTestNoMemCopyable, TestTypesCopyable);
172 //! Typed test fixture for tests requiring a copyable type
173 template<typename T>
174 using HostAllocatorTestCopyable = HostAllocatorTest<T>;
175 TYPED_TEST_CASE(HostAllocatorTestCopyable, TestTypesCopyable);
177 // Note that in GoogleTest typed tests, the use of TestFixture:: and
178 // this-> is sometimes required to get access to things in the fixture
179 // class (or its base classes).
181 // Note also that aspects of this code can be tested even when a GPU
182 // device is not available.
184 TYPED_TEST(HostAllocatorTest, EmptyMemoryAlwaysWorks)
186 typename TestFixture::VectorType v;
189 TYPED_TEST(HostAllocatorTestCopyable, VectorsWithDefaultHostAllocatorAlwaysWorks)
191 typename TestFixture::VectorType input(3), output;
192 output.resizeWithPadding(input.size());
195 // Several tests actually do CUDA transfers. This is not necessary
196 // because the state of page alignment or pinning is not currently
197 // relevant to the success of a CUDA transfer. CUDA checks happen only
198 // during cudaHostRegister and cudaHostUnregister. Such tests are of
199 // value only when this behaviour changes, if ever.
201 TYPED_TEST(HostAllocatorTestCopyable, TransfersWithoutPinningWork)
203 for (const auto& testDevice : getTestHardwareEnvironment()->getTestDeviceList())
205 setActiveDevice(testDevice->deviceInfo());
206 typename TestFixture::VectorType input;
207 fillInput(&input, 1);
208 typename TestFixture::VectorType output;
209 output.resizeWithPadding(input.size());
211 runTest(testDevice->deviceInfo(), makeArrayRef(input), makeArrayRef(output));
215 TYPED_TEST(HostAllocatorTestCopyable, FillInputAlsoWorksAfterCallingReserve)
217 typename TestFixture::VectorType input;
218 input.reserveWithPadding(3);
219 fillInput(&input, 1);
222 TYPED_TEST(HostAllocatorTestNoMem, CreateVector)
224 typename TestFixture::VectorType input1;
225 EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
226 typename TestFixture::VectorType input2({ PinningPolicy::PinnedIfSupported });
227 EXPECT_TRUE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
230 TYPED_TEST(HostAllocatorTestNoMem, MoveAssignment)
232 typename TestFixture::VectorType input1({ PinningPolicy::PinnedIfSupported });
233 input1 = typename TestFixture::VectorType();
234 EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
236 typename TestFixture::VectorType input2;
237 input2 = typename TestFixture::VectorType({ PinningPolicy::PinnedIfSupported });
238 EXPECT_TRUE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
241 TYPED_TEST(HostAllocatorTestNoMem, MoveConstruction)
243 typename TestFixture::VectorType input1;
244 typename TestFixture::VectorType input2(std::move(input1));
245 EXPECT_FALSE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
247 typename TestFixture::VectorType input3({ PinningPolicy::PinnedIfSupported });
248 typename TestFixture::VectorType input4(std::move(input3));
249 EXPECT_TRUE(input4.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
252 TYPED_TEST(HostAllocatorTestNoMemCopyable, CopyAssignment)
254 typename TestFixture::VectorType input1;
255 typename TestFixture::VectorType input2({ PinningPolicy::PinnedIfSupported });
256 input1 = input2;
257 EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
258 EXPECT_TRUE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
259 input2 = input1;
260 EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
261 EXPECT_TRUE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
264 TYPED_TEST(HostAllocatorTestNoMemCopyable, CopyConstruction)
266 typename TestFixture::VectorType input1;
267 typename TestFixture::VectorType input2(input1); //NOLINT(performance-unnecessary-copy-initialization)
268 EXPECT_FALSE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
270 typename TestFixture::VectorType input3({ PinningPolicy::PinnedIfSupported });
271 typename TestFixture::VectorType input4(input3); //NOLINT(performance-unnecessary-copy-initialization)
272 EXPECT_FALSE(input4.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
275 TYPED_TEST(HostAllocatorTestNoMem, Swap)
277 typename TestFixture::VectorType input1;
278 typename TestFixture::VectorType input2({ PinningPolicy::PinnedIfSupported });
279 std::swap(input1, input2);
280 EXPECT_TRUE(input1.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
281 EXPECT_FALSE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
282 std::swap(input2, input1);
283 EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
284 EXPECT_TRUE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
287 TYPED_TEST(HostAllocatorTestNoMem, Comparison)
289 using AllocatorType = typename TestFixture::VectorType::allocator_type;
290 EXPECT_EQ(AllocatorType{}, AllocatorType{});
291 // Should be false for different pinning policy
292 EXPECT_NE(AllocatorType{}, AllocatorType{ PinningPolicy::PinnedIfSupported });
295 #if GMX_GPU_CUDA
297 // Policy suitable for pinning is only supported for a CUDA build
299 TYPED_TEST(HostAllocatorTestCopyable, TransfersWithPinningWorkWithCuda)
301 for (const auto& testDevice : getTestHardwareEnvironment()->getTestDeviceList())
303 setActiveDevice(testDevice->deviceInfo());
304 typename TestFixture::VectorType input;
305 changePinningPolicy(&input, PinningPolicy::PinnedIfSupported);
306 fillInput(&input, 1);
307 typename TestFixture::VectorType output;
308 changePinningPolicy(&output, PinningPolicy::PinnedIfSupported);
309 output.resizeWithPadding(input.size());
311 runTest(testDevice->deviceInfo(), makeArrayRef(input), makeArrayRef(output));
315 //! Helper function for wrapping a call to isHostMemoryPinned.
316 template<typename VectorType>
317 bool isPinned(const VectorType& v)
319 void* data = const_cast<void*>(static_cast<const void*>(v.data()));
320 return isHostMemoryPinned(data);
323 TYPED_TEST(HostAllocatorTestCopyable, ManualPinningOperationsWorkWithCuda)
325 for (const auto& testDevice : getTestHardwareEnvironment()->getTestDeviceList())
327 setActiveDevice(testDevice->deviceInfo());
328 typename TestFixture::VectorType input;
329 changePinningPolicy(&input, PinningPolicy::PinnedIfSupported);
330 EXPECT_TRUE(input.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
331 EXPECT_TRUE(input.empty());
332 fillInput(&input, 1);
333 // realloc and copy).
334 auto oldInputData = input.data();
335 changePinningPolicy(&input, PinningPolicy::CannotBePinned);
336 EXPECT_FALSE(isPinned(input));
337 // These cannot be equal as both had to be allocated at the same
338 // time for the contents to be able to be copied.
339 EXPECT_NE(oldInputData, input.data());
341 // Switching policy to PinnedIfSupported must pin the buffer (via
342 // realloc and copy).
343 oldInputData = input.data();
344 changePinningPolicy(&input, PinningPolicy::PinnedIfSupported);
345 EXPECT_TRUE(isPinned(input));
346 // These cannot be equal as both had to be allocated at the same
347 // time for the contents to be able to be copied.
348 EXPECT_NE(oldInputData, input.data());
352 #endif
354 TYPED_TEST(HostAllocatorTest, StatefulAllocatorUsesMemory)
356 // The HostAllocator has state, so a container using it will be
357 // larger than a normal vector, whose default allocator is
358 // stateless.
359 EXPECT_LT(sizeof(std::vector<typename TestFixture::VectorType::value_type>),
360 sizeof(typename TestFixture::VectorType));
363 TEST(HostAllocatorUntypedTest, Comparison)
365 // Should always be true for the same policy, indpendent of value_type
366 EXPECT_EQ(HostAllocator<float>{}, HostAllocator<double>{});
369 //! Declare allocator types to test.
370 using AllocatorTypesToTest =
371 ::testing::Types<HostAllocator<real>, HostAllocator<int32_t>, HostAllocator<RVec>, HostAllocator<MoveOnly>>;
373 TYPED_TEST_CASE(AllocatorTest, AllocatorTypesToTest);
375 } // namespace test
376 } // namespace gmx
378 // Includes tests common to all allocation policies.
379 #include "gromacs/utility/tests/alignedallocator_impl.h"