Added some more analysisdata unit tests.
authorTeemu Murtola <teemu.murtola@cbr.su.se>
Sat, 29 Oct 2011 18:45:01 +0000 (21:45 +0300)
committerTeemu Murtola <teemu.murtola@gmail.com>
Tue, 21 Feb 2012 20:03:33 +0000 (22:03 +0200)
Added support for multipoint data to the analysis data testing
framework, and added some more tests for AnalysisData using it.

Also fixed a bug (found by the new tests) with frameCount() not being
properly updated for multipoint data.  The count is now kept by
AbstractAnalysisData, so derived classes no longer need to worry about
keeping it up-to-date.

Part of issue #823.

Change-Id: I672fbfb2ec5898b8118616574646f1432c7c0ce1

src/gromacs/analysisdata/abstractdata-impl.h
src/gromacs/analysisdata/abstractdata.cpp
src/gromacs/analysisdata/abstractdata.h
src/gromacs/analysisdata/arraydata.cpp
src/gromacs/analysisdata/arraydata.h
src/gromacs/analysisdata/modules/displacement.cpp
src/gromacs/analysisdata/modules/displacement.h
src/gromacs/analysisdata/tests/analysisdata.cpp
src/gromacs/analysisdata/tests/datatest.cpp
src/gromacs/analysisdata/tests/datatest.h
src/gromacs/analysisdata/tests/mock_module.cpp

index eb768271040ac4fca0b8a32c30b588c91cb0fa0c..31d06252996ac03b37788b2ff6de28c5a85ee602 100644 (file)
@@ -92,6 +92,12 @@ class AbstractAnalysisData::Impl
         real                    _currx;
         //! dx value for the current frame.
         real                    _currdx;
+        /*! \brief
+         * Total number of frames in the data.
+         *
+         * The counter is incremented in notifyFrameStart().
+         */
+        int                     _nframes;
 };
 
 /*! \internal \brief
@@ -141,15 +147,12 @@ class AbstractAnalysisDataStored::Impl
          *
          * \param[in] index  Zero-based index for the frame to query.
          *      Negative value counts backwards from the current frame.
+         * \param[in] nframes  Total number of frames in the data.
          * \returns Index in \a _store corresponding to \p index,
          *      or -1 if not available.
          */
-        int getStoreIndex(int index) const;
+        int getStoreIndex(int index, int nframes) const;
 
-        /*! \brief
-         * Total number of complete frames in the data.
-         */
-        int                     _nframes;
         /*! \brief
          * Number of elements in \a _store.
          *
index 4c09c91dfe63bf1b754715efcb19cfefc688f96f..de562e5c848cc47cc434316f3b8caa89c3767a5d 100644 (file)
@@ -57,7 +57,7 @@ namespace gmx
 
 AbstractAnalysisData::Impl::Impl()
     : _bDataStart(false), _bInData(false), _bInFrame(false),
-      _bAllowMissing(true)
+      _bAllowMissing(true), _nframes(0)
 {
 }
 
@@ -126,6 +126,13 @@ AbstractAnalysisData::~AbstractAnalysisData()
 }
 
 
+int
+AbstractAnalysisData::frameCount() const
+{
+    return _impl->_nframes;
+}
+
+
 bool
 AbstractAnalysisData::getData(int index, real *x, const real **y,
                               const bool **missing) const
@@ -258,6 +265,7 @@ AbstractAnalysisData::notifyFrameStart(real x, real dx) const
     _impl->_bInFrame = true;
     _impl->_currx  = x;
     _impl->_currdx = dx;
+    ++_impl->_nframes;
 
     Impl::ModuleList::const_iterator i;
     for (i = _impl->_modules.begin(); i != _impl->_modules.end(); ++i)
@@ -358,7 +366,7 @@ void AnalysisDataFrame::allocate(int ncol)
  */
 
 AbstractAnalysisDataStored::Impl::Impl()
-    : _nframes(0), _nalloc(0), _bStoreAll(false), _nextind(-1)
+    : _nalloc(0), _bStoreAll(false), _nextind(-1)
 {
 }
 
