Improve AbstractAnalysisData data access interface.
[alexxy/gromacs.git] / src / gromacs / analysisdata / tests / datatest.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 /*! \libinternal \file
32  * \brief
33  * Helper classes for testing classes that derive from AbstractAnalysisData.
34  *
35  * \author Teemu Murtola <teemu.murtola@cbr.su.se>
36  * \inlibraryapi
37  * \ingroup module_analysisdata
38  */
39 #ifndef GMX_ANALYSISDATA_TESTS_DATATEST_H
40 #define GMX_ANALYSISDATA_TESTS_DATATEST_H
41
42 #include <limits>
43 #include <vector>
44
45 #include <gtest/gtest.h>
46
47 #include "gromacs/legacyheaders/types/simple.h"
48 #include "gromacs/fatalerror/gmxassert.h"
49
50 #include "testutils/refdata.h"
51
52 namespace gmx
53 {
54
55 class AbstractAnalysisData;
56 class AnalysisData;
57 class AnalysisDataHandle;
58
59 namespace test
60 {
61
62 class MockAnalysisModule;
63
64 //! Constant to use to signify end of data for AnalysisDataTestInput.
65 const real END_OF_DATA = std::numeric_limits<real>::max();
66 //! Constant to use to signify end of one data frame for AnalysisDataTestInput.
67 const real END_OF_FRAME = std::numeric_limits<real>::min();
68 //! Constant to use to signify end of one multipoint set for AnalysisDataTestInput.
69 const real MPSTOP = -std::numeric_limits<real>::max();
70
71 /*! \libinternal \brief
72  * Represents a single set of points in AnalysisDataTestInputFrame structure.
73  *
74  * If the data is not multipoint, each frame contains exactly one set of
75  * points.  If there is more than one set of points, each of these sets is
76  * passed separately to AnalysisDataHandle::addPoints().
77  *
78  * \ingroup module_analysisdata
79  */
80 class AnalysisDataTestInputPointSet
81 {
82     public:
83         AnalysisDataTestInputPointSet();
84
85         int size() const { return y_.size(); }
86         real y(int i) const { return y_[i]; }
87         const std::vector<real> &yvector() const { return y_; }
88         const real *yptr() const { return &y_[0]; }
89         const real *dyptr() const { return NULL; }
90         const bool *presentptr() const { return NULL; }
91
92     private:
93         std::vector<real>       y_;
94
95         friend class AnalysisDataTestInput;
96 };
97
98 /*! \libinternal \brief
99  * Represents a single frame in AnalysisDataTestInput structure.
100  *
101  * \ingroup module_analysisdata
102  */
103 class AnalysisDataTestInputFrame
104 {
105     public:
106         AnalysisDataTestInputFrame();
107
108         bool isMultipoint() const { return points_.size() > 1; }
109
110         int index() const { return index_; }
111         real x() const { return x_; }
112         real dx() const { return 0.0; }
113
114         int pointSetCount() const { return points_.size(); }
115         const AnalysisDataTestInputPointSet &points(int index = 0) const
116         {
117             GMX_ASSERT(index >= 0 && static_cast<size_t>(index) < points_.size(),
118                        "Point set index out of range");
119             return points_[index];
120         }
121
122     private:
123         int                     index_;
124         real                    x_;
125         std::vector<AnalysisDataTestInputPointSet>  points_;
126
127         friend class AnalysisDataTestInput;
128 };
129
130 /*! \libinternal \brief
131  * Represents static input data for AbstractAnalysisData tests.
132  *
133  * Used to construct structured test input data from a static array of reals,
134  * and then typically used as input to methods in AnalysisDataTestFixture.
135  *
136  * \see AnalysisDataTestFixture
137  *
138  * \ingroup module_analysisdata
139  */
140 class AnalysisDataTestInput
141 {
142     public:
143         /*! \brief
144          * Constructs data representation from a simple array.
145          *
146          * \param[in] data  Array to construct data from.
147          *
148          * The input array should consist of a set of frames, separated by a
149          * END_OF_FRAME marker.  The first value for a frame is the X value,
150          * all following values are Y values.
151          * For multipoint data, one frame can contain several point sets,
152          * separated by MPSTOP markers.  There should be no MPSTOP marker after
153          * the last point set, only an END_OF_FRAME marker.  All point sets are
154          * assumed to start from column zero, but the sets may contain
155          * different number of columns.  For non-multipoint data, all frames
156          * must containt the same number of columns.
157          * The final element in the array (after the last END_OF_FRAME) should
158          * be END_OF_DATA.
159          */
160         explicit AnalysisDataTestInput(const real *data);
161         ~AnalysisDataTestInput();
162
163         int frameCount() const { return frames_.size(); }
164         int columnCount() const { return columnCount_; }
165         bool isMultipoint() const { return bMultipoint_; }
166         const AnalysisDataTestInputFrame &frame(int index) const;
167
168     private:
169         int                     columnCount_;
170         bool                    bMultipoint_;
171         std::vector<AnalysisDataTestInputFrame> frames_;
172 };
173
174 /*! \libinternal \brief
175  * Test fixture for AbstractAnalysisData testing.
176  *
177  * This test fixture is designed to help writing tests for objects that
178  * derive from the AbstractAnalysisData class.  Typical flow in such tests is
179  * that first the test creates the required data objects, then call static
180  * methods in this class to add mock modules (using
181  * AbstractAnalysisData::addModule()) to the data objects to check that they
182  * produce the correct data, and then invokes methods in the data object to
183  * produce the data to be checked.  Static methods are also provided for
184  * pushing data from an AnalysisDataTestInput object to some generic types
185  * derived from AbstractAnalysisData.
186  *
187  * Methods addStaticCheckerModule(), addStaticColumnCheckerModule() and
188  * addStaticStorageCheckerModule() create and add mock modules that check the
189  * data against a given AnalysisDataTestInput instance.
190  *
191  * Method addReferenceCheckerModule() creates and adds a mock module that
192  * checks the output against reference data produced by a previous test
193  * execution (see TestReferenceData).  Two versions are provided, a static
194  * method to be used with any TestReferenceChecker, and a non-static method
195  * that uses the protected \p data_ member.
196  *
197  * presentAllData() and presentDataFrame() are provided to push data from an
198  * AnalysisDataTestInput into an AnalysisData object.  In typical tests, most
199  * checks are done during the these methods, by the added mock modules.
200  * setupArrayData() performs the same function for classes derived from
201  * AbstractAnalysisArrayData.  In that case, the test should separately ensure
202  * that AbstractAnalysisArrayData::valuesReady() gets called.
203  *
204  * \todo
205  * Support for errors and for arbitrary multipoint data.
206  *
207  * \see AnalysisDataTestInput
208  *
209  * \ingroup module_analysisdata
210  */
211 class AnalysisDataTestFixture : public ::testing::Test
212 {
213     public:
214         AnalysisDataTestFixture();
215
216         /*! \brief
217          * Adds all data from AnalysisDataTestInput into an AnalysisData.
218          */
219         static void presentAllData(const AnalysisDataTestInput &input,
220                                    AnalysisData *data);
221         /*! \brief
222          * Adds a single frame from AnalysisDataTestInput into an AnalysisData.
223          */
224         static void presentDataFrame(const AnalysisDataTestInput &input, int row,
225                                      AnalysisDataHandle *handle);
226         /*! \brief
227          * Initializes an array data object from AnalysisDataTestInput.
228          *
229          * \tparam ArrayData  Class derived from AbstractAnalysisArrayData.
230          *
231          * The ArrayData class should expose the setter methods
232          * (setColumnCount(), setRowCount(), allocateValues(), setValue())
233          * publicly or declare the fixture class as a friend.
234          * The X axis in \p data must be configured to match \p input before
235          * calling this method.
236          *
237          * Does not call AbstractAnalysisArrayData::valuesReady().
238          * The test must ensure that this method gets called, otherwise the
239          * mock modules never get called.
240          */
241         template <class ArrayData>
242         static void setupArrayData(const AnalysisDataTestInput &input,
243                                    ArrayData *data);
244
245         /*! \brief
246          * Adds a mock module that verifies output against
247          * AnalysisDataTestInput.
248          *
249          * \param[in]  data     Data to compare against.
250          * \param      source   Data object to verify.
251          *
252          * Creates a mock module that verifies that the
253          * AnalysisDataModuleInterface methods are called correctly by
254          * \p source.  Parameters for the calls are verified against \p data.
255          * Adds the created module to \p source using \p data->addModule().
256          * Any exceptions from the called functions should be caught by the
257          * caller.
258          *
259          * \see AbstractAnalysisData::addModule()
260          */
261         static void addStaticCheckerModule(const AnalysisDataTestInput &data,
262                                            AbstractAnalysisData *source);
263         /*! \brief
264          * Adds a column mock module that verifies output against
265          * AnalysisDataTestInput.
266          *
267          * \param[in]  data     Data to compare against.
268          * \param[in]  firstcol First column to check.
269          * \param[in]  n        Number of columns to check.
270          * \param      source   Data object to verify.
271          *
272          * Creates a mock module that verifies that the
273          * AnalysisDataModuleInterface methods are called correctly by
274          * \p source.  Parameters for the calls are verified against \p data.
275          * Adds the created module to \p source using
276          * \p data->addColumnModule().
277          * Any exceptions from the called functions should be caught by the
278          * caller.
279          *
280          * \see AbstractAnalysisData::addColumnModule()
281          */
282         static void addStaticColumnCheckerModule(const AnalysisDataTestInput &data,
283                                                  int firstcol, int n,
284                                                  AbstractAnalysisData *source);
285         /*! \brief
286          * Adds a mock module that verifies output and storage against
287          * AnalysisDataTestInput.
288          *
289          * \param[in]  data     Data to compare against.
290          * \param[in]  storageCount  Number of previous frames to check
291          *                      (-1 = all).
292          * \param      source   Data object to verify.
293          *
294          * Works like addStaticCheckerModule(), except that in addition, for
295          * each frame, the mock module also checks that previous frames can be
296          * accessed using AbstractAnalysisData::getDataFrame().  In the
297          * AnalysisDataModuleInterface::dataStarted() callback, the mock module
298          * calls AbstractAnalysisData::requestStorage() with \p storageCount as
299          * the parameter.
300          */
301         static void addStaticStorageCheckerModule(const AnalysisDataTestInput &data,
302                                                   int storageCount,
303                                                   AbstractAnalysisData *source);
304         /*! \brief
305          * Adds a mock module that verifies output against reference data.
306          *
307          * \param[in]  checker  Reference data checker to use for comparison.
308          * \param      source   Data object to verify.
309          *
310          * Creates a mock module that verifies that the
311          * AnalysisDataModuleInterface methods are called correctly by
312          * \p source.  Parameters for the calls are verified against reference
313          * data using \p checker.
314          * Adds the created module to \p source using \p data->addModule().
315          * Any exceptions from the called functions should be caught by the
316          * caller.
317          *
318          * \see TestReferenceData
319          */
320         static void addReferenceCheckerModule(const TestReferenceChecker &checker,
321                                               AbstractAnalysisData *source);
322
323         /*! \brief
324          * Adds a mock module that verifies output against reference data.
325          *
326          * \param[in]  id       Identifier for reference data compound to use.
327          * \param      source   Data object to verify.
328          *
329          * Creates a reference checker module using a compound checker with id
330          * \p id at the root level of \p data_.
331          *
332          * See the static overload for other details.
333          */
334         void addReferenceCheckerModule(const char *id,
335                                        AbstractAnalysisData *source);
336
337     protected:
338         gmx::test::TestReferenceData  data_;
339 };
340
341
342 template <class ArrayData>
343 void AnalysisDataTestFixture::setupArrayData(const AnalysisDataTestInput &input,
344                                              ArrayData *data)
345 {
346     GMX_RELEASE_ASSERT(!input.isMultipoint(),
347                        "Array data cannot be initialized from multipoint data");
348     GMX_RELEASE_ASSERT(data->columnCount() == 0 || data->columnCount() == input.columnCount(),
349                        "Mismatching input and target data");
350     GMX_RELEASE_ASSERT(data->rowCount() == 0 || data->rowCount() == input.frameCount(),
351                        "Mismatching input and target data");
352     data->setColumnCount(input.columnCount());
353     data->setRowCount(input.frameCount());
354     data->allocateValues();
355     for (int row = 0; row < input.frameCount(); ++row)
356     {
357         const AnalysisDataTestInputFrame &frame = input.frame(row);
358         EXPECT_FLOAT_EQ(frame.x(), data->xvalue(row));
359         const AnalysisDataTestInputPointSet &points = frame.points();
360         for (int column = 0; column < input.columnCount(); ++column)
361         {
362             data->setValue(row, column, points.y(column));
363         }
364     }
365 }
366
367 } // namespace test
368
369 } // namespace gmx
370
371 #endif