Merge "Remove obsolete pointers in options handling."
[alexxy/gromacs.git] / src / gromacs / analysisdata / tests / datatest.h
1 /*
2  *
3  *                This source code is part of
4  *
5  *                 G   R   O   M   A   C   S
6  *
7  *          GROningen MAchine for Chemical Simulations
8  *
9  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
10  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
11  * Copyright (c) 2001-2009, The GROMACS development team,
12  * check out http://www.gromacs.org for more information.
13
14  * This program is free software; you can redistribute it and/or
15  * modify it under the terms of the GNU General Public License
16  * as published by the Free Software Foundation; either version 2
17  * of the License, or (at your option) any later version.
18  *
19  * If you want to redistribute modifications, please consider that
20  * scientific software is very special. Version control is crucial -
21  * bugs must be traceable. We will be happy to consider code for
22  * inclusion in the official distribution, but derived work must not
23  * be called official GROMACS. Details are found in the README & COPYING
24  * files - if they are missing, get the official version at www.gromacs.org.
25  *
26  * To help us fund GROMACS development, we humbly ask that you cite
27  * the papers on the package - you can find them in the top README file.
28  *
29  * For more info, check our website at http://www.gromacs.org
30  */
31 /*! \libinternal \file
32  * \brief
33  * Helper classes for testing classes that derive from AbstractAnalysisData.
34  *
35  * \author Teemu Murtola <teemu.murtola@cbr.su.se>
36  * \inlibraryapi
37  * \ingroup module_analysisdata
38  */
39 #ifndef GMX_ANALYSISDATA_TESTS_DATATEST_H
40 #define GMX_ANALYSISDATA_TESTS_DATATEST_H
41
42 #include <limits>
43 #include <vector>
44
45 #include <gtest/gtest.h>
46
47 #include "gromacs/legacyheaders/types/simple.h"
48 #include "gromacs/fatalerror/gmxassert.h"
49
50 #include "testutils/refdata.h"
51
52 namespace gmx
53 {
54
55 class AbstractAnalysisData;
56 class AnalysisData;
57 class AnalysisDataHandle;
58
59 namespace test
60 {
61
62 class MockAnalysisModule;
63
64 //! Constant to use to signify end of data for AnalysisDataTestInput.
65 const real END_OF_DATA = std::numeric_limits<real>::max();
66 //! Constant to use to signify end of one data frame for AnalysisDataTestInput.
67 const real END_OF_FRAME = std::numeric_limits<real>::min();
68 //! Constant to use to signify end of one multipoint set for AnalysisDataTestInput.
69 const real MPSTOP = -std::numeric_limits<real>::max();
70
71 /*! \libinternal \brief
72  * Represents a single set of points in AnalysisDataTestInputFrame structure.
73  *
74  * If the data is not multipoint, each frame contains exactly one set of
75  * points.  If there is more than one set of points, each of these sets is
76  * passed separately and AnalysisDataHandle::finishPointSet() called in
77  * between.
78  *
79  * \ingroup module_analysisdata
80  */
81 class AnalysisDataTestInputPointSet
82 {
83     public:
84         //! Returns the number of columns in the point set.
85         int size() const { return y_.size(); }
86         //! Returns the value in column \p i.
87         real y(int i) const { return y_[i]; }
88         //! Returns the error in column \p i.
89         real dy(int i) const { return 0.0; }
90         //! Returns whether the value in column \p i is present.
91         real present(int i) const { return true; }
92         //! Returns a vector of values for all columns.
93         const std::vector<real> &yvector() const { return y_; }
94
95     private:
96         //! Creates an empty point set.
97         AnalysisDataTestInputPointSet();
98
99         std::vector<real>       y_;
100
101         friend class AnalysisDataTestInput;
102 };
103
104 /*! \libinternal \brief
105  * Represents a single frame in AnalysisDataTestInput structure.
106  *
107  * \ingroup module_analysisdata
108  */
109 class AnalysisDataTestInputFrame
110 {
111     public:
112         //! Returns zero-based index for the frame.
113         int index() const { return index_; }
114         //! Returns x coordinate for the frame.
115         real x() const { return x_; }
116         //! Returns error in the x coordinate for the frame.
117         real dx() const { return 0.0; }
118
119         //! Number of individual point sets in the frame.
120         int pointSetCount() const { return points_.size(); }
121         //! Returns a point set object for a given point set.
122         const AnalysisDataTestInputPointSet &points(int index = 0) const
123         {
124             GMX_ASSERT(index >= 0 && static_cast<size_t>(index) < points_.size(),
125                        "Point set index out of range");
126             return points_[index];
127         }
128
129     private:
130         //! Constructs a new frame object with the given values.
131         AnalysisDataTestInputFrame(int index, real x);
132
133         int                     index_;
134         real                    x_;
135         std::vector<AnalysisDataTestInputPointSet>  points_;
136
137         friend class AnalysisDataTestInput;
138 };
139
140 /*! \libinternal \brief
141  * Represents static input data for AbstractAnalysisData tests.
142  *
143  * Used to construct structured test input data from a static array of reals,
144  * and then typically used as input to methods in AnalysisDataTestFixture.
145  *
146  * \see AnalysisDataTestFixture
147  *
148  * \ingroup module_analysisdata
149  */
150 class AnalysisDataTestInput
151 {
152     public:
153         /*! \brief
154          * Constructs data representation from a simple array.
155          *
156          * \param[in] data  Array to construct data from.
157          *
158          * The input array should consist of a set of frames, separated by a
159          * END_OF_FRAME marker.  The first value for a frame is the X value,
160          * all following values are Y values.
161          * For multipoint data, one frame can contain several point sets,
162          * separated by MPSTOP markers.  There should be no MPSTOP marker after
163          * the last point set, only an END_OF_FRAME marker.  All point sets are
164          * assumed to start from column zero, but the sets may contain
165          * different number of columns.  For non-multipoint data, all frames
166          * must containt the same number of columns.
167          * The final element in the array (after the last END_OF_FRAME) should
168          * be END_OF_DATA.
169          */
170         explicit AnalysisDataTestInput(const real *data);
171         ~AnalysisDataTestInput();
172
173         //! Returns the number of frames in the input data.
174         int frameCount() const { return frames_.size(); }
175         //! Returns the number of columns in the input data.
176         int columnCount() const { return columnCount_; }
177         //! Whether the input data is multipoint.
178         bool isMultipoint() const { return bMultipoint_; }
179         //! Returns a frame object for the given input frame.
180         const AnalysisDataTestInputFrame &frame(int index) const;
181
182     private:
183         int                     columnCount_;
184         bool                    bMultipoint_;
185         std::vector<AnalysisDataTestInputFrame> frames_;
186 };
187
188 /*! \libinternal \brief
189  * Test fixture for AbstractAnalysisData testing.
190  *
191  * This test fixture is designed to help writing tests for objects that
192  * derive from the AbstractAnalysisData class.  Typical flow in such tests is
193  * that first the test creates the required data objects, then call static
194  * methods in this class to add mock modules (using
195  * AbstractAnalysisData::addModule()) to the data objects to check that they
196  * produce the correct data, and then invokes methods in the data object to
197  * produce the data to be checked.  Static methods are also provided for
198  * pushing data from an AnalysisDataTestInput object to some generic types
199  * derived from AbstractAnalysisData.
200  *
201  * Methods addStaticCheckerModule(), addStaticColumnCheckerModule() and
202  * addStaticStorageCheckerModule() create and add mock modules that check the
203  * data against a given AnalysisDataTestInput instance.
204  *
205  * Method addReferenceCheckerModule() creates and adds a mock module that
206  * checks the output against reference data produced by a previous test
207  * execution (see TestReferenceData).  Two versions are provided, a static
208  * method to be used with any TestReferenceChecker, and a non-static method
209  * that uses the protected \p data_ member.
210  *
211  * presentAllData() and presentDataFrame() are provided to push data from an
212  * AnalysisDataTestInput into an AnalysisData object.  In typical tests, most
213  * checks are done during the these methods, by the added mock modules.
214  * setupArrayData() performs the same function for classes derived from
215  * AbstractAnalysisArrayData.  In that case, the test should separately ensure
216  * that AbstractAnalysisArrayData::valuesReady() gets called.
217  *
218  * \todo
219  * Support for errors and for arbitrary multipoint data.
220  *
221  * \see AnalysisDataTestInput
222  *
223  * \ingroup module_analysisdata
224  */
225 class AnalysisDataTestFixture : public ::testing::Test
226 {
227     public:
228         AnalysisDataTestFixture();
229
230         /*! \brief
231          * Adds all data from AnalysisDataTestInput into an AnalysisData.
232          */
233         static void presentAllData(const AnalysisDataTestInput &input,
234                                    AnalysisData *data);
235         /*! \brief
236          * Adds a single frame from AnalysisDataTestInput into an AnalysisData.
237          */
238         static void presentDataFrame(const AnalysisDataTestInput &input, int row,
239                                      AnalysisDataHandle handle);
240         /*! \brief
241          * Initializes an array data object from AnalysisDataTestInput.
242          *
243          * \tparam ArrayData  Class derived from AbstractAnalysisArrayData.
244          *
245          * The ArrayData class should expose the setter methods
246          * (setColumnCount(), setRowCount(), allocateValues(), setValue())
247          * publicly or declare the fixture class as a friend.
248          * The X axis in \p data must be configured to match \p input before
249          * calling this method.
250          *
251          * Does not call AbstractAnalysisArrayData::valuesReady().
252          * The test must ensure that this method gets called, otherwise the
253          * mock modules never get called.
254          */
255         template <class ArrayData>
256         static void setupArrayData(const AnalysisDataTestInput &input,
257                                    ArrayData *data);
258
259         /*! \brief
260          * Adds a mock module that verifies output against
261          * AnalysisDataTestInput.
262          *
263          * \param[in]  data     Data to compare against.
264          * \param      source   Data object to verify.
265          *
266          * Creates a mock module that verifies that the
267          * AnalysisDataModuleInterface methods are called correctly by
268          * \p source.  Parameters for the calls are verified against \p data.
269          * Adds the created module to \p source using \p data->addModule().
270          * Any exceptions from the called functions should be caught by the
271          * caller.
272          *
273          * \see AbstractAnalysisData::addModule()
274          */
275         static void addStaticCheckerModule(const AnalysisDataTestInput &data,
276                                            AbstractAnalysisData *source);
277         /*! \brief
278          * Adds a column mock module that verifies output against
279          * AnalysisDataTestInput.
280          *
281          * \param[in]  data     Data to compare against.
282          * \param[in]  firstcol First column to check.
283          * \param[in]  n        Number of columns to check.
284          * \param      source   Data object to verify.
285          *
286          * Creates a mock module that verifies that the
287          * AnalysisDataModuleInterface methods are called correctly by
288          * \p source.  Parameters for the calls are verified against \p data.
289          * Adds the created module to \p source using
290          * \p data->addColumnModule().
291          * Any exceptions from the called functions should be caught by the
292          * caller.
293          *
294          * \see AbstractAnalysisData::addColumnModule()
295          */
296         static void addStaticColumnCheckerModule(const AnalysisDataTestInput &data,
297                                                  int firstcol, int n,
298                                                  AbstractAnalysisData *source);
299         /*! \brief
300          * Adds a mock module that verifies output and storage against
301          * AnalysisDataTestInput.
302          *
303          * \param[in]  data     Data to compare against.
304          * \param[in]  storageCount  Number of previous frames to check
305          *                      (-1 = all).
306          * \param      source   Data object to verify.
307          *
308          * Works like addStaticCheckerModule(), except that in addition, for
309          * each frame, the mock module also checks that previous frames can be
310          * accessed using AbstractAnalysisData::getDataFrame().  In the
311          * AnalysisDataModuleInterface::dataStarted() callback, the mock module
312          * calls AbstractAnalysisData::requestStorage() with \p storageCount as
313          * the parameter.
314          */
315         static void addStaticStorageCheckerModule(const AnalysisDataTestInput &data,
316                                                   int storageCount,
317                                                   AbstractAnalysisData *source);
318         /*! \brief
319          * Adds a mock module that verifies output against reference data.
320          *
321          * \param[in]  checker  Reference data checker to use for comparison.
322          * \param      source   Data object to verify.
323          *
324          * Creates a mock module that verifies that the
325          * AnalysisDataModuleInterface methods are called correctly by
326          * \p source.  Parameters for the calls are verified against reference
327          * data using \p checker.
328          * Adds the created module to \p source using \p data->addModule().
329          * Any exceptions from the called functions should be caught by the
330          * caller.
331          *
332          * \see TestReferenceData
333          */
334         static void addReferenceCheckerModule(const TestReferenceChecker &checker,
335                                               AbstractAnalysisData *source);
336
337         /*! \brief
338          * Adds a mock module that verifies output against reference data.
339          *
340          * \param[in]  id       Identifier for reference data compound to use.
341          * \param      source   Data object to verify.
342          *
343          * Creates a reference checker module using a compound checker with id
344          * \p id at the root level of \p data_.
345          *
346          * See the static overload for other details.
347          */
348         void addReferenceCheckerModule(const char *id,
349                                        AbstractAnalysisData *source);
350
351     protected:
352         /*! \brief
353          * Reference data object used for the reference checker modules.
354          *
355          * Tests can use the data object also for their own purposes if needed.
356          */
357         gmx::test::TestReferenceData  data_;
358 };
359
360
361 template <class ArrayData>
362 void AnalysisDataTestFixture::setupArrayData(const AnalysisDataTestInput &input,
363                                              ArrayData *data)
364 {
365     GMX_RELEASE_ASSERT(!input.isMultipoint(),
366                        "Array data cannot be initialized from multipoint data");
367     GMX_RELEASE_ASSERT(data->columnCount() == 0 || data->columnCount() == input.columnCount(),
368                        "Mismatching input and target data");
369     GMX_RELEASE_ASSERT(data->rowCount() == 0 || data->rowCount() == input.frameCount(),
370                        "Mismatching input and target data");
371     data->setColumnCount(input.columnCount());
372     data->setRowCount(input.frameCount());
373     data->allocateValues();
374     for (int row = 0; row < input.frameCount(); ++row)
375     {
376         const AnalysisDataTestInputFrame &frame = input.frame(row);
377         EXPECT_FLOAT_EQ(frame.x(), data->xvalue(row));
378         const AnalysisDataTestInputPointSet &points = frame.points();
379         for (int column = 0; column < input.columnCount(); ++column)
380         {
381             data->setValue(row, column, points.y(column));
382         }
383     }
384 }
385
386 } // namespace test
387
388 } // namespace gmx
389
390 #endif