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