@@ -374,11 +382,11 @@ AbstractAnalysisDataStored::Impl::~Impl()
 
 
 int
-AbstractAnalysisDataStored::Impl::getStoreIndex(int index) const
+AbstractAnalysisDataStored::Impl::getStoreIndex(int index, int nframes) const
 {
     // Check that the requested index is available.
-    if ((index < 0 && (-index > _nalloc || -index > _nframes))
-        || index >= _nframes || (index >= 0 && index < _nframes - _nalloc))
+    if ((index < 0 && (-index > _nalloc || -index > nframes))
+        || index >= nframes || (index >= 0 && index < nframes - _nalloc))
     {
         return -1;
     }
@@ -391,7 +399,7 @@ AbstractAnalysisDataStored::Impl::getStoreIndex(int index) const
             index += _nalloc;
         }
     }
-    else if (_nframes > _nalloc)
+    else if (nframes > _nalloc)
     {
         index %= _nalloc;
     }
@@ -415,19 +423,12 @@ AbstractAnalysisDataStored::~AbstractAnalysisDataStored()
 }
 
 
-int
-AbstractAnalysisDataStored::frameCount() const
-{
-    return _impl->_nframes;
-}
-
-
 bool
 AbstractAnalysisDataStored::getDataWErr(int index, real *x, real *dx,
                                         const real **y, const real **dy,
                                         const bool **present) const
 {
-    index = _impl->getStoreIndex(index);
+    index = _impl->getStoreIndex(index, frameCount());
     if (index < 0)
     {
         return false;
@@ -578,7 +579,6 @@ AbstractAnalysisDataStored::storeThisFrame(const real *y, const real *dy,
             }
         }
     }
-    ++_impl->_nframes;
 
     // Notify modules of new data.
     notifyPointsAdd(0, ncol, y, dy, present);
index a55af2fe13f144e36004e9bfd23469a4f8305b45..884ab055182311aaf5c9f4510a19dd720d72465a 100644 (file)
@@ -131,7 +131,7 @@ class AbstractAnalysisData
          * produced. If requestStorage() has been successfully called,
          * getData() can be used to access some or all of these frames.
          */
-        virtual int frameCount() const = 0;
+        int frameCount() const;
         /*! \brief
          * Access stored data.
          *
@@ -348,7 +348,6 @@ class AbstractAnalysisDataStored : public AbstractAnalysisData
     public:
         virtual ~AbstractAnalysisDataStored();
 
-        virtual int frameCount() const;
         virtual bool getDataWErr(int index, real *x, real *dx,
                                  const real **y, const real **dy,
                                  const bool **present = 0) const;
index f1f24110e2b3298d73066768e3ebb20167df0a5b..add59fb310270ab67e1bbe6e6de1bf4caaa6b4c6 100644 (file)
@@ -57,13 +57,6 @@ AbstractAnalysisArrayData::~AbstractAnalysisArrayData()
 }
 
 
-int
-AbstractAnalysisArrayData::frameCount() const
-{
-    return _bReady ? _nrows : 0;
-}
-
-
 bool
 AbstractAnalysisArrayData::getDataWErr(int index, real *x, real *dx,
                                        const real **y, const real **dy,
index 5405bdf851db2de2b918638c86fef464ffd95fc9..60eed055c66b1c234c1dacf8c7af84f48329c712 100644 (file)
@@ -64,7 +64,6 @@ class AbstractAnalysisArrayData : public AbstractAnalysisData
     public:
         virtual ~AbstractAnalysisArrayData();
 
-        virtual int frameCount() const;
         virtual bool getDataWErr(int index, real *x, real *dx,
                                  const real **y, const real **dy,
                                  const bool **present = 0) const;
index a799522d77a45c171471db4203cfe7ed4b08baa9..8451e07dfc6084255ad8738ff1f7de992f11fed3 100644 (file)
@@ -106,13 +106,6 @@ AnalysisDataDisplacementModule::setMSDHistogram(AnalysisDataBinAverageModule *hi
 }
 
 
-int
-AnalysisDataDisplacementModule::frameCount() const
-{
-    return _impl->nstored > 1 ? _impl->nstored - 1 : 0;
-}
-
-
 bool
 AnalysisDataDisplacementModule::getDataWErr(int index, real *x, real *dx,
                                             const real **y, const real **dy,
index 9c44dc0996e12123e44281742b73676886b2060a..6a0d9a07722724327e9aa3fdce17d015289464b5 100644 (file)
@@ -79,7 +79,6 @@ class AnalysisDataDisplacementModule : public AbstractAnalysisData,
          */
         void setMSDHistogram(AnalysisDataBinAverageModule *histm);
 
