6f1b57ee446e9694d5fb990897485031cad33f17
[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,2015,2016,2018,2019, 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 __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_.emplace_back(y); }
119     //! Appends a value with an error estimate to this point set.
120     void addValueWithError(real y, real error) { values_.emplace_back(y, error); }
121
122 private:
123     //! Creates an empty point set.
124     AnalysisDataTestInputPointSet(int index, int dataSetIndex, 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, real y1, real y2);
178     //! Adds a point set with given values to this frame.
179     void addPointSetWithValues(int dataSet, int firstColumn, real y1, real y2, real y3);
180     //! Adds a point set with given values to this frame.
181     void addPointSetWithValueAndError(int dataSet, int firstColumn, real y1, real e1);
182
183 private:
184     //! Constructs a new frame object with the given values.
185     AnalysisDataTestInputFrame(int index, real x);
186
187     int                                        index_;
188     real                                       x_;
189     std::vector<AnalysisDataTestInputPointSet> pointSets_;
190
191     //! For constructing new frames.
192     friend class AnalysisDataTestInput;
193 };
194
195 /*! \libinternal \brief
196  * Represents static input data for AbstractAnalysisData tests.
197  *
198  * Used to construct structured test input data for analysis data unit tests.
199  * Typically used as input to methods in AnalysisDataTestFixture.
200  *
201  * \see AnalysisDataTestFixture
202  *
203  * \inlibraryapi
204  * \ingroup module_analysisdata
205  */
206 class AnalysisDataTestInput
207 {
208 public:
209     /*! \brief
210      * Constructs empty input data.
211      *
212      * \param[in] dataSetCount Number of data sets in the data.
213      * \param[in] bMultipoint  Whether the data will be multipoint.
214      *
215      * The column count for each data set must be set with
216      * setColumnCount().
217      */
218     AnalysisDataTestInput(int dataSetCount, bool bMultipoint);
219     ~AnalysisDataTestInput();
220
221     //! Whether the input data is multipoint.
222     bool isMultipoint() const { return bMultipoint_; }
223     //! Returns the number of data sets in the input data.
224     int dataSetCount() const { return columnCounts_.size(); }
225     //! Returns the number of columns in a given data set.
226     int columnCount(int dataSet) const { return columnCounts_[dataSet]; }
227     //! Returns the number of frames in the input data.
228     int frameCount() const { return frames_.size(); }
229     //! Returns a frame object for the given input frame.
230     const AnalysisDataTestInputFrame& frame(int index) const;
231
232     //! Sets the number of columns in a data set.
233     void setColumnCount(int dataSet, int columnCount);
234     //! Appends an empty frame to this data.
235     AnalysisDataTestInputFrame& addFrame(real x);
236     //! Adds a frame with a single point set and the given values.
237     void addFrameWithValues(real x, real y1);
238     //! Adds a frame with a single point set and the given values.
239     void addFrameWithValues(real x, real y1, real y2);
240     //! Adds a frame with a single point set and the given values.
241     void addFrameWithValues(real x, real y1, real y2, real y3);
242     //! Adds a frame with a single point set and the given values.
243     void addFrameWithValueAndError(real x, real y1, real e1);
244
245 private:
246     std::vector<int>                        columnCounts_;
247     bool                                    bMultipoint_;
248     std::vector<AnalysisDataTestInputFrame> frames_;
249 };
250
251 /*! \libinternal \brief
252  * Test fixture for AbstractAnalysisData testing.
253  *
254  * This test fixture is designed to help writing tests for objects that
255  * derive from the AbstractAnalysisData class.  Typical flow in such tests is
256  * that first the test creates the required data objects, then call static
257  * methods in this class to add mock modules (using
258  * AbstractAnalysisData::addModule()) to the data objects to check that they
259  * produce the correct data, and then invokes methods in the data object to
260  * produce the data to be checked.  Static methods are also provided for
261  * pushing data from an AnalysisDataTestInput object to some generic types
262  * derived from AbstractAnalysisData.
263  *
264  * Methods addStaticCheckerModule(), addStaticColumnCheckerModule() and
265  * addStaticStorageCheckerModule() create and add mock modules that check the
266  * data against a given AnalysisDataTestInput instance.
267  *
268  * Method addReferenceCheckerModule() creates and adds a mock module that
269  * checks the output against reference data produced by a previous test
270  * execution (see TestReferenceData).  Two versions are provided, a static
271  * method to be used with any TestReferenceChecker, and a non-static method
272  * that uses the protected \p data_ member.
273  *
274  * presentAllData() and presentDataFrame() are provided to push data from an
275  * AnalysisDataTestInput into an AnalysisData object.  In typical tests, most
276  * checks are done during these methods, by the added mock modules.
277  * setupArrayData() performs the same function for classes derived from
278  * AbstractAnalysisArrayData.  In that case, the test should separately ensure
279  * that AbstractAnalysisArrayData::valuesReady() gets called.
280  *
281  * \todo
282  * Support for arbitrary AnalysisDataValues (errors and missing values).
283  *
284  * \see AnalysisDataTestInput
285  *
286  * \inlibraryapi
287  * \ingroup module_analysisdata
288  */
289 class AnalysisDataTestFixture : public ::testing::Test
290 {
291 public:
292     AnalysisDataTestFixture();
293
294     /*! \brief
295      * Initializes an AnalysisData object from input data.
296      *
297      * Sets the column count and other properties based on the input data.
298      */
299     static void setupDataObject(const AnalysisDataTestInput& input, AnalysisData* data);
300
301     /*! \brief
302      * Adds all data from AnalysisDataTestInput into an AnalysisData.
303      */
304     static void presentAllData(const AnalysisDataTestInput& input, AnalysisData* data);
305     /*! \brief
306      * Adds a single frame from AnalysisDataTestInput into an AnalysisData.
307      */
308     static void presentDataFrame(const AnalysisDataTestInput& input, int row, AnalysisDataHandle handle);
309     /*! \brief
310      * Initializes an array data object from AnalysisDataTestInput.
311      *
312      * \tparam ArrayData  Class derived from AbstractAnalysisArrayData.
313      *
314      * The ArrayData class should expose the setter methods
315      * (setColumnCount(), setRowCount(), allocateValues(), setValue())
316      * publicly or declare the fixture class as a friend.
317      * The X axis in \p data must be configured to match \p input before
318      * calling this method.
319      *
320      * Does not call AbstractAnalysisArrayData::valuesReady().
321      * The test must ensure that this method gets called, otherwise the
322      * mock modules never get called.
323      */
324     template<class ArrayData>
325     static void setupArrayData(const AnalysisDataTestInput& input, ArrayData* data);
326
327     /*! \brief
328      * Adds a mock module that verifies output against
329      * AnalysisDataTestInput.
330      *
331      * \param[in]  data     Data to compare against.
332      * \param      source   Data object to verify.
333      *
334      * Creates a mock module that verifies that the
335      * IAnalysisDataModule methods are called correctly by
336      * \p source.  Parameters for the calls are verified against \p data.
337      * Adds the created module to \p source using \p data->addModule().
338      * Any exceptions from the called functions should be caught by the
339      * caller.
340      *
341      * \see AbstractAnalysisData::addModule()
342      */
343     static void addStaticCheckerModule(const AnalysisDataTestInput& data, AbstractAnalysisData* source);
344     /*! \brief
345      * Adds a mock module that verifies parallel output against
346      * AnalysisDataTestInput.
347      *
348      * \param[in]  data     Data to compare against.
349      * \param      source   Data object to verify.
350      *
351      * Creates a parallel mock module that verifies that the
352      * IAnalysisDataModule methods are called correctly by
353      * \p source.  Parameters for the calls are verified against \p data.
354      * Adds the created module to \p source using \p data->addModule().
355      * Any exceptions from the called functions should be caught by the
356      * caller.
357      *
358      * Differs from addStaticCheckerModule() in that the created mock
359      * module reports that it accepts parallel input data, and accepts and
360      * verifies notification calls following the parallel pattern.
361      *
362      * \see AbstractAnalysisData::addModule()
363      */
364     static void addStaticParallelCheckerModule(const AnalysisDataTestInput& data,
365                                                AbstractAnalysisData*        source);
366     /*! \brief
367      * Adds a column mock module that verifies output against
368      * AnalysisDataTestInput.
369      *
370      * \param[in]  data     Data to compare against.
371      * \param[in]  firstcol First column to check.
372      * \param[in]  n        Number of columns to check.
373      * \param      source   Data object to verify.
374      *
375      * Creates a mock module that verifies that the
376      * IAnalysisDataModule methods are called correctly by
377      * \p source.  Parameters for the calls are verified against \p data.
378      * Adds the created module to \p source using
379      * \p data->addColumnModule().
380      * Any exceptions from the called functions should be caught by the
381      * caller.
382      *
383      * \see AbstractAnalysisData::addColumnModule()
384      */
385     static void addStaticColumnCheckerModule(const AnalysisDataTestInput& data,
386                                              int                          firstcol,
387                                              int                          n,
388                                              AbstractAnalysisData*        source);
389     /*! \brief
390      * Adds a mock module that verifies output and storage against
391      * AnalysisDataTestInput.
392      *
393      * \param[in]  data     Data to compare against.
394      * \param[in]  storageCount  Number of previous frames to check
395      *                      (-1 = all).
396      * \param      source   Data object to verify.
397      *
398      * Works like addStaticCheckerModule(), except that in addition, for
399      * each frame, the mock module also checks that previous frames can be
400      * accessed using AbstractAnalysisData::getDataFrame().  In the
401      * IAnalysisDataModule::dataStarted() callback, the mock module
402      * calls AbstractAnalysisData::requestStorage() with \p storageCount as
403      * the parameter.
404      */
405     static void addStaticStorageCheckerModule(const AnalysisDataTestInput& data,
406                                               int                          storageCount,
407                                               AbstractAnalysisData*        source);
408     /*! \brief
409      * Adds a mock module that verifies output against reference data.
410      *
411      * \param[in]  checker   Reference data checker to use for comparison.
412      * \param[in]  id        Identifier for reference data compound to use.
413      * \param      source    Data object to verify.
414      * \param[in]  tolerance Tolerance to use for comparison.
415      *
416      * Creates a mock module that verifies that the
417      * IAnalysisDataModule methods are called correctly by
418      * \p source.  Parameters for the calls are verified against reference
419      * data using a child compound \p id of \p checker.
420      * Adds the created module to \p source using \p data->addModule().
421      * Any exceptions from the called functions should be caught by the
422      * caller.
423      *
424      * \see TestReferenceData
425      */
426     static void addReferenceCheckerModule(const TestReferenceChecker&   checker,
427                                           const char*                   id,
428                                           AbstractAnalysisData*         source,
429                                           const FloatingPointTolerance& tolerance);
430
431     /*! \brief
432      * Adds a mock module that verifies output against reference data.
433      *
434      * \param[in]  id       Identifier for reference data compound to use.
435      * \param      source   Data object to verify.
436      *
437      * Creates a reference checker module using a compound checker with id
438      * \p id at the root level of \p data_.
439      *
440      * See the static overload for other details.
441      */
442     void addReferenceCheckerModule(const char* id, AbstractAnalysisData* source);
443
444 private:
445     /*! \brief
446      * Reference data object used for the reference checker modules.
447      *
448      * Tests can use the data object also for their own purposes if needed.
449      */
450     gmx::test::TestReferenceData data_;
451 };
452
453
454 template<class ArrayData>
455 void AnalysisDataTestFixture::setupArrayData(const AnalysisDataTestInput& input, ArrayData* data)
456 {
457     GMX_RELEASE_ASSERT(!input.isMultipoint(),
458                        "Array data cannot be initialized from multipoint data");
459     GMX_RELEASE_ASSERT(input.dataSetCount() == 1,
460                        "Array data cannot be initialized from multiple data sets");
461     GMX_RELEASE_ASSERT(data->columnCount() == 0 || data->columnCount() == input.columnCount(0),
462                        "Mismatching input and target data");
463     GMX_RELEASE_ASSERT(data->rowCount() == 0 || data->rowCount() == input.frameCount(),
464                        "Mismatching input and target data");
465     data->setColumnCount(input.columnCount(0));
466     data->setRowCount(input.frameCount());
467     data->allocateValues();
468     for (int row = 0; row < input.frameCount(); ++row)
469     {
470         const AnalysisDataTestInputFrame& frame = input.frame(row);
471         EXPECT_FLOAT_EQ(frame.x(), data->xvalue(row));
472         GMX_RELEASE_ASSERT(frame.pointSetCount() == 1,
473                            "Multiple point sets not supported by array data");
474         const AnalysisDataTestInputPointSet& points = frame.pointSet(0);
475         for (int column = 0; column < points.size(); ++column)
476         {
477             data->value(row, column + points.firstColumn()) = points.value(column);
478         }
479     }
480 }
481
482 } // namespace test
483
484 } // namespace gmx
485
486 #endif