Enable uncrustify for OpenCL
[gromacs.git] / src / gromacs / analysisdata / arraydata.h
blobf255ae70f612e34f2e1dd7829fbccbdcb97854c0
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2010,2011,2012,2013,2014, 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::AbstractAnalysisArrayData and gmx::AnalysisArrayData.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \inpublicapi
41 * \ingroup module_analysisdata
43 #ifndef GMX_ANALYSISDATA_ARRAYDATA_H
44 #define GMX_ANALYSISDATA_ARRAYDATA_H
46 #include <vector>
48 #include "gromacs/analysisdata/abstractdata.h"
49 #include "gromacs/analysisdata/dataframe.h"
50 #include "gromacs/utility/gmxassert.h"
52 namespace gmx
55 /*! \brief
56 * Abstract base class for data objects that present in-memory data.
58 * This class implements a subclass of AbstractAnalysisData that presents an
59 * in-memory array through the AbstractAnalysisData interface. Subclasses
60 * should initialize the in-memory array through the provided protected member
61 * functions. This class provides public accessor methods for read access to
62 * the data.
64 * Public accessor methods in this class do not throw, but assert if data is
65 * accessed before it is available.
67 * \todo
68 * Add support for multiple data sets.
70 * \inlibraryapi
71 * \ingroup module_analysisdata
73 class AbstractAnalysisArrayData : public AbstractAnalysisData
75 public:
76 virtual ~AbstractAnalysisArrayData();
78 virtual int frameCount() const
80 return bReady_ ? rowCount_ : 0;
83 /*! \brief
84 * Returns the number of rows in the data array.
86 * This function is identical to frameCount(), except that frameCount()
87 * returns 0 before valuesReady() has been called.
89 int rowCount() const { return rowCount_; }
90 //! Returns true if values have been allocated.
91 bool isAllocated() const { return !value_.empty(); }
92 //! Returns the x value of the first frame.
93 real xstart() const { return xvalue_[0]; }
94 //! Returns the step between frame x values.
95 real xstep() const
97 GMX_ASSERT(bUniformX_, "Accessing x step for non-uniform data");
98 return xstep_;
100 //! Returns the x value of a row.
101 real xvalue(int row) const
103 GMX_ASSERT(row >= 0 && row < rowCount(), "Row index out of range");
104 return xvalue_[row];
106 //! Returns a given array element.
107 const AnalysisDataValue &value(int row, int col) const
109 GMX_ASSERT(row >= 0 && row < rowCount(), "Row index out of range");
110 GMX_ASSERT(col >= 0 && col < columnCount(), "Column index out of range");
111 GMX_ASSERT(isAllocated(), "Data array not allocated");
112 return value_[row * columnCount() + col];
115 protected:
116 /*! \brief
117 * Initializes an empty array data object.
119 * \throws std::bad_alloc if out of memory.
121 AbstractAnalysisArrayData();
123 /*! \brief
124 * Sets the number of columns in the data array.
126 * \param[in] ncols Number of columns in the data.
128 * Cannot be called after allocateValues().
130 * See AbstractAnalysisData::setColumnCount() for exception behavior.
132 void setColumnCount(int ncols);
133 /*! \brief
134 * Sets the number of rows in the data array.
136 * \param[in] rowCount Number of rows in the data.
138 * Cannot be called after allocateValues().
140 * Does not throw.
142 void setRowCount(int rowCount);
143 /*! \brief
144 * Allocates memory for the values.
146 * \throws std::bad_alloc if memory allocation fails.
148 * setColumnCount() and setRowCount() must have been called.
150 * Strong exception safety guarantee.
152 void allocateValues();
153 /*! \brief
154 * Sets the values reported as x values for frames.
156 * \param[in] start x value for the first frame.
157 * \param[in] step Step between x values of successive frames.
159 * Must not be called after valuesReady().
160 * Any values set with setXAxisValue() are overwritten.
162 * Does not throw.
164 void setXAxis(real start, real step);
165 /*! \brief
166 * Sets a single value reported as x value for frames.
168 * \param[in] row Row/frame for which to set the value.
169 * \param[in] value x value for the frame specified by \p row.
171 * Must not be called after valuesReady().
173 * Does not throw.
175 void setXAxisValue(int row, real value);
176 //! Returns a reference to a given array element.
177 AnalysisDataValue &value(int row, int col)
179 GMX_ASSERT(row >= 0 && row < rowCount(), "Row index out of range");
180 GMX_ASSERT(col >= 0 && col < columnCount(), "Column index out of range");
181 GMX_ASSERT(isAllocated(), "Data array not allocated");
182 return value_[row * columnCount() + col];
184 /*! \brief
185 * Notifies modules of the data.
187 * \throws unspecified Any exception thrown by attached data modules
188 * in data notification methods.
190 * This function should be called once the values in the array
191 * have been initialized. The values should not be changed after this
192 * function has been called.
194 void valuesReady();
196 /*! \brief
197 * Copies the contents into a new object.
199 * \param[in] src Object to copy data from.
200 * \param[in,out] dest Empty array data object to copy data to.
201 * \throws std::bad_alloc if memory allocation for \p dest fails.
203 * \p dest should not have previous contents.
205 static void copyContents(const AbstractAnalysisArrayData *src,
206 AbstractAnalysisArrayData *dest);
208 private:
209 virtual AnalysisDataFrameRef tryGetDataFrameInternal(int index) const;
210 virtual bool requestStorageInternal(int nframes);
212 int rowCount_;
213 AnalysisDataPointSetInfo pointSetInfo_;
214 std::vector<AnalysisDataValue> value_;
215 std::vector<real> xvalue_;
216 real xstep_;
217 bool bUniformX_;
218 bool bReady_;
220 // Copy and assign disallowed by base.
223 /*! \brief
224 * Simple in-memory data array.
226 * This class is a simple alternative to AnalysisData for in-memory data arrays
227 * that are constructed in-place.
229 * Public accessor methods in this class do not throw, but assert if data is
230 * accessed before it is available.
232 * \if libapi
233 * This class exposes the protected functions of AbstractAnalysisArrayData for
234 * users.
235 * \endif
237 * \inpublicapi
238 * \ingroup module_analysisdata
240 class AnalysisArrayData : public AbstractAnalysisArrayData
242 public:
243 /*! \brief
244 * Initializes an empty array data object.
246 * \throws std::bad_alloc if out of memory.
248 AnalysisArrayData() {}
250 // TODO: These statements cause Doxygen to generate confusing
251 // documentation.
252 using AbstractAnalysisArrayData::setColumnCount;
253 using AbstractAnalysisArrayData::setRowCount;
254 using AbstractAnalysisArrayData::allocateValues;
255 using AbstractAnalysisArrayData::setXAxis;
256 using AbstractAnalysisArrayData::setXAxisValue;
257 using AbstractAnalysisArrayData::value;
258 using AbstractAnalysisArrayData::valuesReady;
260 // Copy and assign disallowed by base.
263 } // namespace gmx
265 #endif