Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / mdspan / tests / mdspan.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2018,2019, 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 /*
36  * This file is a modified version of original work of Sandia Corporation.
37  * In the spirit of the original code, this particular file can be distributed
38  * on the terms of Sandia Corporation.
39  */
40 /*
41  *                         Kokkos v. 2.0
42  *               Copyright (2014) Sandia Corporation
43  *
44  * Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
45  * the U.S. Government retains certain rights in this software.
46  *
47  * Kokkos is licensed under 3-clause BSD terms of use:
48  *
49  * Redistribution and use in source and binary forms, with or without
50  * modification, are permitted provided that the following conditions are
51  * met:
52  *
53  * 1. Redistributions of source code must retain the above copyright
54  * notice, this list of conditions and the following disclaimer.
55  *
56  * 2. Redistributions in binary form must reproduce the above copyright
57  * notice, this list of conditions and the following disclaimer in the
58  * documentation and/or other materials provided with the distribution.
59  *
60  * 3. Neither the name of the Corporation nor the names of the
61  * contributors may be used to endorse or promote products derived from
62  * this software without specific prior written permission.
63  *
64  * THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
65  * EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
66  * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
67  * PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
68  * CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
69  * EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
70  * PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
71  * PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
72  * LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
73  * NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
74  * SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
75  *
76  * Questions? Contact Christian R. Trott (crtrott@sandia.gov)
77  */
78 /*! \internal \file
79  * \brief Testing gmx::basic_mdspan.
80  *
81  * \author Christian Trott <crtrott@sandia.gov>
82  * \author Carter Edwards <hedwards@nvidia.com>
83  * \author David Hollman <dshollm@sandia.gov>
84  * \author Christian Blau <cblau@gwdg.de>
85  */
86 #include "gmxpre.h"
87
88 #include "gromacs/mdspan/mdspan.h"
89
90 #include <cstdio>
91
92 #include <gtest/gtest.h>
93
94 namespace gmx
95 {
96
97 namespace
98 {
99
100
101 // Test basic_mdspan with a mixture of dynamic and static extents, as well as extent of one.
102 // The dynamic extents will be set to 4 and 2 repsectively so that we'll tests
103 // a multidimensional array of extent 5 x 4 x 3 x 2 x 1, i.e. 120 elements
104
105 //! View on int data with mixed static and dynamic extents
106 using mdspan_int =
107         basic_mdspan<int, extents<5, dynamic_extent, 3, dynamic_extent, 1>, layout_right, accessor_basic<int>>;
108
109 TEST(MdSpanTest, MdSpanWrapsBasicMdSpanCorrectly)
110 {
111     // Check that mdspan wraps basic_mdspan as expected
112     ::testing::StaticAssertTypeEq<mdspan_int, mdspan<int, 5, dynamic_extent, 3, dynamic_extent, 1>>();
113 }
114
115 //! View on float data with mixed static and dynamic extents
116 using mdspan_float =
117         basic_mdspan<float, extents<5, dynamic_extent, 3, dynamic_extent, 1>, layout_right, accessor_basic<float>>;
118 //! Types to be tested
119 using MdSpanTypes = ::testing::Types<mdspan_int, mdspan_float>;
120
121 template<class ElementType, class Mapping>
122 struct fill_raw_data
123 {
124     static void fill(ElementType* p, Mapping m)
125     {
126         typename Mapping::extents_type e = m.extents();
127         for (ptrdiff_t i0 = 0; i0 < e.extent(0); i0++)
128         {
129             for (ptrdiff_t i1 = 0; i1 < e.extent(1); i1++)
130             {
131                 for (ptrdiff_t i2 = 0; i2 < e.extent(2); i2++)
132                 {
133                     for (ptrdiff_t i3 = 0; i3 < e.extent(3); i3++)
134                     {
135                         for (ptrdiff_t i4 = 0; i4 < e.extent(4); i4++)
136                         {
137                             p[i0 * m.stride(0) + i1 * m.stride(1) + i2 * m.stride(2)
138                               + i3 * m.stride(3) + i4 * m.stride(4)] =
139                                     ElementType(i0 * 10000 + i1 * 1000 + i2 * 100 + i3 * 10 + i4);
140                         }
141                     }
142                 }
143             }
144         }
145     }
146 };
147 template<class MDSPAN>
148 struct MdSpanTest : public ::testing::Test
149 {
150     using mdspan_type   = MDSPAN;
151     using element_type  = typename mdspan_type::element_type;
152     using extents_type  = typename mdspan_type::extents_type;
153     using mapping_type  = typename mdspan_type::mapping_type;
154     using accessor_type = typename mdspan_type::accessor_type;
155     using pointer_type  = typename mdspan_type::pointer;
156
157     mdspan_type my_mdspan_extents;
158     mdspan_type my_mdspan_array;
159     mdspan_type my_mdspan_mapping;
160     mdspan_type my_mdspan_map_acc;
161     mdspan_type my_mdspan_copy;
162
163     std::vector<element_type> rawData;
164
165     template<class... ED>
166     void SetUp(ED... e)
167     {
168         mapping_type  map{ extents_type(e...) };
169         accessor_type acc;
170         rawData.resize(map.required_span_size());
171         fill_raw_data<element_type, mapping_type>::fill(rawData.data(), map);
172         pointer_type p(rawData.data());
173
174         my_mdspan_array   = mdspan_type(p, std::array<ptrdiff_t, sizeof...(ED)>({ { e... } }));
175         my_mdspan_mapping = mdspan_type(p, map);
176         my_mdspan_map_acc = mdspan_type(p, map, acc);
177         my_mdspan_extents = mdspan_type(p, e...);
178         my_mdspan_copy    = my_mdspan_mapping;
179     }
180
181     void check_rank(ptrdiff_t r)
182     {
183         EXPECT_EQ(my_mdspan_mapping.rank(), r);
184         EXPECT_EQ(my_mdspan_map_acc.rank(), r);
185         EXPECT_EQ(my_mdspan_extents.rank(), r);
186         EXPECT_EQ(my_mdspan_copy.rank(), r);
187     }
188     void check_rank_dynamic(ptrdiff_t r)
189     {
190         EXPECT_EQ(my_mdspan_mapping.rank_dynamic(), r);
191         EXPECT_EQ(my_mdspan_map_acc.rank_dynamic(), r);
192         EXPECT_EQ(my_mdspan_extents.rank_dynamic(), r);
193         EXPECT_EQ(my_mdspan_copy.rank_dynamic(), r);
194     }
195     template<class... E>
196     void check_extents(E... e)
197     {
198         std::array<ptrdiff_t, extents_type::rank()> a{ { e... } };
199         for (size_t r = 0; r < extents_type::rank(); r++)
200         {
201             EXPECT_EQ(my_mdspan_mapping.extent(r), a[r]);
202             EXPECT_EQ(my_mdspan_map_acc.extent(r), a[r]);
203             EXPECT_EQ(my_mdspan_extents.extent(r), a[r]);
204             EXPECT_EQ(my_mdspan_copy.extent(r), a[r]);
205         }
206     }
207     template<class... E>
208     void check_strides(E... e)
209     {
210         std::array<ptrdiff_t, extents_type::rank()> a{ { e... } };
211         for (size_t r = 0; r < extents_type::rank(); r++)
212         {
213             EXPECT_EQ(my_mdspan_mapping.stride(r), a[r]);
214             EXPECT_EQ(my_mdspan_map_acc.stride(r), a[r]);
215             EXPECT_EQ(my_mdspan_extents.stride(r), a[r]);
216             EXPECT_EQ(my_mdspan_copy.stride(r), a[r]);
217         }
218     }
219
220     void check_properties_internal(mdspan_type my_mdspan,
221                                    bool        always_unique,
222                                    bool        always_contiguous,
223                                    bool        always_strided,
224                                    bool        unique,
225                                    bool        contiguous,
226                                    bool        strided)
227     {
228         EXPECT_EQ(my_mdspan.is_always_unique(), always_unique);
229         EXPECT_EQ(my_mdspan.is_always_contiguous(), always_contiguous);
230         EXPECT_EQ(my_mdspan.is_always_strided(), always_strided);
231         EXPECT_EQ(my_mdspan.is_unique(), unique);
232         EXPECT_EQ(my_mdspan.is_contiguous(), contiguous);
233         EXPECT_EQ(my_mdspan.is_strided(), strided);
234     }
235
236     void check_properties(bool always_unique,
237                           bool always_contiguous,
238                           bool always_strided,
239                           bool unique,
240                           bool contiguous,
241                           bool strided)
242     {
243         check_properties_internal(my_mdspan_mapping, always_unique, always_contiguous,
244                                   always_strided, unique, contiguous, strided);
245         check_properties_internal(my_mdspan_map_acc, always_unique, always_contiguous,
246                                   always_strided, unique, contiguous, strided);
247         check_properties_internal(my_mdspan_extents, always_unique, always_contiguous,
248                                   always_strided, unique, contiguous, strided);
249         check_properties_internal(my_mdspan_copy, always_unique, always_contiguous, always_strided,
250                                   unique, contiguous, strided);
251     }
252
253     void check_operator()
254     {
255         extents_type e = my_mdspan_mapping.extents();
256         for (ptrdiff_t i0 = 0; i0 < e.extent(0); i0++)
257         {
258             for (ptrdiff_t i1 = 0; i1 < e.extent(1); i1++)
259             {
260                 for (ptrdiff_t i2 = 0; i2 < e.extent(2); i2++)
261                 {
262                     for (ptrdiff_t i3 = 0; i3 < e.extent(3); i3++)
263                     {
264                         for (ptrdiff_t i4 = 0; i4 < e.extent(4); i4++)
265                         {
266                             element_type value = i0 * 10000 + i1 * 1000 + i2 * 100 + i3 * 10 + i4;
267                             EXPECT_EQ(my_mdspan_mapping(i0, i1, i2, i3, i4), value);
268                             EXPECT_EQ(my_mdspan_map_acc(i0, i1, i2, i3, i4), value);
269                             EXPECT_EQ(my_mdspan_extents(i0, i1, i2, i3, i4), value);
270                             EXPECT_EQ(my_mdspan_copy(i0, i1, i2, i3, i4), value);
271                         }
272                     }
273                 }
274             }
275         }
276     }
277 };
278
279 TYPED_TEST_CASE(MdSpanTest, MdSpanTypes);
280
281 TYPED_TEST(MdSpanTest, Rank)
282 {
283     this->SetUp(4, 2);
284     this->check_rank(5);
285 }
286
287 TYPED_TEST(MdSpanTest, DynamicRank)
288 {
289     this->SetUp(4, 2);
290     this->check_rank_dynamic(2);
291 }
292
293 TYPED_TEST(MdSpanTest, Extents)
294 {
295     this->SetUp(4, 2);
296     this->check_extents(5, 4, 3, 2, 1);
297 }
298
299 TYPED_TEST(MdSpanTest, Strides)
300 {
301     this->SetUp(4, 2);
302     this->check_strides(24, 6, 2, 1, 1);
303 }
304
305 TYPED_TEST(MdSpanTest, Properties)
306 {
307     this->SetUp(4, 2);
308     const bool always_unique     = true;
309     const bool always_contiguous = true;
310     const bool always_strided    = true;
311     const bool unique            = true;
312     const bool contiguous        = true;
313     const bool strided           = true;
314     this->check_properties(always_unique, always_contiguous, always_strided, unique, contiguous, strided);
315 }
316
317 TYPED_TEST(MdSpanTest, Operator)
318 {
319     this->SetUp(4, 2);
320     this->check_operator();
321 }
322
323 } // namespace
324
325 } // namespace gmx