689987bec79b65217b228af613bbfe8642628c1c
[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, 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 <type_traits>
46 #include <vector>
47
48 #include <gtest/gtest.h>
49
50 #include "gromacs/math/vectypes.h"
51 #include "gromacs/utility/real.h"
52
53 #include "devicetransfers.h"
54 #include "gputest.h"
55
56 namespace gmx
57 {
58
59 namespace
60 {
61
62 //! The types used in testing.
63 typedef ::testing::Types<int, real, RVec> TestTypes;
64
65 //! Typed test fixture
66 template <typename T>
67 class HostAllocatorTest : public test::GpuTest
68 {
69     public:
70         //! Convenience type
71         using ValueType = T;
72         //! Convenience type
73         using AllocatorType = HostAllocator<T>;
74         //! Convenience type
75         using VectorType = std::vector<ValueType, AllocatorType>;
76         //! Convenience type
77         using ViewType = ArrayRef<ValueType>;
78         //! Convenience type
79         using ConstViewType = ArrayRef<const ValueType>;
80         //! Prepare contents of a VectorType.
81         void fillInput(VectorType *input) const;
82         //! Compares input and output vectors.
83         void compareVectors(ConstViewType input,
84                             ConstViewType output) const;
85         //! Do some transfers and test the results.
86         void runTest(ConstViewType input, ViewType output) const;
87 };
88
89 // Already documented
90 template <typename T>
91 void HostAllocatorTest<T>::fillInput(VectorType *input) const
92 {
93     input->push_back(1);
94     input->push_back(2);
95     input->push_back(3);
96 }
97
98 //! Initialization specialization for RVec
99 template <>
100 void HostAllocatorTest<RVec>::fillInput(VectorType *input) const
101 {
102     input->push_back({1, 2, 3});
103 }
104
105 // Already documented
106 template <typename T>
107 void HostAllocatorTest<T>::compareVectors(ConstViewType input,
108                                           ConstViewType output) const
109 {
110     for (size_t i = 0; i != input.size(); ++i)
111     {
112         EXPECT_EQ(input[i], output[i]) << "for index " << i;
113     }
114 }
115
116 //! Comparison specialization for RVec
117 template <>
118 void HostAllocatorTest<RVec>::compareVectors(ConstViewType input,
119                                              ConstViewType output) const
120 {
121     for (size_t i = 0; i != input.size(); ++i)
122     {
123         EXPECT_EQ(input[i][XX], output[i][XX]) << "for index " << i;
124         EXPECT_EQ(input[i][YY], output[i][YY]) << "for index " << i;
125         EXPECT_EQ(input[i][ZZ], output[i][ZZ]) << "for index " << i;
126     }
127 }
128
129 /*! \brief Convenience function to transform a view into one with base
130  * type of (non-const) char.
131  *
132  * This transformation is useful for using containers with C APIs
133  * where the function signature is not declared const even where the
134  * semantics of the usage actually are const.
135  *
136  * \param[in]    data   The data pointer.
137  * \param[in]    size   The size of the data pointer (in T).
138  * \tparam       T      The base type of the container
139  * */
140 template <typename T>
141 ArrayRef<char> charArrayRefFromArray(T *data, size_t size)
142 {
143     // Make a type like T, but without its possible const qualifier.
144     using NonConstT = typename std::remove_const<T>::type;
145     return arrayRefFromArray<char>(reinterpret_cast<char *>(const_cast<NonConstT *>(data)), size * sizeof(T));
146 }
147
148 template <typename T>
149 void HostAllocatorTest<T>::runTest(ConstViewType input, ViewType output) const
150 {
151     // We can't do a test that does a transfer unless we have a
152     // compatible device.
153     if (!this->haveValidGpus())
154     {
155         return;
156     }
157
158     // Convert the views of input and output to flat non-const chars,
159     // so that there's no templating when we call doDeviceTransfers.
160     auto inputRef  = charArrayRefFromArray(input.data(), input.size());
161     auto outputRef = charArrayRefFromArray(output.data(), output.size());
162
163     doDeviceTransfers(*this->gpuInfo_, inputRef, outputRef);
164     this->compareVectors(input, output);
165 }
166
167 TYPED_TEST_CASE(HostAllocatorTest, TestTypes);
168
169 // Note that in GoogleTest typed tests, the use of TestFixture:: and
170 // this-> is sometimes required to get access to things in the fixture
171 // class (or its base classes).
172
173 // Note also that aspects of this code can be tested even when a GPU
174 // device is not available.
175
176 TYPED_TEST(HostAllocatorTest, EmptyMemoryAlwaysWorks)
177 {
178     typename TestFixture::VectorType v;
179 }
180
181 TYPED_TEST(HostAllocatorTest, TransfersUsingDefaultHostAllocatorWork)
182 {
183     typename TestFixture::VectorType input = {{1, 2, 3}}, output;
184     output.resize(input.size());
185
186     this->runTest(input, output);
187 }
188
189 TYPED_TEST(HostAllocatorTest, TransfersUsingNormalCpuHostAllocatorWork)
190 {
191     // Make an allocator with a 'normal CPU' allocation policy. This
192     // might be slower than another policy, but still works.
193     using AllocatorType       = typename TestFixture::AllocatorType;
194     using AllocatorPolicyType = typename AllocatorType::allocation_policy;
195     AllocatorPolicyType              policy(AllocatorPolicyType::Impl::AllocateAligned);
196     AllocatorType                    allocator(policy);
197
198     typename TestFixture::VectorType input(allocator);
199     this->fillInput(&input);
200     typename TestFixture::VectorType output(allocator);
201     output.resize(input.size());
202
203     this->runTest(input, output);
204 }
205
206 TYPED_TEST(HostAllocatorTest, TransfersUsingGpuHostAllocatorWork)
207 {
208     // Make an allocator with a 'for GPU' allocation policy. This
209     // should be more efficient, but we can't test that.
210     using AllocatorType       = typename TestFixture::AllocatorType;
211     using AllocatorPolicyType = typename AllocatorType::allocation_policy;
212     AllocatorPolicyType              policy(AllocatorPolicyType::Impl::AllocateForGpu);
213     AllocatorType                    allocator(policy);
214
215     typename TestFixture::VectorType input(allocator);
216     this->fillInput(&input);
217     typename TestFixture::VectorType output(allocator);
218     output.resize(input.size());
219
220     this->runTest(input, output);
221 }
222
223 TYPED_TEST(HostAllocatorTest, StatefulAllocatorUsesMemory)
224 {
225     // The HostAllocator has state, so a container using it will be
226     // larger than a normal vector, whose default allocator is
227     // stateless.
228     EXPECT_LT(sizeof(std::vector<typename TestFixture::ValueType>),
229               sizeof(typename TestFixture::VectorType));
230 }
231
232 } // namespace
233 } // namespace