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