Support for storing multipoint analysis data.
authorTeemu Murtola <teemu.murtola@gmail.com>
Mon, 10 Jun 2013 18:24:23 +0000 (21:24 +0300)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Wed, 17 Jul 2013 18:49:27 +0000 (20:49 +0200)
Main changes:
- The storage frames now store a vector of "point set info" objects,
  which partition the vector of values into point sets, and annotate
  those sets with extra information (currently, the index of the first
  column).
- The implementation of the frame access wrappers in dataframe.* are
  adapted similarly, and methods are added to allow access to multiple
  point sets within a single frame.
- Added tests for multipoint storage.
- Other changes are adapting to these interface changes.

Still doesn't work perfectly; full support for addColumnModule()
requires reworking the logic in dataproxy.cpp.

Part of #869.

Change-Id: I21141bad83b3f2bb72ece9146aaec02792e6cd0e

12 files changed:
src/gromacs/analysisdata/abstractdata.cpp
src/gromacs/analysisdata/abstractdata.h
src/gromacs/analysisdata/analysisdata.cpp
src/gromacs/analysisdata/analysisdata.h
src/gromacs/analysisdata/arraydata.cpp
src/gromacs/analysisdata/arraydata.h
src/gromacs/analysisdata/dataframe.cpp
src/gromacs/analysisdata/dataframe.h
src/gromacs/analysisdata/datastorage.cpp
src/gromacs/analysisdata/datastorage.h
src/gromacs/analysisdata/tests/analysisdata.cpp
src/testutils/mock_datamodule.cpp

index 11a478fe93fb3c6dbedf773868db388de2be1f2d..2ef54e21f33423ee7629c41be2f8d4eb4c21db2e 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -134,7 +134,10 @@ AbstractAnalysisData::Impl::presentData(AbstractAnalysisData        *data,
             GMX_THROW(APIError("Missing data not supported by a module"));
         }
         module->frameStarted(frame.header());
-        module->pointsAdded(frame.points());
+        for (int j = 0; j < frame.pointSetCount(); ++j)
+        {
+            module->pointsAdded(frame.pointSet(j));
+        }
         module->frameFinished(frame.header());
     }
     if (!bInData_)
index 84bee683de7dafb0422125258ab2a70c35fd584c..ea0ff3041a52e228ffffb75ed06bd31ad173c5b5 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -239,6 +239,11 @@ class AbstractAnalysisData
          * \param     module  Module to add.
          * \throws    APIError in same situations as addModule().
          *
+         * \todo
+         * This method doesn't currently work in all cases with multipoint
+         * data.  In particular, if the added module requests storage and uses
+         * getDataFrame(), it will behave unpredictably (most likely asserts).
+         *
          * \see addModule()
          */
         void addColumnModule(int col, int span, AnalysisDataModulePointer module);
index cf6d903261450231b777cc1e001ba4c42d1b979d..2ec211576956197d42a54203e8e5aacbf1c32937 100644 (file)
@@ -140,7 +140,6 @@ AnalysisData::setMultipoint(bool bMultipoint)
     GMX_RELEASE_ASSERT(impl_->handles_.empty(),
                        "Cannot change data type after creating handles");
     AbstractAnalysisData::setMultipoint(bMultipoint);
-    impl_->storage_.setMultipoint(bMultipoint);
 }
 
 
@@ -155,10 +154,6 @@ AnalysisData::startData(const AnalysisDataParallelOptions &opt)
         impl_->storage_.setParallelOptions(opt);
         impl_->storage_.startDataStorage(this);
     }
-    else if (isMultipoint())
-    {
-        GMX_THROW(NotImplementedError("Parallelism not supported for multipoint data"));
-    }
 
     Impl::HandlePointer handle(new internal::AnalysisDataHandleImpl(this));
     impl_->handles_.push_back(move(handle));
index c5ca4dac221e90fe85998b0076b5a1e24d6d1d22..bb5f570e0cd878305d5289064cba18c5a2918d57 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -68,9 +68,6 @@ class AnalysisDataParallelOptions;
  * expect them.
  *
  * \todo
- * Currently, multiple handles with multipoint data are not implemented.
- *
- * \todo
  * Parallel implementation is not complete.
  *
  * \if internal
index 2f4a16fbe6f1e473a62ddae266caabbacda5e0ef..4e5b5f3cdc6786f5af4971adb0b574c2d00d6cd6 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -51,7 +51,8 @@ namespace gmx
 {
 
 AbstractAnalysisArrayData::AbstractAnalysisArrayData()
-    : rowCount_(0), xstart_(0.0), xstep_(1.0), bReady_(false)
+    : rowCount_(0), pointSetInfo_(0, 0, 0), xstart_(0.0), xstep_(1.0),
+      bReady_(false)
 {
 }
 
@@ -71,7 +72,8 @@ AbstractAnalysisArrayData::tryGetDataFrameInternal(int index) const
         = value_.begin() + index * columnCount();
     return AnalysisDataFrameRef(
             AnalysisDataFrameHeader(index, xvalue(index), 0.0),
-            AnalysisDataValuesRef(begin, begin + columnCount()));
+            AnalysisDataValuesRef(begin, begin + columnCount()),
+            AnalysisDataPointSetInfosRef(&pointSetInfo_, 1));
 }
 
 
@@ -88,6 +90,7 @@ AbstractAnalysisArrayData::setColumnCount(int ncols)
     GMX_RELEASE_ASSERT(!isAllocated(),
                        "Cannot change column count after data has been allocated");
     AbstractAnalysisData::setColumnCount(ncols);
+    pointSetInfo_ = AnalysisDataPointSetInfo(0, ncols, 0);
 }
 
 
