Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / analysisdata / tests / datatest.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 /*! \libinternal \file
36  * \brief
37  * Helper classes for testing classes that derive from AbstractAnalysisData.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \inlibraryapi
41  * \ingroup module_analysisdata
42  */
43 #ifndef GMX_ANALYSISDATA_TESTS_DATATEST_H
44 #define GMX_ANALYSISDATA_TESTS_DATATEST_H
45
46 #include <vector>
47
48 #include <gtest/gtest.h>
49
50 #include "gromacs/analysisdata/dataframe.h"
51 #include "gromacs/utility/gmxassert.h"
52 #include "gromacs/utility/real.h"
53
54 #include "testutils/refdata.h"
55
56 namespace gmx
57 {
58
59 class AbstractAnalysisData;
60 class AnalysisData;
61 class AnalysisDataHandle;
62
63 namespace test
64 {
65
66 class FloatingPointTolerance;
67
68 /*! \libinternal \brief
69  * Represents a single set of points in AnalysisDataTestInputFrame structure.
70  *
71  * If the data is not multipoint, each frame contains exactly one set of
72  * points.  If there is more than one set of points, each of these sets is
73  * passed separately and AnalysisDataHandle::finishPointSet() called in
74  * between.
75  *
76  * \inlibraryapi
77  * \ingroup module_analysisdata
78  */
79 class AnalysisDataTestInputPointSet
80 {
81     public:
82         //! Returns zero-based index of this point set in its frame.
83         int index() const { return index_; }
84         //! Returns zero-based index of the data set of this point set.
85         int dataSetIndex() const { return dataSetIndex_; }
86         //! Returns zero-based index of the first column in this point set.
87         int firstColumn() const { return firstColumn_; }
88         //! Returns zero-based index of the last column in this point set.
89         int lastColumn() const { return firstColumn_ + size() - 1; }
90         //! Returns the number of columns in the point set.
91         int size() const { return values_.size(); }
92         //! Returns the value in column \p i.
93         real y(int i) const { return values_[i].y; }
94         //! Returns whether the error is present for column \p i.
95         bool hasError(int i) const { return values_[i].bError; }
96         //! Returns the error in column \p i.
97         real error(int i) const { return values_[i].error; }
98         //! Returns whether the value in column \p i is present.
99         bool present(int /*i*/) const { return true; }
100         //! Returns an AnalysisDataValue for column \p i.
101         AnalysisDataValue value(int i) const
102         {
103             AnalysisDataValue result;
104             result.setValue(values_[i].y);
105             if (values_[i].bError)
106             {
107                 result.setError(values_[i].error);
108             }
109             return result;
110         }
111
112         //! Appends a value to this point set.
113         void addValue(real y) { values_.push_back(Value(y)); }
114         //! Appends a value with an error estimate to this point set.
115         void addValueWithError(real y, real error)
116         {
117             values_.push_back(Value(y, error));
118         }
119
120     private:
121         //! Creates an empty point set.
122         AnalysisDataTestInputPointSet(int index, int dataSetIndex,
123                                       int firstColumn);
124
125         struct Value
126         {
127             Value() : y(0.0), error(0.0), bError(false) {}
128             explicit Value(real y) : y(y), error(0.0), bError(false) {}
129             Value(real y, real error) : y(y), error(error), bError(true) {}
130
131             real                y;
132             real                error;
133             bool                bError;
134         };
135
136         int                     index_;
137         int                     dataSetIndex_;
138         int                     firstColumn_;
139         std::vector<Value>      values_;
140
141         //! For constructing new point sets.
142         friend class AnalysisDataTestInputFrame;
143 };
144
145 /*! \libinternal \brief
146  * Represents a single frame in AnalysisDataTestInput structure.
147  *
148  * \inlibraryapi
149  * \ingroup module_analysisdata
150  */
151 class AnalysisDataTestInputFrame
152 {
153     public:
154         //! Returns zero-based index for the frame.
155         int index() const { return index_; }
156         //! Returns x coordinate for the frame.
157         real x() const { return x_; }
158         //! Returns error in the x coordinate for the frame.
159         real dx() const { return 0.0; }
160
161         //! Number of individual point sets in the frame.
162         int pointSetCount() const { return pointSets_.size(); }
163         //! Returns a point set object for a given point set.
164         const AnalysisDataTestInputPointSet &pointSet(int index) const
165         {
166             GMX_ASSERT(index >= 0 && static_cast<size_t>(index) < pointSets_.size(),
167                        "Point set index out of range");
168             return pointSets_[index];
169         }
170
171         //! Appends an empty point set to this frame.
172         AnalysisDataTestInputPointSet &addPointSet(int dataSet, int firstColumn);
173         //! Adds a point set with given values to this frame.
174         void addPointSetWithValues(int dataSet, int firstColumn, real y1);
175         //! Adds a point set with given values to this frame.
176         void addPointSetWithValues(int dataSet, int firstColumn,
177                                    real y1, real y2);
178         //! Adds a point set with given values to this frame.
179         void addPointSetWithValues(int dataSet, int firstColumn,
180                                    real y1, real y2, real y3);
181         //! Adds a point set with given values to this frame.
182         void addPointSetWithValueAndError(int dataSet, int firstColumn,
183                                           real y1, real e1);
184
185     private:
186         //! Constructs a new frame object with the given values.
187         AnalysisDataTestInputFrame(int index, real x);
188
189         int                                         index_;
190         real                                        x_;
191         std::vector<AnalysisDataTestInputPointSet>  pointSets_;
192
193         //! For constructing new frames.
194         friend class AnalysisDataTestInput;
195 };
196
197 /*! \libinternal \brief
198  * Represents static input data for AbstractAnalysisData tests.
199  *
200  * Used to construct structured test input data for analysis data unit tests.
201  * Typically used as input to methods in AnalysisDataTestFixture.
202  *
203  * \see AnalysisDataTestFixture
204  *
205  * \inlibraryapi
206  * \ingroup module_analysisdata
207  */
208 class AnalysisDataTestInput
209 {
210     public:
211         /*! \brief
212          * Constructs empty input data.
213          *
214          * \param[in] dataSetCount Number of data sets in the data.
215          * \param[in] bMultipoint  Whether the data will be multipoint.
216          *
217          * The column count for each data set must be set with
218          * setColumnCount().
219          */
220         AnalysisDataTestInput(int dataSetCount, bool bMultipoint);
221         ~AnalysisDataTestInput();
222
223         //! Whether the input data is multipoint.
224         bool isMultipoint() const { return bMultipoint_; }
225         //! Returns the number of data sets in the input data.
226         int dataSetCount() const { return columnCounts_.size(); }
227         //! Returns the number of columns in a given data set.
228         int columnCount(int dataSet) const { return columnCounts_[dataSet]; }
229         //! Returns the number of frames in the input data.
230         int frameCount() const { return frames_.size(); }
231         //! Returns a frame object for the given input frame.
232         const AnalysisDataTestInputFrame &frame(int index) const;
233
234         //! Sets the number of columns in a data set.
235         void setColumnCount(int dataSet, int columnCount);
236         //! Appends an empty frame to this data.
237         AnalysisDataTestInputFrame &addFrame(real x);
238         //! Adds a frame with a single point set and the given values.
239         void addFrameWithValues(real x, real y1);
240         //! Adds a frame with a single point set and the given values.
241         void addFrameWithValues(real x, real y1, real y2);
242         //! Adds a frame with a single point set and the given values.
243         void addFrameWithValues(real x, real y1, real y2, real y3);
244         //! Adds a frame with a single point set and the given values.
245         void addFrameWithValueAndError(real x, real y1, real e1);
246
247     private:
248         std::vector<int>                        columnCounts_;
249         bool                                    bMultipoint_;
250         std::vector<AnalysisDataTestInputFrame> frames_;
251 };
252
253 /*! \libinternal \brief
254  * Test fixture for AbstractAnalysisData testing.
255  *
256  * This test fixture is designed to help writing tests for objects that
257  * derive from the AbstractAnalysisData class.  Typical flow in such tests is
258  * that first the test creates the required data objects, then call static
259  * methods in this class to add mock modules (using
260  * AbstractAnalysisData::addModule()) to the data objects to check that they
261  * produce the correct data, and then invokes methods in the data object to
262  * produce the data to be checked.  Static methods are also provided for
263  * pushing data from an AnalysisDataTestInput object to some generic types
264  * derived from AbstractAnalysisData.
265  *
266  * Methods addStaticCheckerModule(), addStaticColumnCheckerModule() and
267  * addStaticStorageCheckerModule() create and add mock modules that check the
268  * data against a given AnalysisDataTestInput instance.
269  *
270  * Method addReferenceCheckerModule() creates and adds a mock module that
271  * checks the output against reference data produced by a previous test
272  * execution (see TestReferenceData).  Two versions are provided, a static
273  * method to be used with any TestReferenceChecker, and a non-static method
274  * that uses the protected \p data_ member.
275  *
276  * presentAllData() and presentDataFrame() are provided to push data from an
277  * AnalysisDataTestInput into an AnalysisData object.  In typical tests, most
278  * checks are done during these methods, by the added mock modules.
279  * setupArrayData() performs the same function for classes derived from
280  * AbstractAnalysisArrayData.  In that case, the test should separately ensure
281  * that AbstractAnalysisArrayData::valuesReady() gets called.
282  *
283  * \todo
284  * Support for arbitrary AnalysisDataValues (errors and missing values).
285  *
286  * \see AnalysisDataTestInput
287  *
288  * \inlibraryapi
289  * \ingroup module_analysisdata
290  */
291 class AnalysisDataTestFixture : public ::testing::Test
292 {
293     public:
294         AnalysisDataTestFixture();
295
296         /*! \brief
297          * Initializes an AnalysisData object from input data.
298          *
299          * Sets the column count and other properties based on the input data.
300          */
301         static void setupDataObject(const AnalysisDataTestInput &input,
302                                     AnalysisData                *data);
303
304         /*! \brief
305          * Adds all data from AnalysisDataTestInput into an AnalysisData.
306          */
307         static void presentAllData(const AnalysisDataTestInput &input,
308                                    AnalysisData                *data);
309         /*! \brief
310          * Adds a single frame from AnalysisDataTestInput into an AnalysisData.
311          */
312         static void presentDataFrame(const AnalysisDataTestInput &input, int row,
313                                      AnalysisDataHandle handle);
314         /*! \brief
315          * Initializes an array data object from AnalysisDataTestInput.
316          *
317          * \tparam ArrayData  Class derived from AbstractAnalysisArrayData.
318          *
319          * The ArrayData class should expose the setter methods
320          * (setColumnCount(), setRowCount(), allocateValues(), setValue())
321          * publicly or declare the fixture class as a friend.
322          * The X axis in \p data must be configured to match \p input before
323          * calling this method.
324          *
325          * Does not call AbstractAnalysisArrayData::valuesReady().
326          * The test must ensure that this method gets called, otherwise the
327          * mock modules never get called.
328          */
329         template <class ArrayData>
330         static void setupArrayData(const AnalysisDataTestInput &input,
331                                    ArrayData                   *data);
332
333         /*! \brief
334          * Adds a mock module that verifies output against
335          * AnalysisDataTestInput.
336          *
337          * \param[in]  data     Data to compare against.
338          * \param      source   Data object to verify.
339          *
340          * Creates a mock module that verifies that the
341          * AnalysisDataModuleInterface methods are called correctly by
342          * \p source.  Parameters for the calls are verified against \p data.
343          * Adds the created module to \p source using \p data->addModule().
344          * Any exceptions from the called functions should be caught by the
345          * caller.
346          *
347          * \see AbstractAnalysisData::addModule()
348          */
349         static void addStaticCheckerModule(const AnalysisDataTestInput &data,
350                                            AbstractAnalysisData        *source);
351         /*! \brief
352          * Adds a mock module that verifies parallel output against
353          * AnalysisDataTestInput.
354          *
355          * \param[in]  data     Data to compare against.
356          * \param      source   Data object to verify.
357          *
358          * Creates a parallel mock module that verifies that the
359          * AnalysisDataModuleInterface methods are called correctly by
360          * \p source.  Parameters for the calls are verified against \p data.
361          * Adds the created module to \p source using \p data->addModule().
362          * Any exceptions from the called functions should be caught by the
363          * caller.
364          *
365          * Differs from addStaticCheckerModule() in that the created mock
366          * module reports that it accepts parallel input data, and accepts and
367          * verifies notification calls following the parallel pattern.
368          *
369          * \see AbstractAnalysisData::addModule()
370          */
371         static void addStaticParallelCheckerModule(
372             const AnalysisDataTestInput &data,
373             AbstractAnalysisData        *source);
374         /*! \brief
375          * Adds a column mock module that verifies output against
376          * AnalysisDataTestInput.
377          *
378          * \param[in]  data     Data to compare against.
379          * \param[in]  firstcol First column to check.
380          * \param[in]  n        Number of columns to check.
381          * \param      source   Data object to verify.
382          *
383          * Creates a mock module that verifies that the
384          * AnalysisDataModuleInterface methods are called correctly by
385          * \p source.  Parameters for the calls are verified against \p data.
386          * Adds the created module to \p source using
387          * \p data->addColumnModule().
388          * Any exceptions from the called functions should be caught by the
389          * caller.
390          *
391          * \see AbstractAnalysisData::addColumnModule()
392          */
393         static void addStaticColumnCheckerModule(const AnalysisDataTestInput &data,
394                                                  int firstcol, int n,
395                                                  AbstractAnalysisData *source);
396         /*! \brief
397          * Adds a mock module that verifies output and storage against
398          * AnalysisDataTestInput.
399          *
400          * \param[in]  data     Data to compare against.
401          * \param[in]  storageCount  Number of previous frames to check
402          *                      (-1 = all).
403          * \param      source   Data object to verify.
404          *
405          * Works like addStaticCheckerModule(), except that in addition, for
406          * each frame, the mock module also checks that previous frames can be
407          * accessed using AbstractAnalysisData::getDataFrame().  In the
408          * AnalysisDataModuleInterface::dataStarted() callback, the mock module
409          * calls AbstractAnalysisData::requestStorage() with \p storageCount as
410          * the parameter.
411          */
412         static void addStaticStorageCheckerModule(const AnalysisDataTestInput &data,
413                                                   int                          storageCount,
414                                                   AbstractAnalysisData        *source);
415         /*! \brief
416          * Adds a mock module that verifies output against reference data.
417          *
418          * \param[in]  checker   Reference data checker to use for comparison.
419          * \param[in]  id        Identifier for reference data compound to use.
420          * \param      source    Data object to verify.
421          * \param[in]  tolerance Tolerance to use for comparison.
422          *
423          * Creates a mock module that verifies that the
424          * AnalysisDataModuleInterface methods are called correctly by
425          * \p source.  Parameters for the calls are verified against reference
426          * data using a child compound \p id of \p checker.
427          * Adds the created module to \p source using \p data->addModule().
428          * Any exceptions from the called functions should be caught by the
429          * caller.
430          *
431          * \see TestReferenceData
432          */
433         static void addReferenceCheckerModule(TestReferenceChecker          checker,
434                                               const char                   *id,
435                                               AbstractAnalysisData         *source,
436                                               const FloatingPointTolerance &tolerance);
437
438         /*! \brief
439          * Adds a mock module that verifies output against reference data.
440          *
441          * \param[in]  id       Identifier for reference data compound to use.
442          * \param      source   Data object to verify.
443          *
444          * Creates a reference checker module using a compound checker with id
445          * \p id at the root level of \p data_.
446          *
447          * See the static overload for other details.
448          */
449         void addReferenceCheckerModule(const char           *id,
450                                        AbstractAnalysisData *source);
451
452     protected:
453         /*! \brief
454          * Reference data object used for the reference checker modules.
455          *
456          * Tests can use the data object also for their own purposes if needed.
457          */
458         gmx::test::TestReferenceData  data_;
459 };
460
461
462 template <class ArrayData>
463 void AnalysisDataTestFixture::setupArrayData(const AnalysisDataTestInput &input,
464                                              ArrayData                   *data)
465 {
466     GMX_RELEASE_ASSERT(!input.isMultipoint(),
467                        "Array data cannot be initialized from multipoint data");
468     GMX_RELEASE_ASSERT(input.dataSetCount() == 1,
469                        "Array data cannot be initialized from multiple data sets");
470     GMX_RELEASE_ASSERT(data->columnCount() == 0 || data->columnCount() == input.columnCount(0),
471                        "Mismatching input and target data");
472     GMX_RELEASE_ASSERT(data->rowCount() == 0 || data->rowCount() == input.frameCount(),
473                        "Mismatching input and target data");
474     data->setColumnCount(input.columnCount(0));
475     data->setRowCount(input.frameCount());
476     data->allocateValues();
477     for (int row = 0; row < input.frameCount(); ++row)
478     {
479         const AnalysisDataTestInputFrame    &frame = input.frame(row);
480         EXPECT_FLOAT_EQ(frame.x(), data->xvalue(row));
481         GMX_RELEASE_ASSERT(frame.pointSetCount() == 1,
482                            "Multiple point sets not supported by array data");
483         const AnalysisDataTestInputPointSet &points = frame.pointSet(0);
484         for (int column = 0; column < points.size(); ++column)
485         {
486             data->value(row, column + points.firstColumn()) = points.value(column);
487         }
488     }
489 }
490
491 } // namespace test
492
493 } // namespace gmx
494
495 #endif