Add replacements for pbc enumerations
[gromacs.git] / src / gromacs / math / paddedvector.h
blobfa66983af828a899ea0eb6530d05bd42a44bb621
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,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.
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 /*! \file
36 * \brief
37 * Declares gmx::PaddedRVecVector
39 * \author Mark Abraham <mark.j.abraham@gmail.com>
40 * \inpublicapi
41 * \ingroup module_math
43 #ifndef GMX_MATH_PADDEDVECTOR_H
44 #define GMX_MATH_PADDEDVECTOR_H
46 #include <algorithm>
47 #include <vector>
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"
55 namespace gmx
58 namespace detail
61 /*! \brief Traits classes for handling padding for types used with PaddedVector
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.
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.
74 template<typename T>
75 struct PaddingTraits {};
77 template<>
78 struct PaddingTraits<int32_t>
80 using SimdBaseType = int32_t;
81 static constexpr int maxSimdWidthOfBaseType = 16;
84 template<>
85 struct PaddingTraits<float>
87 using SimdBaseType = float;
88 static constexpr int maxSimdWidthOfBaseType = GMX_FLOAT_MAX_SIMD_WIDTH;
91 template<>
92 struct PaddingTraits<double>
94 using SimdBaseType = double;
95 static constexpr int maxSimdWidthOfBaseType = GMX_DOUBLE_MAX_SIMD_WIDTH;
98 template<>
99 struct PaddingTraits < BasicVector < float>>
101 using SimdBaseType = float;
102 static constexpr int maxSimdWidthOfBaseType = GMX_FLOAT_MAX_SIMD_WIDTH;
105 template<>
106 struct PaddingTraits < BasicVector < double>>
108 using SimdBaseType = double;
109 static constexpr int maxSimdWidthOfBaseType = GMX_DOUBLE_MAX_SIMD_WIDTH;
112 /*! \brief Returns the allocation size for PaddedVector that contains
113 * \c numElements elements plus padding for SIMD operations.
115 * \param[in] numElements The number of T elements for which data will be stored.
116 * \returns The number of T elements that must be allocated
117 * (ie >= numElements).
119 template <typename T>
120 index computePaddedSize(index numElements)
122 // We don't need padding if there is no access.
123 if (numElements == 0)
125 return 0;
128 // We sometimes load a whole extra element when doing 4-wide SIMD
129 // operations (which might e.g. be an RVec) so we need to pad for
130 // that.
131 index simdScatterAccessSize = numElements + 1;
133 // For SIMD updates based on RVec, we might load starting from
134 // the last RVec element, so that sets the minimum extent of the
135 // padding. That extent must take the initialized allocation up to
136 // the SIMD width of the base type multiplied by the width of T in
137 // that base type. But since storage_ contains RVec, we only have
138 // to tell it the number of elements, which means to round up to
139 // the next SIMD width.
141 // We don't want a dependence on the SIMD module for the actual
142 // SIMD width of the base type, so we use maximum for the base
143 // type via the traits. A little extra padding won't really hurt.
144 constexpr int maxSimdWidth = PaddingTraits<T>::maxSimdWidthOfBaseType;
145 index simdFlatAccessSize = (numElements + (maxSimdWidth-1)) / maxSimdWidth * maxSimdWidth;
147 return std::max(simdScatterAccessSize, simdFlatAccessSize);
150 //! Helper function to insert padding elements for most T.
151 template <typename T, typename AllocatorType>
152 inline void insertPaddingElements(std::vector<T, AllocatorType> *v,
153 index newPaddedSize)
155 // Ensure the padding region is initialized to zero. There is no
156 // way to insert a number of default-initialized elements. So we
157 // have to provide a value for those elements, which anyway suits
158 // this use case.
159 v->insert(v->end(), newPaddedSize - v->size(), 0);
162 //! Specialization of helper function to insert padding elements, used for BasicVector<T>.
163 template <typename T, typename AllocatorType>
164 inline void insertPaddingElements(std::vector<BasicVector<T>, AllocatorType> *v,
165 index newPaddedSize)
167 // Ensure the padding region is initialized to zero.
168 v->insert(v->end(), newPaddedSize - v->size(), BasicVector<T>(0, 0, 0));
171 } // namespace detail
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.
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.
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().
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.
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.
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.
211 template <typename T, typename Allocator = Allocator < T, AlignedAllocationPolicy > >
212 class PaddedVector
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 //! \}
230 PaddedVector() :
231 storage_(),
232 unpaddedEnd_(begin())
234 /*! \brief Constructor that specifes the initial size. */
235 explicit PaddedVector(size_type count,
236 const allocator_type &allocator = Allocator()) :
237 storage_(count, allocator),
238 unpaddedEnd_(begin() + count)
240 // The count elements have been default inserted, and now
241 // the padding elements are added
242 resizeWithPadding(count);
244 /*! \brief Constructor that specifes the initial size and an element to copy. */
245 explicit PaddedVector(size_type count,
246 value_type const &v,
247 const allocator_type &allocator = Allocator()) :
248 storage_(count, v, allocator),
249 unpaddedEnd_(begin() + count)
251 // The count elements have been default inserted, and now
252 // the padding elements are added
253 resizeWithPadding(count);
255 //! Default constructor with allocator
256 explicit PaddedVector(allocator_type const &allocator) :
257 storage_(allocator),
258 unpaddedEnd_(begin())
260 //! Copy constructor
261 PaddedVector(PaddedVector const &o) :
262 storage_(o.storage_),
263 unpaddedEnd_(begin() + o.size())
265 //! Move constructor
266 PaddedVector(PaddedVector &&o) noexcept :
267 storage_(std::move(o.storage_)),
268 unpaddedEnd_(std::move(o.unpaddedEnd_))
270 unpaddedEnd_ = begin();
272 //! Move constructor using \c alloc for the new vector.
273 PaddedVector(PaddedVector &&o, const Allocator &alloc) noexcept :
274 storage_(std::move(alloc)),
275 unpaddedEnd_(begin())
277 auto unpaddedSize = o.size();
278 if (alloc == o.storage_.get_allocator())
280 storage_ = std::move(o.storage_);
282 else
284 // If the allocator compares differently, we must
285 // reallocate and copy.
286 resizeWithPadding(unpaddedSize);
287 std::copy(o.begin(), o.end(), storage_.begin());
289 unpaddedEnd_ = begin() + unpaddedSize;
291 //! Construct from an initializer list
292 PaddedVector(std::initializer_list<value_type> const &il) :
293 storage_(il),
294 unpaddedEnd_(storage_.end())
296 // We can't choose the padding until we know the size of
297 // the normal vector, so we have to make the storage_ and
298 // then resize it.
299 resizeWithPadding(storage_.size());
301 //! Reserve storage for the container to contain newExtent elements, plus the required padding.
302 void reserveWithPadding(const size_type newExtent)
304 auto unpaddedSize = end() - begin();
305 /* v.reserve(13) should allocate enough memory so that
306 v.resize(13) does not reallocate. This means that the
307 new extent should be large enough for the padded
308 storage for a vector whose size is newExtent. */
309 auto newPaddedExtent = detail::computePaddedSize<T>(newExtent);
310 storage_.reserve(newPaddedExtent);
311 unpaddedEnd_ = begin() + unpaddedSize;
313 //! Resize the container to contain newSize elements, plus the required padding.
314 void resizeWithPadding(const size_type newSize)
316 // When the contained type is e.g. a scalar, then the
317 // default initialization behaviour is to zero all
318 // elements, which is OK, but we have to make sure that it
319 // happens for the elements in the padded region when the
320 // vector is shrinking.
321 auto newPaddedSize = detail::computePaddedSize<T>(newSize);
322 // Make sure there is room for padding if we need to grow.
323 storage_.reserve(newPaddedSize);
324 // Make the unpadded size correct, with any additional
325 // elements initialized by the default constructor. It is
326 // particularly important to destruct former elements when
327 // newSize is smaller than the old size.
328 storage_.resize(newSize);
329 // Ensure the padding region is zeroed if required.
330 detail::insertPaddingElements(&storage_, newPaddedSize);
331 unpaddedEnd_ = begin() + newSize;
333 //! Return the size of the view without the padding.
334 size_type size() const { return end() - begin(); }
335 //! Return the container size including the padding.
336 size_type paddedSize() const { return storage_.size(); }
337 //! Return whether the storage is empty.
338 bool empty() const { return storage_.empty(); }
339 //! Swap two PaddedVectors
340 void swap(PaddedVector &x)
342 std::swap(storage_, x.storage_);
343 std::swap(unpaddedEnd_, x.unpaddedEnd_);
345 //! Clear the vector, ie. set size to zero and remove padding.
346 void clear()
348 storage_.clear();
349 unpaddedEnd_ = begin();
351 //! Iterator getters refer to a view without padding.
352 //! \{
353 pointer data() noexcept { return storage_.data(); }
354 const_pointer data() const noexcept { return storage_.data(); }
356 iterator begin() { return storage_.begin(); }
357 iterator end() { return iterator(unpaddedEnd_); }
359 const_iterator cbegin() { return const_iterator(begin()); }
360 const_iterator cend() { return const_iterator(unpaddedEnd_); }
362 const_iterator begin() const { return storage_.begin(); }
363 const_iterator end() const { return const_iterator(unpaddedEnd_); }
365 const_iterator cbegin() const { return const_iterator(begin()); }
366 const_iterator cend() const { return const_iterator(unpaddedEnd_); }
367 //! \}
368 // TODO should these do bounds checking for the unpadded range? In debug mode?
369 //! Indexing operator.
370 reference operator[](int i) { return storage_[i]; }
371 //! Indexing operator as const.
372 const_reference operator[](int i) const { return storage_[i]; }
373 //! Returns an ArrayRef of elements that includes the padding region, e.g. for use in SIMD code.
374 ArrayRefWithPadding<T> arrayRefWithPadding()
376 return ArrayRefWithPadding<T>(data(), data()+size(), data()+paddedSize());
378 //! Returns an ArrayRef of const elements that includes the padding region, e.g. for use in SIMD code.
379 ArrayRefWithPadding<const T> constArrayRefWithPadding() const
381 return ArrayRefWithPadding<const T>(data(), data()+size(), data()+paddedSize());
383 //! Returns an rvec * pointer for containers of RVec, for use with legacy code.
384 template <typename AlsoT = T,
385 typename = typename std::enable_if<std::is_same<AlsoT, RVec>::value> >
386 rvec *rvec_array()
388 return as_rvec_array(data());
390 //! Returns a const rvec * pointer for containers of RVec, for use with legacy code.
391 template <typename AlsoT = T,
392 typename = typename std::enable_if<std::is_same<AlsoT, RVec>::value> >
393 const rvec *rvec_array() const
395 return as_rvec_array(data());
397 //! Copy assignment operator
398 PaddedVector &operator=(PaddedVector const &o)
400 if (&o != this)
402 storage_ = o.storage_;
403 unpaddedEnd_ = begin() + o.size();
405 return *this;
407 //! Move assignment operator
408 PaddedVector &operator=(PaddedVector &&o) noexcept
410 if (&o != this)
412 auto oSize = o.size();
413 storage_ = std::move(o.storage_);
414 unpaddedEnd_ = begin() + oSize;
415 o.unpaddedEnd_ = o.begin();
417 return *this;
419 //! Getter for the allocator
420 allocator_type
421 get_allocator() const
423 return storage_.get_allocator();
426 private:
427 storage_type storage_;
428 iterator unpaddedEnd_;
431 } // namespace gmx
433 // TODO These are hacks to avoid littering gmx:: all over code that is
434 // almost all destined to move into the gmx namespace at some point.
435 // An alternative would be about 20 files with using statements.
436 using gmx::PaddedVector; //NOLINT(google-global-names-in-headers)
438 #endif