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