Add unit tests for ArrayRef
authorMark Abraham <mark.j.abraham@gmail.com>
Wed, 18 Mar 2015 08:52:15 +0000 (08:52 +0000)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Tue, 16 Jun 2015 09:06:45 +0000 (11:06 +0200)
This fixes mdrun signalling on BG/Q. The templated constructor for
ArrayRef does not work with xlc 12 on BG/Q with array fields of
structs unless the base type has size equal to char. Presumably this
is another example of the way that attributes of struct fields just
get thrown away by compilers.

Changed name of local variable to sig, just in case "signal" clashes
with a preprocessor symbol somewhere...

Some changes to the arrayref.h implementation to make it easier to
write the tests as type-parameterized, by making the factory functions
also appear as members of the corresponding classes, with the same
name for both ArrayRef and ConstArrayRef. While there, remove some
documentation duplication.

Fixes #1701

Change-Id: I6894706b224dc5f3db7893503371107f1ff324d2

src/gromacs/mdlib/mdrun_signalling.cpp
src/gromacs/mdlib/mdrun_signalling.h
src/gromacs/utility/arrayref.h
src/gromacs/utility/tests/CMakeLists.txt
src/gromacs/utility/tests/arrayref.cpp [new file with mode: 0644]

index c76c4f7980c93cebb0fbe3315c05dd5ae925be24..1714753f6bbc74077ba6c7f435d8d66cd60aaf02 100644 (file)
@@ -88,10 +88,11 @@ prepareSignalBuffer(struct gmx_signalling_t *gs)
 {
     if (gs)
     {
-        gmx::ArrayRef<int>  signal(gs->sig);
+        gmx::ArrayRef<char> sig(gs->sig);
         gmx::ArrayRef<real> temp(gs->mpiBuffer);
 
-        std::copy(signal.begin(), signal.end(), temp.begin());
+        std::copy(sig.begin(), sig.end(), temp.begin());
+
         return temp;
     }
     else
@@ -130,9 +131,9 @@ handleSignals(struct gmx_signalling_t  *gs,
             /* Set the communicated signal only when it is non-zero,
              * since signals might not be processed at each MD step.
              */
-            int gsi = (gs->mpiBuffer[i] >= 0.0 ?
-                       (int)(gs->mpiBuffer[i] + 0.5) :
-                       (int)(gs->mpiBuffer[i] - 0.5));
+            char gsi = (gs->mpiBuffer[i] >= 0.0 ?
+                        static_cast<char>(gs->mpiBuffer[i] + 0.5) :
+                        static_cast<char>(gs->mpiBuffer[i] - 0.5));
             if (gsi != 0)
             {
                 gs->set[i] = gsi;
index 2fdde5d1bf4056bc18f9bf49b39efeea3751eaa0..8f3c0340daacebe737c1cb4c9a3ca5f5cbd2c145 100644 (file)
@@ -72,11 +72,15 @@ enum {
     eglsCHKPT, eglsSTOPCOND, eglsRESETCOUNTERS, eglsNR
 };
 
-/*! \internal \brief Object used by mdrun ranks to signal to each other */
+/*! \internal
+ * \brief Object used by mdrun ranks to signal to each other
+ *
+ * Note that xlc on BG/Q requires sig to be of size char (see unit tests
+ * of ArrayRef for details). */
 struct gmx_signalling_t {
     int  nstms;             /**< The frequency for inter-simulation communication */
-    int  sig[eglsNR];       /**< The signal set by this rank in do_md */
-    int  set[eglsNR];       /**< The communicated signal, equal for all ranks once communication has occurred */
+    char sig[eglsNR];       /**< The signal set by this rank in do_md */
+    char set[eglsNR];       /**< The communicated signal, equal for all ranks once communication has occurred */
     real mpiBuffer[eglsNR]; /**< Buffer for communication */
 };
 
index d59a3548577305800b7836d5f6a7d1acc1ca58f5..ccf6e5e9d04202aed0ac46cf265c7b1dec1b6b2a 100644 (file)
@@ -125,6 +125,51 @@ class ArrayRef
         //! Standard reverse iterator.
         typedef std::reverse_iterator<const_iterator> const_reverse_iterator;
 
+        /*! \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.
+         */
+        static ArrayRef<value_type>
+        fromPointers(value_type *begin, value_type *end)
+        {
+            return ArrayRef<value_type>(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.
+         */
+        static ArrayRef<value_type>
+        fromArray(value_type *begin, size_t size)
+        {
+            return ArrayRef<value_type>(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.
+         */
+        static ArrayRef<value_type>
+        fromVector(typename std::vector<value_type>::iterator begin,
+                   typename std::vector<value_type>::iterator end)
+        {
+            value_type *p_begin = (begin != end) ? &*begin : NULL;
+            value_type *p_end   = p_begin + (end-begin);
+            return ArrayRef<value_type>(p_begin, p_end);
+        }
+
         /*! \brief
          * Constructs an empty reference.
          */
@@ -186,6 +231,11 @@ class ArrayRef
          *
          * This constructor is not explicit to allow directly passing
          * a C array to a function that takes an ArrayRef parameter.
+         *
+         * xlc on BG/Q compiles wrong code if the C array is a struct
+         * field, unless value_type is char or unsigned char. There's
+         * no good way to assert on this before C++11 (which that
+         * compiler will never support).
          */
         template <size_t count>
         ArrayRef(value_type (&array)[count])
@@ -272,60 +322,6 @@ class ArrayRef
 };
 
 
-/*! \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
@@ -381,6 +377,28 @@ class ConstArrayRef
         //! Standard reverse iterator.
         typedef std::reverse_iterator<const_iterator> const_reverse_iterator;
 
+        //! \copydoc ArrayRef::fromPointers()
+        static ConstArrayRef<value_type>
+        fromPointers(const value_type *begin, const value_type *end)
+        {
+            return ConstArrayRef<value_type>(begin, end);
+        }
+        //! \copydoc ArrayRef::fromArray()
+        static ConstArrayRef<value_type>
+        fromArray(const value_type *begin, size_t size)
+        {
+            return ConstArrayRef<value_type>(begin, begin+size);
+        }
+        //! \copydoc ArrayRef::fromVector()
+        static ConstArrayRef<value_type>
+        fromVector(typename std::vector<value_type>::const_iterator begin,
+                   typename std::vector<value_type>::const_iterator end)
+        {
+            const value_type *p_begin = (begin != end) ? &*begin : NULL;
+            const value_type *p_end   = p_begin + (end-begin);
+            return ConstArrayRef<value_type>(p_begin, p_end);
+        }
+
         /*! \brief
          * Constructs an empty reference.
          */
@@ -437,10 +455,16 @@ class ConstArrayRef
          * pointer).  It constructs a reference to the whole array, without
          * a need to pass the number of elements explicitly.  The compiler
          * must be able to deduce the array size.
+         *
          * Passed array must remain valid for the lifetime of this object.
          *
          * This constructor is not explicit to allow directly passing
          * a C array to a function that takes a ConstArrayRef parameter.
+         *
+         * xlc on BG/Q compiles wrong code if the C array is a struct
+         * field, unless value_type is char or unsigned char. There's
+         * no good way to assert on this before C++11 (which that
+         * compiler will never support).
          */
         template <size_t count>
         ConstArrayRef(const value_type (&array)[count])
@@ -502,57 +526,51 @@ class ConstArrayRef
 };
 
 
-/*! \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
- */
+//! \copydoc ArrayRef::fromPointers()
+//! \related ArrayRef
 template <typename T>
-ConstArrayRef<T> constArrayRefFromPointers(const T * begin, const T * end)
+ArrayRef<T> arrayRefFromPointers(T *begin, T *end)
 {
-    return ConstArrayRef<T>(begin, end);
+    return ArrayRef<T>::fromPointers(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
- */
+//! \copydoc ArrayRef::fromArray()
+//! \related ArrayRef
 template <typename T>
-ConstArrayRef<T> constArrayRefFromArray(const T * begin, size_t size)
+ArrayRef<T> arrayRefFromArray(T *begin, size_t size)
 {
-    return constArrayRefFromPointers<T>(begin, begin+size);
+    return ArrayRef<T>::fromArray(begin, size);
+}
+//! \copydoc ArrayRef::fromVector()
+//! \related ArrayRef
+template <typename T>
+ArrayRef<T> arrayRefFromVector(typename std::vector<T>::iterator begin,
+                               typename std::vector<T>::iterator end)
+{
+    return ArrayRef<T>::fromVector(begin, end);
 }
 
-/*! \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
- */
+
+//! \copydoc ConstArrayRef::fromPointers()
+//! \related ConstArrayRef
+template <typename T>
+ConstArrayRef<T> constArrayRefFromPointers(const T *begin, const T *end)
+{
+    return ConstArrayRef<T>::fromPointers(begin, end);
+}
+//! \copydoc ConstArrayRef::fromArray()
+//! \related ConstArrayRef
+template <typename T>
+ConstArrayRef<T> constArrayRefFromArray(const T *begin, size_t size)
+{
+    return ConstArrayRef<T>::fromArray(begin, size);
+}
+//! \copydoc ConstArrayRef::fromVector()
+//! \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);
+    return ConstArrayRef<T>::fromVector(begin, end);
 }
 
 /*! \brief
index 9024e56ea608d12acd0320e8af4243ba0b88c705..a002c762c3f3d91d2d001eca20826f516356bf27 100644 (file)
@@ -1,7 +1,7 @@
 #
 # This file is part of the GROMACS molecular simulation package.
 #
-# Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
+# Copyright (c) 2012,2013,2014,2015, 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.
@@ -33,5 +33,6 @@
 # the research papers on the package. Check out http://www.gromacs.org.
 
 gmx_add_unit_test(UtilityUnitTests utility-test
+                  arrayref.cpp
                   bitmask32.cpp bitmask64.cpp bitmask128.cpp
                   stringutil.cpp)
diff --git a/src/gromacs/utility/tests/arrayref.cpp b/src/gromacs/utility/tests/arrayref.cpp
new file mode 100644 (file)
index 0000000..d5f00bc
--- /dev/null
@@ -0,0 +1,266 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2015, 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.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \internal \file
+ * \brief Tests for gmx::ArrayRef and gmx::ConstArrayRef.
+ *
+ * \author Mark Abraham <mark.j.abraham@gmail.com>
+ * \ingroup module_utility
+ */
+#include "gmxpre.h"
+
+#include "gromacs/utility/arrayref.h"
+
+#include "config.h"
+
+#include <vector>
+
+#include <gtest/gtest.h>
+
+#include "gromacs/utility/basedefinitions.h"
+#include "gromacs/utility/real.h"
+
+namespace gmx
+{
+
+namespace
+{
+
+TEST(EmptyArrayRefTest, IsEmpty)
+{
+    EmptyArrayRef  emptyArray = {};
+    ArrayRef<real> empty(emptyArray);
+
+    EXPECT_EQ(0U, empty.size());
+    EXPECT_TRUE(empty.empty());
+}
+
+TEST(EmptyConstArrayRefTest, IsEmpty)
+{
+    EmptyArrayRef       emptyArray = {};
+    ConstArrayRef<real> empty(emptyArray);
+
+    EXPECT_EQ(0U, empty.size());
+    EXPECT_TRUE(empty.empty());
+}
+
+#ifdef GTEST_HAS_TYPED_TEST
+
+//! Define the types that end up being available as TypeParam in the test cases for both kinds of ArrayRef
+typedef ::testing::Types<
+        ArrayRef<char>,
+        ArrayRef<unsigned char>,
+        ArrayRef<int>,
+        ArrayRef<unsigned int>,
+        ArrayRef<long>,
+        ArrayRef<unsigned long>,
+        ArrayRef<gmx_int64_t>,
+        ArrayRef<gmx_uint64_t>,
+        ArrayRef<float>,
+        ArrayRef<double>,
+        ConstArrayRef<char>,
+        ConstArrayRef<unsigned char>,
+        ConstArrayRef<int>,
+        ConstArrayRef<unsigned int>,
+        ConstArrayRef<long>,
+        ConstArrayRef<unsigned long>,
+        ConstArrayRef<gmx_int64_t>,
+        ConstArrayRef<gmx_uint64_t>,
+        ConstArrayRef<float>,
+        ConstArrayRef<double>
+        > ArrayRefTypes;
+
+/*! \brief Permit all the tests to run on all kinds of ArrayRefs
+ *
+ * The main objective is to verify that all the different kinds of
+ * construction lead to the expected result. */
+template <typename TypeParam>
+class ArrayRefTest : public ::testing::Test
+{
+    public:
+        typedef TypeParam ArrayRefType;
+        typedef typename ArrayRefType::value_type ValueType;
+
+        /*! \brief Run the same tests all the time
+         *
+         * Note that test cases must call this->runTests(), because
+         * that's how the derived-class templates that implement
+         * type-parameterized tests actually work. */
+        void runTests(ValueType     *a,
+                      size_t         aSize,
+                      ValueType     *aData,
+                      ArrayRefType  &arrayRef)
+        {
+            ASSERT_EQ(aSize, arrayRef.size());
+            ASSERT_FALSE(arrayRef.empty());
+            EXPECT_EQ(aData, arrayRef.data());
+            EXPECT_EQ(a[0], arrayRef.front());
+            EXPECT_EQ(a[aSize-1], arrayRef.back());
+            for (size_t i = 0; i != aSize; ++i)
+            {
+                EXPECT_EQ(a[i], arrayRef[i]);
+            }
+        }
+};
+
+TYPED_TEST_CASE(ArrayRefTest, ArrayRefTypes);
+
+/* Welcome back to the past. While you can declare a static array[] of
+   templated type in a class, in C++98, you have to define it outside
+   the class, and when you do, the compiler knows the declaration is
+   incomplete and can't match the types to actual functions. So,
+   declaring locals is the only choice available, so we need macros to
+   avoid duplication. Lovely. */
+#define DEFINE_ARRAY(a, aSize)                                  \
+    typename TestFixture::ValueType (a)[] = {                   \
+        static_cast<typename TestFixture::ValueType>(1.2),      \
+        static_cast<typename TestFixture::ValueType>(2.4),      \
+        static_cast<typename TestFixture::ValueType>(3.1)       \
+    };                                                          \
+    size_t (aSize) = sizeof((a)) / sizeof(typename TestFixture::ValueType);
+
+TYPED_TEST(ArrayRefTest, MakeWithAssignmentWorks)
+{
+    DEFINE_ARRAY(a, aSize);
+    typename TestFixture::ArrayRefType arrayRef = a;
+    this->runTests(a, aSize, a, arrayRef);
+}
+
+TYPED_TEST(ArrayRefTest, ConstructWithTemplateConstructorWorks)
+{
+    DEFINE_ARRAY(a, aSize);
+    typename TestFixture::ArrayRefType arrayRef(a);
+    this->runTests(a, aSize, a, arrayRef);
+}
+
+TYPED_TEST(ArrayRefTest, ConstructFromPointersWorks)
+{
+    DEFINE_ARRAY(a, aSize);
+    typename TestFixture::ArrayRefType arrayRef(a, a + aSize);
+    this->runTests(a, aSize, a, arrayRef);
+}
+
+TYPED_TEST(ArrayRefTest, MakeFromPointersWorks)
+{
+    DEFINE_ARRAY(a, aSize);
+    typename TestFixture::ArrayRefType arrayRef
+        = TestFixture::ArrayRefType::fromPointers(a, a + aSize);
+    this->runTests(a, aSize, a, arrayRef);
+}
+
+TYPED_TEST(ArrayRefTest, MakeFromArrayWorks)
+{
+    DEFINE_ARRAY(a, aSize);
+    typename TestFixture::ArrayRefType arrayRef
+        = TestFixture::ArrayRefType::fromArray(a, aSize);
+    this->runTests(a, aSize, a, arrayRef);
+}
+
+TYPED_TEST(ArrayRefTest, ConstructFromVectorWorks)
+{
+    DEFINE_ARRAY(a, aSize);
+    std::vector<typename TestFixture::ValueType> v(a, a + aSize);
+    typename TestFixture::ArrayRefType           arrayRef(v);
+    this->runTests(a, aSize, v.data(), arrayRef);
+}
+
+TYPED_TEST(ArrayRefTest, MakeFromVectorWorks)
+{
+    DEFINE_ARRAY(a, aSize);
+    std::vector<typename TestFixture::ValueType> v(a, a + aSize);
+    typename TestFixture::ArrayRefType arrayRef
+        = TestFixture::ArrayRefType::fromVector(v.begin(), v.end());
+    this->runTests(a, aSize, v.data(), arrayRef);
+}
+
+//! Helper struct for the case actually used in mdrun signalling
+template <typename T>
+struct Helper
+{
+    public:
+        T   a[3];
+        int size;
+};
+
+/*! \brief Test of the case actually used in mdrun signalling
+ *
+ * There, we take a non-const struct-field array of static length and
+ * make an ArrayRef to it using the template constructor that is
+ * supposed to infer the length from the static size. But on xlc on
+ * BlueGene/Q, if the base type is not char (or unsigned char), the
+ * generated code ends up with an ArrayRef of zero size, so everything
+ * breaks. Presumably the default code path accidentally works for
+ * char.
+ *
+ * Fortunately, all current uses of that constructor have a base type
+ * of char, so there's no big problem. Using a pointer-based
+ * constructor does work, if there's ever a problem (and that is
+ * tested above). */
+
+TYPED_TEST(ArrayRefTest, ConstructFromStructFieldWithTemplateConstructorWorks)
+{
+    DEFINE_ARRAY(a, aSize);
+    Helper<typename TestFixture::ValueType> h;
+    h.size = aSize;
+    for (int i = 0; i != h.size; ++i)
+    {
+        h.a[i] = a[i];
+    }
+    typename TestFixture::ArrayRefType arrayRef(h.a);
+#if defined(GMX_TARGET_BGQ) && defined(__xlC__)
+    if (sizeof(typename TestFixture::ValueType) == sizeof(char))
+#endif
+    {
+        this->runTests(h.a, h.size, h.a, arrayRef);
+    }
+}
+
+#else   // GTEST_HAS_TYPED_TEST
+
+/* A dummy test that at least signals that something is missing if one runs the
+ * unit test executable itself.
+ */
+TEST(DISABLED_ArrayRefTest, GenericTests)
+{
+    ADD_FAILURE()
+    << "Tests for generic ArrayRef functionality require support for "
+    << "Google Test typed tests, which was not available when the tests "
+    << "were compiled.";
+}
+
+#endif // GTEST_HAS_TYPED_TEST
+
+}      // namespace
+
+}      // namespace