128-bit AVX2 SIMD for AMD Ryzen
[gromacs.git] / src / gromacs / analysisdata / framelocaldata.h
blob99fb08c6ac52066c7d658f7301aa99ca3d2bab99
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2014,2015, 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 /*! \internal \file
36 * \brief
37 * Defines gmx::AnalysisDataFrameLocalData and supporting types.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_analysisdata
42 #ifndef GMX_ANALYSISDATA_FRAMELOCALDATA_H
43 #define GMX_ANALYSISDATA_FRAMELOCALDATA_H
45 #include <algorithm>
46 #include <numeric>
47 #include <vector>
49 #include "gromacs/analysisdata/paralleloptions.h"
50 #include "gromacs/utility/arrayref.h"
51 #include "gromacs/utility/gmxassert.h"
53 namespace gmx
56 //! \addtogroup module_analysisdata
57 //! \{
59 /*! \internal
60 * \brief
61 * Handle to a single data set within frame-local data array.
63 * Methods in this class do not throw.
65 * \see AnalysisDataFrameLocalData
67 template<typename ValueType>
68 class AnalysisDataFrameLocalDataSetHandle
70 public:
71 //! Constructs a handle from an array of values.
72 explicit AnalysisDataFrameLocalDataSetHandle(ArrayRef<ValueType> values)
73 : values_(values)
77 //! Clears all values in the data set.
78 void clear()
80 std::fill(values_.begin(), values_.end(), ValueType());
83 //! Accesses a single value in the data set.
84 ValueType &value(int column)
86 GMX_ASSERT(column >= 0 && column < static_cast<int>(values_.size()),
87 "Invalid column index");
88 return values_[column];
91 private:
92 ArrayRef<ValueType> values_;
95 /*! \internal
96 * \brief
97 * Handle to a single frame data within frame-local data array.
99 * Methods in this class do not throw.
101 * \see AnalysisDataFrameLocalData
103 template<typename ValueType>
104 class AnalysisDataFrameLocalDataHandle
106 public:
107 //! Shorthand for the internal array of values.
108 typedef std::vector<ValueType> ValueArray;
109 //! Shorthand for a handle to a single data set.
110 typedef AnalysisDataFrameLocalDataSetHandle<ValueType> DataSetHandle;
112 //! Constructs a handle from specified frame data.
113 AnalysisDataFrameLocalDataHandle(const std::vector<int> *dataSetIndices,
114 ValueArray *values)
115 : dataSetIndices_(dataSetIndices), values_(values)
119 //! Returns the number of data sets in the array.
120 int dataSetCount() const
122 return dataSetIndices_->size() - 1;
124 //! Clears all values in the frame.
125 void clear()
127 std::fill(values_->begin(), values_->end(), ValueType());
130 //! Returns a handle for a single data set.
131 DataSetHandle dataSet(int dataSet)
133 GMX_ASSERT(dataSet >= 0 && dataSet < dataSetCount(),
134 "Invalid data set index");
135 const int firstIndex = (*dataSetIndices_)[dataSet];
136 const int lastIndex = (*dataSetIndices_)[dataSet + 1];
137 typename ValueArray::iterator begin = values_->begin() + firstIndex;
138 typename ValueArray::iterator end = values_->begin() + lastIndex;
139 return DataSetHandle(arrayRefFromVector<ValueType>(begin, end));
141 //! Accesses a single value in the frame.
142 ValueType &value(int dataSet, int column)
144 GMX_ASSERT(dataSet >= 0 && dataSet < dataSetCount(),
145 "Invalid data set index");
146 const int firstIndex = (*dataSetIndices_)[dataSet];
147 GMX_ASSERT(column >= 0
148 && column < (*dataSetIndices_)[dataSet+1] - firstIndex,
149 "Invalid column index");
150 return (*values_)[firstIndex + column];
153 private:
154 const std::vector<int> *dataSetIndices_;
155 ValueArray *values_;
158 /*! \internal \brief
159 * Container for an array of frame-local values that supports parallel data
160 * processing.
162 * \tparam ValueType Type of values to store.
164 * This class provides a convenient interface to create an array of frame-local
165 * data for use in analysis data modules that support parallel processing.
166 * The object is initialized by setting the desired dimensionality with
167 * setDataSetCount() and setColumnCount(), followed by a call to init(),
168 * typically in IAnalysisDataModule::parallelDataStarted(),
170 * After initialization, frameData() can be used to access the data for a given
171 * frame, independently from other frames. This works if the assumptions about
172 * parallelism hold: if `N` is the parallelization factor given for init() with
173 * AnalysisDataParallelOptions::parallelizationFactor(), then frame `i+N` must
174 * not be accessed before all processing for frame `i` is finished.
175 * Technically, the data for different frames is kept in a ring buffer of size
176 * `N`.
178 * The data for a frame is not cleared after it is reused for a new frame (but
179 * is initially cleared). This allows using the data for accumulating values
180 * over all frames in a lock-free manner.
182 * frameDataSet() is provided for convenience when only a single data set
183 * needs to be accessed (typically in IAnalysisDataModule::pointsAdded()).
185 * Methods in this class do not throw except where indicated.
187 * \see AnalysisDataFrameLocalData
189 template<typename ValueType>
190 class AnalysisDataFrameLocalData
192 public:
193 //! Shorthand for the internal array of values for a frame.
194 typedef std::vector<ValueType> ValueArray;
195 //! Shorthand for a handle to a single frame.
196 typedef AnalysisDataFrameLocalDataHandle<ValueType> FrameHandle;
197 //! Shorthand for a handle to a single data set.
198 typedef AnalysisDataFrameLocalDataSetHandle<ValueType> DataSetHandle;
200 //! Constructs an empty container with a single data set.
201 AnalysisDataFrameLocalData()
203 dataSetColumns_.resize(2);
206 //! Whether init() has been called.
207 bool isInitialized() const { return !values_.empty(); }
208 /*! \brief
209 * Returns number of independent data frames in this object.
211 * This supports looping over all the frame arrays to, e.g., sum them
212 * up at the end in accumulation scenarios.
214 int frameCount() const { return values_.size(); }
216 /*! \brief
217 * Sets the number of data sets stored for each frame.
219 * \throws std::bad_alloc if out of memory.
221 * If not called, there is a single data set in the object.
222 * Cannot be called after init().
224 void setDataSetCount(int dataSetCount)
226 GMX_RELEASE_ASSERT(!isInitialized(),
227 "Cannot change value count after init()");
228 GMX_RELEASE_ASSERT(dataSetCount >= 0,
229 "Invalid data set count");
230 dataSetColumns_.resize(dataSetCount + 1);
232 /*! \brief
233 * Sets the number of columns stored for a data set.
235 * Must be called for each data set that needs to have values,
236 * otherwise there will be zero columns for that data set.
237 * Cannot be called after init().
239 void setColumnCount(int dataSet, int columnCount)
241 GMX_RELEASE_ASSERT(!isInitialized(),
242 "Cannot change value count after init()");
243 GMX_RELEASE_ASSERT(dataSet >= 0 && dataSet < static_cast<int>(dataSetColumns_.size()) - 1,
244 "Invalid data set index");
245 GMX_RELEASE_ASSERT(columnCount >= 0,
246 "Invalid column count");
247 dataSetColumns_[dataSet + 1] = columnCount;
250 /*! \brief
251 * Initializes the storage to support specified parallelism.
253 * \throws std::bad_alloc if out of memory.
255 void init(const AnalysisDataParallelOptions &opt)
257 GMX_RELEASE_ASSERT(!isInitialized(), "init() called multiple times");
258 std::partial_sum(dataSetColumns_.begin(), dataSetColumns_.end(),
259 dataSetColumns_.begin());
260 values_.resize(opt.parallelizationFactor());
261 typename std::vector<ValueArray>::iterator i;
262 for (i = values_.begin(); i != values_.end(); ++i)
264 i->resize(dataSetColumns_.back());
268 //! Returns a handle to access data for a frame.
269 FrameHandle frameData(int frameIndex)
271 GMX_ASSERT(frameIndex >= 0, "Invalid frame index");
272 GMX_ASSERT(isInitialized(), "Cannot access data before init()");
273 return FrameHandle(&dataSetColumns_,
274 &values_[frameIndex % values_.size()]);
276 //! Returns a handle to access a single data set within a frame.
277 DataSetHandle frameDataSet(int frameIndex, int dataSet)
279 return frameData(frameIndex).dataSet(dataSet);
282 private:
283 /*! \brief
284 * Index to find data sets within a per-frame array in `values_`.
286 * The first entry is always zero, followed by one entry for each data
287 * set. Before init(), the data set entries hold the numbers set with
288 * setColumnCount(). After init(), the data set entries hold the
289 * indices of the first column for that data set in the per-frame
290 * arrays in `values_`.
292 std::vector<int> dataSetColumns_;
293 /*! \brief
294 * Data array for each frame.
296 * This is a ring buffer whose size is specified by the desired
297 * parallelism level. For each frame, there is a single array of
298 * values, where the individual data sets are indexed with
299 * `dataSetColumns_`.
301 std::vector<ValueArray> values_;
304 //! \}
306 } // namespace gmx
308 #endif