-        virtual int frameCount() const;
         virtual bool getDataWErr(int index, real *x, real *dx,
                                  const real **y, const real **dy,
                                  const bool **present = 0) const;
index f7c66fb4d7ccb13b3e709a3d0cc3a9edb590faf9..42c958b97279d7e1fe613fdf5f1f99e1d41f9588 100644 (file)
@@ -121,6 +121,7 @@ typedef gmx::test::AnalysisDataTestFixture AnalysisDataTest;
 // Input data for the tests below.
 using gmx::test::END_OF_DATA;
 using gmx::test::END_OF_FRAME;
+using gmx::test::MPSTOP;
 static const real inputdata[] = {
     1.0,  0.0, 1.0, 2.0, END_OF_FRAME,
     2.0,  1.0, 1.0, 1.0, END_OF_FRAME,
@@ -176,10 +177,11 @@ TEST_F(AnalysisDataTest, CallsModuleCorrectlyWithIndividualPoints)
     for (int row = 0; row < input.frameCount(); ++row)
     {
         const gmx::test::AnalysisDataTestInputFrame &frame = input.frame(row);
+        const gmx::test::AnalysisDataTestInputPointSet &points = frame.points();
         ASSERT_NO_THROW(handle->startFrame(row, frame.x(), frame.dx()));
         for (int column = 0; column < input.columnCount(); ++column)
         {
-            ASSERT_NO_THROW(handle->addPoint(column, frame.yptr()[column]));
+            ASSERT_NO_THROW(handle->addPoint(column, points.y(column)));
         }
         ASSERT_NO_THROW(handle->finishFrame());
         EXPECT_EQ(row + 1, data.frameCount());
@@ -253,4 +255,43 @@ TEST_F(AnalysisDataTest, LimitedStorageWorks)
     ASSERT_NO_THROW(presentAllData(input, &data));
 }
 
+// Input data for the tests below.
+static const real multipointinputdata[] = {
+    1.0,  0.0, 1.0, 2.0, MPSTOP, 1.1, 2.1, 1.1, MPSTOP, 2.2, 1.2, 0.2, END_OF_FRAME,
+    2.0,  1.0, 1.0, 1.0, MPSTOP, 2.1, 1.1, 0.1, MPSTOP, 1.2, 0.2, 1.2, END_OF_FRAME,
+    3.0,  2.0, 0.0, 0.0, MPSTOP, 3.1, 2.1, 1.1, MPSTOP, 0.2, 2.2, 1.2, END_OF_FRAME,
+    END_OF_DATA
+};
+
+/*
+ * Tests that multipoint data is forwarded correctly to modules using two
+ * independent modules.
+ */
+TEST_F(AnalysisDataTest, MultipointCallsModuleCorrectly)
+{
+    gmx::test::AnalysisDataTestInput input(multipointinputdata);
+    gmx::AnalysisData data;
+    data.setColumns(input.columnCount(), true);
+
+    ASSERT_NO_THROW(addStaticCheckerModule(input, &data));
+    ASSERT_NO_THROW(addStaticCheckerModule(input, &data));
+    ASSERT_NO_THROW(presentAllData(input, &data));
+}
+
+/*
+ * Tests that multipoint data is forwarded correctly to modules that are added
+ * using addColumnModule().
+ * Uses two independent modules.
+ */
+TEST_F(AnalysisDataTest, MultipointCallsColumnModuleCorrectly)
+{
+    gmx::test::AnalysisDataTestInput input(multipointinputdata);
+    gmx::AnalysisData data;
+    data.setColumns(input.columnCount(), true);
+
+    ASSERT_NO_THROW(addStaticColumnCheckerModule(input, 0, 2, &data));
+    ASSERT_NO_THROW(addStaticColumnCheckerModule(input, 2, 1, &data));
+    ASSERT_NO_THROW(presentAllData(input, &data));
+}
+
 } // namespace
