SYCL: Avoid using no_init read accessor in rocFFT
[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 by the GROMACS development team.
5  * Copyright (c) 2016,2018,2019,2020, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \libinternal \file
37  * \brief
38  * Helper classes for testing classes that derive from AbstractAnalysisData.
39  *
40  * \author Teemu Murtola <teemu.murtola@gmail.com>
41  * \inlibraryapi
42  * \ingroup module_analysisdata
43  */
44 #ifndef GMX_ANALYSISDATA_TESTS_DATATEST_H
45 #define GMX_ANALYSISDATA_TESTS_DATATEST_H
46
47 #include <vector>
48
49 #include <gtest/gtest.h>
50
51 #include "gromacs/analysisdata/dataframe.h"
52 #include "gromacs/utility/gmxassert.h"
53 #include "gromacs/utility/real.h"
54
55 #include "testutils/refdata.h"
56
57 // currently the bug manifests itself only in AbstractAnalysisData testing
58 #if defined __PATHSCALE__
59 #    define 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     static bool present(int /*i*/) { 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_.emplace_back(y); }
120     //! Appends a value with an error estimate to this point set.
121     void addValueWithError(real y, real error) { values_.emplace_back(y, error); }
122
123 private:
124     //! Creates an empty point set.
125     AnalysisDataTestInputPointSet(int index, int dataSetIndex, int firstColumn);
126
127     struct Value
128     {
129         Value() : y(0.0), error(0.0), bError(false) {}
130         explicit Value(real y) : y(y), error(0.0), bError(false) {}
131         Value(real y, real error) : y(y), error(error), bError(true) {}
132
133         real y;
134         real error;
135         bool bError;
136     };
137
138     int                index_;
139     int                dataSetIndex_;
140     int                firstColumn_;
141     std::vector<Value> values_;
142
143     //! For constructing new point sets.
144     friend class AnalysisDataTestInputFrame;
145 };
146
147 /*! \libinternal \brief
148  * Represents a single frame in AnalysisDataTestInput structure.
149  *
150  * \inlibraryapi
151  * \ingroup module_analysisdata
152  */
153 class AnalysisDataTestInputFrame
154 {
155 public:
156     //! Returns zero-based index for the frame.
157     int index() const { return index_; }
158     //! Returns x coordinate for the frame.
159     real x() const { return x_; }
160     //! Returns error in the x coordinate for the frame.
161     static real dx() { return 0.0; }
162
163     //! Number of individual point sets in the frame.
164     int pointSetCount() const { return pointSets_.size(); }
165     //! Returns a point set object for a given point set.
166     const AnalysisDataTestInputPointSet& pointSet(int index) const
167     {
168         GMX_ASSERT(index >= 0 && static_cast<size_t>(index) < pointSets_.size(),
169                    "Point set index out of range");
170         return pointSets_[index];
171     }
172
173     //! Appends an empty point set to this frame.
174     AnalysisDataTestInputPointSet& addPointSet(int dataSet, int firstColumn);
175     //! Adds a point set with given values to this frame.
176     void addPointSetWithValues(int dataSet, int firstColumn, real y1);
177     //! Adds a point set with given values to this frame.
178     void addPointSetWithValues(int dataSet, int firstColumn, real y1, real y2);
179     //! Adds a point set with given values to this frame.
180     void addPointSetWithValues(int dataSet, int firstColumn, real y1, real y2, real y3);
181     //! Adds a point set with given values to this frame.
182     void addPointSetWithValueAndError(int dataSet, int firstColumn, 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, AnalysisData* data);
301
302     /*! \brief
303      * Adds all data from AnalysisDataTestInput into an AnalysisData.
304      */
305     static void presentAllData(const AnalysisDataTestInput& input, AnalysisData* data);
306     /*! \brief
307      * Adds a single frame from AnalysisDataTestInput into an AnalysisData.
308      */
309     static void presentDataFrame(const AnalysisDataTestInput& input, int row, AnalysisDataHandle handle);
310     /*! \brief
311      * Initializes an array data object from AnalysisDataTestInput.
312      *
313      * \tparam ArrayData  Class derived from AbstractAnalysisArrayData.
314      *
315      * The ArrayData class should expose the setter methods
316      * (setColumnCount(), setRowCount(), allocateValues(), setValue())
317      * publicly or declare the fixture class as a friend.
318      * The X axis in \p data must be configured to match \p input before
319      * calling this method.
320      *
321      * Does not call AbstractAnalysisArrayData::valuesReady().
322      * The test must ensure that this method gets called, otherwise the
323      * mock modules never get called.
324      */
325     template<class ArrayData>
326     static void setupArrayData(const AnalysisDataTestInput& input, ArrayData* data);
327
328     /*! \brief
329      * Adds a mock module that verifies output against
330      * AnalysisDataTestInput.
331      *
332      * \param[in]  data     Data to compare against.
333      * \param      source   Data object to verify.
334      *
335      * Creates a mock module that verifies that the
336      * IAnalysisDataModule methods are called correctly by
337      * \p source.  Parameters for the calls are verified against \p data.
338      * Adds the created module to \p source using \p data->addModule().
339      * Any exceptions from the called functions should be caught by the
340      * caller.
341      *
342      * \see AbstractAnalysisData::addModule()
343      */
344     static void addStaticCheckerModule(const AnalysisDataTestInput& data, AbstractAnalysisData* source);
345     /*! \brief
346      * Adds a mock module that verifies parallel output against
347      * AnalysisDataTestInput.
348      *
349      * \param[in]  data     Data to compare against.
350      * \param      source   Data object to verify.
351      *
352      * Creates a parallel mock module that verifies that the
353      * IAnalysisDataModule methods are called correctly by
354      * \p source.  Parameters for the calls are verified against \p data.
355      * Adds the created module to \p source using \p data->addModule().
356      * Any exceptions from the called functions should be caught by the
357      * caller.
358      *
359      * Differs from addStaticCheckerModule() in that the created mock
360      * module reports that it accepts parallel input data, and accepts and
361      * verifies notification calls following the parallel pattern.
362      *
363      * \see AbstractAnalysisData::addModule()
364      */
365     static void addStaticParallelCheckerModule(const AnalysisDataTestInput& data,
366                                                AbstractAnalysisData*        source);
367     /*! \brief
368      * Adds a column mock module that verifies output against
369      * AnalysisDataTestInput.
370      *
371      * \param[in]  data     Data to compare against.
372      * \param[in]  firstcol First column to check.
373      * \param[in]  n        Number of columns to check.
374      * \param      source   Data object to verify.
375      *
376      * Creates a mock module that verifies that the
377      * IAnalysisDataModule methods are called correctly by
378      * \p source.  Parameters for the calls are verified against \p data.
379      * Adds the created module to \p source using
380      * \p data->addColumnModule().
381      * Any exceptions from the called functions should be caught by the
382      * caller.
383      *
384      * \see AbstractAnalysisData::addColumnModule()
385      */
386     static void addStaticColumnCheckerModule(const AnalysisDataTestInput& data,
387                                              int                          firstcol,
388                                              int                          n,
389                                              AbstractAnalysisData*        source);
390     /*! \brief
391      * Adds a mock module that verifies output and storage against
392      * AnalysisDataTestInput.
393      *
394      * \param[in]  data     Data to compare against.
395      * \param[in]  storageCount  Number of previous frames to check
396      *                      (-1 = all).
397      * \param      source   Data object to verify.
398      *
399      * Works like addStaticCheckerModule(), except that in addition, for
400      * each frame, the mock module also checks that previous frames can be
401      * accessed using AbstractAnalysisData::getDataFrame().  In the
402      * IAnalysisDataModule::dataStarted() callback, the mock module
403      * calls AbstractAnalysisData::requestStorage() with \p storageCount as
404      * the parameter.
405      */
406     static void addStaticStorageCheckerModule(const AnalysisDataTestInput& data,
407                                               int                          storageCount,
408                                               AbstractAnalysisData*        source);
409     /*! \brief
410      * Adds a mock module that verifies output against reference data.
411      *
412      * \param[in]  checker   Reference data checker to use for comparison.
413      * \param[in]  id        Identifier for reference data compound to use.
414      * \param      source    Data object to verify.
415      * \param[in]  tolerance Tolerance to use for comparison.
416      *
417      * Creates a mock module that verifies that the
418      * IAnalysisDataModule methods are called correctly by
419      * \p source.  Parameters for the calls are verified against reference
420      * data using a child compound \p id of \p checker.
421      * Adds the created module to \p source using \p data->addModule().
422      * Any exceptions from the called functions should be caught by the
423      * caller.
424      *
425      * \see TestReferenceData
426      */
427     static void addReferenceCheckerModule(const TestReferenceChecker&   checker,
428                                           const char*                   id,
429                                           AbstractAnalysisData*         source,
430                                           const FloatingPointTolerance& tolerance);
431
432     /*! \brief
433      * Adds a mock module that verifies output against reference data.
434      *
435      * \param[in]  id       Identifier for reference data compound to use.
436      * \param      source   Data object to verify.
437      *
438      * Creates a reference checker module using a compound checker with id
439      * \p id at the root level of \p data_.
440      *
441      * See the static overload for other details.
442      */
443     void addReferenceCheckerModule(const char* id, AbstractAnalysisData* source);
444
445 private:
446     /*! \brief
447      * Reference data object used for the reference checker modules.
448      *
449      * Tests can use the data object also for their own purposes if needed.
450      */
451     gmx::test::TestReferenceData data_;
452 };
453
454
455 template<class ArrayData>
456 void AnalysisDataTestFixture::setupArrayData(const AnalysisDataTestInput& input, ArrayData* data)
457 {
458     GMX_RELEASE_ASSERT(!input.isMultipoint(),
459                        "Array data cannot be initialized from multipoint data");
460     GMX_RELEASE_ASSERT(input.dataSetCount() == 1,
461                        "Array data cannot be initialized from multiple data sets");
462     GMX_RELEASE_ASSERT(data->columnCount() == 0 || data->columnCount() == input.columnCount(0),
463                        "Mismatching input and target data");
464     GMX_RELEASE_ASSERT(data->rowCount() == 0 || data->rowCount() == input.frameCount(),
465                        "Mismatching input and target data");
466     data->setColumnCount(input.columnCount(0));
467     data->setRowCount(input.frameCount());
468     data->allocateValues();
469     for (int row = 0; row < input.frameCount(); ++row)
470     {
471         const AnalysisDataTestInputFrame& frame = input.frame(row);
472         EXPECT_FLOAT_EQ(frame.x(), data->xvalue(row));
473         GMX_RELEASE_ASSERT(frame.pointSetCount() == 1,
474                            "Multiple point sets not supported by array data");
475         const AnalysisDataTestInputPointSet& points = frame.pointSet(0);
476         for (int column = 0; column < points.size(); ++column)
477         {
478             data->value(row, column + points.firstColumn()) = points.value(column);
479         }
480     }
481 }
482
483 } // namespace test
484
485 } // namespace gmx
486
487 #endif