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