Tidy: modernize-use-nullptr
[alexxy/gromacs.git] / src / gromacs / analysisdata / tests / mock_datamodule.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2011,2012,2013,2014,2015,2017, 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 /*! \internal \file
36  * \brief
37  * Implements classes in mock_datamodule.h.
38  *
39  * \author Teemu Murtola <teemu.murtola@gmail.com>
40  * \ingroup module_analysisdata
41  */
42 #include "gmxpre.h"
43
44 #include "mock_datamodule.h"
45
46 #include <gmock/gmock.h>
47 #include <gtest/gtest.h>
48
49 #include "gromacs/analysisdata/analysisdata.h"
50 #include "gromacs/analysisdata/dataframe.h"
51 #include "gromacs/utility/gmxassert.h"
52 #include "gromacs/utility/stringutil.h"
53
54 #include "gromacs/analysisdata/tests/datatest.h"
55 #include "testutils/refdata.h"
56 #include "testutils/testasserts.h"
57
58 namespace gmx
59 {
60 namespace test
61 {
62
63 /********************************************************************
64  * MockAnalysisDataModule::Impl
65  */
66
67 /*! \internal \brief
68  * Private implementation class for gmx::test::MockAnalysisDataModule.
69  *
70  * \ingroup module_analysisdata
71  */
72 class MockAnalysisDataModule::Impl
73 {
74     public:
75         //! Initializes a mock object with the given flags.
76         explicit Impl(int flags);
77
78         /*! \brief
79          * Callback used to initialize reference data checks
80          *
81          * Called in response to dataStarted().
82          * Records the source data for later use (for access to data properties).
83          */
84         void startReferenceData(AbstractAnalysisData *data);
85         /*! \brief
86          * Callback used to check frame start against reference data.
87          *
88          * Called to check parameters and order of calls to frameStarted().
89          * In addition to reference data checks, this method checks statically
90          * that the new frame matches \a frameIndex_.
91          */
92         void startReferenceFrame(const AnalysisDataFrameHeader &header);
93         /*! \brief
94          * Callback used to check frame points against reference data.
95          *
96          * Called to check parameters and order of calls to pointsAdded().
97          */
98         void checkReferencePoints(const AnalysisDataPointSetRef &points);
99         /*! \brief
100          * Callback used to check frame finish against reference data.
101          *
102          * Called to check parameters and order of calls to frameFinished().
103          */
104         void finishReferenceFrame(const AnalysisDataFrameHeader &header);
105         /*! \brief
106          * Callback used to check serial frame finish with reference data.
107          *
108          * Called to check parameters and order of calls to
109          * frameFinishedSerial().
110          * \a frameIndex_ is incremented here.
111          */
112         void finishReferenceFrameSerial(int frameIndex);
113
114         /*! \brief
115          * Reference data checker to use for checking frames.
116          *
117          * Must be initialized if startReferenceFrame() is called.
118          */
119         TestReferenceChecker         rootChecker_;
120         /*! \brief
121          * Reference data checker to use to check the current frame.
122          *
123          * Initialized between startReferenceFrame() and finishReferenceFrame()
124          * calls.
125          */
126         TestReferenceChecker         frameChecker_;
127         //! Source data.
128         const AbstractAnalysisData  *source_;
129         //! Flags that will be returned by the mock module.
130         int                          flags_;
131         //! Index of the current/next frame.
132         int                          frameIndex_;
133 };
134
135 namespace
136 {
137
138 /*! \internal \brief
139  * Checks a single AnalysisDataValue.
140  *
141  * \ingroup module_analysisdata
142  */
143 void checkReferenceDataPoint(TestReferenceChecker    *checker,
144                              const AnalysisDataValue &value)
145 {
146     TestReferenceChecker compound(checker->checkCompound("DataValue", nullptr));
147     compound.checkReal(value.value(), "Value");
148     if (compound.checkPresent(value.hasError(), "Error"))
149     {
150         compound.checkReal(value.error(), "Error");
151     }
152     if (compound.checkPresent(!value.isPresent(), "Present"))
153     {
154         compound.checkBoolean(value.isPresent(), "Present");
155     }
156 }
157
158 }       // namespace
159
160 MockAnalysisDataModule::Impl::Impl(int flags)
161     : source_(nullptr), flags_(flags), frameIndex_(0)
162 {
163 }
164
165
166 void
167 MockAnalysisDataModule::Impl::startReferenceData(AbstractAnalysisData *data)
168 {
169     source_ = data;
170 }
171
172
173 void
174 MockAnalysisDataModule::Impl::startReferenceFrame(
175         const AnalysisDataFrameHeader &header)
176 {
177     GMX_RELEASE_ASSERT(rootChecker_,
178                        "Root checker not set, but reference data used");
179     EXPECT_FALSE(frameChecker_.isValid());
180     EXPECT_EQ(frameIndex_, header.index());
181     frameChecker_ = rootChecker_.checkCompound("DataFrame",
182                                                formatString("Frame%d", frameIndex_).c_str());
183     frameChecker_.checkReal(header.x(), "X");
184 }
185
186
187 void
188 MockAnalysisDataModule::Impl::checkReferencePoints(
189         const AnalysisDataPointSetRef &points)
190 {
191     EXPECT_TRUE(frameChecker_.isValid());
192     if (frameChecker_)
193     {
194         TestReferenceChecker checker(
195                 frameChecker_.checkCompound("DataValues", nullptr));
196         checker.checkInteger(points.columnCount(), "Count");
197         if (checker.checkPresent(source_->dataSetCount() > 1, "DataSet"))
198         {
199             checker.checkInteger(points.dataSetIndex(), "DataSet");
200         }
201         const int  sourceColumnCount = source_->columnCount(points.dataSetIndex());
202         const bool bAllColumns       = (points.firstColumn() == 0
203                                         && points.columnCount() == sourceColumnCount);
204         if (checker.checkPresent(!bAllColumns, "FirstColumn"))
205         {
206             checker.checkInteger(points.firstColumn(), "FirstColumn");
207             checker.checkInteger(points.lastColumn(),  "LastColumn");
208         }
209
210         AnalysisDataValuesRef::const_iterator value;
211         for (value = points.values().begin(); value != points.values().end(); ++value)
212         {
213             checkReferenceDataPoint(&checker, *value);
214         }
215     }
216 }
217
218
219 void
220 MockAnalysisDataModule::Impl::finishReferenceFrame(
221         const AnalysisDataFrameHeader &header)
222 {
223     EXPECT_TRUE(frameChecker_.isValid());
224     EXPECT_EQ(frameIndex_, header.index());
225     frameChecker_ = TestReferenceChecker();
226 }
227
228
229 void
230 MockAnalysisDataModule::Impl::finishReferenceFrameSerial(int frameIndex)
231 {
232     EXPECT_FALSE(frameChecker_.isValid());
233     EXPECT_EQ(frameIndex_, frameIndex);
234     ++frameIndex_;
235 }
236
237
238 /********************************************************************
239  * MockAnalysisDataModule
240  */
241
242 namespace
243 {
244
245 /*! \brief
246  * Helper function for checking the data frame header against static data.
247  *
248  * \param[in] header    Frame header to check.
249  * \param[in] refFrame  Data to check against.
250  */
251 void checkHeader(const AnalysisDataFrameHeader    &header,
252                  const AnalysisDataTestInputFrame &refFrame)
253 {
254     EXPECT_EQ(refFrame.index(), header.index());
255     EXPECT_FLOAT_EQ(refFrame.x(), header.x());
256     EXPECT_FLOAT_EQ(refFrame.dx(), header.dx());
257 }
258
259 /*! \brief
260  * Helper function for checking a point set against static data.
261  *
262  * \param[in] points       Point set to check.
263  * \param[in] refPoints    Data to check against.
264  * \param[in] columnOffset Offset of first column of \p points in \p refPoints.
265  */
266 void checkPoints(const AnalysisDataPointSetRef       &points,
267                  const AnalysisDataTestInputPointSet &refPoints,
268                  int                                  columnOffset)
269 {
270     for (int i = 0; i < points.columnCount(); ++i)
271     {
272         const int column = points.firstColumn() - refPoints.firstColumn() + i + columnOffset;
273         EXPECT_FLOAT_EQ(refPoints.y(column),
274                         points.y(i))
275         << "  Column: " << i+1 << " / " << points.columnCount()
276         << " (+" << points.firstColumn() << ")\n"
277         << "Ref. col: " << column+1 << " / " << refPoints.size()
278         << " (+" << refPoints.firstColumn() << ", offs " << columnOffset << ")";
279     }
280 }
281
282 /*! \brief
283  * Helper function for checking a full frame against static data.
284  *
285  * \param[in] frame     Frame to check.
286  * \param[in] refFrame  Data to check against.
287  */
288 void checkFrame(const AnalysisDataFrameRef       &frame,
289                 const AnalysisDataTestInputFrame &refFrame)
290 {
291     checkHeader(frame.header(), refFrame);
292     ASSERT_EQ(refFrame.pointSetCount(), frame.pointSetCount());
293     for (int i = 0; i < frame.pointSetCount(); ++i)
294     {
295         const AnalysisDataPointSetRef       &points    = frame.pointSet(i);
296         const AnalysisDataTestInputPointSet &refPoints = refFrame.pointSet(i);
297         EXPECT_EQ(refPoints.firstColumn(), points.firstColumn());
298         checkPoints(points, refPoints, 0);
299     }
300 }
301
302 /*! \brief
303  * Functor for checking data frame header against static test input data.
304  *
305  * This functor is designed to be invoked as a handled for
306  * IAnalysisDataModule::frameStarted().
307  */
308 class StaticDataFrameHeaderChecker
309 {
310     public:
311         /*! \brief
312          * Constructs a checker against a given input data frame.
313          *
314          * \param[in] frame Frame to check against.
315          *
316          * \p frame must exist for the lifetime of this object.
317          */
318         StaticDataFrameHeaderChecker(const AnalysisDataTestInputFrame *frame)
319             : frame_(frame)
320         {
321         }
322
323         //! Function call operator for the functor.
324         void operator()(const AnalysisDataFrameHeader &header) const
325         {
326             SCOPED_TRACE(formatString("Frame %d", frame_->index()));
327             checkHeader(header, *frame_);
328         }
329
330     private:
331         const AnalysisDataTestInputFrame *frame_;
332 };
333
334 /*! \brief
335  * Functor for checking data frame points against static test input data.
336  *
337  * This functor is designed to be invoked as a handled for
338  * IAnalysisDataModule::pointsAdded().
339  */
340 class StaticDataPointsChecker
341 {
342     public:
343         /*! \brief
344          * Constructs a checker against a given input data frame and point set.
345          *
346          * \param[in] frame    Frame to check against.
347          * \param[in] points   Point set in \p frame to check against.
348          * \param[in] firstcol Expected first column.
349          * \param[in] n        Expected number of columns.
350          *
351          * \p firstcol and \p n are used to create a checker that only expects
352          * to be called for a subset of columns.
353          * \p frame and \p points must exist for the lifetime of this object.
354          */
355         StaticDataPointsChecker(const AnalysisDataTestInputFrame *frame,
356                                 const AnalysisDataTestInputPointSet *points,
357                                 int firstcol, int n)
358             : frame_(frame), points_(points), firstcol_(firstcol), n_(n)
359         {
360         }
361
362         //! Function call operator for the functor.
363         void operator()(const AnalysisDataPointSetRef &points) const
364         {
365             SCOPED_TRACE(formatString("Frame %d, point set %d",
366                                       frame_->index(), points_->index()));
367             EXPECT_EQ(points_->dataSetIndex(), points.dataSetIndex());
368             const int expectedFirstColumn
369                 = std::max(0, points_->firstColumn() - firstcol_);
370             const int expectedLastColumn
371                 = std::min(n_ - 1, points_->lastColumn() - firstcol_);
372             EXPECT_EQ(expectedFirstColumn, points.firstColumn());
373             EXPECT_EQ(expectedLastColumn,  points.lastColumn());
374             checkHeader(points.header(), *frame_);
375             checkPoints(points, *points_, firstcol_);
376         }
377
378     private:
379         const AnalysisDataTestInputFrame    *frame_;
380         const AnalysisDataTestInputPointSet *points_;
381         int                                  firstcol_;
382         int                                  n_;
383 };
384
385 /*! \brief
386  * Functor for requesting data storage.
387  *
388  * This functor is designed to be invoked as a handled for
389  * IAnalysisDataModule::dataStarted().
390  */
391 class DataStorageRequester
392 {
393     public:
394         /*! \brief
395          * Constructs a functor that requests the given amount of storage.
396          *
397          * \param[in] count  Number of frames of storage to request, or
398          *      -1 for all frames.
399          *
400          * \see AbstractAnalysisData::requestStorage()
401          */
402         explicit DataStorageRequester(int count) : count_(count) {}
403
404         //! Function call operator for the functor.
405         void operator()(AbstractAnalysisData *data) const
406         {
407             EXPECT_TRUE(data->requestStorage(count_));
408         }
409
410     private:
411         int                     count_;
412 };
413
414 /*! \brief
415  * Functor for checking data frame points and storage against static test input
416  * data.
417  *
418  * This functor is designed to be invoked as a handled for
419  * IAnalysisDataModule::pointsAdded().
420  */
421 class StaticDataPointsStorageChecker
422 {
423     public:
424         /*! \brief
425          * Constructs a checker for a given frame.
426          *
427          * \param[in] source     Data object that is being checked.
428          * \param[in] data       Test input data to check against.
429          * \param[in] frameIndex Frame index for which this functor expects
430          *      to be called.
431          * \param[in] pointSetIndex Point set for which this functor expects
432          *      to be called.
433          * \param[in] storageCount How many past frames should be checked for
434          *      storage (-1 = check all frames).
435          *
436          * This checker works as StaticDataPointsChecker, but additionally
437          * checks that previous frames can be accessed using access methods
438          * in AbstractAnalysisData and that correct data is returned.
439          *
440          * \p source and \p data must exist for the lifetime of this object.
441          */
442         StaticDataPointsStorageChecker(AbstractAnalysisData        *source,
443                                        const AnalysisDataTestInput *data,
444                                        int frameIndex, int pointSetIndex,
445                                        int storageCount)
446             : source_(source), data_(data),
447               frameIndex_(frameIndex), pointSetIndex_(pointSetIndex),
448               storageCount_(storageCount)
449         {
450         }
451
452         //! Function call operator for the functor.
453         void operator()(const AnalysisDataPointSetRef &points) const
454         {
455             SCOPED_TRACE(formatString("Frame %d", frameIndex_));
456             const AnalysisDataTestInputFrame    &refFrame  = data_->frame(frameIndex_);
457             const AnalysisDataTestInputPointSet &refPoints = refFrame.pointSet(pointSetIndex_);
458             EXPECT_EQ(refPoints.firstColumn(), points.firstColumn());
459             EXPECT_EQ(refPoints.size(), points.columnCount());
460             checkHeader(points.header(), refFrame);
461             checkPoints(points, refPoints, 0);
462             for (int past = 1;
463                  (storageCount_ < 0 || past <= storageCount_) && past <= frameIndex_;
464                  ++past)
465             {
466                 int   index = frameIndex_ - past;
467                 SCOPED_TRACE(formatString("Checking storage of frame %d", index));
468                 ASSERT_NO_THROW_GMX({
469                                         AnalysisDataFrameRef frame = source_->getDataFrame(index);
470                                         ASSERT_TRUE(frame.isValid());
471                                         checkFrame(frame, data_->frame(index));
472                                     });
473             }
474         }
475
476     private:
477         AbstractAnalysisData        *source_;
478         const AnalysisDataTestInput *data_;
479         int                          frameIndex_;
480         int                          pointSetIndex_;
481         int                          storageCount_;
482 };
483
484 /*! \brief
485  * Sets the mock object expectation to mimick AnalysisDataModuleSerial.
486  *
487  * Makes MockAnalysisDataModule::parallelDataStarted() behave as if the mock
488  * object was an AnalysisDataModuleSerial object: forward the call to
489  * MockAnalysisDataModule::dataStarted() and return false.
490  */
491 void setSerialExpectationForParallelDataStarted(MockAnalysisDataModule *mock)
492 {
493     using ::testing::_;
494     using ::testing::AtMost;
495     using ::testing::DoAll;
496     using ::testing::Invoke;
497     using ::testing::Return;
498     using ::testing::WithArg;
499     EXPECT_CALL(*mock, parallelDataStarted(_, _))
500         .Times(AtMost(1))
501         .WillOnce(DoAll(WithArg<0>(Invoke(mock, &MockAnalysisDataModule::dataStarted)),
502                         Return(false)));
503 }
504
505 }       // anonymous namespace
506
507
508 MockAnalysisDataModule::MockAnalysisDataModule(int flags)
509     : impl_(new Impl(flags))
510 {
511 }
512
513
514 MockAnalysisDataModule::~MockAnalysisDataModule()
515 {
516 }
517
518
519 int MockAnalysisDataModule::flags() const
520 {
521     return impl_->flags_;
522 }
523
524
525 void
526 MockAnalysisDataModule::setupStaticCheck(const AnalysisDataTestInput &data,
527                                          AbstractAnalysisData        *source,
528                                          bool                         bParallel)
529 {
530     impl_->flags_ |= efAllowMulticolumn | efAllowMultipoint | efAllowMultipleDataSets;
531
532     using ::testing::_;
533     using ::testing::Invoke;
534     using ::testing::Property;
535     using ::testing::Return;
536
537     if (bParallel)
538     {
539         ::testing::Expectation init =
540             EXPECT_CALL(*this, parallelDataStarted(source, _))
541                 .WillOnce(Return(true));
542         ::testing::ExpectationSet framesFinished;
543         ::testing::Expectation    prevFinish;
544         for (int row = 0; row < data.frameCount(); ++row)
545         {
546             ::testing::InSequence frameSequence;
547             const AnalysisDataTestInputFrame &frame = data.frame(row);
548             EXPECT_CALL(*this, frameStarted(Property(&AnalysisDataFrameHeader::index, row)))
549                 .After(init)
550                 .WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
551             for (int ps = 0; ps < frame.pointSetCount(); ++ps)
552             {
553                 const AnalysisDataTestInputPointSet &points = frame.pointSet(ps);
554                 StaticDataPointsChecker              checker(&frame, &points, 0,
555                                                              data.columnCount(points.dataSetIndex()));
556                 EXPECT_CALL(*this, pointsAdded(Property(&AnalysisDataPointSetRef::frameIndex, row)))
557                     .WillOnce(Invoke(checker));
558             }
559             EXPECT_CALL(*this, frameFinished(Property(&AnalysisDataFrameHeader::index, row)))
560                 .WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
561             ::testing::Expectation finish;
562             if (row > 0)
563             {
564                 finish = EXPECT_CALL(*this, frameFinishedSerial(row))
565                         .After(prevFinish);
566             }
567             else
568             {
569                 finish = EXPECT_CALL(*this, frameFinishedSerial(row));
570             }
571             framesFinished += finish;
572             prevFinish      = finish;
573         }
574         EXPECT_CALL(*this, dataFinished())
575             .After(framesFinished);
576     }
577     else
578     {
579         ::testing::InSequence dummy;
580         setSerialExpectationForParallelDataStarted(this);
581         EXPECT_CALL(*this, dataStarted(source));
582         for (int row = 0; row < data.frameCount(); ++row)
583         {
584             const AnalysisDataTestInputFrame &frame = data.frame(row);
585             EXPECT_CALL(*this, frameStarted(_))
586                 .WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
587             for (int ps = 0; ps < frame.pointSetCount(); ++ps)
588             {
589                 const AnalysisDataTestInputPointSet &points = frame.pointSet(ps);
590                 StaticDataPointsChecker              checker(&frame, &points, 0,
591                                                              data.columnCount(points.dataSetIndex()));
592                 EXPECT_CALL(*this, pointsAdded(_)).WillOnce(Invoke(checker));
593             }
594             EXPECT_CALL(*this, frameFinished(_))
595                 .WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
596             EXPECT_CALL(*this, frameFinishedSerial(row));
597         }
598         EXPECT_CALL(*this, dataFinished());
599     }
600 }
601
602
603 void
604 MockAnalysisDataModule::setupStaticColumnCheck(
605         const AnalysisDataTestInput &data,
606         int firstcol, int n, AbstractAnalysisData * /*source*/)
607 {
608     impl_->flags_ |= efAllowMulticolumn | efAllowMultipoint | efAllowMultipleDataSets;
609
610     ::testing::InSequence dummy;
611     using ::testing::_;
612     using ::testing::Invoke;
613
614     setSerialExpectationForParallelDataStarted(this);
615     EXPECT_CALL(*this, dataStarted(_));
616     for (int row = 0; row < data.frameCount(); ++row)
617     {
618         const AnalysisDataTestInputFrame &frame = data.frame(row);
619         EXPECT_CALL(*this, frameStarted(_))
620             .WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
621         for (int ps = 0; ps < frame.pointSetCount(); ++ps)
622         {
623             const AnalysisDataTestInputPointSet &points = frame.pointSet(ps);
624             if (points.lastColumn() >= firstcol
625                 && points.firstColumn() <= firstcol + n - 1)
626             {
627                 EXPECT_CALL(*this, pointsAdded(_))
628                     .WillOnce(Invoke(StaticDataPointsChecker(&frame, &points, firstcol, n)));
629             }
630         }
631         EXPECT_CALL(*this, frameFinished(_))
632             .WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
633         EXPECT_CALL(*this, frameFinishedSerial(row));
634     }
635     EXPECT_CALL(*this, dataFinished());
636 }
637
638
639 void
640 MockAnalysisDataModule::setupStaticStorageCheck(
641         const AnalysisDataTestInput &data,
642         int storageCount, AbstractAnalysisData *source)
643 {
644     GMX_RELEASE_ASSERT(data.isMultipoint() == source->isMultipoint(),
645                        "Mismatching multipoint properties");
646     impl_->flags_ |= efAllowMulticolumn | efAllowMultipoint | efAllowMultipleDataSets;
647
648     ::testing::InSequence dummy;
649     using ::testing::_;
650     using ::testing::Invoke;
651
652     setSerialExpectationForParallelDataStarted(this);
653     EXPECT_CALL(*this, dataStarted(source))
654         .WillOnce(Invoke(DataStorageRequester(storageCount)));
655     for (int row = 0; row < data.frameCount(); ++row)
656     {
657         const AnalysisDataTestInputFrame &frame = data.frame(row);
658         EXPECT_CALL(*this, frameStarted(_))
659             .WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
660         for (int pointSet = 0; pointSet < frame.pointSetCount(); ++pointSet)
661         {
662             StaticDataPointsStorageChecker checker(source, &data, row, pointSet,
663                                                    storageCount);
664             EXPECT_CALL(*this, pointsAdded(_)).WillOnce(Invoke(checker));
665         }
666         EXPECT_CALL(*this, frameFinished(_))
667             .WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
668         EXPECT_CALL(*this, frameFinishedSerial(row));
669     }
670     EXPECT_CALL(*this, dataFinished());
671 }
672
673
674 void
675 MockAnalysisDataModule::setupReferenceCheck(const TestReferenceChecker &checker,
676                                             AbstractAnalysisData       *source)
677 {
678     impl_->flags_ |= efAllowMulticolumn | efAllowMultipoint | efAllowMissing
679         | efAllowMultipleDataSets;
680
681     impl_->rootChecker_ = TestReferenceChecker(checker);
682     // Google Mock does not support checking the order fully, because
683     // the number of frames is not known.
684     // Order of the frameStarted(), pointsAdded() and frameFinished()
685     // calls is checked using Google Test assertions in the invoked methods.
686     using ::testing::_;
687     using ::testing::AnyNumber;
688     using ::testing::Expectation;
689     using ::testing::Invoke;
690
691     setSerialExpectationForParallelDataStarted(this);
692     Expectation dataStart = EXPECT_CALL(*this, dataStarted(source))
693             .WillOnce(Invoke(impl_.get(), &Impl::startReferenceData));
694     Expectation frameStart = EXPECT_CALL(*this, frameStarted(_))
695             .After(dataStart)
696             .WillRepeatedly(Invoke(impl_.get(), &Impl::startReferenceFrame));
697     Expectation pointsAdd = EXPECT_CALL(*this, pointsAdded(_))
698             .After(dataStart)
699             .WillRepeatedly(Invoke(impl_.get(), &Impl::checkReferencePoints));
700     Expectation frameFinish = EXPECT_CALL(*this, frameFinished(_))
701             .After(dataStart)
702             .WillRepeatedly(Invoke(impl_.get(), &Impl::finishReferenceFrame));
703     Expectation frameFinishSerial = EXPECT_CALL(*this, frameFinishedSerial(_))
704             .After(dataStart)
705             .WillRepeatedly(Invoke(impl_.get(), &Impl::finishReferenceFrameSerial));
706     EXPECT_CALL(*this, dataFinished())
707         .After(frameStart, pointsAdd, frameFinish, frameFinishSerial);
708 }
709
710 } // namespace test
711 } // namespace gmx