index 8a202824ffdf847fbb5163f665b15cc18e0eebef..589232f71b603d3a9d14b8820fa17038b7f13262 100644 (file)
@@ -57,73 +57,70 @@ namespace test
 {
 
 /********************************************************************
- * AnalysisDataTestInputFrame
+ * AnalysisDataTestInputPointSet
  */
 
-AnalysisDataTestInputFrame::AnalysisDataTestInputFrame()
+AnalysisDataTestInputPointSet::AnalysisDataTestInputPointSet()
 {
 }
 
 
 /********************************************************************
- * AnalysisDataTestInput
+ * AnalysisDataTestInputFrame
  */
 
-namespace
+AnalysisDataTestInputFrame::AnalysisDataTestInputFrame()
 {
-    void checkValidDataItem(real data)
-    {
-        GMX_RELEASE_ASSERT(data != END_OF_DATA && data != END_OF_FRAME,
-                           "Inconsistent data array");
-    }
 }
 
+
+/********************************************************************
+ * AnalysisDataTestInput
+ */
+
 AnalysisDataTestInput::AnalysisDataTestInput(const real *data)
-    : columnCount_(0)
+    : columnCount_(0), bMultipoint_(false)
 {
-    // Count rows and columns.
-    int rows = 0, columns = -1;
+    size_t columns = 0;
     const real *dataptr = data;
-    for (int i = 0; dataptr[i] != END_OF_DATA; ++i)
+
+    while (*dataptr != END_OF_DATA)
     {
-        if (dataptr[i] == END_OF_FRAME)
+        frames_.push_back(AnalysisDataTestInputFrame());
+        AnalysisDataTestInputFrame &frame = frames_.back();
+        GMX_RELEASE_ASSERT(*dataptr != END_OF_FRAME && *dataptr != MPSTOP,
+                           "Empty data frame");
+        frame.index_ = frames_.size() - 1;
+        frame.x_ = *dataptr;
+        while (*dataptr != END_OF_FRAME)
         {
-            GMX_RELEASE_ASSERT(i > 0, "Empty data frame");
-            if (columns == -1)
+            ++dataptr;
+            frame.points_.push_back(AnalysisDataTestInputPointSet());
+            AnalysisDataTestInputPointSet &points = frame.points_.back();
+            while (*dataptr != MPSTOP && *dataptr != END_OF_FRAME)
+            {
+                GMX_RELEASE_ASSERT(*dataptr != END_OF_DATA,
+                                   "Premature end of data marker");
+                points.y_.push_back(*dataptr);
+                ++dataptr;
+            }
+            size_t frameColumns = points.y_.size();
+            GMX_RELEASE_ASSERT(frameColumns > 0U, "Empty data point set");
+            if (*dataptr == MPSTOP)
             {
-                columns = i - 1;
+                bMultipoint_ = true;
             }
-            GMX_RELEASE_ASSERT(columns == i - 1,
+            GMX_RELEASE_ASSERT(!(!bMultipoint_ && columns > 0U && columns != frameColumns),
                                "Different frames have different number of columns");
-            ++rows;
-            dataptr += i + 1;
-            i = -1;
+            if (columns < frameColumns)
+            {
+                columns = frameColumns;
+            }
         }
-    }
-    GMX_RELEASE_ASSERT(rows > 0, "Empty data");
-    columnCount_ = columns;
-
-    // Store the data.
-    frames_.resize(rows);
-    dataptr = data;
-    for (int r = 0; r < rows; ++r)
-    {
-        AnalysisDataTestInputFrame &frame = frames_[r];
-        checkValidDataItem(*dataptr);
-        frame.index_ = r;
-        frame.x_ = *dataptr;
         ++dataptr;
-        frame.y_.resize(columns);
-        for (int c = 0; c < columns; ++c)
-        {
-            checkValidDataItem(dataptr[c]);
-            frame.y_[c] = dataptr[c];
-        }
-        GMX_RELEASE_ASSERT(dataptr[columns] == END_OF_FRAME,
-                           "Inconsistent data array");
-        dataptr += columns + 1;
     }
-    GMX_RELEASE_ASSERT(*dataptr == END_OF_DATA, "Inconsistent data array");
+    GMX_RELEASE_ASSERT(frames_.size() > 0U, "Empty data");
+    columnCount_ = columns;
 }
 
 
@@ -168,8 +165,12 @@ void AnalysisDataTestFixture::presentDataFrame(const AnalysisDataTestInput &inpu
 {
     const AnalysisDataTestInputFrame &frame = input.frame(row);
     handle->startFrame(row, frame.x(), frame.dx());
-    handle->addPoints(0, input.columnCount(), frame.yptr(), frame.dyptr(),
-                      frame.presentptr());
+    for (int i = 0; i < frame.pointSetCount(); ++i)
+    {
+        const AnalysisDataTestInputPointSet &points = frame.points(i);
+        handle->addPoints(0, points.size(), points.yptr(), points.dyptr(),
+                          points.presentptr());
+    }
     handle->finishFrame();
 }
 
@@ -179,7 +180,8 @@ AnalysisDataTestFixture::addStaticCheckerModule(const AnalysisDataTestInput &dat
                                                 AbstractAnalysisData *source)
 {
     std::auto_ptr<MockAnalysisModule> module(
-            new MockAnalysisModule(gmx::AnalysisDataModuleInterface::efAllowMulticolumn));
+            new MockAnalysisModule(gmx::AnalysisDataModuleInterface::efAllowMulticolumn |
+                                   gmx::AnalysisDataModuleInterface::efAllowMultipoint));
     module->setupStaticCheck(data, source);
     source->addModule(module.release());
 }
@@ -191,7 +193,8 @@ AnalysisDataTestFixture::addStaticColumnCheckerModule(const AnalysisDataTestInpu
                                                       AbstractAnalysisData *source)
 {
     std::auto_ptr<MockAnalysisModule> module(
-            new MockAnalysisModule(gmx::AnalysisDataModuleInterface::efAllowMulticolumn));
+            new MockAnalysisModule(gmx::AnalysisDataModuleInterface::efAllowMulticolumn |
+                                   gmx::AnalysisDataModuleInterface::efAllowMultipoint));
     module->setupStaticColumnCheck(data, firstcol, n, source);
     source->addColumnModule(firstcol, n, module.release());
 }
index f0daad61bd46770587b7a714ddcc7ad5475ca51a..7081ce442d162e408d91b36c2fe560310ff5e051 100644 (file)
@@ -45,6 +45,7 @@
 #include <gtest/gtest.h>
 
 #include "gromacs/legacyheaders/types/simple.h"
+#include "gromacs/fatalerror/gmxassert.h"
 
 #include "testutils/refdata.h"
 
@@ -64,6 +65,35 @@ class MockAnalysisModule;
 const real END_OF_DATA = std::numeric_limits<real>::max();
 //! Constant to use to signify end of one data frame for AnalysisDataTestInput.
 const real END_OF_FRAME = std::numeric_limits<real>::min();
+//! Constant to use to signify end of one multipoint set for AnalysisDataTestInput.
+const real MPSTOP = -std::numeric_limits<real>::max();
+
+/*! \libinternal \brief
+ * Represents a single set of points in AnalysisDataTestInputFrame structure.
+ *
+ * If the data is not multipoint, each frame contains exactly one set of
+ * points.  If there is more than one set of points, each of these sets is
+ * passed separately to AnalysisDataHandle::addPoints().
+ *
+ * \ingroup module_analysisdata
+ */
+class AnalysisDataTestInputPointSet
+{
+    public:
+        AnalysisDataTestInputPointSet();
+
+        int size() const { return y_.size(); }
+        real y(int i) const { return y_[i]; }
+        const std::vector<real> &yvector() const { return y_; }
+        const real *yptr() const { return &y_[0]; }
+        const real *dyptr() const { return NULL; }
+        const bool *presentptr() const { return NULL; }
+
+    private:
+        std::vector<real>       y_;
+
+        friend class AnalysisDataTestInput;
+};
 
 /*! \libinternal \brief
  * Represents a single frame in AnalysisDataTestInput structure.
@@ -75,19 +105,24 @@ class AnalysisDataTestInputFrame
     public:
         AnalysisDataTestInputFrame();
 
+        bool isMultipoint() const { return points_.size() > 1; }
+
         int index() const { return index_; }
         real x() const { return x_; }
         real dx() const { return 0.0; }
-        real y(int i) const { return y_[i]; }
-        const std::vector<real> &yvector() const { return y_; }
-        const real *yptr() const { return &y_[0]; }
-        const real *dyptr() const { return NULL; }
-        const bool *presentptr() const { return NULL; }
+
+        int pointSetCount() const { return points_.size(); }
+        const AnalysisDataTestInputPointSet &points(int index = 0) const
+        {
+            GMX_ASSERT(index >= 0 && static_cast<size_t>(index) < points_.size(),
+                       "Point set index out of range");
+            return points_[index];
+        }
 
     private:
         int                     index_;
         real                    x_;
-        std::vector<real>       y_;
+        std::vector<AnalysisDataTestInputPointSet>  points_;
 
         friend class AnalysisDataTestInput;
 };
@@ -113,6 +148,12 @@ class AnalysisDataTestInput
          * The input array should consist of a set of frames, separated by a
          * END_OF_FRAME marker.  The first value for a frame is the X value,
          * all following values are Y values.
+         * For multipoint data, one frame can contain several point sets,
+         * separated by MPSTOP markers.  There should be no MPSTOP marker after
+         * the last point set, only an END_OF_FRAME marker.  All point sets are
+         * assumed to start from column zero, but the sets may contain
+         * different number of columns.  For non-multipoint data, all frames
+         * must containt the same number of columns.
          * The final element in the array (after the last END_OF_FRAME) should
          * be END_OF_DATA.
          */
@@ -121,10 +162,12 @@ class AnalysisDataTestInput
 
         int frameCount() const { return frames_.size(); }
         int columnCount() const { return columnCount_; }
+        bool isMultipoint() const { return bMultipoint_; }
         const AnalysisDataTestInputFrame &frame(int index) const;
 
     private:
         int                     columnCount_;
+        bool                    bMultipoint_;
         std::vector<AnalysisDataTestInputFrame> frames_;
 };
 