@@ -141,9 +144,11 @@ AbstractAnalysisArrayData::valuesReady()
     {
         AnalysisDataFrameHeader header(i, xvalue(i), 0);
         notifyFrameStart(header);
-        notifyPointsAdd(AnalysisDataPointSetRef(header, 0,
-                                                AnalysisDataValuesRef(valueIter,
-                                                                      valueIter + columnCount())));
+        notifyPointsAdd(
+                AnalysisDataPointSetRef(
+                        header, pointSetInfo_,
+                        AnalysisDataValuesRef(valueIter,
+                                              valueIter + columnCount())));
         notifyFrameFinish(header);
     }
     notifyDataFinish();
index df30162ac71ecbad62adb2be6cfd80ec3dfedc42..be0a8785700bae3f25aed5ef1186c91d38c03ecb 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -203,6 +203,7 @@ class AbstractAnalysisArrayData : public AbstractAnalysisData
         virtual bool requestStorageInternal(int nframes);
 
         int                            rowCount_;
+        AnalysisDataPointSetInfo       pointSetInfo_;
         std::vector<AnalysisDataValue> value_;
         real                           xstart_;
         real                           xstep_;
index 543ae428310b58d7e1ed349a335c9b028e9ef20c..51663a71d7f61dc74eea1c6e018e88d8080d792d 100644 (file)
@@ -68,13 +68,15 @@ AnalysisDataFrameHeader::AnalysisDataFrameHeader(int index, real x, real dx)
  */
 
 AnalysisDataPointSetRef::AnalysisDataPointSetRef(
-        const AnalysisDataFrameHeader &header, int firstColumn,
-        const AnalysisDataValuesRef &values)
-    : header_(header), firstColumn_(firstColumn), values_(values)
+        const AnalysisDataFrameHeader  &header,
+        const AnalysisDataPointSetInfo &pointSetInfo,
+        const AnalysisDataValuesRef    &values)
+    : header_(header), firstColumn_(pointSetInfo.firstColumn()),
+      values_(&*values.begin() + pointSetInfo.valueOffset(),
+              pointSetInfo.valueCount())
 {
     GMX_ASSERT(header_.isValid(),
                "Invalid point set reference should not be constructed");
-    GMX_ASSERT(firstColumn >= 0, "Invalid first column");
 }
 
 
@@ -147,28 +149,40 @@ AnalysisDataFrameRef::AnalysisDataFrameRef()
 
 
 AnalysisDataFrameRef::AnalysisDataFrameRef(
-        const AnalysisDataFrameHeader &header,
-        const AnalysisDataValuesRef   &values)
-    : header_(header), values_(values)
+        const AnalysisDataFrameHeader      &header,
+        const AnalysisDataValuesRef        &values,
+        const AnalysisDataPointSetInfosRef &pointSets)
+    : header_(header), values_(values), pointSets_(pointSets)
 {
+    GMX_ASSERT(!pointSets_.empty(), "There must always be a point set");
 }
 
 
 AnalysisDataFrameRef::AnalysisDataFrameRef(
-        const AnalysisDataFrameHeader        &header,
-        const std::vector<AnalysisDataValue> &values)
-    : header_(header), values_(values.begin(), values.end())
+        const AnalysisDataFrameHeader               &header,
+        const std::vector<AnalysisDataValue>        &values,
+        const std::vector<AnalysisDataPointSetInfo> &pointSets)
+    : header_(header), values_(values.begin(), values.end()),
+      pointSets_(pointSets.begin(), pointSets.end())
 {
+    GMX_ASSERT(!pointSets_.empty(), "There must always be a point set");
 }
 
 
 AnalysisDataFrameRef::AnalysisDataFrameRef(
         const AnalysisDataFrameRef &frame, int firstColumn, int columnCount)
-    : header_(frame.header()), values_(&frame.values_[firstColumn], columnCount)
+    : header_(frame.header()),
+      values_(&frame.values_[firstColumn], columnCount),
+      pointSets_(frame.pointSets_)
 {
+    // FIXME: This doesn't produce a valid internal state, although it does
+    // work in some cases. The point sets cannot be correctly managed here, but
+    // need to be handles by the data proxy class.
     GMX_ASSERT(firstColumn >= 0, "Invalid first column");
     GMX_ASSERT(columnCount >= 0, "Invalid column count");
-    GMX_ASSERT(firstColumn + columnCount <= frame.columnCount(),
+    GMX_ASSERT(pointSets_.size() == 1U,
+               "Subsets of frames only supported with simple data");
+    GMX_ASSERT(firstColumn + columnCount <= static_cast<int>(values_.size()),
                "Invalid last column");
 }
 
index 02cf100267e519c1e305942017d983d3b82cb7fb..1205fe98e0c00100206c41bd41c87ececa8b5e9e 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2012, by the GROMACS development team, led by
+ * Copyright (c) 2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -264,6 +264,67 @@ class AnalysisDataFrameHeader
 };
 
 
