Implement common analysis data value object.
[alexxy/gromacs.git] / src / gromacs / analysisdata / arraydata.h
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11  * Copyright (c) 2001-2009, The GROMACS development team,
12  * check out http://www.gromacs.org for more information.
13
14  * This program is free software; you can redistribute it and/or
15  * modify it under the terms of the GNU General Public License
16  * as published by the Free Software Foundation; either version 2
17  * of the License, or (at your option) any later version.
18  *
19  * If you want to redistribute modifications, please consider that
20  * scientific software is very special. Version control is crucial -
21  * bugs must be traceable. We will be happy to consider code for
22  * inclusion in the official distribution, but derived work must not
23  * be called official GROMACS. Details are found in the README & COPYING
24  * files - if they are missing, get the official version at www.gromacs.org.
25  *
26  * To help us fund GROMACS development, we humbly ask that you cite
27  * the papers on the package - you can find them in the top README file.
28  *
29  * For more info, check our website at http://www.gromacs.org
30  */
31 /*! \file
32  * \brief
33  * Declares gmx::AbstractAnalysisArrayData and gmx::AnalysisArrayData.
34  *
35  * \author Teemu Murtola <teemu.murtola@cbr.su.se>
36  * \inpublicapi
37  * \ingroup module_analysisdata
38  */
39 #ifndef GMX_ANALYSISDATA_ARRAYDATA_H
40 #define GMX_ANALYSISDATA_ARRAYDATA_H
41
42 #include <vector>
43
44 #include "../fatalerror/gmxassert.h"
45
46 #include "abstractdata.h"
47 #include "dataframe.h"
48
49 namespace gmx
50 {
51
52 /*! \brief
53  * Abstract base class for data objects that present in-memory data.
54  *
55  * This class implements a subclass of AbstractAnalysisData that presents an
56  * in-memory array through the AbstractAnalysisData interface.  Subclasses
57  * should initialize the in-memory array through the provided protected member
58  * functions.
59  *
60  * \todo
61  * Add methods to take full advantage of AnalysisDataValue features.
62  *
63  * \inlibraryapi
64  * \ingroup module_analysisdata
65  */
66 class AbstractAnalysisArrayData : public AbstractAnalysisData
67 {
68     public:
69         virtual ~AbstractAnalysisArrayData();
70
71         /*! \brief
72          * Returns the number of rows in the data array.
73          *
74          * This function is identical to frameCount(), except that frameCount()
75          * returns 0 before valuesReady() has been called.
76          */
77         int rowCount() const { return _nrows; }
78         //! Returns true if values have been allocated.
79         bool isAllocated() const { return !_value.empty(); }
80         //! Returns the x value of the first frame.
81         real xstart() const { return _xstart; }
82         //! Returns the step between frame x values.
83         real xstep() const { return _xstep; }
84         //! Returns the x value of a row.
85         real xvalue(int row) const
86         {
87             GMX_ASSERT(row >= 0 && row < rowCount(), "Row index out of range");
88             return xstart() + row * xstep();
89         }
90         //! Returns a given array element.
91         real value(int row, int col) const
92         {
93             GMX_ASSERT(row >= 0 && row < rowCount(), "Row index out of range");
94             GMX_ASSERT(col >= 0 && col < columnCount(), "Column index out of range");
95             GMX_ASSERT(isAllocated(), "Data array not allocated");
96             return _value[row * columnCount() + col].value();
97         }
98
99     protected:
100         AbstractAnalysisArrayData();
101
102         /*! \brief
103          * Sets the number of columns in the data array.
104          *
105          * Cannot be called after allocateValues().
106          */
107         void setColumnCount(int ncols);
108         /*! \brief
109          * Sets the number of rows in the data array.
110          *
111          * Cannot be called after allocateValues().
112          */
113         void setRowCount(int nrows);
114         /*! \brief
115          * Allocates memory for the values.
116          *
117          * setColumnCount() and setRowCount() must have been called.
118          */
119         void allocateValues();
120         /*! \brief
121          * Sets the values reported as x values for frames.
122          */
123         void setXAxis(real start, real step);
124         //! Returns a reference to a given array element.
125         real &value(int row, int col)
126         {
127             GMX_ASSERT(row >= 0 && row < rowCount(), "Row index out of range");
128             GMX_ASSERT(col >= 0 && col < columnCount(), "Column index out of range");
129             GMX_ASSERT(isAllocated(), "Data array not allocated");
130             return _value[row * columnCount() + col].value();
131         }
132         /*! \brief
133          * Sets the value of an element in the array.
134          */
135         void setValue(int row, int col, real val)
136         {
137             value(row, col) = val;
138         }
139         /*! \brief
140          * Notifies modules of the data.
141          *
142          * This function should be called once the values in the array
143          * have been initialized. The values should not be changed after this
144          * function has been called.
145          */
146         void valuesReady();
147
148         /*! \brief
149          * Copies the contents into a new object.
150          *
151          * \param[in]     src  Object to copy data from.
152          * \param[in,out] dest Empty array data object to copy data to.
153          *
154          * \p dest should not have previous contents.
155          */
156         static void copyContents(const AbstractAnalysisArrayData *src,
157                                  AbstractAnalysisArrayData *dest);
158
159     private:
160         virtual AnalysisDataFrameRef tryGetDataFrameInternal(int index) const;
161         virtual bool requestStorageInternal(int nframes);
162
163         int                  _nrows;
164         std::vector<AnalysisDataValue> _value;
165         real                 _xstart;
166         real                 _xstep;
167         bool                 _bReady;
168
169         // Copy and assign disallowed by base.
170 };
171
172 /*! \brief
173  * Simple in-memory data array.
174  *
175  * This class simply exposes the protected functions of
176  * AbstractAnalysisArrayData to allow construction of simple in-memory data
177  * arrays for input into data modules.
178  *
179  * \inpublicapi
180  * \ingroup module_analysisdata
181  */
182 class AnalysisArrayData : public AbstractAnalysisArrayData
183 {
184     public:
185         AnalysisArrayData() {}
186
187         // TODO: These statements cause Doxygen to generate confusing
188         // documentation.
189         using AbstractAnalysisArrayData::setColumnCount;
190         using AbstractAnalysisArrayData::setRowCount;
191         using AbstractAnalysisArrayData::allocateValues;
192         using AbstractAnalysisArrayData::setXAxis;
193         using AbstractAnalysisArrayData::value;
194         using AbstractAnalysisArrayData::setValue;
195         using AbstractAnalysisArrayData::valuesReady;
196
197         // Copy and assign disallowed by base.
198 };
199
200 } // namespace gmx
201
202 #endif