Fix portability issue with ArrayRef initializer overloading
authorErik Lindahl <erik@kth.se>
Wed, 30 Jul 2014 12:46:57 +0000 (14:46 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Thu, 7 Aug 2014 17:59:16 +0000 (19:59 +0200)
Some compilers (in particular the Fujitsu compilers that are
derived from an earlier version of Clang) will not allow
overloading functions where one version uses a pointer to a
type, and the other a vector iterator to the same type - likely
because they have implemented iterators with pointers.
Regardless of what the C++ standard says, we need this working
on K computer, so this patch replaces the overloaded
initializers with non-member functions that create an ArrayRef
or ConstArrayRef either from pointers, an array,
or iterators.

Change-Id: I4c4e327c869920cc08e3f955e88cb3a5b28c7e87

src/gromacs/analysisdata/arraydata.cpp
src/gromacs/analysisdata/dataframe.cpp
src/gromacs/analysisdata/datastorage.cpp
src/gromacs/commandline/pargs.cpp
src/gromacs/options/filenameoption.cpp
src/gromacs/selection/nbsearch.cpp
src/gromacs/selection/selection.h
src/gromacs/utility/arrayref.h

index c14fece0aec9b2fb6f5f9cc898f377f83823d622..55dc204df5c45c4da35c4c85fc1452ed600fdf83 100644 (file)
@@ -74,8 +74,8 @@ AbstractAnalysisArrayData::tryGetDataFrameInternal(int index) const
         = value_.begin() + index * columnCount();
     return AnalysisDataFrameRef(
             AnalysisDataFrameHeader(index, xvalue(index), 0.0),
-            AnalysisDataValuesRef(begin, begin + columnCount()),
-            AnalysisDataPointSetInfosRef(&pointSetInfo_, 1));
+            constArrayRefFromVector<AnalysisDataValue>(begin, begin + columnCount()),
+            constArrayRefFromArray(&pointSetInfo_, 1));
 }
 
 
@@ -184,8 +184,8 @@ AbstractAnalysisArrayData::valuesReady()
         modules.notifyPointsAdd(
                 AnalysisDataPointSetRef(
                         header, pointSetInfo_,
-                        AnalysisDataValuesRef(valueIter,
-                                              valueIter + columnCount())));
+                        constArrayRefFromVector<AnalysisDataValue>(valueIter,
+                                                                   valueIter + columnCount())));
         modules.notifyFrameFinish(header);
     }
     modules.notifyDataFinish();