@@ -156,7 +199,7 @@ class AnalysisDataTestInput
  * checks are done during the these methods, by the added mock modules.
  *
  * \todo
- * Support for errors and for multipoint data.
+ * Support for errors and for arbitrary multipoint data.
  *
  * \see AnalysisDataTestInput
  *
index 5e8c6c8fccebd3db2ab05d5d2848c3d8b771893e..507e2b9004e6c150d4ae0159455c075d07acc38c 100644 (file)
@@ -84,23 +84,35 @@ MockAnalysisModule::Impl::checkReferencePoints(real x, real dx, int firstcol, in
 namespace
 {
 
-void checkFrame(real x, real dx, int firstcol, int n, const real *y,
-                const AnalysisDataTestInputFrame &frame)
+void checkFrame(real x, real dx, const AnalysisDataTestInputFrame &frame)
 {
     EXPECT_FLOAT_EQ(frame.x(), x);
     EXPECT_FLOAT_EQ(frame.dx(), dx);
+}
+
+void checkPoints(int firstcol, int n, const real *y,
+                 const AnalysisDataTestInputPointSet &points)
+{
     for (int i = 0; i < n; ++i)
     {
-        EXPECT_FLOAT_EQ(frame.y(firstcol + i), y[i]);
+        EXPECT_FLOAT_EQ(points.y(firstcol + i), y[i]);
     }
 }
 
+void checkFrame(real x, real dx, int firstcol, int n, const real *y,
+                const AnalysisDataTestInputFrame &frame)
+{
+    checkFrame(x, dx, frame);
+    checkPoints(firstcol, n, y, frame.points());
+}
+
 class StaticDataPointsChecker
 {
     public:
         StaticDataPointsChecker(const AnalysisDataTestInputFrame *frame,
+                                const AnalysisDataTestInputPointSet *points,
                                 int firstcol, int n)
-            : frame_(frame), firstcol_(firstcol), n_(n)
+            : frame_(frame), points_(points), firstcol_(firstcol), n_(n)
         {
         }
 
@@ -111,11 +123,13 @@ class StaticDataPointsChecker
             SCOPED_TRACE(formatString("Frame %d", frame_->index()));
             EXPECT_EQ(0, firstcol);
             EXPECT_EQ(n_, n);
-            checkFrame(x, dx, firstcol_ + firstcol, n, y, *frame_);
+            checkFrame(x, dx, *frame_);
+            checkPoints(firstcol_ + firstcol, n, y, *points_);
         }
 
     private:
         const AnalysisDataTestInputFrame *frame_;
+        const AnalysisDataTestInputPointSet *points_;
         int                     firstcol_;
         int                     n_;
 };
@@ -221,9 +235,14 @@ MockAnalysisModule::setupStaticCheck(const AnalysisDataTestInput &data,
     {
         const AnalysisDataTestInputFrame &frame = data.frame(row);
         EXPECT_CALL(*this, frameStarted(frame.x(), frame.dx()));
-        EXPECT_CALL(*this, pointsAdded(frame.x(), frame.dx(), 0,
-                                       data.columnCount(), _, _, _))
-            .WillOnce(Invoke(StaticDataPointsChecker(&frame, 0, data.columnCount())));
+        for (int ps = 0; ps < frame.pointSetCount(); ++ps)
+        {
+            const AnalysisDataTestInputPointSet &points = frame.points(ps);
+            EXPECT_CALL(*this, pointsAdded(frame.x(), frame.dx(), 0,
+                                           points.size(), _, _, _))
+                .WillOnce(Invoke(StaticDataPointsChecker(&frame, &points, 0,
+                                                         data.columnCount())));
+        }
         EXPECT_CALL(*this, frameFinished());
     }
     EXPECT_CALL(*this, dataFinished());
@@ -249,8 +268,12 @@ MockAnalysisModule::setupStaticColumnCheck(const AnalysisDataTestInput &data,
     {
         const AnalysisDataTestInputFrame &frame = data.frame(row);
         EXPECT_CALL(*this, frameStarted(frame.x(), frame.dx()));
-        EXPECT_CALL(*this, pointsAdded(frame.x(), frame.dx(), 0, n, _, _, _))
-            .WillOnce(Invoke(StaticDataPointsChecker(&frame, firstcol, n)));
+        for (int ps = 0; ps < frame.pointSetCount(); ++ps)
+        {
+            const AnalysisDataTestInputPointSet &points = frame.points(ps);
+            EXPECT_CALL(*this, pointsAdded(frame.x(), frame.dx(), 0, n, _, _, _))
+                .WillOnce(Invoke(StaticDataPointsChecker(&frame, &points, firstcol, n)));
+        }
         EXPECT_CALL(*this, frameFinished());
     }
     EXPECT_CALL(*this, dataFinished());
@@ -264,6 +287,8 @@ MockAnalysisModule::setupStaticStorageCheck(const AnalysisDataTestInput &data,
 {
     GMX_RELEASE_ASSERT(data.columnCount() == source->columnCount(),
                        "Mismatching data column count");
+    GMX_RELEASE_ASSERT(!data.isMultipoint() && !source->isMultipoint(),
+                       "Storage testing not implemented for multipoint data");
 
     ::testing::InSequence dummy;
     using ::testing::_;