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