Fix random doxygen warnings and typos
[alexxy/gromacs.git] / src / gromacs / math / paddedvector.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,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 /*! \file
36  * \brief
37  * Declares gmx::PaddedRVecVector
38  *
39  * \author Mark Abraham <mark.j.abraham@gmail.com>
40  * \inpublicapi
41  * \ingroup module_math
42  */
43 #ifndef GMX_MATH_PADDEDVECTOR_H
44 #define GMX_MATH_PADDEDVECTOR_H
45
46 #include <algorithm>
47 #include <vector>
48
49 #include "gromacs/math/arrayrefwithpadding.h"
50 #include "gromacs/math/vectypes.h"
51 #include "gromacs/utility/alignedallocator.h"
52 #include "gromacs/utility/allocator.h"
53 #include "gromacs/utility/arrayref.h"
54
55 namespace gmx
56 {
57
58 namespace detail
59 {
60
61 /*! \brief Traits classes for handling padding for types used with PaddedVector
62  *
63  * Only the base types of the SIMD module are supported for
64  * PaddedVector, because the purpose of the padding is to permit
65  * SIMD-width operations from the SIMD module.
66  *
67  * \todo Consider explicitly tying these types to the SimdTrait
68  * types. That would require depending on the SIMD module, or
69  * extracting the traits from it. This would also permit
70  * maxSimdWidthOfBaseType to be set more efficiently, e.g. as a
71  * metaprogramming max over the maximum width from different
72  * implementations.
73  */
74 template<typename T>
75 struct PaddingTraits
76 {
77 };
78
79 template<>
80 struct PaddingTraits<int32_t>
81 {
82     using SimdBaseType                          = int32_t;
83     static constexpr int maxSimdWidthOfBaseType = 16;
84 };
85
86 template<>
87 struct PaddingTraits<float>
88 {
89     using SimdBaseType                          = float;
90     static constexpr int maxSimdWidthOfBaseType = GMX_FLOAT_MAX_SIMD_WIDTH;
91 };
92
93 template<>
94 struct PaddingTraits<double>
95 {
96     using SimdBaseType                          = double;
97     static constexpr int maxSimdWidthOfBaseType = GMX_DOUBLE_MAX_SIMD_WIDTH;
98 };
99
100 template<>
101 struct PaddingTraits<BasicVector<float>>
102 {
103     using SimdBaseType                          = float;
104     static constexpr int maxSimdWidthOfBaseType = GMX_FLOAT_MAX_SIMD_WIDTH;
105 };
106
107 template<>
108 struct PaddingTraits<BasicVector<double>>
109 {
110     using SimdBaseType                          = double;
111     static constexpr int maxSimdWidthOfBaseType = GMX_DOUBLE_MAX_SIMD_WIDTH;
112 };
113
114 /*! \brief Returns the allocation size for PaddedVector that contains
115  * \c numElements elements plus padding for SIMD operations.
116  *
117  * \param[in] numElements  The number of T elements for which data will be stored.
118  * \returns                The number of T elements that must be allocated
119  *                         (ie >= numElements).
120  */
121 template<typename T>
122 index computePaddedSize(index numElements)
123 {
124     // We don't need padding if there is no access.
125     if (numElements == 0)
126     {
127         return 0;
128     }
129
130     // We sometimes load a whole extra element when doing 4-wide SIMD
131     // operations (which might e.g. be an RVec) so we need to pad for
132     // that.
133     index simdScatterAccessSize = numElements + 1;
134
135     // For SIMD updates based on RVec, we might load starting from
136     // the last RVec element, so that sets the minimum extent of the
137     // padding. That extent must take the initialized allocation up to
138     // the SIMD width of the base type multiplied by the width of T in
139     // that base type. But since storage_ contains RVec, we only have
140     // to tell it the number of elements, which means to round up to
141     // the next SIMD width.
142     //
143     // We don't want a dependence on the SIMD module for the actual
144     // SIMD width of the base type, so we use maximum for the base
145     // type via the traits. A little extra padding won't really hurt.
146     constexpr int maxSimdWidth = PaddingTraits<T>::maxSimdWidthOfBaseType;
147     index simdFlatAccessSize   = (numElements + (maxSimdWidth - 1)) / maxSimdWidth * maxSimdWidth;
148
149     return std::max(simdScatterAccessSize, simdFlatAccessSize);
150 }
151
152 //! Helper function to insert padding elements for most T.
153 template<typename T, typename AllocatorType>
154 inline void insertPaddingElements(std::vector<T, AllocatorType>* v, index newPaddedSize)
155 {
156     // Ensure the padding region is initialized to zero. There is no
157     // way to insert a number of default-initialized elements. So we
158     // have to provide a value for those elements, which anyway suits
159     // this use case.
160     v->insert(v->end(), newPaddedSize - v->size(), 0);
161 }
162
163 //! Specialization of helper function to insert padding elements, used for BasicVector<T>.
164 template<typename T, typename AllocatorType>
165 inline void insertPaddingElements(std::vector<BasicVector<T>, AllocatorType>* v, index newPaddedSize)
166 {
167     // Ensure the padding region is initialized to zero.
168     v->insert(v->end(), newPaddedSize - v->size(), BasicVector<T>(0, 0, 0));
169 }
170
171 } // namespace detail
172
173 /*! \brief PaddedVector is a container of elements in contiguous
174  * storage that allocates extra memory for safe SIMD-style loads for
175  * operations used in GROMACS.
176  *
177  * \tparam T the type of objects within the container
178  * \tparam Allocator the allocator used. Can be any standard-compliant
179  * allocator, such gmx::Allocator used for alignment and/or pinning.
180  *
181  * The interface resembles std::vector. However, access
182  * intended to include padded elements must be via ArrayRef objects
183  * explicitly created to view those elements. Most other aspects of
184  * this vector refer to the unpadded view, e.g. iterators, data(),
185  * size().
186  *
187  * The underlying storage is allocated with extra elements, properly
188  * initialized, that ensure that any operations accessing the any
189  * non-additional element that operate on memory equivalent to a full
190  * SIMD lane do so on allocated memory that has been initialized, so
191  * that memory traps will not occur, and arithmetic operations will
192  * not cause e.g. floating-point exceptions so long as the values in
193  * the padded elements are properly managed.
194  *
195  * Proper initialization is tricker than it would first appear, since
196  * we intend this container to be used with scalar and class types
197  * (e.g. RVec). Resize and construction operations use "default
198  * insertion" which leads to zero initialization for the former, and
199  * calling the default constructor for the latter. BasicVector has a
200  * default constructor that leaves the elements uninitialized, which
201  * is particularly risky for elements only present as padding. Thus
202  * the implementation specifically initializes the padded elements to
203  * zero, which makes no difference to the scalar template
204  * instantiations, and makes the BasicVector ones safer to use.
205  *
206  * Because the allocator can be configured, the memory allocation can
207  * have other attributes such as SIMD alignment or being pinned to
208  * physical memory for efficient transfers. The default allocator
209  * ensures alignment, but std::allocator also works.
210  */
211 template<typename T, typename Allocator = Allocator<T, AlignedAllocationPolicy>>
212 class PaddedVector
213 {
214 public:
215     //! Standard helper types
216     //! \{
217     using value_type      = T;
218     using allocator_type  = Allocator;
219     using size_type       = index;
220     using reference       = value_type&;
221     using const_reference = const value_type&;
222     using storage_type    = std::vector<T, allocator_type>;
223     using pointer         = typename storage_type::pointer;
224     using const_pointer   = typename storage_type::const_pointer;
225     using iterator        = typename storage_type::iterator;
226     using const_iterator  = typename storage_type::const_iterator;
227     using difference_type = typename storage_type::iterator::difference_type;
228     //! \}
229
230     PaddedVector() : storage_(), unpaddedEnd_(begin()) {}
231     /*! \brief Constructor that specifies the initial size. */
232     explicit PaddedVector(size_type count, const allocator_type& allocator = Allocator()) :
233         storage_(count, allocator),
234         unpaddedEnd_(begin() + count)
235     {
236         // The count elements have been default inserted, and now
237         // the padding elements are added
238         resizeWithPadding(count);
239     }
240     /*! \brief Constructor that specifies the initial size and an element to copy. */
241     explicit PaddedVector(size_type count, value_type const& v, const allocator_type& allocator = Allocator()) :
242         storage_(count, v, allocator),
243         unpaddedEnd_(begin() + count)
244     {
245         // The count elements have been default inserted, and now
246         // the padding elements are added
247         resizeWithPadding(count);
248     }
249     //! Default constructor with allocator
250     explicit PaddedVector(allocator_type const& allocator) :
251         storage_(allocator),
252         unpaddedEnd_(begin())
253     {
254     }
255     //! Copy constructor
256     PaddedVector(PaddedVector const& o) : storage_(o.storage_), unpaddedEnd_(begin() + o.size()) {}
257     //! Move constructor
258     PaddedVector(PaddedVector&& o) noexcept :
259         storage_(std::move(o.storage_)),
260         unpaddedEnd_(std::move(o.unpaddedEnd_))
261     {
262         unpaddedEnd_ = begin();
263     }
264     //! Move constructor using \c alloc for the new vector.
265     PaddedVector(PaddedVector&& o, const Allocator& alloc) noexcept :
266         storage_(std::move(alloc)),
267         unpaddedEnd_(begin())
268     {
269         auto unpaddedSize = o.size();
270         if (alloc == o.storage_.get_allocator())
271         {
272             storage_ = std::move(o.storage_);
273         }
274         else
275         {
276             // If the allocator compares differently, we must
277             // reallocate and copy.
278             resizeWithPadding(unpaddedSize);
279             std::copy(o.begin(), o.end(), storage_.begin());
280         }
281         unpaddedEnd_ = begin() + unpaddedSize;
282     }
283     //! Construct from an initializer list
284     PaddedVector(std::initializer_list<value_type> const& il) :
285         storage_(il),
286         unpaddedEnd_(storage_.end())
287     {
288         // We can't choose the padding until we know the size of
289         // the normal vector, so we have to make the storage_ and
290         // then resize it.
291         resizeWithPadding(storage_.size());
292     }
293     //! Reserve storage for the container to contain newExtent elements, plus the required padding.
294     void reserveWithPadding(const size_type newExtent)
295     {
296         auto unpaddedSize = end() - begin();
297         /* v.reserve(13) should allocate enough memory so that
298            v.resize(13) does not reallocate. This means that the
299            new extent should be large enough for the padded
300            storage for a vector whose size is newExtent. */
301         auto newPaddedExtent = detail::computePaddedSize<T>(newExtent);
302         storage_.reserve(newPaddedExtent);
303         unpaddedEnd_ = begin() + unpaddedSize;
304     }
305     //! Resize the container to contain newSize elements, plus the required padding.
306     void resizeWithPadding(const size_type newSize)
307     {
308         // When the contained type is e.g. a scalar, then the
309         // default initialization behaviour is to zero all
310         // elements, which is OK, but we have to make sure that it
311         // happens for the elements in the padded region when the
312         // vector is shrinking.
313         auto newPaddedSize = detail::computePaddedSize<T>(newSize);
314         // Make sure there is room for padding if we need to grow.
315         storage_.reserve(newPaddedSize);
316         // Make the unpadded size correct, with any additional
317         // elements initialized by the default constructor. It is
318         // particularly important to destruct former elements when
319         // newSize is smaller than the old size.
320         storage_.resize(newSize);
321         // Ensure the padding region is zeroed if required.
322         detail::insertPaddingElements(&storage_, newPaddedSize);
323         unpaddedEnd_ = begin() + newSize;
324     }
325     //! Return the size of the view without the padding.
326     size_type size() const { return end() - begin(); }
327     //! Return the container size including the padding.
328     size_type paddedSize() const { return storage_.size(); }
329     //! Return whether the storage is empty.
330     bool empty() const { return storage_.empty(); }
331     //! Swap two PaddedVectors
332     void swap(PaddedVector& x)
333     {
334         std::swap(storage_, x.storage_);
335         std::swap(unpaddedEnd_, x.unpaddedEnd_);
336     }
337     //! Clear the vector, ie. set size to zero and remove padding.
338     void clear()
339     {
340         storage_.clear();
341         unpaddedEnd_ = begin();
342     }
343     //! Iterator getters refer to a view without padding.
344     //! \{
345     pointer       data() noexcept { return storage_.data(); }
346     const_pointer data() const noexcept { return storage_.data(); }
347
348     iterator begin() { return storage_.begin(); }
349     iterator end() { return iterator(unpaddedEnd_); }
350
351     const_iterator cbegin() { return const_iterator(begin()); }
352     const_iterator cend() { return const_iterator(unpaddedEnd_); }
353
354     const_iterator begin() const { return storage_.begin(); }
355     const_iterator end() const { return const_iterator(unpaddedEnd_); }
356
357     const_iterator cbegin() const { return const_iterator(begin()); }
358     const_iterator cend() const { return const_iterator(unpaddedEnd_); }
359     //! \}
360     // TODO should these do bounds checking for the unpadded range? In debug mode?
361     //! Indexing operator.
362     reference operator[](int i) { return storage_[i]; }
363     //! Indexing operator as const.
364     const_reference operator[](int i) const { return storage_[i]; }
365     //! Returns an ArrayRef of elements that includes the padding region, e.g. for use in SIMD code.
366     ArrayRefWithPadding<T> arrayRefWithPadding()
367     {
368         return ArrayRefWithPadding<T>(data(), data() + size(), data() + paddedSize());
369     }
370     //! Returns an ArrayRef of const elements that includes the padding region, e.g. for use in SIMD code.
371     ArrayRefWithPadding<const T> constArrayRefWithPadding() const
372     {
373         return ArrayRefWithPadding<const T>(data(), data() + size(), data() + paddedSize());
374     }
375     /*! \brief Returns an rvec * pointer for containers of RVec, for use with legacy code.
376      *
377      * \todo Use std::is_same_v when CUDA 11 is a requirement.
378      */
379     template<typename AlsoT = T, typename = typename std::enable_if<std::is_same<AlsoT, RVec>::value>>
380     rvec* rvec_array()
381     {
382         return as_rvec_array(data());
383     }
384     /*! \brief Returns a const rvec * pointer for containers of RVec, for use with legacy code.
385      *
386      * \todo Use std::is_same_v when CUDA 11 is a requirement.
387      */
388     template<typename AlsoT = T, typename = typename std::enable_if<std::is_same<AlsoT, RVec>::value>>
389     const rvec* rvec_array() const
390     {
391         return as_rvec_array(data());
392     }
393     //! Copy assignment operator
394     PaddedVector& operator=(PaddedVector const& o)
395     {
396         if (&o != this)
397         {
398             storage_     = o.storage_;
399             unpaddedEnd_ = begin() + o.size();
400         }
401         return *this;
402     }
403     //! Move assignment operator
404     PaddedVector& operator=(PaddedVector&& o) noexcept
405     {
406         if (&o != this)
407         {
408             auto oSize     = o.size();
409             storage_       = std::move(o.storage_);
410             unpaddedEnd_   = begin() + oSize;
411             o.unpaddedEnd_ = o.begin();
412         }
413         return *this;
414     }
415     //! Getter for the allocator
416     allocator_type get_allocator() const { return storage_.get_allocator(); }
417
418 private:
419     storage_type storage_;
420     iterator     unpaddedEnd_;
421 };
422
423 } // namespace gmx
424
425 // TODO These are hacks to avoid littering gmx:: all over code that is
426 // almost all destined to move into the gmx namespace at some point.
427 // An alternative would be about 20 files with using statements.
428 using gmx::PaddedVector; //NOLINT(google-global-names-in-headers)
429
430 #endif