Access the device status directly, remove the getter
[alexxy/gromacs.git] / src / gromacs / gpu_utils / tests / hostallocator.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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 /*! \internal \file
36  * \brief
37  * Tests for GPU host allocator.
38  *
39  * \author Mark Abraham <mark.j.abraham@gmail.com>
40  */
41 #include "gmxpre.h"
42
43 #include "gromacs/gpu_utils/hostallocator.h"
44
45 #include "config.h"
46
47 #include <type_traits>
48 #include <vector>
49
50 #include <gtest/gtest.h>
51
52 #include "gromacs/gpu_utils/gpu_utils.h"
53 #include "gromacs/math/vectypes.h"
54 #include "gromacs/utility/arrayref.h"
55 #include "gromacs/utility/real.h"
56
57 #include "gromacs/math/tests/testarrayrefs.h"
58
59 #include "devicetransfers.h"
60 #include "gputest.h"
61
62 namespace gmx
63 {
64
65 namespace test
66 {
67
68 /*! \internal \brief Typed test fixture for infrastructure for
69  * host-side memory used for GPU transfers. */
70 template<typename T>
71 class HostMemoryTest : public test::GpuTest
72 {
73 public:
74     //! Convenience type
75     using ValueType = T;
76 };
77
78 /*! \brief Convenience function to transform a view into one with base
79  * type of (non-const) char.
80  *
81  * This transformation is useful for using containers with C APIs
82  * where the function signature is not declared const even where the
83  * semantics of the usage actually are const.
84  *
85  * \param[in]    data   The data pointer.
86  * \param[in]    size   The size of the data pointer (in T).
87  * \tparam       T      The base type of the container
88  * */
89 template<typename T>
90 ArrayRef<char> charArrayRefFromArray(T* data, size_t size)
91 {
92     // Make a type like T, but without its possible const qualifier.
93     using NonConstT = std::remove_const_t<T>;
94     return arrayRefFromArray<char>(reinterpret_cast<char*>(const_cast<NonConstT*>(data)),
95                                    size * sizeof(T));
96 }
97
98 //! Does a device transfer of \c input to the device in \c gpuInfo, and back to \c output.
99 template<typename T>
100 void runTest(const DeviceInformation& deviceInfo, ArrayRef<T> input, ArrayRef<T> output)
101 {
102     // Convert the views of input and output to flat non-const chars,
103     // so that there's no templating when we call doDeviceTransfers.
104     auto inputRef  = charArrayRefFromArray(input.data(), input.size());
105     auto outputRef = charArrayRefFromArray(output.data(), output.size());
106
107     ASSERT_EQ(inputRef.size(), outputRef.size());
108
109     doDeviceTransfers(deviceInfo, inputRef, outputRef);
110     compareViews(input, output);
111 }
112
113 struct MoveOnly
114 {
115     MoveOnly(real x = 0) : x(x) {}
116     MoveOnly(const MoveOnly&) = delete;
117     MoveOnly(MoveOnly&&)      = default;
118     MoveOnly& operator=(const MoveOnly&) = delete;
119     MoveOnly& operator=(MoveOnly&&) = default;
120     bool      operator==(const MoveOnly& o) const { return x == o.x; }
121     real      operator*=(int scaleFactor) { return x *= scaleFactor; }
122     real      x;
123 };
124
125 } // namespace test
126
127 namespace detail
128 {
129
130 template<>
131 struct PaddingTraits<test::MoveOnly>
132 {
133     using SimdBaseType                          = real;
134     static constexpr int maxSimdWidthOfBaseType = GMX_REAL_MAX_SIMD_WIDTH;
135 };
136
137 } // namespace detail
138
139 namespace test
140 {
141
142 //! The types used in testing of all operations.
143 typedef ::testing::Types<int32_t, real, RVec, test::MoveOnly> TestTypes;
144
145 //! Typed test fixture
146 template<typename T>
147 struct HostAllocatorTest : HostMemoryTest<T>
148 {
149     using VectorType = PaddedHostVector<T>; //!< PaddedHostVector of type tested
150 };
151 TYPED_TEST_CASE(HostAllocatorTest, TestTypes);
152
153 //! Typed test fixture (no mem/gpu initializtion - much faster)
154 template<typename T>
155 struct HostAllocatorTestNoMem : ::testing::Test
156 {
157     using VectorType = PaddedHostVector<T>; //!< PaddedHostVector of type tested
158 };
159 TYPED_TEST_CASE(HostAllocatorTestNoMem, TestTypes);
160
161 //! Typed test fixture for tests requiring a copyable type
162 template<typename T>
163 struct HostAllocatorTestNoMemCopyable : HostAllocatorTestNoMem<T>
164 {
165 };
166 //! The types used in testing minus move only types
167 using TestTypesCopyable = ::testing::Types<int32_t, real, RVec>;
168
169 TYPED_TEST_CASE(HostAllocatorTestNoMemCopyable, TestTypesCopyable);
170
171 //! Typed test fixture for tests requiring a copyable type
172 template<typename T>
173 using HostAllocatorTestCopyable = HostAllocatorTest<T>;
174 TYPED_TEST_CASE(HostAllocatorTestCopyable, TestTypesCopyable);
175
176 // Note that in GoogleTest typed tests, the use of TestFixture:: and
177 // this-> is sometimes required to get access to things in the fixture
178 // class (or its base classes).
179
180 // Note also that aspects of this code can be tested even when a GPU
181 // device is not available.
182
183 TYPED_TEST(HostAllocatorTest, EmptyMemoryAlwaysWorks)
184 {
185     typename TestFixture::VectorType v;
186 }
187
188 TYPED_TEST(HostAllocatorTestCopyable, VectorsWithDefaultHostAllocatorAlwaysWorks)
189 {
190     typename TestFixture::VectorType input(3), output;
191     output.resizeWithPadding(input.size());
192 }
193
194 // Several tests actually do CUDA transfers. This is not necessary
195 // because the state of page alignment or pinning is not currently
196 // relevant to the success of a CUDA transfer. CUDA checks happen only
197 // during cudaHostRegister and cudaHostUnregister. Such tests are of
198 // value only when this behaviour changes, if ever.
199
200 TYPED_TEST(HostAllocatorTestCopyable, TransfersWithoutPinningWork)
201 {
202     for (const DeviceInformation& compatibleDeviceInfo : getCompatibleDevices(this->deviceInfoList_))
203     {
204         typename TestFixture::VectorType input;
205         fillInput(&input, 1);
206         typename TestFixture::VectorType output;
207         output.resizeWithPadding(input.size());
208
209         runTest(compatibleDeviceInfo, makeArrayRef(input), makeArrayRef(output));
210     }
211 }
212
213 TYPED_TEST(HostAllocatorTestCopyable, FillInputAlsoWorksAfterCallingReserve)
214 {
215     typename TestFixture::VectorType input;
216     input.reserveWithPadding(3);
217     fillInput(&input, 1);
218 }
219
220 TYPED_TEST(HostAllocatorTestNoMem, CreateVector)
221 {
222     typename TestFixture::VectorType input1;
223     EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
224     typename TestFixture::VectorType input2({ PinningPolicy::PinnedIfSupported });
225     EXPECT_TRUE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
226 }
227
228 TYPED_TEST(HostAllocatorTestNoMem, MoveAssignment)
229 {
230     typename TestFixture::VectorType input1({ PinningPolicy::PinnedIfSupported });
231     input1 = typename TestFixture::VectorType();
232     EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
233
234     typename TestFixture::VectorType input2;
235     input2 = typename TestFixture::VectorType({ PinningPolicy::PinnedIfSupported });
236     EXPECT_TRUE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
237 }
238
239 TYPED_TEST(HostAllocatorTestNoMem, MoveConstruction)
240 {
241     typename TestFixture::VectorType input1;
242     typename TestFixture::VectorType input2(std::move(input1));
243     EXPECT_FALSE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
244
245     typename TestFixture::VectorType input3({ PinningPolicy::PinnedIfSupported });
246     typename TestFixture::VectorType input4(std::move(input3));
247     EXPECT_TRUE(input4.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
248 }
249
250 TYPED_TEST(HostAllocatorTestNoMemCopyable, CopyAssignment)
251 {
252     typename TestFixture::VectorType input1;
253     typename TestFixture::VectorType input2({ PinningPolicy::PinnedIfSupported });
254     input1 = input2;
255     EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
256     EXPECT_TRUE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
257     input2 = input1;
258     EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
259     EXPECT_TRUE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
260 }
261
262 TYPED_TEST(HostAllocatorTestNoMemCopyable, CopyConstruction)
263 {
264     typename TestFixture::VectorType input1;
265     typename TestFixture::VectorType input2(input1); //NOLINT(performance-unnecessary-copy-initialization)
266     EXPECT_FALSE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
267
268     typename TestFixture::VectorType input3({ PinningPolicy::PinnedIfSupported });
269     typename TestFixture::VectorType input4(input3); //NOLINT(performance-unnecessary-copy-initialization)
270     EXPECT_FALSE(input4.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
271 }
272
273 TYPED_TEST(HostAllocatorTestNoMem, Swap)
274 {
275     typename TestFixture::VectorType input1;
276     typename TestFixture::VectorType input2({ PinningPolicy::PinnedIfSupported });
277     std::swap(input1, input2);
278     EXPECT_TRUE(input1.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
279     EXPECT_FALSE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
280     std::swap(input2, input1);
281     EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
282     EXPECT_TRUE(input2.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
283 }
284
285 TYPED_TEST(HostAllocatorTestNoMem, Comparison)
286 {
287     using AllocatorType = typename TestFixture::VectorType::allocator_type;
288     EXPECT_EQ(AllocatorType{}, AllocatorType{});
289     // Should be false for different pinning policy
290     EXPECT_NE(AllocatorType{}, AllocatorType{ PinningPolicy::PinnedIfSupported });
291 }
292
293 #if GMX_GPU_CUDA
294
295 // Policy suitable for pinning is only supported for a CUDA build
296
297 TYPED_TEST(HostAllocatorTestCopyable, TransfersWithPinningWorkWithCuda)
298 {
299     for (auto& deviceInfo : this->deviceInfoList_)
300     {
301         typename TestFixture::VectorType input;
302         changePinningPolicy(&input, PinningPolicy::PinnedIfSupported);
303         fillInput(&input, 1);
304         typename TestFixture::VectorType output;
305         changePinningPolicy(&output, PinningPolicy::PinnedIfSupported);
306         output.resizeWithPadding(input.size());
307
308         runTest(*deviceInfo, makeArrayRef(input), makeArrayRef(output));
309     }
310 }
311
312 //! Helper function for wrapping a call to isHostMemoryPinned.
313 template<typename VectorType>
314 bool isPinned(const VectorType& v)
315 {
316     void* data = const_cast<void*>(static_cast<const void*>(v.data()));
317     return isHostMemoryPinned(data);
318 }
319
320 TYPED_TEST(HostAllocatorTestCopyable, ManualPinningOperationsWorkWithCuda)
321 {
322     if (!canComputeOnDevice())
323     {
324         return;
325     }
326
327     typename TestFixture::VectorType input;
328     changePinningPolicy(&input, PinningPolicy::PinnedIfSupported);
329     EXPECT_TRUE(input.get_allocator().pinningPolicy() == PinningPolicy::PinnedIfSupported);
330     EXPECT_EQ(0, input.size());
331     EXPECT_EQ(0, input.paddedSize());
332     EXPECT_TRUE(input.empty());
333     EXPECT_FALSE(isPinned(input)) << "should not be pinned before allocation";
334
335     // Fill some contents, which will be pinned because of the policy.
336     fillInput(&input, 1);
337     EXPECT_TRUE(isPinned(input)) << "should be pinned after allocation";
338
339     // Switching policy to CannotBePinned must unpin the buffer (via
340     // realloc and copy).
341     auto oldInputData = input.data();
342     changePinningPolicy(&input, PinningPolicy::CannotBePinned);
343     EXPECT_FALSE(isPinned(input)) << "should not be pinned after changing policy to CannotBePinned";
344     // These cannot be equal as both had to be allocated at the same
345     // time for the contents to be able to be copied.
346     EXPECT_NE(oldInputData, input.data());
347
348     // Switching policy to PinnedIfSupported must pin the buffer (via
349     // realloc and copy).
350     oldInputData = input.data();
351     changePinningPolicy(&input, PinningPolicy::PinnedIfSupported);
352     EXPECT_TRUE(isPinned(input)) << "should be pinned after changing policy to PinnedIfSupported";
353     // These cannot be equal as both had to be allocated at the same
354     // time for the contents to be able to be copied.
355     EXPECT_NE(oldInputData, input.data());
356 }
357
358 #endif
359
360 TYPED_TEST(HostAllocatorTest, StatefulAllocatorUsesMemory)
361 {
362     // The HostAllocator has state, so a container using it will be
363     // larger than a normal vector, whose default allocator is
364     // stateless.
365     EXPECT_LT(sizeof(std::vector<typename TestFixture::VectorType::value_type>),
366               sizeof(typename TestFixture::VectorType));
367 }
368
369 TEST(HostAllocatorUntypedTest, Comparison)
370 {
371     // Should always be true for the same policy, indpendent of value_type
372     EXPECT_EQ(HostAllocator<float>{}, HostAllocator<double>{});
373 }
374
375 //! Declare allocator types to test.
376 using AllocatorTypesToTest =
377         ::testing::Types<HostAllocator<real>, HostAllocator<int32_t>, HostAllocator<RVec>, HostAllocator<MoveOnly>>;
378
379 TYPED_TEST_CASE(AllocatorTest, AllocatorTypesToTest);
380
381 } // namespace test
382 } // namespace gmx
383
384 // Includes tests common to all allocation policies.
385 #include "gromacs/utility/tests/alignedallocator_impl.h"