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