Apply clang-format to source tree
[alexxy/gromacs.git] / src / gromacs / math / arrayrefwithpadding.h
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 /*! \file
36  * \brief Declares gmx::ArrayRefWithPadding that refers to memory whose
37  * size includes padding for SIMD operations.
38  *
39  * \author Mark Abraham <mark.j.abraham@gmail.com>
40  *
41  * \inpublicapi
42  * \ingroup module_math
43  */
44 #ifndef GMX_MATH_ARRAYREFWITHPADDING_H
45 #define GMX_MATH_ARRAYREFWITHPADDING_H
46
47 #include <cstddef>
48
49 #include "gromacs/utility/arrayref.h"
50 #include "gromacs/utility/basedefinitions.h"
51 #include "gromacs/utility/gmxassert.h"
52
53 namespace gmx
54 {
55
56 /*! \brief Interface to a C array of T (or part of a std container of
57  * T), that includes padding that is suitable for the kinds of SIMD
58  * operations GROMACS uses.
59  *
60  * High level code e.g. in the force calculation routines need to hold
61  * a non-owning view of memory and be able to create ArrayRef objects
62  * that view padded or unpadded memory, suitable for the various
63  * component routines. This class implements that non-owning view,
64  * and the only available functionality refers to its size, and the
65  * methods to create such ArrayRef objects.
66  *
67  * \copydoc ArrayRef
68  * \inlibraryapi
69  * \ingroup module_math
70  */
71 template<typename T>
72 class ArrayRefWithPadding
73 {
74 public:
75     //! Type of values stored in the reference.
76     using value_type = T;
77     //! Type for representing size of the reference.
78     using size_type = index;
79     //! Const pointer to an element.
80     using const_pointer = const T*;
81     //! Const iterator type to an element.
82     using const_iterator = const T*;
83     //! Pointer to an element.
84     using pointer = T*;
85     //! Iterator type to an element.
86     using iterator = T*;
87
88     /*! \brief
89      * Constructs an empty reference.
90      */
91     ArrayRefWithPadding() : begin_(nullptr), end_(nullptr), paddedEnd_(nullptr) {}
92     /*! \brief
93      * Constructs a reference to a particular range.
94      *
95      * \param[in] begin        Pointer to the beginning of a range.
96      * \param[in] end          Pointer to the end of a range without padding.
97      * \param[in] paddedEnd    Pointer to the end of a range including padding.
98      *
99      * Passed pointers must remain valid for the lifetime of this object.
100      */
101     ArrayRefWithPadding(pointer begin, pointer end, pointer paddedEnd) :
102         begin_(begin),
103         end_(end),
104         paddedEnd_(paddedEnd)
105     {
106         GMX_ASSERT(end >= begin, "Invalid range");
107         GMX_ASSERT(paddedEnd >= end, "Invalid range");
108     }
109     //! Copy constructor
110     ArrayRefWithPadding(const ArrayRefWithPadding& o) :
111         begin_(o.begin_),
112         end_(o.end_),
113         paddedEnd_(o.paddedEnd_)
114     {
115     }
116     //! Move constructor
117     ArrayRefWithPadding(ArrayRefWithPadding&& o) noexcept :
118         begin_(std::move(o.begin_)),
119         end_(std::move(o.end_)),
120         paddedEnd_(std::move(o.paddedEnd_))
121     {
122     }
123     //! Convenience overload constructor to make an ArrayRefWithPadding<const T> from a non-const one.
124     template<typename U, typename = std::enable_if_t<std::is_same<value_type, const typename std::remove_reference_t<U>::value_type>::value>>
125     ArrayRefWithPadding(U&& o)
126     {
127         auto constArrayRefWithPadding = o.constArrayRefWithPadding();
128         this->swap(constArrayRefWithPadding);
129     }
130     //! Copy assignment operator
131     ArrayRefWithPadding& operator=(ArrayRefWithPadding const& o)
132     {
133         if (&o != this)
134         {
135             begin_     = o.begin_;
136             end_       = o.end_;
137             paddedEnd_ = o.paddedEnd_;
138         }
139         return *this;
140     }
141     //! Move assignment operator
142     ArrayRefWithPadding& operator=(ArrayRefWithPadding&& o) noexcept
143     {
144         if (&o != this)
145         {
146             begin_     = std::move(o.begin_);
147             end_       = std::move(o.end_);
148             paddedEnd_ = std::move(o.paddedEnd_);
149         }
150         return *this;
151     }
152
153     //! Return the size of the view, i.e with the padding.
154     size_type size() const { return paddedEnd_ - begin_; }
155     //! Whether the reference refers to no memory.
156     bool empty() const { return begin_ == end_; }
157
158     //! Returns an ArrayRef of elements that does not include the padding region.
159     ArrayRef<T> unpaddedArrayRef() { return { begin_, end_ }; }
160     //! Returns an ArrayRef of const elements that does not include the padding region.
161     ArrayRef<const T> unpaddedConstArrayRef() const { return { begin_, end_ }; }
162     //! Returns an ArrayRef of elements that does include the padding region.
163     ArrayRef<T> paddedArrayRef() { return { begin_, paddedEnd_ }; }
164     //! Returns an ArrayRef of const elements that does include the padding region.
165     ArrayRef<const T> paddedConstArrayRef() const { return { begin_, paddedEnd_ }; }
166     //! Returns an identical ArrayRefWithPadding that refers to const elements.
167     ArrayRefWithPadding<const T> constArrayRefWithPadding() const
168     {
169         return { begin_, end_, paddedEnd_ };
170     }
171     /*! \brief
172      * Swaps referenced memory with the other object.
173      *
174      * The actual memory areas are not modified, only the references are
175      * swapped.
176      */
177     void swap(ArrayRefWithPadding<T>& other)
178     {
179         std::swap(begin_, other.begin_);
180         std::swap(end_, other.end_);
181         std::swap(paddedEnd_, other.paddedEnd_);
182     }
183
184 private:
185     pointer begin_;
186     pointer end_;
187     pointer paddedEnd_;
188 };
189
190 } // namespace gmx
191
192 #endif