+/*! \cond libinternal */
+/*! \libinternal \brief
+ * Value type for internal indexing of point sets.
+ *
+ * This class contains the necessary data to split an array of
+ * AnalysisDataValue objects into point sets.  It is always specified in the
+ * context of an array of AnalysisDataValues: the point set specified by this
+ * class contains valueCount() values, starting from the array index
+ * valueOffset().
+ * The value at location valueOffset() corresponds to column firstColumn().
+ * It is not necessary for code using the analysis data framework to know of
+ * this class, but it is declared in a public header to allow using it in other
+ * types.
+ *
+ * Default copy constructor and assignment operator are used and work as
+ * intended.
+ * Typically new objects of this type are only constructed internally by the
+ * library and in classes that are derived from AbstractAnalysisData.
+ *
+ * Methods in this class do not throw, but may contain asserts for incorrect
+ * usage.
+ *
+ * Note that it is not possible to change the contents of an initialized
+ * object, except by assigning a new object to replace it completely.
+ *
+ * \inlibraryapi
+ * \ingroup module_analysisdata
+ */
+class AnalysisDataPointSetInfo
+{
+    public:
+        //! Construct point set data object with the given values.
+        AnalysisDataPointSetInfo(int valueOffset, int valueCount,
+                                 int firstColumn)
+            : valueOffset_(valueOffset), valueCount_(valueCount),
+              firstColumn_(firstColumn)
+        {
+            GMX_ASSERT(valueOffset >= 0, "Negative value offsets are invalid");
+            GMX_ASSERT(valueCount  >= 0, "Negative value counts are invalid");
+            GMX_ASSERT(firstColumn >= 0, "Negative column indices are invalid");
+        }
+
+        //! Returns the offset of the first value in the referenced value array.
+        int valueOffset() const { return valueOffset_; }
+        //! Returns the number of values in this point set.
+        int valueCount() const { return valueCount_; }
+        //! Returns the index of the first column in this point set.
+        int firstColumn() const { return firstColumn_; }
+
+    private:
+        int                     valueOffset_;
+        int                     valueCount_;
+        int                     firstColumn_;
+};
+
+//! Shorthand for reference to an array of point set data objects.
+typedef ConstArrayRef<AnalysisDataPointSetInfo> AnalysisDataPointSetInfosRef;
+
+//! \endcond
+
+
 /*! \brief
  * Value type wrapper for non-mutable access to a set of data column values.
  *
@@ -290,16 +351,16 @@ class AnalysisDataPointSetRef
         /*! \brief
          * Constructs a point set reference from given values.
          *
-         * \param[in] header      Header for the frame.
-         * \param[in] firstColumn Zero-based index of the first column.
-         *     Must be >= 0.
-         * \param[in] values      Values for each column.
+         * \param[in] header       Header for the frame.
+         * \param[in] pointSetInfo Information about the point set.
+         * \param[in] values       Values for each column.
          *
-         * The first element in \p values should correspond to \p firstColumn.
+         * The first element of the point set should be found from \p values
+         * using the offset in \p pointSetInfo.
          */
-        AnalysisDataPointSetRef(const AnalysisDataFrameHeader &header,
-                                int                            firstColumn,
-                                const AnalysisDataValuesRef   &values);
+        AnalysisDataPointSetRef(const AnalysisDataFrameHeader  &header,
+                                const AnalysisDataPointSetInfo &pointSetInfo,
+                                const AnalysisDataValuesRef    &values);
         /*! \brief
          * Constructs a point set reference from given values.
          *
@@ -445,9 +506,6 @@ class AnalysisDataPointSetRef
  * Note that it is not possible to change the contents of an initialized
  * object, except by assigning a new object to replace it completely.
  *
- * \todo
- * Support for multipoint data.
- *
  * \inpublicapi
  * \ingroup module_analysisdata
  */
@@ -466,17 +524,21 @@ class AnalysisDataFrameRef
          *
          * \param[in] header      Header for the frame.
          * \param[in] values      Values for each column.
+         * \param[in] pointSets   Point set data.
          */
-        AnalysisDataFrameRef(const AnalysisDataFrameHeader &header,
-                             const AnalysisDataValuesRef   &values);
+        AnalysisDataFrameRef(const AnalysisDataFrameHeader      &header,
+                             const AnalysisDataValuesRef        &values,
+                             const AnalysisDataPointSetInfosRef &pointSets);
         /*! \brief
          * Constructs a frame reference from given values.
          *
          * \param[in] header      Header for the frame.
          * \param[in] values      Values for each column.
+         * \param[in] pointSets   Point set data.
          */
-        AnalysisDataFrameRef(const AnalysisDataFrameHeader        &header,
-                             const std::vector<AnalysisDataValue> &values);
+        AnalysisDataFrameRef(const AnalysisDataFrameHeader               &header,
+                             const std::vector<AnalysisDataValue>        &values,
+                             const std::vector<AnalysisDataPointSetInfo> &pointSets);
         /*! \brief
          * Constructs a frame reference to a subset of columns.
          *
@@ -523,63 +585,54 @@ class AnalysisDataFrameRef
             return header().dx();
         }
         /*! \brief
-         * Returns point set reference to the column values of this frame.
-         *
-         * Should not be called for invalid frames.
-         */
-        AnalysisDataPointSetRef points() const
-        {
-            GMX_ASSERT(isValid(), "Invalid data frame accessed");
-            return AnalysisDataPointSetRef(header_, 0, values_);
-        }
-        /*! \brief
-         * Returns number of columns in this frame.
+         * Returns the number of point sets for this frame.
          *
          * Returns zero for an invalid frame.
          */
-        int columnCount() const
+        int pointSetCount() const
         {
-            return values_.size();
+            return pointSets_.size();
         }
         /*! \brief
-         * Returns reference container for all column values.
+         * Returns point set reference for a given point set.
+         *
+         * Should not be called for invalid frames.
          */