index 12ff2c03c7b9ad33f378c8e8a364def0047a8caf..cfa09e86a3811889cfe6adbabfe52bd4f235e492 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2012,2013, by the GROMACS development team, led by
+ * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -74,8 +74,8 @@ AnalysisDataPointSetRef::AnalysisDataPointSetRef(
     : header_(header),
       dataSetIndex_(pointSetInfo.dataSetIndex()),
       firstColumn_(pointSetInfo.firstColumn()),
-      values_(&*values.begin() + pointSetInfo.valueOffset(),
-              pointSetInfo.valueCount())
+      values_(constArrayRefFromArray(&*values.begin() + pointSetInfo.valueOffset(),
+                                     pointSetInfo.valueCount()))
 {
     GMX_ASSERT(header_.isValid(),
                "Invalid point set reference should not be constructed");
@@ -86,7 +86,7 @@ AnalysisDataPointSetRef::AnalysisDataPointSetRef(
         const AnalysisDataFrameHeader        &header,
         const std::vector<AnalysisDataValue> &values)
     : header_(header), dataSetIndex_(0), firstColumn_(0),
-      values_(values.begin(), values.end())
+      values_(constArrayRefFromVector<AnalysisDataValue>(values.begin(), values.end()))
 {
     GMX_ASSERT(header_.isValid(),
                "Invalid point set reference should not be constructed");
@@ -166,8 +166,8 @@ AnalysisDataFrameRef::AnalysisDataFrameRef(
         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())
+    : header_(header), values_(constArrayRefFromVector<AnalysisDataValue>(values.begin(), values.end())),
+      pointSets_(constArrayRefFromVector<AnalysisDataPointSetInfo>(pointSets.begin(), pointSets.end()))
 {
     GMX_ASSERT(!pointSets_.empty(), "There must always be a point set");
 }
@@ -176,7 +176,7 @@ AnalysisDataFrameRef::AnalysisDataFrameRef(
 AnalysisDataFrameRef::AnalysisDataFrameRef(
         const AnalysisDataFrameRef &frame, int firstColumn, int columnCount)
     : header_(frame.header()),
-      values_(&frame.values_[firstColumn], columnCount),
+      values_(constArrayRefFromArray(&frame.values_[firstColumn], columnCount)),
       pointSets_(frame.pointSets_)
 {
     // FIXME: This doesn't produce a valid internal state, although it does
index ec80d9b4c8c8552d1e94cb254ab32d23031f00a5..57dfea17c7bbd6564a1d5e11edaf3bfd52ed6415 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2012,2013, by the GROMACS development team, led by
+ * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -635,7 +635,7 @@ AnalysisDataStorageFrameData::addPointSet(int dataSetIndex, int firstColumn,
     AnalysisDataPointSetInfo pointSetInfo(0, valueCount,
                                           dataSetIndex, firstColumn);
     AnalysisDataPointSetRef  pointSet(header(), pointSetInfo,
-                                      AnalysisDataValuesRef(begin, end));
+                                      constArrayRefFromVector<AnalysisDataValue>(begin, end));
     storageImpl().modules_->notifyParallelPointsAdd(pointSet);
     if (storageImpl().shouldNotifyImmediately())
     {
@@ -684,7 +684,7 @@ AnalysisDataStorageFrameData::pointSet(int index) const
                "Invalid point set index");
     return AnalysisDataPointSetRef(
             header_, pointSets_[index],
-            AnalysisDataValuesRef(values_.begin(), values_.end()));
+            constArrayRefFromVector<AnalysisDataValue>(values_.begin(), values_.end()));
 }
 
 }   // namespace internal
index e7236740c4f1d7666f01d81bf4e1719e75d461b9..944b3089b5ea935dd24de645142239c42ff38dcd 100644 (file)
@@ -798,7 +798,7 @@ gmx_bool parse_common_args(int *argc, char *argv[], unsigned long Flags,
         if (context != NULL && !(FF(PCA_QUIET)))
         {
             gmx::Options options(NULL, NULL);
-            options.setDescription(gmx::ConstArrayRef<const char *>(desc, ndesc));
+            options.setDescription(gmx::constArrayRefFromArray(desc, ndesc));
             for (i = 0; i < nfile; i++)
             {
                 gmx::filenmToOptions(&options, &fnm[i]);
@@ -810,7 +810,7 @@ gmx_bool parse_common_args(int *argc, char *argv[], unsigned long Flags,
             gmx::CommandLineHelpWriter(options)
                 .setShowDescriptions(true)
                 .setTimeUnitString(output_env_get_time_unit(*oenv))
-                .setKnownIssues(gmx::ConstArrayRef<const char *>(bugs, nbugs))
+                .setKnownIssues(gmx::constArrayRefFromArray(bugs, nbugs))
                 .writeHelp(*context);
         }
     }
index 4bc2810ba960958a41638e2cd376e55bc84bf62b..f5b58f0665236e1ff7253f1bc251f4cc924d6e09 100644 (file)
@@ -371,7 +371,7 @@ ConstArrayRef<const char *> FileNameOptionStorage::extensions() const
     const FileTypeRegistry &registry    = FileTypeRegistry::instance();
     const FileTypeHandler  &typeHandler = registry.handlerForType(filetype_, legacyType_);
     const ExtensionList    &extensions  = typeHandler.extensions();
-    return ConstArrayRef<const char *>(extensions.begin(), extensions.end());
+    return constArrayRefFromVector<const char *>(extensions.begin(), extensions.end());
 }
 
 /********************************************************************
index 06ecbe62a37e417ebd5118ecff94f79711fd84c8..24a66ddbda7538ff88f2e9f9e03c72ed9356bc5f 100644 (file)
@@ -606,14 +606,14 @@ void AnalysisNeighborhoodPairSearchImpl::startSearch(
 {
     if (positions.index_ < 0)
     {
-        testPositions_ = ConstArrayRef<rvec>(positions.x_, positions.count_);
+        testPositions_ = constArrayRefFromArray<rvec>(positions.x_, positions.count_);
         reset(0);
     }
     else
     {
         // Somewhat of a hack: setup the array such that only the last position
         // will be used.
-        testPositions_ = ConstArrayRef<rvec>(positions.x_, positions.index_ + 1);
+        testPositions_ = constArrayRefFromArray<rvec>(positions.x_, positions.index_ + 1);
         reset(positions.index_);
     }
 }
index a09167d5a53875dc5d5ba89e1b2352c232046767..2f36c1ecf216453bbe9b575356eab376569c6fb3 100644 (file)
@@ -333,8 +333,8 @@ class Selection
         //! Returns atom indices of all atoms in the selection.
         ConstArrayRef<int> atomIndices() const
         {
-            return ConstArrayRef<int>(sel_->rawPositions_.m.mapb.a,
-                                      sel_->rawPositions_.m.mapb.nra);
+            return constArrayRefFromArray(sel_->rawPositions_.m.mapb.a,
+                                          sel_->rawPositions_.m.mapb.nra);
         }
         //! Number of positions in the selection.
         int posCount() const { return data().posCount(); }
@@ -343,7 +343,7 @@ class Selection
         //! Returns coordinates for this selection as a continuous array.
         ConstArrayRef<rvec> coordinates() const
         {
-            return ConstArrayRef<rvec>(data().rawPositions_.x, posCount());
+            return constArrayRefFromArray(data().rawPositions_.x, posCount());
         }
         //! Returns whether velocities are available for this selection.
         bool hasVelocities() const { return data().rawPositions_.v != NULL; }
@@ -355,7 +355,7 @@ class Selection
         ConstArrayRef<rvec> velocities() const
         {
             GMX_ASSERT(hasVelocities(), "Velocities accessed, but unavailable");
-            return ConstArrayRef<rvec>(data().rawPositions_.v, posCount());
+            return constArrayRefFromArray(data().rawPositions_.v, posCount());
         }
         //! Returns whether forces are available for this selection.
         bool hasForces() const { return sel_->rawPositions_.f != NULL; }
@@ -367,7 +367,7 @@ class Selection
         ConstArrayRef<rvec> forces() const
         {
             GMX_ASSERT(hasForces(), "Forces accessed, but unavailable");
-            return ConstArrayRef<rvec>(data().rawPositions_.f, posCount());
+            return constArrayRefFromArray(data().rawPositions_.f, posCount());
         }
         //! Returns masses for this selection as a continuous array.
         ConstArrayRef<real> masses() const
@@ -377,8 +377,8 @@ class Selection
             // (and thus the masses and charges are fixed).
             GMX_ASSERT(data().posMass_.size() >= static_cast<size_t>(posCount()),
                        "Internal inconsistency");
-            return ConstArrayRef<real>(data().posMass_.begin(),
-                                       data().posMass_.begin() + posCount());
+            return constArrayRefFromVector<real>(data().posMass_.begin(),
+                                                 data().posMass_.begin() + posCount());
         }
         //! Returns charges for this selection as a continuous array.
         ConstArrayRef<real> charges() const
@@ -388,8 +388,8 @@ class Selection
             // (and thus the masses and charges are fixed).
             GMX_ASSERT(data().posCharge_.size() >= static_cast<size_t>(posCount()),
                        "Internal inconsistency");
-            return ConstArrayRef<real>(data().posCharge_.begin(),
-                                       data().posCharge_.begin() + posCount());
+            return constArrayRefFromVector<real>(data().posCharge_.begin(),
+                                                 data().posCharge_.begin() + posCount());
         }
         /*! \brief
          * Returns reference IDs for this selection as a continuous array.
@@ -398,7 +398,7 @@ class Selection
          */
         ConstArrayRef<int> refIds() const
         {
-            return ConstArrayRef<int>(data().rawPositions_.m.refid, posCount());
+            return constArrayRefFromArray(data().rawPositions_.m.refid, posCount());
         }
         /*! \brief
          * Returns mapped IDs for this selection as a continuous array.
@@ -407,7 +407,7 @@ class Selection
          */
         ConstArrayRef<int> mappedIds() const
         {
-            return ConstArrayRef<int>(data().rawPositions_.m.mapid, posCount());
+            return constArrayRefFromArray(data().rawPositions_.m.mapid, posCount());
         }
 
         //! Returns whether the covered fraction can change between frames.
@@ -654,7 +654,7 @@ class SelectionPosition
                 return ConstArrayRef<int>();
             }
             const int first = sel_->rawPositions_.m.mapb.index[i_];
-            return ConstArrayRef<int>(&atoms[first], atomCount());
+            return constArrayRefFromArray(&atoms[first], atomCount());
         }
         /*! \brief
          * Returns whether this position is selected in the current frame.
index c02f120ab0570bb1c86deb033fd4dd3f72cae886..982c712b7d4128cdd72027290082d6d6e29fc467 100644 (file)
@@ -144,41 +144,15 @@ class ArrayRef
          * \param[in] end    Pointer to the end of a range.
          *
          * Passed pointers must remain valid for the lifetime of this object.
+         *
+         * \note For clarity, use the non-member function arrayRefFromPointers
+         * instead.
          */
         ArrayRef(pointer begin, pointer end)
             : begin_(begin), end_(end)
         {
             GMX_ASSERT(end >= begin, "Invalid range");
         }
-        /*! \brief
-         * Constructs a reference to a particular range in a std::vector.
-         *
-         * \param[in] begin  Iterator to the beginning of a range.
-         * \param[in] end    Iterator to the end of a range.
-         *
-         * The referenced vector must remain valid and not be reallocated for
-         * the lifetime of this object.
-         */
-        ArrayRef(typename std::vector<value_type>::iterator begin,
-                 typename std::vector<value_type>::iterator end)
-            : begin_((begin != end) ? &*begin : NULL),
-              end_(begin_+(end-begin))
-        {
-            GMX_ASSERT(end >= begin, "Invalid range");
-        }
-        /*! \brief
-         * Constructs a reference to an array.
-         *
-         * \param[in] begin  Pointer to the beginning of the array.
-         *      May be NULL if \p size is zero.
-         * \param[in] size   Number of elements in the array.
-         *
-         * Passed pointer must remain valid for the lifetime of this object.
-         */
-        ArrayRef(pointer begin, size_type size)
-            : begin_(begin), end_(begin + size)
-        {
-        }
         //! \cond
         // Doxygen 1.8.5 doesn't parse the declaration correctly...
         /*! \brief
@@ -281,6 +255,62 @@ class ArrayRef
         pointer           end_;
 };
 
+
+/*! \brief
+ * Constructs a reference to a particular range from two pointers.
+ *
+ * \param[in] begin  Pointer to the beginning of a range.
+ * \param[in] end    Pointer to the end of a range.
+ *
+ * Passed pointers must remain valid for the lifetime of this object.
+ *
+ * \related ArrayRef
+ */
+template <typename T>
+ArrayRef<T> arrayRefFromPointers(T * begin, T * end)
+{
+    return ArrayRef<T>(begin, end);
+}
+
+/*! \brief
+ * Constructs a reference to an array
+ *
+ * \param[in] begin  Pointer to the beginning of the array.
+ *                   May be NULL if \p size is zero.
+ * \param[in] size   Number of elements in the array.
+ *
+ * Passed pointer must remain valid for the lifetime of this object.
+ *
+ * \related ArrayRef
+ */
+template <typename T>
+ArrayRef<T> arrayRefFromArray(T * begin, size_t size)
+{
+    return arrayRefFromPointers<T>(begin, begin+size);
+}
+
+/*! \brief
+ * Constructs a reference to a particular range in a std::vector.
+ *
+ * \param[in] begin  Iterator to the beginning of a range.
+ * \param[in] end    Iterator to the end of a range.
+ *
+ * The referenced vector must remain valid and not be reallocated for
+ * the lifetime of this object.
+ *
+ * \related ArrayRef
+ */
+template <typename T>
+ArrayRef<T> arrayRefFromVector(typename std::vector<T>::iterator begin,
+                               typename std::vector<T>::iterator end)
+{
+    T * p_begin = (begin != end) ? &*begin : NULL;
+    T * p_end   = p_begin + (end-begin);
+    return arrayRefFromPointers<T>(p_begin, p_end);
+}
+
+
+
 /*! \brief
  * STL-like container for non-mutable interface to a C array (or part of a
  * std::vector).
@@ -354,41 +384,15 @@ class ConstArrayRef
          * \param[in] end    Pointer to the end of a range.
          *
          * Passed pointers must remain valid for the lifetime of this object.
+         *
+         * \note For clarity, use the non-member function constArrayRefFromPointers
+         * instead.
          */
         ConstArrayRef(const_pointer begin, const_pointer end)
             : begin_(begin), end_(end)
         {
             GMX_ASSERT(end >= begin, "Invalid range");
         }
-        /*! \brief
-         * Constructs a reference to a particular range in a std::vector.
-         *
-         * \param[in] begin  Iterator to the beginning of a range.
-         * \param[in] end    Iterator to the end of a range.
-         *
-         * The referenced vector must remain valid and not be reallocated for
-         * the lifetime of this object.
-         */
-        ConstArrayRef(typename std::vector<value_type>::const_iterator begin,
-                      typename std::vector<value_type>::const_iterator end)
-            : begin_((begin != end) ? &*begin : NULL),
-              end_(begin_+(end-begin))
-        {
-            GMX_ASSERT(end >= begin, "Invalid range");
-        }
-        /*! \brief
-         * Constructs a reference to an array.
-         *
-         * \param[in] begin  Pointer to the beginning of the array.
-         *      May be NULL if \p size is zero.
-         * \param[in] size   Number of elements in the array.
-         *
-         * Passed pointer must remain valid for the lifetime of this object.
-         */
-        ConstArrayRef(const_pointer begin, size_type size)
-            : begin_(begin), end_(begin + size)
-        {
-        }
         //! \cond
         // Doxygen 1.8.5 doesn't parse the declaration correctly...
         /*! \brief
@@ -465,6 +469,60 @@ class ConstArrayRef
         const_pointer           end_;
 };
 
+
+/*! \brief
+ * Constructs a reference to a particular range from two pointers.
+ *
+ * \param[in] begin  Pointer to the beginning of a range.
+ * \param[in] end    Pointer to the end of a range.
+ *
+ * Passed pointers must remain valid for the lifetime of this object.
+ *
+ * \related ConstArrayRef
+ */
+template <typename T>
+ConstArrayRef<T> constArrayRefFromPointers(const T * begin, const T * end)
+{
+    return ConstArrayRef<T>(begin, end);
+}
+
+/*! \brief
+ * Constructs a reference to an array.
+ *
+ * \param[in] begin  Pointer to the beginning of the array.
+ *                   May be NULL if \p size is zero.
+ * \param[in] size   Number of elements in the array.
+ *
+ * Passed pointer must remain valid for the lifetime of this object.
+ *
+ * \related ConstArrayRef
+ */
+template <typename T>
+ConstArrayRef<T> constArrayRefFromArray(const T * begin, size_t size)
+{
+    return constArrayRefFromPointers<T>(begin, begin+size);
+}
+
+/*! \brief
+ * Constructs a reference to a particular range in a std::vector.
+ *
+ * \param[in] begin  Iterator to the beginning of a range.
+ * \param[in] end    Iterator to the end of a range.
+ *
+ * The referenced vector must remain valid and not be reallocated for
+ * the lifetime of this object.
+ *
+ * \related ConstArrayRef
+ */
+template <typename T>
+ConstArrayRef<T> constArrayRefFromVector(typename std::vector<T>::const_iterator begin,
+                                         typename std::vector<T>::const_iterator end)
+{
+    const T * p_begin = (begin != end) ? &*begin : NULL;
+    const T * p_end   = p_begin + (end-begin);
+    return constArrayRefFromPointers<T>(p_begin, p_end);
+}
+
 /*! \brief
  * Simple swap method for ArrayRef objects.
  *