Merge branch release-2018
[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, 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 "devicetransfers.h"
58 #include "gputest.h"
59
60 namespace gmx
61 {
62
63 namespace
64 {
65
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
70 {
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;
86 };
87
88 // Already documented
89 template <typename T> template <typename VectorType>
90 void HostMemoryTest<T>::fillInput(VectorType *input) const
91 {
92     input->resize(3);
93     (*input)[0] = 1;
94     (*input)[1] = 2;
95     (*input)[2] = 3;
96 }
97
98 //! Initialization specialization for RVec
99 template <> template <typename VectorType>
100 void HostMemoryTest<RVec>::fillInput(VectorType *input) const
101 {
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};
107 }
108
109 // Already documented
110 template <typename T>
111 void HostMemoryTest<T>::compareVectors(ConstViewType input,
112                                        ConstViewType output) const
113 {
114     for (size_t i = 0; i != input.size(); ++i)
115     {
116         EXPECT_EQ(input[i], output[i]) << "for index " << i;
117     }
118 }
119
120 //! Comparison specialization for RVec
121 template <>
122 void HostMemoryTest<RVec>::compareVectors(ConstViewType input,
123                                           ConstViewType output) const
124 {
125     for (size_t i = 0; i != input.size(); ++i)
126     {
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;
130     }
131 }
132
133 /*! \brief Convenience function to transform a view into one with base
134  * type of (non-const) char.
135  *
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.
139  *
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)
146 {
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));
150 }
151
152 template <typename T>
153 void HostMemoryTest<T>::runTest(ConstViewType input, ViewType output) const
154 {
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());
159
160     doDeviceTransfers(*this->gpuInfo_, inputRef, outputRef);
161     this->compareVectors(input, output);
162 }
163
164 //! The types used in testing.
165 typedef ::testing::Types<int, real, RVec> TestTypes;
166
167 //! Typed test fixture
168 template <typename T>
169 class HostAllocatorTest : public HostMemoryTest<T>
170 {
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>;
178 };
179
180 TYPED_TEST_CASE(HostAllocatorTest, TestTypes);
181
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).
185
186 // Note also that aspects of this code can be tested even when a GPU
187 // device is not available.
188
189 TYPED_TEST(HostAllocatorTest, EmptyMemoryAlwaysWorks)
190 {
191     typename TestFixture::VectorType v;
192 }
193
194 TYPED_TEST(HostAllocatorTest, VectorsWithDefaultHostAllocatorAlwaysWorks)
195 {
196     typename TestFixture::VectorType input = {{1, 2, 3}}, output;
197     output.resize(input.size());
198 }
199
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.
205
206 TYPED_TEST(HostAllocatorTest, TransfersWithoutPinningWork)
207 {
208     typename TestFixture::VectorType input;
209     this->fillInput(&input);
210     typename TestFixture::VectorType output;
211     output.resize(input.size());
212
213     this->runTest(input, output);
214 }
215
216 TYPED_TEST(HostAllocatorTest, FillInputAlsoWorksAfterCallingReserve)
217 {
218     typename TestFixture::VectorType input;
219     input.reserve(3);
220     this->fillInput(&input);
221 }
222
223 #if GMX_GPU == GMX_GPU_CUDA
224
225 // Policy suitable for pinning is only supported for a CUDA build
226
227 TYPED_TEST(HostAllocatorTest, TransfersWithPinningWorkWithCuda)
228 {
229     if (!this->haveValidGpus())
230     {
231         return;
232     }
233
234     typename TestFixture::VectorType input;
235     changePinningPolicy(&input, PinningPolicy::CanBePinned);
236     this->fillInput(&input);
237     typename TestFixture::VectorType output;
238     changePinningPolicy(&output, PinningPolicy::CanBePinned);
239     output.resize(input.size());
240
241     this->runTest(input, output);
242 }
243
244 //! Helper function for wrapping a call to isHostMemoryPinned.
245 template <typename VectorType>
246 bool isPinned(const VectorType &v)
247 {
248     void *data = const_cast<void *>(static_cast<const void *>(v.data()));
249     return isHostMemoryPinned(data);
250 }
251
252 TYPED_TEST(HostAllocatorTest, ManualPinningOperationsWorkWithCuda)
253 {
254     if (!this->haveValidGpus())
255     {
256         return;
257     }
258
259     typename TestFixture::VectorType input;
260     changePinningPolicy(&input, PinningPolicy::CanBePinned);
261     EXPECT_FALSE(isPinned(input));
262
263     // Unpin before allocation is fine, but does nothing.
264     input.get_allocator().getPolicy().unpin();
265     EXPECT_FALSE(isPinned(input));
266
267     // Pin with no contents is fine, but does nothing.
268     input.get_allocator().getPolicy().pin();
269     EXPECT_FALSE(isPinned(input));
270
271     // Fill some contents, which will be pinned because of the policy.
272     this->fillInput(&input);
273     EXPECT_TRUE(isPinned(input));
274
275     // Unpin after pin is fine.
276     input.get_allocator().getPolicy().unpin();
277     EXPECT_FALSE(isPinned(input));
278
279     // Repeated unpin should be a no-op.
280     input.get_allocator().getPolicy().unpin();
281
282     // Pin after unpin is fine.
283     input.get_allocator().getPolicy().pin();
284     EXPECT_TRUE(isPinned(input));
285
286     // Repeated pin should be a no-op, and still pinned.
287     input.get_allocator().getPolicy().pin();
288     EXPECT_TRUE(isPinned(input));
289
290     // Switching policy to CannotBePinned must unpin the buffer (via
291     // realloc and copy).
292     auto oldInputData = input.data();
293     changePinningPolicy(&input, PinningPolicy::CannotBePinned);
294     EXPECT_FALSE(isPinned(input));
295     // These cannot be equal as both had to be allocated at the same
296     // time for the contents to be able to be copied.
297     EXPECT_NE(oldInputData, input.data());
298
299     // Switching policy to CanBePinned must pin the buffer (via
300     // realloc and copy).
301     oldInputData = input.data();
302     changePinningPolicy(&input, PinningPolicy::CanBePinned);
303     EXPECT_TRUE(isPinned(input));
304     // These cannot be equal as both had to be allocated at the same
305     // time for the contents to be able to be copied.
306     EXPECT_NE(oldInputData, input.data());
307 }
308
309 #else
310
311 TYPED_TEST(HostAllocatorTest, ChangingPinningPolicyRequiresCuda)
312 {
313     typename TestFixture::VectorType input;
314     EXPECT_DEATH(changePinningPolicy(&input, PinningPolicy::CanBePinned),
315                  ".*A suitable build of GROMACS.* is required.*");
316 }
317
318 TYPED_TEST(HostAllocatorTest, ManualPinningOperationsWorkEvenWithoutCuda)
319 {
320     typename TestFixture::VectorType input;
321
322     // Since the buffer can't be pinned and isn't pinned, and the
323     // calling code can't be unhappy about this, these are OK.
324     input.get_allocator().getPolicy().pin();
325     input.get_allocator().getPolicy().unpin();
326 }
327
328 #endif
329
330 TYPED_TEST(HostAllocatorTest, StatefulAllocatorUsesMemory)
331 {
332     // The HostAllocator has state, so a container using it will be
333     // larger than a normal vector, whose default allocator is
334     // stateless.
335     EXPECT_LT(sizeof(std::vector<typename TestFixture::ValueType>),
336               sizeof(typename TestFixture::VectorType));
337 }
338
339 //! Declare allocator types to test.
340 using AllocatorTypesToTest = ::testing::Types<HostAllocator<real>,
341                                               HostAllocator<int>,
342                                               HostAllocator<RVec>
343                                               >;
344
345 TYPED_TEST_CASE(AllocatorTest, AllocatorTypesToTest);
346
347 } // namespace
348 } // namespace
349
350 // Includes tests common to all allocation policies.
351 #include "gromacs/utility/tests/alignedallocator-impl.h"