-        const AnalysisDataValuesRef &values() const
+        AnalysisDataPointSetRef pointSet(int index) const
         {
-            return values_;
+            GMX_ASSERT(isValid(), "Invalid data frame accessed");
+            GMX_ASSERT(index >= 0 && index < pointSetCount(),
+                       "Out of range data access");
+            return AnalysisDataPointSetRef(header_, pointSets_[index], values_);
         }
         /*! \brief
-         * Convenience method for accessing a column value.
+         * Convenience method for accessing a column value in simple data.
          *
          * \copydetails AnalysisDataPointSetRef::y()
          */
         real y(int i) const
         {
-            GMX_ASSERT(isValid(), "Invalid data frame accessed");
-            GMX_ASSERT(i >= 0 && i < columnCount(), "Out of range data access");
-            return values_[i].value();
+            return singleColumnValue(i).value();
         }
         /*! \brief
-         * Convenience method for accessing error for a column value.
+         * Convenience method for accessing error for a column value in simple
+         * data.
          *
          * \copydetails AnalysisDataPointSetRef::dy()
          */
         real dy(int i) const
         {
-            GMX_ASSERT(isValid(), "Invalid data frame accessed");
-            GMX_ASSERT(i >= 0 && i < columnCount(), "Out of range data access");
-            return values_[i].error();
+            return singleColumnValue(i).error();
         }
         /*! \brief
-         * Convenience method for accessing present status for a column.
+         * Convenience method for accessing present status for a column in
+         * simple data.
          *
          * \copydetails AnalysisDataPointSetRef::present()
          */
         bool present(int i) const
         {
-            GMX_ASSERT(isValid(), "Invalid data frame accessed");
-            GMX_ASSERT(i >= 0 && i < columnCount(), "Out of range data access");
-            return values_[i].isPresent();
+            return singleColumnValue(i).isPresent();
         }
         /*! \brief
          * Returns true if all points in this frame are present.
@@ -587,8 +640,20 @@ class AnalysisDataFrameRef
         bool allPresent() const;
 
     private:
-        AnalysisDataFrameHeader header_;
-        AnalysisDataValuesRef   values_;
+        //! Helper method for accessing single columns in simple data.
+        const AnalysisDataValue &singleColumnValue(int i) const
+        {
+            GMX_ASSERT(isValid(), "Invalid data frame accessed");
+            GMX_ASSERT(pointSets_.size() == 1U && pointSets_[0].firstColumn() == 0,
+                       "Convenience method not available for multiple point sets");
+            GMX_ASSERT(i >= 0 && i < static_cast<int>(values_.size()),
+                       "Out of range data access");
+            return values_[i];
+        }
+
+        AnalysisDataFrameHeader      header_;
+        AnalysisDataValuesRef        values_;
+        AnalysisDataPointSetInfosRef pointSets_;
 };
 
 } // namespace gmx
index 0c53e61a3bbe10fafba22bef7e9a79bb328aa2c0..1b80d95430284930e39c0122d346400b84ad63a4 100644 (file)
@@ -41,6 +41,8 @@
  */
 #include "datastorage.h"
 
+#include <algorithm>
+#include <iterator>
 #include <limits>
 #include <vector>
 
@@ -165,6 +167,18 @@ class AnalysisDataStorage::Impl
          */
         FrameBuilderPointer getFrameBuilder();
 
+        /*! \brief
+         * Returns whether notifications should be immediately fired.
+         *
+         * This is used to optimize multipoint handling for non-parallel cases,
+         * where it is not necessary to store even a single frame.
+         *
+         * Does not throw.
+         */
+        bool shouldNotifyImmediately() const
+        {
+            return isMultipoint() && storageLimit_ == 0 && pendingLimit_ == 1;
+        }
         /*! \brief
          * Calls notification method in \a data_.
          *
@@ -190,14 +204,6 @@ class AnalysisDataStorage::Impl
 
         //! Data object to use for notification calls.
         AbstractAnalysisData   *data_;
-        /*! \brief
-         * Whether the storage has been set to allow multipoint.
-         *
-         * Should be possible to remove once full support for multipoint data
-         * has been implemented;  isMultipoint() can simply return
-         * \c data_->isMultipoint() in that case.
-         */
-        bool                    bMultipoint_;
         /*! \brief
          * Number of past frames that need to be stored.
          *
@@ -321,6 +327,8 @@ class AnalysisDataStorageFrameData
         const AnalysisDataFrameHeader &header() const { return header_; }
         //! Returns zero-based index of the frame.
         int frameIndex() const { return header().index(); }
+        //! Returns the number of point sets for the frame.
+        int pointSetCount() const { return pointSets_.size(); }
 
         //! Clears the frame for reusing as a new frame.
         void clearFrame(int newIndex);
@@ -338,21 +346,25 @@ class AnalysisDataStorageFrameData
             GMX_ASSERT(builder_, "Accessing builder for not-in-progress frame");
             return *builder_;
         }
+        /*! \brief
+         * Adds a new point set to this frame.
+         */
+        void addPointSet(int firstColumn, ValueIterator begin, ValueIterator end);
         /*! \brief
          * Finalizes the frame during AnalysisDataStorage::finishFrame().
          *
          * \returns The builder object used by the frame, for reusing it for
          *      other frames.
          */
-        AnalysisDataFrameBuilderPointer finishFrame();
+        AnalysisDataFrameBuilderPointer finishFrame(bool bMultipoint);
 
         //! Returns frame reference to this frame.
         AnalysisDataFrameRef frameReference() const
         {
-            return AnalysisDataFrameRef(header_, values_);
+            return AnalysisDataFrameRef(header_, values_, pointSets_);
         }
-        //! Returns point set reference to currently set values.
-        AnalysisDataPointSetRef currentPoints() const;
+        //! Returns point set reference to a given point set.
+        AnalysisDataPointSetRef pointSet(int index) const;
 
     private:
         //! Storage object that contains this frame.
@@ -361,6 +373,8 @@ class AnalysisDataStorageFrameData
         AnalysisDataFrameHeader                 header_;
         //! Values for the frame.
         std::vector<AnalysisDataValue>          values_;
