Support pinning in HostAllocator
[gromacs.git] / src / gromacs / gpu_utils / tests / hostallocator.cpp
blob809cce2b5a593292c82126e9493dc83f9943de32
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2017, 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/math/vectypes.h"
54 #include "gromacs/utility/arrayref.h"
55 #include "gromacs/utility/real.h"
57 #include "devicetransfers.h"
58 #include "gputest.h"
60 namespace gmx
63 namespace
66 /*! \internal \brief Typed test fixture for infrastructure for
67 * host-side memory used for GPU transfers. */
68 template <typename T>
69 class HostMemoryTest : public test::GpuTest
71 public:
72 //! Convenience type
73 using ValueType = T;
74 //! Convenience type
75 using ViewType = ArrayRef<ValueType>;
76 //! Convenience type
77 using ConstViewType = ArrayRef<const ValueType>;
78 //! Prepare contents of a VectorType.
79 template <typename VectorType>
80 void fillInput(VectorType *input) const;
81 //! Compares input and output vectors.
82 void compareVectors(ConstViewType input,
83 ConstViewType output) const;
84 //! Do some transfers and test the results.
85 void runTest(ConstViewType input, ViewType output) const;
88 // Already documented
89 template <typename T> template <typename VectorType>
90 void HostMemoryTest<T>::fillInput(VectorType *input) const
92 input->resize(3);
93 (*input)[0] = 1;
94 (*input)[1] = 2;
95 (*input)[2] = 3;
98 //! Initialization specialization for RVec
99 template <> template <typename VectorType>
100 void HostMemoryTest<RVec>::fillInput(VectorType *input) const
102 input->reserve(3);
103 input->resize(3);
104 (*input)[0] = {1, 2, 3};
105 (*input)[1] = {4, 5, 6};
106 (*input)[2] = {7, 8, 9};
109 // Already documented
110 template <typename T>
111 void HostMemoryTest<T>::compareVectors(ConstViewType input,
112 ConstViewType output) const
114 for (size_t i = 0; i != input.size(); ++i)
116 EXPECT_EQ(input[i], output[i]) << "for index " << i;
120 //! Comparison specialization for RVec
121 template <>
122 void HostMemoryTest<RVec>::compareVectors(ConstViewType input,
123 ConstViewType output) const
125 for (size_t i = 0; i != input.size(); ++i)
127 EXPECT_EQ(input[i][XX], output[i][XX]) << "for index " << i;
128 EXPECT_EQ(input[i][YY], output[i][YY]) << "for index " << i;
129 EXPECT_EQ(input[i][ZZ], output[i][ZZ]) << "for index " << i;
133 /*! \brief Convenience function to transform a view into one with base
134 * type of (non-const) char.
136 * This transformation is useful for using containers with C APIs
137 * where the function signature is not declared const even where the
138 * semantics of the usage actually are const.
140 * \param[in] data The data pointer.
141 * \param[in] size The size of the data pointer (in T).
142 * \tparam T The base type of the container
143 * */
144 template <typename T>
145 ArrayRef<char> charArrayRefFromArray(T *data, size_t size)
147 // Make a type like T, but without its possible const qualifier.
148 using NonConstT = typename std::remove_const<T>::type;
149 return arrayRefFromArray<char>(reinterpret_cast<char *>(const_cast<NonConstT *>(data)), size * sizeof(T));
152 template <typename T>
153 void HostMemoryTest<T>::runTest(ConstViewType input, ViewType output) const
155 // Convert the views of input and output to flat non-const chars,
156 // so that there's no templating when we call doDeviceTransfers.
157 auto inputRef = charArrayRefFromArray(input.data(), input.size());
158 auto outputRef = charArrayRefFromArray(output.data(), output.size());
160 doDeviceTransfers(*this->gpuInfo_, inputRef, outputRef);
161 this->compareVectors(input, output);
164 //! The types used in testing.
165 typedef ::testing::Types<int, real, RVec> TestTypes;
167 //! Typed test fixture
168 template <typename T>
169 class HostAllocatorTest : public HostMemoryTest<T>
171 public:
172 //! Convenience type
173 using ValueType = T;
174 //! Convenience type
175 using AllocatorType = HostAllocator<T>;
176 //! Convenience type
177 using VectorType = std::vector<ValueType, AllocatorType>;
180 TYPED_TEST_CASE(HostAllocatorTest, TestTypes);
182 // Note that in GoogleTest typed tests, the use of TestFixture:: and
183 // this-> is sometimes required to get access to things in the fixture
184 // class (or its base classes).
186 // Note also that aspects of this code can be tested even when a GPU
187 // device is not available.
189 TYPED_TEST(HostAllocatorTest, EmptyMemoryAlwaysWorks)
191 typename TestFixture::VectorType v;
194 TYPED_TEST(HostAllocatorTest, VectorsWithDefaultHostAllocatorAlwaysWorks)
196 typename TestFixture::VectorType input = {{1, 2, 3}}, output;
197 output.resize(input.size());
200 // Several tests actually do CUDA transfers. This is not necessary
201 // because the state of page alignment or pinning is not currently
202 // relevant to the success of a CUDA transfer. CUDA checks happen only
203 // during cudaHostRegister and cudaHostUnregister. Such tests are of
204 // value only when this behaviour changes, if ever.
206 TYPED_TEST(HostAllocatorTest, TransfersWithoutPinningWork)
208 typename TestFixture::VectorType input;
209 this->fillInput(&input);
210 typename TestFixture::VectorType output;
211 output.resize(input.size());
213 this->runTest(input, output);
216 TYPED_TEST(HostAllocatorTest, FillInputAlsoWorksAfterCallingReserve)
218 typename TestFixture::VectorType input;
219 input.reserve(3);
220 this->fillInput(&input);
223 #if GMX_GPU == GMX_GPU_CUDA
225 // Policy suitable for pinning is only supported for a CUDA build
227 TYPED_TEST(HostAllocatorTest, TransfersWithPinningWorkWithCuda)
229 typename TestFixture::VectorType input;
230 changePinningPolicy(&input, PinningPolicy::CanBePinned);
231 this->fillInput(&input);
232 typename TestFixture::VectorType output;
233 changePinningPolicy(&output, PinningPolicy::CanBePinned);
234 output.resize(input.size());
236 this->runTest(input, output);
239 //! Helper function for wrapping a call to isHostMemoryPinned.
240 template <typename VectorType>
241 bool isPinned(const VectorType &v)
243 void *data = const_cast<void *>(static_cast<const void *>(v.data()));
244 return isHostMemoryPinned(data);
247 TYPED_TEST(HostAllocatorTest, ManualPinningOperationsWorkWithCuda)
249 typename TestFixture::VectorType input;
250 changePinningPolicy(&input, PinningPolicy::CanBePinned);
251 EXPECT_FALSE(isPinned(input));
253 // Unpin before allocation is fine, but does nothing.
254 input.get_allocator().getPolicy().unpin();
255 EXPECT_FALSE(isPinned(input));
257 // Pin with no contents is fine, but does nothing.
258 input.get_allocator().getPolicy().pin();
259 EXPECT_FALSE(isPinned(input));
261 // Fill some contents, which will be pinned because of the policy.
262 this->fillInput(&input);
263 EXPECT_TRUE(isPinned(input));
265 // Unpin after pin is fine.
266 input.get_allocator().getPolicy().unpin();
267 EXPECT_FALSE(isPinned(input));
269 // Repeated unpin should be a no-op.
270 input.get_allocator().getPolicy().unpin();
272 // Pin after unpin is fine.
273 input.get_allocator().getPolicy().pin();
274 EXPECT_TRUE(isPinned(input));
276 // Repeated pin should be a no-op, and still pinned.
277 input.get_allocator().getPolicy().pin();
278 EXPECT_TRUE(isPinned(input));
280 // Switching policy to CannotBePinned must unpin the buffer (via
281 // realloc and copy).
282 auto oldInputData = input.data();
283 changePinningPolicy(&input, PinningPolicy::CannotBePinned);
284 EXPECT_FALSE(isPinned(input));
285 // These cannot be equal as both had to be allocated at the same
286 // time for the contents to be able to be copied.
287 EXPECT_NE(oldInputData, input.data());
289 // Switching policy to CanBePinned must pin the buffer (via
290 // realloc and copy).
291 oldInputData = input.data();
292 changePinningPolicy(&input, PinningPolicy::CanBePinned);
293 EXPECT_TRUE(isPinned(input));
294 // These cannot be equal as both had to be allocated at the same
295 // time for the contents to be able to be copied.
296 EXPECT_NE(oldInputData, input.data());
299 #else
301 TYPED_TEST(HostAllocatorTest, ChangingPinningPolicyRequiresCuda)
303 typename TestFixture::VectorType input;
304 EXPECT_DEATH(changePinningPolicy(&input, PinningPolicy::CanBePinned),
305 ".*A suitable build of GROMACS.* is required.*");
308 TYPED_TEST(HostAllocatorTest, ManualPinningOperationsWorkEvenWithoutCuda)
310 typename TestFixture::VectorType input;
312 // Since the buffer can't be pinned and isn't pinned, and the
313 // calling code can't be unhappy about this, these are OK.
314 input.get_allocator().getPolicy().pin();
315 input.get_allocator().getPolicy().unpin();
318 #endif
320 TYPED_TEST(HostAllocatorTest, StatefulAllocatorUsesMemory)
322 // The HostAllocator has state, so a container using it will be
323 // larger than a normal vector, whose default allocator is
324 // stateless.
325 EXPECT_LT(sizeof(std::vector<typename TestFixture::ValueType>),
326 sizeof(typename TestFixture::VectorType));
329 //! Declare allocator types to test.
330 using AllocatorTypesToTest = ::testing::Types<HostAllocator<real>,
331 HostAllocator<int>,
332 HostAllocator<RVec>
335 TYPED_TEST_CASE(AllocatorTest, AllocatorTypesToTest);
337 } // namespace
338 } // namespace
340 // Includes tests common to all allocation policies.
341 #include "gromacs/utility/tests/alignedallocator-impl.h"