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