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