Added unit tests for array data and minor improvements.
[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 <cstddef>
43
44 #include <vector>
45
46 #include "../fatalerror/gmxassert.h"
47
48 #include "abstractdata.h"
49
50 namespace gmx
51 {
52
53 /*! \brief
54  * Abstract base class for data objects that present in-memory data.
55  *
56  * This class implements a subclass of AbstractAnalysisData that presents an
57  * in-memory array through the AbstractAnalysisData interface.  Subclasses
58  * should initialize the in-memory array through the provided protected member
59  * functions.
60  *
61  * \inlibraryapi
62  * \ingroup module_analysisdata
63  */
64 class AbstractAnalysisArrayData : public AbstractAnalysisData
65 {
66     public:
67         virtual ~AbstractAnalysisArrayData();
68
69         virtual bool getDataWErr(int index, real *x, real *dx,
70                                  const real **y, const real **dy,
71                                  const bool **present = 0) const;
72         virtual bool requestStorage(int nframes = -1);
73
74         /*! \brief
75          * Returns the number of rows in the data array.
76          *
77          * This function is identical to frameCount(), except that frameCount()
78          * returns 0 before valuesReady() has been called.
79          */
80         int rowCount() const { return _nrows; }
81         //! Returns true if values have been allocated.
82         bool isAllocated() const { return !_value.empty(); }
83         //! Returns the x value of the first frame.
84         real xstart() const { return _xstart; }
85         //! Returns the step between frame x values.
86         real xstep() const { return _xstep; }
87         //! Returns the x value of a row.
88         real xvalue(int row) const
89         {
90             GMX_ASSERT(row >= 0 && row < rowCount(), "Row index out of range");
91             return xstart() + row * xstep();
92         }
93         //! Returns a given array element.
94         real value(int row, int col) const
95         {
96             GMX_ASSERT(row >= 0 && row < rowCount(), "Row index out of range");
97             GMX_ASSERT(col >= 0 && col < columnCount(), "Column index out of range");
98             GMX_ASSERT(isAllocated(), "Data array not allocated");
99             return _value[row * columnCount() + col];
100         }
101
102     protected:
103         AbstractAnalysisArrayData();
104
105         /*! \brief
106          * Sets the number of columns in the data array.
107          *
108          * Cannot be called after allocateValues().
109          */
110         void setColumnCount(int ncols);
111         /*! \brief
112          * Sets the number of rows in the data array.
113          *
114          * Cannot be called after allocateValues().
115          */
116         void setRowCount(int nrows);
117         /*! \brief
118          * Allocates memory for the values.
119          *
120          * setColumnCount() and setRowCount() must have been called.
121          */
122         void allocateValues();
123         /*! \brief
124          * Sets the values reported as x values for frames.
125          */
126         void setXAxis(real start, real step);
127         //! Returns a reference to a given array element.
128         real &value(int row, int col)
129         {
130             GMX_ASSERT(row >= 0 && row < rowCount(), "Row index out of range");
131             GMX_ASSERT(col >= 0 && col < columnCount(), "Column index out of range");
132             GMX_ASSERT(isAllocated(), "Data array not allocated");
133             return _value[row * columnCount() + col];
134         }
135         /*! \brief
136          * Sets the value of an element in the array.
137          */
138         void setValue(int row, int col, real val)
139         {
140             value(row, col) = val;
141         }
142         /*! \brief
143          * Notifies modules of the data.
144          *
145          * This function should be called once the values in the array
146          * have been initialized. The values should not be changed after this
147          * function has been called.
148          */
149         void valuesReady();
150
151         /*! \brief
152          * Copies the contents into a new object.
153          *
154          * \param[in]     src  Object to copy data from.
155          * \param[in,out] dest Empty array data object to copy data to.
156          *
157          * \p dest should not have previous contents.
158          */
159         static void copyContents(const AbstractAnalysisArrayData *src,
160                                  AbstractAnalysisArrayData *dest);
161
162     private:
163         int                  _nrows;
164         std::vector<real>    _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