+        //! Information about each point set in the frame.
+        std::vector<AnalysisDataPointSetInfo>   pointSets_;
         /*! \brief
          * Builder object for the frame.
          *
@@ -372,8 +386,6 @@ class AnalysisDataStorageFrameData
         Status                                  status_;
 
         GMX_DISALLOW_COPY_AND_ASSIGN(AnalysisDataStorageFrameData);
-
-        friend class gmx::AnalysisDataStorageFrame;
 };
 
 }   // namespace internal
@@ -383,7 +395,7 @@ class AnalysisDataStorageFrameData
  */
 
 AnalysisDataStorage::Impl::Impl()
-    : data_(NULL), bMultipoint_(false),
+    : data_(NULL),
       storageLimit_(0), pendingLimit_(1), firstFrameLocation_(0), nextIndex_(0)
 {
 }
@@ -400,7 +412,8 @@ AnalysisDataStorage::Impl::columnCount() const
 bool
 AnalysisDataStorage::Impl::isMultipoint() const
 {
-    return bMultipoint_;
+    GMX_ASSERT(data_ != NULL, "isMultipoint() called too early");
+    return data_->isMultipoint();
 }
 
 
@@ -517,7 +530,10 @@ AnalysisDataStorage::Impl::notifyNextFrames(size_t firstLocation)
         if (!storedFrame.isNotified())
         {
             data_->notifyFrameStart(storedFrame.header());
-            data_->notifyPointsAdd(storedFrame.currentPoints());
+            for (int j = 0; j < storedFrame.pointSetCount(); ++j)
+            {
+                data_->notifyPointsAdd(storedFrame.pointSet(j));
+            }
             data_->notifyFrameFinish(storedFrame.header());
             storedFrame.markNotified();
             if (storedFrame.frameIndex() >= storageLimit_)
@@ -546,10 +562,9 @@ AnalysisDataStorage::Impl::finishFrame(int index)
                        "finishFrame() called twice for the same frame");
     GMX_RELEASE_ASSERT(storedFrame.frameIndex() == index,
                        "Inconsistent internal frame indexing");
-    builders_.push_back(storedFrame.finishFrame());
-    if (isMultipoint())
+    builders_.push_back(storedFrame.finishFrame(isMultipoint()));
+    if (shouldNotifyImmediately())
     {
-        // TODO: Check that the last point set has been finished
         data_->notifyFrameFinish(storedFrame.header());
         if (storedFrame.frameIndex() >= storageLimit_)
         {
@@ -575,6 +590,15 @@ AnalysisDataStorageFrameData::AnalysisDataStorageFrameData(
         int                        index)
     : storageImpl_(*storageImpl), header_(index, 0.0, 0.0), status_(eMissing)
 {
+    GMX_RELEASE_ASSERT(storageImpl->data_ != NULL,
+                       "Storage frame constructed before data started");
+    // With non-multipoint data, the point set structure is static,
+    // so initialize it only once here.
+    if (!baseData().isMultipoint())
+    {
+        int columnCount = baseData().columnCount();
+        pointSets_.push_back(AnalysisDataPointSetInfo(0, columnCount, 0));
+    }
 }
 
 
@@ -585,6 +609,10 @@ AnalysisDataStorageFrameData::clearFrame(int newIndex)
     status_ = eMissing;
     header_ = AnalysisDataFrameHeader(newIndex, 0.0, 0.0);
     values_.clear();
+    if (baseData().isMultipoint())
+    {
+        pointSets_.clear();
+    }
 }
 
 
@@ -600,12 +628,43 @@ AnalysisDataStorageFrameData::startFrame(
 }
 
 
+void
+AnalysisDataStorageFrameData::addPointSet(int firstColumn,
+                                          ValueIterator begin, ValueIterator end)
+{
+    const int valueCount  = end - begin;
+    if (storageImpl().shouldNotifyImmediately())
+    {
+        AnalysisDataPointSetInfo pointSetInfo(0, valueCount, firstColumn);
+        storageImpl().notifyPointSet(
+                AnalysisDataPointSetRef(header(), pointSetInfo,
+                                        AnalysisDataValuesRef(begin, end)));
+    }
+    else
+    {
+        pointSets_.push_back(
+                AnalysisDataPointSetInfo(values_.size(), valueCount, firstColumn));
+        std::copy(begin, end, std::back_inserter(values_));
+    }
+}
+
+
 AnalysisDataFrameBuilderPointer
-AnalysisDataStorageFrameData::finishFrame()
+AnalysisDataStorageFrameData::finishFrame(bool bMultipoint)
 {
     status_ = eFinished;
-    values_ = builder_->values_;
-    builder_->clearValues();
+    if (!bMultipoint)
+    {
+        GMX_RELEASE_ASSERT(pointSets_.size() == 1U,
+                           "Point sets created for non-multipoint data");
+        values_ = builder_->values_;
+        builder_->clearValues();
+    }
+    else
+    {
+        GMX_RELEASE_ASSERT(!builder_->bPointSetInProgress_,
+                           "Unfinished point set");
+    }
     AnalysisDataFrameBuilderPointer builder(move(builder_));
     builder_.reset();
     return move(builder);
@@ -613,21 +672,13 @@ AnalysisDataStorageFrameData::finishFrame()
 
 
 AnalysisDataPointSetRef
-AnalysisDataStorageFrameData::currentPoints() const
+AnalysisDataStorageFrameData::pointSet(int index) const
 {
-    std::vector<AnalysisDataValue>::const_iterator begin = values_.begin();
-    std::vector<AnalysisDataValue>::const_iterator end   = values_.end();
-    while (begin != end && !begin->isSet())
-    {
-        ++begin;
-    }
-    while (end != begin && !(end-1)->isSet())
-    {
-        --end;
-    }
-    int firstColumn = (begin != end) ? begin - values_.begin() : 0;
-    return AnalysisDataPointSetRef(header_, firstColumn,
-                                   AnalysisDataValuesRef(begin, end));
+    GMX_ASSERT(index >= 0 && index < pointSetCount(),
+               "Invalid point set index");
+    return AnalysisDataPointSetRef(
+            header_, pointSets_[index],
+            AnalysisDataValuesRef(values_.begin(), values_.end()));
 }
 
 }   // namespace internal
@@ -637,7 +688,7 @@ AnalysisDataStorageFrameData::currentPoints() const
  */
 
 AnalysisDataStorageFrame::AnalysisDataStorageFrame(int columnCount)
-    : data_(NULL), values_(columnCount)
+    : data_(NULL), values_(columnCount), bPointSetInProgress_(false)
 {
 }
 
@@ -650,11 +701,15 @@ AnalysisDataStorageFrame::~AnalysisDataStorageFrame()
 void
 AnalysisDataStorageFrame::clearValues()
 {
-    std::vector<AnalysisDataValue>::iterator i;
-    for (i = values_.begin(); i != values_.end(); ++i)
+    if (bPointSetInProgress_)
     {
-        i->clear();
+        std::vector<AnalysisDataValue>::iterator i;
+        for (i = values_.begin(); i != values_.end(); ++i)
+        {
+            i->clear();
+        }
     }
+    bPointSetInProgress_ = false;
 }
 
 
@@ -664,8 +719,21 @@ AnalysisDataStorageFrame::finishPointSet()
     GMX_RELEASE_ASSERT(data_ != NULL, "Invalid frame accessed");
     GMX_RELEASE_ASSERT(data_->baseData().isMultipoint(),
                        "Should not be called for non-multipoint data");
-    data_->values_ = values_;
-    data_->storageImpl().notifyPointSet(data_->currentPoints());
+    if (bPointSetInProgress_)
+    {
+        std::vector<AnalysisDataValue>::const_iterator begin = values_.begin();
+        std::vector<AnalysisDataValue>::const_iterator end   = values_.end();
+        while (begin != end && !begin->isSet())
+        {
+            ++begin;
+        }
+        while (end != begin && !(end-1)->isSet())
+        {
+            --end;
+        }
+        int firstColumn = (begin != end) ? begin - values_.begin() : 0;
+        data_->addPointSet(firstColumn, begin, end);
+    }
     clearValues();
 }
 
