Apply clang-format to source tree
[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,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 /*! \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, const AnalysisDataValue& value)
144 {
145     TestReferenceChecker compound(checker->checkCompound("DataValue", nullptr));
146     compound.checkReal(value.value(), "Value");
147     if (compound.checkPresent(value.hasError(), "Error"))
148     {
149         compound.checkReal(value.error(), "Error");
150     }
151     if (compound.checkPresent(!value.isPresent(), "Present"))
152     {
153         compound.checkBoolean(value.isPresent(), "Present");
154     }
155 }
156
157 } // namespace
158
159 MockAnalysisDataModule::Impl::Impl(int flags) : source_(nullptr), flags_(flags), frameIndex_(0) {}
160
161
162 void MockAnalysisDataModule::Impl::startReferenceData(AbstractAnalysisData* data)
163 {
164     source_ = data;
165 }
166
167
168 void MockAnalysisDataModule::Impl::startReferenceFrame(const AnalysisDataFrameHeader& header)
169 {
170     GMX_RELEASE_ASSERT(rootChecker_, "Root checker not set, but reference data used");
171     EXPECT_FALSE(frameChecker_.isValid());
172     EXPECT_EQ(frameIndex_, header.index());
173     frameChecker_ =
174             rootChecker_.checkCompound("DataFrame", formatString("Frame%d", frameIndex_).c_str());
175     frameChecker_.checkReal(header.x(), "X");
176 }
177
178
179 void MockAnalysisDataModule::Impl::checkReferencePoints(const AnalysisDataPointSetRef& points)
180 {
181     EXPECT_TRUE(frameChecker_.isValid());
182     if (frameChecker_)
183     {
184         TestReferenceChecker checker(frameChecker_.checkCompound("DataValues", nullptr));
185         checker.checkInteger(points.columnCount(), "Count");
186         if (checker.checkPresent(source_->dataSetCount() > 1, "DataSet"))
187         {
188             checker.checkInteger(points.dataSetIndex(), "DataSet");
189         }
190         const int sourceColumnCount = source_->columnCount(points.dataSetIndex());
191         const bool bAllColumns = (points.firstColumn() == 0 && points.columnCount() == sourceColumnCount);
192         if (checker.checkPresent(!bAllColumns, "FirstColumn"))
193         {
194             checker.checkInteger(points.firstColumn(), "FirstColumn");
195             checker.checkInteger(points.lastColumn(), "LastColumn");
196         }
197
198         AnalysisDataValuesRef::const_iterator value;
199         for (value = points.values().begin(); value != points.values().end(); ++value)
200         {
201             checkReferenceDataPoint(&checker, *value);
202         }
203     }
204 }
205
206
207 void MockAnalysisDataModule::Impl::finishReferenceFrame(const AnalysisDataFrameHeader& header)
208 {
209     EXPECT_TRUE(frameChecker_.isValid());
210     EXPECT_EQ(frameIndex_, header.index());
211     frameChecker_ = TestReferenceChecker();
212 }
213
214
215 void MockAnalysisDataModule::Impl::finishReferenceFrameSerial(int frameIndex)
216 {
217     EXPECT_FALSE(frameChecker_.isValid());
218     EXPECT_EQ(frameIndex_, frameIndex);
219     ++frameIndex_;
220 }
221
222
223 /********************************************************************
224  * MockAnalysisDataModule
225  */
226
227 namespace
228 {
229
230 /*! \brief
231  * Helper function for checking the data frame header against static data.
232  *
233  * \param[in] header    Frame header to check.
234  * \param[in] refFrame  Data to check against.
235  */
236 void checkHeader(const AnalysisDataFrameHeader& header, const AnalysisDataTestInputFrame& refFrame)
237 {
238     EXPECT_EQ(refFrame.index(), header.index());
239     EXPECT_FLOAT_EQ(refFrame.x(), header.x());
240     EXPECT_FLOAT_EQ(refFrame.dx(), header.dx());
241 }
242
243 /*! \brief
244  * Helper function for checking a point set against static data.
245  *
246  * \param[in] points       Point set to check.
247  * \param[in] refPoints    Data to check against.
248  * \param[in] columnOffset Offset of first column of \p points in \p refPoints.
249  */
250 void checkPoints(const AnalysisDataPointSetRef&       points,
251                  const AnalysisDataTestInputPointSet& refPoints,
252                  int                                  columnOffset)
253 {
254     for (int i = 0; i < points.columnCount(); ++i)
255     {
256         const int column = points.firstColumn() - refPoints.firstColumn() + i + columnOffset;
257         EXPECT_FLOAT_EQ(refPoints.y(column), points.y(i))
258                 << "  Column: " << i + 1 << " / " << points.columnCount() << " (+"
259                 << points.firstColumn() << ")\n"
260                 << "Ref. col: " << column + 1 << " / " << refPoints.size() << " (+"
261                 << refPoints.firstColumn() << ", offs " << columnOffset << ")";
262     }
263 }
264
265 /*! \brief
266  * Helper function for checking a full frame against static data.
267  *
268  * \param[in] frame     Frame to check.
269  * \param[in] refFrame  Data to check against.
270  */
271 void checkFrame(const AnalysisDataFrameRef& frame, const AnalysisDataTestInputFrame& refFrame)
272 {
273     checkHeader(frame.header(), refFrame);
274     ASSERT_EQ(refFrame.pointSetCount(), frame.pointSetCount());
275     for (int i = 0; i < frame.pointSetCount(); ++i)
276     {
277         const AnalysisDataPointSetRef&       points    = frame.pointSet(i);
278         const AnalysisDataTestInputPointSet& refPoints = refFrame.pointSet(i);
279         EXPECT_EQ(refPoints.firstColumn(), points.firstColumn());
280         checkPoints(points, refPoints, 0);
281     }
282 }
283
284 /*! \brief
285  * Functor for checking data frame header against static test input data.
286  *
287  * This functor is designed to be invoked as a handled for
288  * IAnalysisDataModule::frameStarted().
289  */
290 class StaticDataFrameHeaderChecker
291 {
292 public:
293     /*! \brief
294      * Constructs a checker against a given input data frame.
295      *
296      * \param[in] frame Frame to check against.
297      *
298      * \p frame must exist for the lifetime of this object.
299      */
300     StaticDataFrameHeaderChecker(const AnalysisDataTestInputFrame* frame) : frame_(frame) {}
301
302     //! Function call operator for the functor.
303     void operator()(const AnalysisDataFrameHeader& header) const
304     {
305         SCOPED_TRACE(formatString("Frame %d", frame_->index()));
306         checkHeader(header, *frame_);
307     }
308
309 private:
310     const AnalysisDataTestInputFrame* frame_;
311 };
312
313 /*! \brief
314  * Functor for checking data frame points against static test input data.
315  *
316  * This functor is designed to be invoked as a handled for
317  * IAnalysisDataModule::pointsAdded().
318  */
319 class StaticDataPointsChecker
320 {
321 public:
322     /*! \brief
323      * Constructs a checker against a given input data frame and point set.
324      *
325      * \param[in] frame    Frame to check against.
326      * \param[in] points   Point set in \p frame to check against.
327      * \param[in] firstcol Expected first column.
328      * \param[in] n        Expected number of columns.
329      *
330      * \p firstcol and \p n are used to create a checker that only expects
331      * to be called for a subset of columns.
332      * \p frame and \p points must exist for the lifetime of this object.
333      */
334     StaticDataPointsChecker(const AnalysisDataTestInputFrame*    frame,
335                             const AnalysisDataTestInputPointSet* points,
336                             int                                  firstcol,
337                             int                                  n) :
338         frame_(frame),
339         points_(points),
340         firstcol_(firstcol),
341         n_(n)
342     {
343     }
344
345     //! Function call operator for the functor.
346     void operator()(const AnalysisDataPointSetRef& points) const
347     {
348         SCOPED_TRACE(formatString("Frame %d, point set %d", frame_->index(), points_->index()));
349         EXPECT_EQ(points_->dataSetIndex(), points.dataSetIndex());
350         const int expectedFirstColumn = std::max(0, points_->firstColumn() - firstcol_);
351         const int expectedLastColumn  = std::min(n_ - 1, points_->lastColumn() - firstcol_);
352         EXPECT_EQ(expectedFirstColumn, points.firstColumn());
353         EXPECT_EQ(expectedLastColumn, points.lastColumn());
354         checkHeader(points.header(), *frame_);
355         checkPoints(points, *points_, firstcol_);
356     }
357
358 private:
359     const AnalysisDataTestInputFrame*    frame_;
360     const AnalysisDataTestInputPointSet* points_;
361     int                                  firstcol_;
362     int                                  n_;
363 };
364
365 /*! \brief
366  * Functor for requesting data storage.
367  *
368  * This functor is designed to be invoked as a handled for
369  * IAnalysisDataModule::dataStarted().
370  */
371 class DataStorageRequester
372 {
373 public:
374     /*! \brief
375      * Constructs a functor that requests the given amount of storage.
376      *
377      * \param[in] count  Number of frames of storage to request, or
378      *      -1 for all frames.
379      *
380      * \see AbstractAnalysisData::requestStorage()
381      */
382     explicit DataStorageRequester(int count) : count_(count) {}
383
384     //! Function call operator for the functor.
385     void operator()(AbstractAnalysisData* data) const { EXPECT_TRUE(data->requestStorage(count_)); }
386
387 private:
388     int count_;
389 };
390
391 /*! \brief
392  * Functor for checking data frame points and storage against static test input
393  * data.
394  *
395  * This functor is designed to be invoked as a handled for
396  * IAnalysisDataModule::pointsAdded().
397  */
398 class StaticDataPointsStorageChecker
399 {
400 public:
401     /*! \brief
402      * Constructs a checker for a given frame.
403      *
404      * \param[in] source     Data object that is being checked.
405      * \param[in] data       Test input data to check against.
406      * \param[in] frameIndex Frame index for which this functor expects
407      *      to be called.
408      * \param[in] pointSetIndex Point set for which this functor expects
409      *      to be called.
410      * \param[in] storageCount How many past frames should be checked for
411      *      storage (-1 = check all frames).
412      *
413      * This checker works as StaticDataPointsChecker, but additionally
414      * checks that previous frames can be accessed using access methods
415      * in AbstractAnalysisData and that correct data is returned.
416      *
417      * \p source and \p data must exist for the lifetime of this object.
418      */
419     StaticDataPointsStorageChecker(AbstractAnalysisData*        source,
420                                    const AnalysisDataTestInput* data,
421                                    int                          frameIndex,
422                                    int                          pointSetIndex,
423                                    int                          storageCount) :
424         source_(source),
425         data_(data),
426         frameIndex_(frameIndex),
427         pointSetIndex_(pointSetIndex),
428         storageCount_(storageCount)
429     {
430     }
431
432     //! Function call operator for the functor.
433     void operator()(const AnalysisDataPointSetRef& points) const
434     {
435         SCOPED_TRACE(formatString("Frame %d", frameIndex_));
436         const AnalysisDataTestInputFrame&    refFrame  = data_->frame(frameIndex_);
437         const AnalysisDataTestInputPointSet& refPoints = refFrame.pointSet(pointSetIndex_);
438         EXPECT_EQ(refPoints.firstColumn(), points.firstColumn());
439         EXPECT_EQ(refPoints.size(), points.columnCount());
440         checkHeader(points.header(), refFrame);
441         checkPoints(points, refPoints, 0);
442         for (int past = 1; (storageCount_ < 0 || past <= storageCount_) && past <= frameIndex_; ++past)
443         {
444             int index = frameIndex_ - past;
445             SCOPED_TRACE(formatString("Checking storage of frame %d", index));
446             ASSERT_NO_THROW_GMX({
447                 AnalysisDataFrameRef frame = source_->getDataFrame(index);
448                 ASSERT_TRUE(frame.isValid());
449                 checkFrame(frame, data_->frame(index));
450             });
451         }
452     }
453
454 private:
455     AbstractAnalysisData*        source_;
456     const AnalysisDataTestInput* data_;
457     int                          frameIndex_;
458     int                          pointSetIndex_;
459     int                          storageCount_;
460 };
461
462 /*! \brief
463  * Sets the mock object expectation to mimick AnalysisDataModuleSerial.
464  *
465  * Makes MockAnalysisDataModule::parallelDataStarted() behave as if the mock
466  * object was an AnalysisDataModuleSerial object: forward the call to
467  * MockAnalysisDataModule::dataStarted() and return false.
468  */
469 void setSerialExpectationForParallelDataStarted(MockAnalysisDataModule* mock)
470 {
471     using ::testing::_;
472     using ::testing::AtMost;
473     using ::testing::DoAll;
474     using ::testing::Invoke;
475     using ::testing::Return;
476     using ::testing::WithArg;
477     EXPECT_CALL(*mock, parallelDataStarted(_, _))
478             .Times(AtMost(1))
479             .WillOnce(DoAll(WithArg<0>(Invoke(mock, &MockAnalysisDataModule::dataStarted)), Return(false)));
480 }
481
482 } // anonymous namespace
483
484
485 MockAnalysisDataModule::MockAnalysisDataModule(int flags) : impl_(new Impl(flags)) {}
486
487
488 MockAnalysisDataModule::~MockAnalysisDataModule() {}
489
490
491 int MockAnalysisDataModule::flags() const
492 {
493     return impl_->flags_;
494 }
495
496
497 void MockAnalysisDataModule::setupStaticCheck(const AnalysisDataTestInput& data,
498                                               AbstractAnalysisData*        source,
499                                               bool                         bParallel)
500 {
501     impl_->flags_ |= efAllowMulticolumn | efAllowMultipoint | efAllowMultipleDataSets;
502
503     using ::testing::_;
504     using ::testing::Invoke;
505     using ::testing::Property;
506     using ::testing::Return;
507
508     if (bParallel)
509     {
510 #ifndef __clang_analyzer__
511         ::testing::Expectation init =
512                 EXPECT_CALL(*this, parallelDataStarted(source, _)).WillOnce(Return(true));
513         ::testing::ExpectationSet framesFinished;
514         ::testing::Expectation    prevFinish;
515         for (int row = 0; row < data.frameCount(); ++row)
516         {
517             ::testing::InSequence             frameSequence;
518             const AnalysisDataTestInputFrame& frame = data.frame(row);
519             EXPECT_CALL(*this, frameStarted(Property(&AnalysisDataFrameHeader::index, row)))
520                     .After(init)
521                     .WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
522             for (int ps = 0; ps < frame.pointSetCount(); ++ps)
523             {
524                 const AnalysisDataTestInputPointSet& points = frame.pointSet(ps);
525                 StaticDataPointsChecker              checker(&frame, &points, 0,
526                                                 data.columnCount(points.dataSetIndex()));
527                 EXPECT_CALL(*this, pointsAdded(Property(&AnalysisDataPointSetRef::frameIndex, row)))
528                         .WillOnce(Invoke(checker));
529             }
530             EXPECT_CALL(*this, frameFinished(Property(&AnalysisDataFrameHeader::index, row)))
531                     .WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
532             ::testing::Expectation finish;
533             if (row > 0)
534             {
535                 finish = EXPECT_CALL(*this, frameFinishedSerial(row)).After(prevFinish);
536             }
537             else
538             {
539                 finish = EXPECT_CALL(*this, frameFinishedSerial(row));
540             }
541             framesFinished += finish;
542             prevFinish = finish;
543         }
544         EXPECT_CALL(*this, dataFinished()).After(framesFinished);
545 #endif
546     }
547     else
548     {
549         ::testing::InSequence dummy;
550         setSerialExpectationForParallelDataStarted(this);
551         EXPECT_CALL(*this, dataStarted(source));
552         for (int row = 0; row < data.frameCount(); ++row)
553         {
554             const AnalysisDataTestInputFrame& frame = data.frame(row);
555             EXPECT_CALL(*this, frameStarted(_)).WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
556             for (int ps = 0; ps < frame.pointSetCount(); ++ps)
557             {
558                 const AnalysisDataTestInputPointSet& points = frame.pointSet(ps);
559                 StaticDataPointsChecker              checker(&frame, &points, 0,
560                                                 data.columnCount(points.dataSetIndex()));
561                 EXPECT_CALL(*this, pointsAdded(_)).WillOnce(Invoke(checker));
562             }
563             EXPECT_CALL(*this, frameFinished(_)).WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
564             EXPECT_CALL(*this, frameFinishedSerial(row));
565         }
566         EXPECT_CALL(*this, dataFinished());
567     }
568 }
569
570
571 void MockAnalysisDataModule::setupStaticColumnCheck(const AnalysisDataTestInput& data,
572                                                     int                          firstcol,
573                                                     int                          n,
574                                                     AbstractAnalysisData* /*source*/)
575 {
576     impl_->flags_ |= efAllowMulticolumn | efAllowMultipoint | efAllowMultipleDataSets;
577
578     ::testing::InSequence dummy;
579     using ::testing::_;
580     using ::testing::Invoke;
581
582     setSerialExpectationForParallelDataStarted(this);
583     EXPECT_CALL(*this, dataStarted(_));
584     for (int row = 0; row < data.frameCount(); ++row)
585     {
586         const AnalysisDataTestInputFrame& frame = data.frame(row);
587         EXPECT_CALL(*this, frameStarted(_)).WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
588         for (int ps = 0; ps < frame.pointSetCount(); ++ps)
589         {
590             const AnalysisDataTestInputPointSet& points = frame.pointSet(ps);
591             if (points.lastColumn() >= firstcol && points.firstColumn() <= firstcol + n - 1)
592             {
593                 EXPECT_CALL(*this, pointsAdded(_))
594                         .WillOnce(Invoke(StaticDataPointsChecker(&frame, &points, firstcol, n)));
595             }
596         }
597         EXPECT_CALL(*this, frameFinished(_)).WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
598         EXPECT_CALL(*this, frameFinishedSerial(row));
599     }
600     EXPECT_CALL(*this, dataFinished());
601 }
602
603
604 void MockAnalysisDataModule::setupStaticStorageCheck(const AnalysisDataTestInput& data,
605                                                      int                          storageCount,
606                                                      AbstractAnalysisData*        source)
607 {
608     GMX_RELEASE_ASSERT(data.isMultipoint() == source->isMultipoint(),
609                        "Mismatching multipoint properties");
610     impl_->flags_ |= efAllowMulticolumn | efAllowMultipoint | efAllowMultipleDataSets;
611
612     ::testing::InSequence dummy;
613     using ::testing::_;
614     using ::testing::Invoke;
615
616     setSerialExpectationForParallelDataStarted(this);
617     EXPECT_CALL(*this, dataStarted(source)).WillOnce(Invoke(DataStorageRequester(storageCount)));
618     for (int row = 0; row < data.frameCount(); ++row)
619     {
620         const AnalysisDataTestInputFrame& frame = data.frame(row);
621         EXPECT_CALL(*this, frameStarted(_)).WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
622         for (int pointSet = 0; pointSet < frame.pointSetCount(); ++pointSet)
623         {
624             StaticDataPointsStorageChecker checker(source, &data, row, pointSet, storageCount);
625             EXPECT_CALL(*this, pointsAdded(_)).WillOnce(Invoke(checker));
626         }
627         EXPECT_CALL(*this, frameFinished(_)).WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
628         EXPECT_CALL(*this, frameFinishedSerial(row));
629     }
630     EXPECT_CALL(*this, dataFinished());
631 }
632
633
634 void MockAnalysisDataModule::setupReferenceCheck(const TestReferenceChecker& checker,
635                                                  AbstractAnalysisData*       source)
636 {
637     impl_->flags_ |= efAllowMulticolumn | efAllowMultipoint | efAllowMissing | efAllowMultipleDataSets;
638
639     impl_->rootChecker_ = TestReferenceChecker(checker);
640     // Google Mock does not support checking the order fully, because
641     // the number of frames is not known.
642     // Order of the frameStarted(), pointsAdded() and frameFinished()
643     // calls is checked using Google Test assertions in the invoked methods.
644     using ::testing::_;
645     using ::testing::AnyNumber;
646     using ::testing::Expectation;
647     using ::testing::Invoke;
648
649     setSerialExpectationForParallelDataStarted(this);
650     Expectation dataStart =
651             EXPECT_CALL(*this, dataStarted(source)).WillOnce(Invoke(impl_.get(), &Impl::startReferenceData));
652     Expectation frameStart = EXPECT_CALL(*this, frameStarted(_))
653                                      .After(dataStart)
654                                      .WillRepeatedly(Invoke(impl_.get(), &Impl::startReferenceFrame));
655     Expectation pointsAdd = EXPECT_CALL(*this, pointsAdded(_))
656                                     .After(dataStart)
657                                     .WillRepeatedly(Invoke(impl_.get(), &Impl::checkReferencePoints));
658     Expectation frameFinish = EXPECT_CALL(*this, frameFinished(_))
659                                       .After(dataStart)
660                                       .WillRepeatedly(Invoke(impl_.get(), &Impl::finishReferenceFrame));
661     Expectation frameFinishSerial =
662             EXPECT_CALL(*this, frameFinishedSerial(_))
663                     .After(dataStart)
664                     .WillRepeatedly(Invoke(impl_.get(), &Impl::finishReferenceFrameSerial));
665     EXPECT_CALL(*this, dataFinished()).After(frameStart, pointsAdd, frameFinish, frameFinishSerial);
666 }
667
668 } // namespace test
669 } // namespace gmx