@@ -693,17 +761,6 @@ AnalysisDataStorage::~AnalysisDataStorage()
 }
 
 
-void
-AnalysisDataStorage::setMultipoint(bool bMultipoint)
-{
-    if (bMultipoint && impl_->storageLimit_ > 0)
-    {
-        GMX_THROW(APIError("Storage of multipoint data not supported"));
-    }
-    impl_->bMultipoint_ = bMultipoint;
-}
-
-
 void
 AnalysisDataStorage::setParallelOptions(const AnalysisDataParallelOptions &opt)
 {
@@ -714,10 +771,6 @@ AnalysisDataStorage::setParallelOptions(const AnalysisDataParallelOptions &opt)
 AnalysisDataFrameRef
 AnalysisDataStorage::tryGetDataFrame(int index) const
 {
-    if (impl_->isMultipoint())
-    {
-        return AnalysisDataFrameRef();
-    }
     int storageIndex = impl_->computeStorageLocation(index);
     if (storageIndex == -1)
     {
@@ -735,11 +788,6 @@ AnalysisDataStorage::tryGetDataFrame(int index) const
 bool
 AnalysisDataStorage::requestStorage(int nframes)
 {
-    if (impl_->isMultipoint())
-    {
-        return false;
-    }
-
     // Handle the case when everything needs to be stored.
     if (nframes == -1)
     {
@@ -761,7 +809,6 @@ AnalysisDataStorage::startDataStorage(AbstractAnalysisData *data)
 {
     // Data needs to be set before calling extendBuffer()
     impl_->data_ = data;
-    setMultipoint(data->isMultipoint());
     if (!impl_->storeAll())
     {
         impl_->extendBuffer(impl_->storageLimit_ + impl_->pendingLimit_ + 1);
@@ -797,7 +844,7 @@ AnalysisDataStorage::startFrame(const AnalysisDataFrameHeader &header)
     GMX_RELEASE_ASSERT(storedFrame->frameIndex() == header.index(),
                        "Inconsistent internal frame indexing");
     storedFrame->startFrame(header, impl_->getFrameBuilder());
-    if (impl_->isMultipoint())
+    if (impl_->shouldNotifyImmediately())
     {
         impl_->data_->notifyFrameStart(header);
     }
index 5f1cb6957389632435b62a58087274e87934e124..650a613b0fd11f0485958aec2782c727e3643ba4 100644 (file)
@@ -107,6 +107,7 @@ class AnalysisDataStorageFrame
             GMX_ASSERT(column >= 0 && column < columnCount(),
                        "Invalid column index");
             values_[column].setValue(value, bPresent);
+            bPointSetInProgress_ = true;
         }
         /*! \brief
          * Sets value for a column.
@@ -126,6 +127,7 @@ class AnalysisDataStorageFrame
             GMX_ASSERT(column >= 0 && column < columnCount(),
                        "Invalid column index");
             values_[column].setValue(value, error, bPresent);
+            bPointSetInProgress_ = true;
         }
         /*! \brief
          * Access value for a column.
@@ -203,6 +205,9 @@ class AnalysisDataStorageFrame
         //! Values for the currently in-progress point set.
         std::vector<AnalysisDataValue>          values_;
 
+        //! Whether any values have been set in the current point set.
+        bool                                    bPointSetInProgress_;
+
         //! Needed for access to the constructor.
         friend class AnalysisDataStorage;
         //! Needed for managing the frame the object points to.
@@ -229,11 +234,6 @@ class AnalysisDataStorageFrame
  * AbstractAnalysisData::notifyFrameFinish() appropriately.
  *
  * \todo
- * Full support for multipoint data.
- * Currently, multipoint data is only supported in serial pass-through mode
- * without any storage.
- *
- * \todo
  * Proper multi-threaded implementation.
  *
  * \inlibraryapi
@@ -246,18 +246,6 @@ class AnalysisDataStorage
         AnalysisDataStorage();
         ~AnalysisDataStorage();
 
-        /*! \brief
-         * Sets whether the data will be multipoint.
-         *
-         * \exception APIError if storage has been requested.
-         *
-         * It is not mandatory to call this method (the multipoint property is
-         * automatically detected in startDataStorage()), but currently calling
-         * it will produce better diagnostic messages.
-         * When full support for multipoint data has been implemented, this
-         * method can be removed.
-         */
-        void setMultipoint(bool bMultipoint);
         /*! \brief
          * Set parallelization options for the storage.
          *
@@ -290,9 +278,6 @@ class AnalysisDataStorage
          * AbstractAnalysisData::requestStorageInternal() can be directly
          * forwarded to this method.  See that method for more documentation.
          *
-         * Currently, multipoint data cannot be stored, but all other storage
-         * request will be honored.
-         *
          * \see AbstractAnalysisData::requestStorageInternal()
          */
         bool requestStorage(int nframes);
@@ -302,8 +287,6 @@ class AnalysisDataStorage
          *
          * \param  data  AbstractAnalysisData object containing this storage.
          * \exception std::bad_alloc if storage allocation fails.
-         * \exception APIError if storage has been requested and \p data is
-         *      multipoint.
          *
          * Lifetime of \p data must exceed the lifetime of the storage object
          * (typically, the storage object will be a member in \p data).
index 5143dd5251c3be7bece7abfeae3b30d2cd7692bc..008dce536e342675380b8ba92891817350cf8ee3 100644 (file)
@@ -211,6 +211,11 @@ class AnalysisDataTest : public AnalysisDataTestFixture
             AnalysisDataTestFixture::addStaticColumnCheckerModule(
                     input_, firstColumn, columnCount, &data_);
         }
+        void addStaticStorageCheckerModule(int storageCount)
+        {
+            AnalysisDataTestFixture::addStaticStorageCheckerModule(
+                    input_, storageCount, &data_);
+        }
         void presentAllData()
         {
             AnalysisDataTestFixture::presentAllData(input_, &data_);
@@ -264,18 +269,18 @@ TYPED_TEST(AnalysisDataCommonTest, CallsColumnModuleCorrectly)
  * Tests that data is forwarded correctly (in frame order) to modules when the
  * data is added through multiple handles in non-increasing order.
  */
-TEST_F(AnalysisDataSimpleTest, CallsModuleCorrectlyWithOutOfOrderFrames)
+TYPED_TEST(AnalysisDataCommonTest, CallsModuleCorrectlyWithOutOfOrderFrames)
 {
-    ASSERT_NO_THROW_GMX(addStaticCheckerModule());
-    ASSERT_NO_THROW_GMX(addStaticColumnCheckerModule(1, 2));
+    ASSERT_NO_THROW_GMX(AnalysisDataTest::addStaticCheckerModule());
+    ASSERT_NO_THROW_GMX(AnalysisDataTest::addStaticColumnCheckerModule(1, 2));
     gmx::AnalysisDataHandle          handle1;
     gmx::AnalysisDataHandle          handle2;
     gmx::AnalysisDataParallelOptions options(2);
-    ASSERT_NO_THROW_GMX(handle1 = data_.startData(options));
-    ASSERT_NO_THROW_GMX(handle2 = data_.startData(options));
-    ASSERT_NO_THROW_GMX(presentDataFrame(input_, 1, handle1));
-    ASSERT_NO_THROW_GMX(presentDataFrame(input_, 0, handle2));
-    ASSERT_NO_THROW_GMX(presentDataFrame(input_, 2, handle1));
+    ASSERT_NO_THROW_GMX(handle1 = this->data_.startData(options));
+    ASSERT_NO_THROW_GMX(handle2 = this->data_.startData(options));
+    ASSERT_NO_THROW_GMX(AnalysisDataTest::presentDataFrame(this->input_, 1, handle1));
+    ASSERT_NO_THROW_GMX(AnalysisDataTest::presentDataFrame(this->input_, 0, handle2));
+    ASSERT_NO_THROW_GMX(AnalysisDataTest::presentDataFrame(this->input_, 2, handle1));
     ASSERT_NO_THROW_GMX(handle1.finishData());
     ASSERT_NO_THROW_GMX(handle2.finishData());
 }
@@ -284,32 +289,32 @@ TEST_F(AnalysisDataSimpleTest, CallsModuleCorrectlyWithOutOfOrderFrames)
  * Tests that data can be accessed correctly from a module that requests
  * storage using AbstractAnalysisData::requestStorage() with parameter -1.
  */
-TEST_F(AnalysisDataSimpleTest, FullStorageWorks)
+TYPED_TEST(AnalysisDataCommonTest, FullStorageWorks)
 {
-    ASSERT_NO_THROW_GMX(addStaticStorageCheckerModule(input_, -1, &data_));
-    ASSERT_NO_THROW_GMX(presentAllData());
+    ASSERT_NO_THROW_GMX(AnalysisDataTest::addStaticStorageCheckerModule(-1));
+    ASSERT_NO_THROW_GMX(AnalysisDataTest::presentAllData());
 }
 
 /*
  * Tests that a data module can be added to an AnalysisData object after data
  * has been added if all data is still available in storage.
  */
-TEST_F(AnalysisDataSimpleTest, CanAddModuleAfterStoredData)
+TYPED_TEST(AnalysisDataCommonTest, CanAddModuleAfterStoredData)
 {
-    ASSERT_TRUE(data_.requestStorage(-1));
+    ASSERT_TRUE(this->data_.requestStorage(-1));
 
-    ASSERT_NO_THROW_GMX(presentAllData());
-    ASSERT_NO_THROW_GMX(addStaticCheckerModule());
+    ASSERT_NO_THROW_GMX(AnalysisDataTest::presentAllData());
+    ASSERT_NO_THROW_GMX(AnalysisDataTest::addStaticCheckerModule());
 }
 
 /*
  * Tests that data can be accessed correctly from a module that requests
  * storage using AbstractAnalysisData::requestStorage() only for one frame.
  */
-TEST_F(AnalysisDataSimpleTest, LimitedStorageWorks)
+TYPED_TEST(AnalysisDataCommonTest, LimitedStorageWorks)
 {
-    ASSERT_NO_THROW_GMX(addStaticStorageCheckerModule(input_, 1, &data_));
-    ASSERT_NO_THROW_GMX(presentAllData());
+    ASSERT_NO_THROW_GMX(AnalysisDataTest::addStaticStorageCheckerModule(1));
+    ASSERT_NO_THROW_GMX(AnalysisDataTest::presentAllData());
 }
 
 } // namespace
index ea9f57cf97344eb0757c2e1afaec2522c1892b5b..4be423ac1b24f03ade918c0751ce94a8b46447a6 100644 (file)
@@ -267,7 +267,14 @@ void checkFrame(const AnalysisDataFrameRef       &frame,
                 const AnalysisDataTestInputFrame &refFrame)
 {
     checkHeader(frame.header(), refFrame);
-    checkPoints(frame.points(), refFrame.pointSet(0), 0);
+    ASSERT_EQ(refFrame.pointSetCount(), frame.pointSetCount());
+    for (int i = 0; i < frame.pointSetCount(); ++i)
+    {
+        const AnalysisDataPointSetRef       &points    = frame.pointSet(i);
+        const AnalysisDataTestInputPointSet &refPoints = refFrame.pointSet(i);
+        EXPECT_EQ(refPoints.firstColumn(), points.firstColumn());
+        checkPoints(points, refPoints, 0);
+    }
 }
 
 /*! \internal \brief
@@ -398,6 +405,8 @@ class StaticDataPointsStorageChecker
          * \param[in] data       Test input data to check against.
          * \param[in] frameIndex Frame index for which this functor expects
          *      to be called.
+         * \param[in] pointSetIndex Point set for which this functor expects
+         *      to be called.
          * \param[in] storageCount How many past frames should be checked for
          *      storage (-1 = check all frames).
          *
@@ -407,11 +416,13 @@ class StaticDataPointsStorageChecker
          *
          * \p source and \p data must exist for the lifetime of this object.
          */
-        StaticDataPointsStorageChecker(AbstractAnalysisData *source,
+        StaticDataPointsStorageChecker(AbstractAnalysisData        *source,
                                        const AnalysisDataTestInput *data,
-                                       int frameIndex, int storageCount)
+                                       int frameIndex, int pointSetIndex,
+                                       int storageCount)
             : source_(source), data_(data),
-              frameIndex_(frameIndex), storageCount_(storageCount)
+              frameIndex_(frameIndex), pointSetIndex_(pointSetIndex),
+              storageCount_(storageCount)
         {
         }
 
@@ -420,9 +431,9 @@ class StaticDataPointsStorageChecker
         {
             SCOPED_TRACE(formatString("Frame %d", frameIndex_));
             const AnalysisDataTestInputFrame    &refFrame  = data_->frame(frameIndex_);
-            const AnalysisDataTestInputPointSet &refPoints = refFrame.pointSet(0);
+            const AnalysisDataTestInputPointSet &refPoints = refFrame.pointSet(pointSetIndex_);
             EXPECT_EQ(refPoints.firstColumn(), points.firstColumn());
-            EXPECT_EQ(data_->columnCount(), points.columnCount());
+            EXPECT_EQ(refPoints.size(), points.columnCount());
             checkHeader(points.header(), refFrame);
             checkPoints(points, refPoints, 0);
             for (int past = 1;
@@ -443,6 +454,7 @@ class StaticDataPointsStorageChecker
         AbstractAnalysisData        *source_;
         const AnalysisDataTestInput *data_;
         int                          frameIndex_;
+        int                          pointSetIndex_;
         int                          storageCount_;
 };
 
@@ -543,9 +555,9 @@ MockAnalysisDataModule::setupStaticStorageCheck(
 {
     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");
-    impl_->flags_ |= efAllowMulticolumn;
+    GMX_RELEASE_ASSERT(data.isMultipoint() == source->isMultipoint(),
+                       "Mismatching multipoint properties");
+    impl_->flags_ |= efAllowMulticolumn | efAllowMultipoint;
 
     ::testing::InSequence dummy;
     using ::testing::_;
@@ -558,9 +570,12 @@ MockAnalysisDataModule::setupStaticStorageCheck(
         const AnalysisDataTestInputFrame &frame = data.frame(row);
         EXPECT_CALL(*this, frameStarted(_))
             .WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
-        EXPECT_CALL(*this, pointsAdded(_))
-            .WillOnce(Invoke(StaticDataPointsStorageChecker(source, &data, row,
-                                                            storageCount)));
+        for (int pointSet = 0; pointSet < frame.pointSetCount(); ++pointSet)
+        {
+            StaticDataPointsStorageChecker checker(source, &data, row, pointSet,
+                                                   storageCount);
+            EXPECT_CALL(*this, pointsAdded(_)).WillOnce(Invoke(checker));
+        }
         EXPECT_CALL(*this, frameFinished(_))
             .WillOnce(Invoke(StaticDataFrameHeaderChecker(&frame)));
     }