Remove meta information from host allocator policy
authorRoland Schulz <roland.schulz@intel.com>
Tue, 14 Aug 2018 22:23:06 +0000 (15:23 -0700)
committerRoland Schulz <roland.schulz@intel.com>
Wed, 5 Sep 2018 23:17:17 +0000 (16:17 -0700)
Also prevent HostAllocator to be copied/propagated.
A HostVector shouldn't be copied without the calling code
having a choice on the pinning policy.
This commit prevents copy construction. Move construction and
copy and move assignment is still allowed.
For move/swap nothing is changed.
For copy assignment the pinning policy of the destination
container is used.

Fixes #2612

Change-Id: I15034e4db2434ea46eda84e70e3f8e6072660d20

cmake/gmxCFlags.cmake
src/gromacs/gpu_utils/hostallocator.cpp
src/gromacs/gpu_utils/hostallocator.h
src/gromacs/gpu_utils/tests/hostallocator.cpp
src/gromacs/utility/allocator.h

index 752392771a32a741629bc5794e95f46df4951930..29623e2befbd227bce417f52a46c380eec022530 100644 (file)
@@ -207,10 +207,11 @@ GMX_TEST_CFLAG(CFLAGS_WARN "/W3 /wd177 /wd411 /wd593 /wd981 /wd1418 /wd1419 /wd1
                     GMX_TEST_CXXFLAG(CXXFLAGS_WARN_OLD_GPU "-wd7 -wd82 -wd193 -wd3346" GMXC_CXXFLAGS)
                 endif()
 #All but the following warnings are identical for the C-compiler (see above)
+# 304: access control not specified
 # 383: value copied to temporary, reference to temporary used
 # 444: destructor for base class ".." is not virtual
 #2282: unrecognized GCC pragma
-                GMX_TEST_CXXFLAG(CXXFLAGS_WARN "-w3 -wd177 -wd383 -wd411 -wd444 -wd981 -wd1418 -wd1572 -wd1599 -wd2259 -wd2547 -wd3280 -wd11074 -wd11076 -wd2282" GMXC_CXXFLAGS)
+                GMX_TEST_CXXFLAG(CXXFLAGS_WARN "-w3 -wd177 -wd304 -wd383 -wd411 -wd444 -wd981 -wd1418 -wd1572 -wd1599 -wd2259 -wd2547 -wd3280 -wd11074 -wd11076 -wd2282" GMXC_CXXFLAGS)
             endif()
             GMX_TEST_CXXFLAG(CXXFLAGS_OPT "-ip -funroll-all-loops -alias-const -ansi-alias -no-prec-div -fimf-domain-exclusion=14 -qoverride-limits" GMXC_CXXFLAGS_RELEASE)
             GMX_TEST_CXXFLAG(CXXFLAGS_DEBUG "-O0" GMXC_CXXFLAGS_DEBUG)
@@ -221,7 +222,7 @@ GMX_TEST_CFLAG(CFLAGS_WARN "/W3 /wd177 /wd411 /wd593 /wd981 /wd1418 /wd1419 /wd1
             endif()
             if (GMX_COMPILER_WARNINGS)
 #809: exception specification for virtual function X is incompatible with that of overridden function
-                GMX_TEST_CXXFLAG(CXXFLAGS_WARN "/W3 /wd177 /wd383 /wd411 /wd444 /wd809 /wd981 /wd1418 /wd1572 /wd1599 /wd1786 /wd2259 /wd2547 /wd3280 /wd11074 /wd11076 /wd2282" GMXC_CXXFLAGS)
+                GMX_TEST_CXXFLAG(CXXFLAGS_WARN "/W3 /wd177 /wd304 /wd383 /wd411 /wd444 /wd809 /wd981 /wd1418 /wd1572 /wd1599 /wd1786 /wd2259 /wd2547 /wd3280 /wd11074 /wd11076 /wd2282" GMXC_CXXFLAGS)
             endif()
             GMX_TEST_CXXFLAG(CXXFLAGS_OPT "/Qip" GMXC_CXXFLAGS_RELEASE)
         endif()
index b983ab0fd496e881f20c77eb7b7ae26757e39fa4..1553aa930a24b0b2f41ce1162b493c7aa0a49d15 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2017, by the GROMACS development team, led by
+ * Copyright (c) 2017,2018, 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.
@@ -48,6 +48,7 @@
 
 #include <memory>
 
+#include "gromacs/gpu_utils/gpu_utils.h"
 #include "gromacs/gpu_utils/pinning.h"
 #include "gromacs/utility/alignedallocator.h"
 #include "gromacs/utility/gmxassert.h"
 namespace gmx
 {
 
-//! Private implementation class.
-class HostAllocationPolicy::Impl
-{
-    public:
-        /*! \brief Pointer to the last unfreed allocation, or nullptr
-         * if no allocation exists.
-         *
-         * Note that during e.g. std::vector.resize() a call to its
-         * allocator's allocate() function precedes the call to its
-         * allocator's deallocate() function for freeing the old
-         * buffer after the data has been copied from it. So in
-         * general, pointer_ will not match the argument received by
-         * free(). */
-        void         *pointer_ = nullptr;
-        //! Number of bytes in the last unfreed allocation.
-        std::size_t   numBytes_ = 0;
-        //! The pointer to any storage that has been pinned, or nullptr if none has been pinned.
-        void         *pinnedPointer_ = nullptr;
-        //! Whether this object is in mode where new allocations will be pinned by default.
-        PinningPolicy pinningPolicy_ = PinningPolicy::CannotBePinned;
-};
-
-HostAllocationPolicy::HostAllocationPolicy() : impl_(std::make_shared<Impl>())
+HostAllocationPolicy::HostAllocationPolicy(PinningPolicy pinningPolicy)
+    : pinningPolicy_(pinningPolicy)
 {
+    if (GMX_GPU != GMX_GPU_CUDA)
+    {
+        GMX_RELEASE_ASSERT(pinningPolicy == PinningPolicy::CannotBePinned,
+                           "A suitable build of GROMACS (e.g. with CUDA) is required for a "
+                           "HostAllocationPolicy to be set to a mode that produces pinning.");
+    }
 }
 
 std::size_t HostAllocationPolicy::alignment()
 {
-    return (impl_->pinningPolicy_ == PinningPolicy::CanBePinned ?
+    return (pinningPolicy_ == PinningPolicy::CanBePinned ?
             PageAlignedAllocationPolicy::alignment() :
             AlignedAllocationPolicy::alignment());
 }
+
 void *HostAllocationPolicy::malloc(std::size_t bytes) const noexcept
 {
-    // A container could have a pinned allocation that is being
-    // extended, in which case we must un-pin while we still know the
-    // old pinned vector, and which also ensures we don't pin two
-    // buffers at the same time. If there's no allocation, or it isn't
-    // pinned, then attempting to unpin it is OK, too.
-    unpin();
-    impl_->pointer_ = (impl_->pinningPolicy_ == PinningPolicy::CanBePinned ?
-                       PageAlignedAllocationPolicy::malloc(bytes) :
-                       AlignedAllocationPolicy::malloc(bytes));
-
-    if (impl_->pointer_ != nullptr)
+    if (pinningPolicy_ == PinningPolicy::CanBePinned)
+    {
+        void *p = PageAlignedAllocationPolicy::malloc(bytes);
+        pin(p, bytes);
+        return p;
+    }
+    else
     {
-        impl_->numBytes_ = bytes;
+        return AlignedAllocationPolicy::malloc(bytes);
     }
-    pin();
-    return impl_->pointer_;
 }
 
 void HostAllocationPolicy::free(void *buffer) const noexcept
 {
-    unpin();
     if (buffer == nullptr)
     {
         // Nothing to do
         return;
     }
-    if (impl_->pinningPolicy_ == PinningPolicy::CanBePinned)
+    unpin(buffer);
+    if (pinningPolicy_ == PinningPolicy::CanBePinned)
     {
         PageAlignedAllocationPolicy::free(buffer);
     }
@@ -124,65 +105,29 @@ void HostAllocationPolicy::free(void *buffer) const noexcept
     {
         AlignedAllocationPolicy::free(buffer);
     }
-    impl_->pointer_  = nullptr;
-    impl_->numBytes_ = 0;
 }
 
-PinningPolicy HostAllocationPolicy::pinningPolicy() const
+void HostAllocationPolicy::pin(void gmx_unused *p, size_t gmx_unused n) const noexcept
 {
-    return impl_->pinningPolicy_;
-}
-
-void HostAllocationPolicy::setPinningPolicy(PinningPolicy pinningPolicy)
-{
-    if (GMX_GPU != GMX_GPU_CUDA)
-    {
-        GMX_RELEASE_ASSERT(pinningPolicy == PinningPolicy::CannotBePinned,
-                           "A suitable build of GROMACS (e.g. with CUDA) is required for a "
-                           "HostAllocationPolicy to be set to a mode that produces pinning.");
-    }
-    impl_->pinningPolicy_ = pinningPolicy;
-}
-
-void HostAllocationPolicy::pin() const noexcept
-{
-    if (impl_->pinningPolicy_ == PinningPolicy::CannotBePinned ||
-        impl_->pointer_ == nullptr ||
-        impl_->pinnedPointer_ != nullptr)
+#if GMX_GPU == GMX_GPU_CUDA
+    // I believe this if statement isn't required for calls from malloc. But it
+    // is required for the unit tests. Which might not make much sense to test
+    // cases which can't actually happen for an allocator.
+    if (p == nullptr || n == 0 || isHostMemoryPinned(p))
     {
-        // Do nothing if we're not in pinning mode, or the allocation
-        // is empty, or it is already pinned.
         return;
     }
-#if GMX_GPU == GMX_GPU_CUDA
-    pinBuffer(impl_->pointer_, impl_->numBytes_);
-    impl_->pinnedPointer_ = impl_->pointer_;
-#else
-    const char *errorMessage = "Could not register the host memory for pinning.";
-
-    GMX_RELEASE_ASSERT(impl_->pinningPolicy_ == PinningPolicy::CannotBePinned,
-                       formatString("%s This build configuration must only have pinning policy "
-                                    "that leads to no pinning.", errorMessage).c_str());
+    pinBuffer(p, n);
 #endif
 }
 
-void HostAllocationPolicy::unpin() const noexcept
+void HostAllocationPolicy::unpin(void gmx_unused *p) const noexcept
 {
-    if (impl_->pinnedPointer_ == nullptr)
+#if GMX_GPU == GMX_GPU_CUDA
+    if (isHostMemoryPinned(p))
     {
-        return;
+        unpinBuffer(p);
     }
-
-#if GMX_GPU == GMX_GPU_CUDA
-    // Note that if the caller deactivated pinning mode, we still want
-    // to be able to unpin if the allocation is still pinned.
-
-    unpinBuffer(impl_->pointer_);
-    impl_->pinnedPointer_ = nullptr;
-#else
-    GMX_RELEASE_ASSERT(impl_->pinnedPointer_ == nullptr,
-                       "Since the build configuration does not support pinning, then "
-                       "the pinned pointer must be nullptr.");
 #endif
 }
 
index 3de80177597c0d61605c9833354bea7a36ce0d84..eeccd27a9ee2e28c49b266ea9d6ae7b3e19242ef 100644 (file)
@@ -134,8 +134,8 @@ using HostVector = std::vector<T, HostAllocator<T> >;
 class HostAllocationPolicy
 {
     public:
-        //! Default constructor.
-        HostAllocationPolicy();
+        //! Constructor
+        HostAllocationPolicy(PinningPolicy policy = PinningPolicy::CannotBePinned);
         /*! \brief Return the alignment size currently used by the active pinning policy. */
         std::size_t alignment();
         /*! \brief Allocate and perhaps pin page-aligned memory suitable for
@@ -178,7 +178,7 @@ class HostAllocationPolicy
          *
          * Does not throw.
          */
-        void pin() const noexcept;
+        void pin(void* p, size_t n) const noexcept;
         /*! \brief Unpin the allocation, if appropriate.
          *
          * Regardless of the allocation policy, unpin the memory if
@@ -186,49 +186,28 @@ class HostAllocationPolicy
          *
          * Does not throw.
          */
-        void unpin() const noexcept;
+        void unpin(void* p) const noexcept;
         /*! \brief Return the current pinning policy (which is semi-independent
          * of whether the buffer is actually pinned).
          *
          * Does not throw.
          */
-        PinningPolicy pinningPolicy() const;
-        //! Specify an allocator trait so that the stateful allocator should propagate.
-        using propagate_on_container_copy_assignment = std::true_type;
-        //! Specify an allocator trait so that the stateful allocator should propagate.
+        PinningPolicy pinningPolicy() const { return pinningPolicy_; }
+        //! Don't propagate for copy
+        using propagate_on_container_copy_assignment = std::false_type;
+        //! Propagate for move
         using propagate_on_container_move_assignment = std::true_type;
-        //! Specify an allocator trait so that the stateful allocator should propagate.
+        //! Propagate for move
         using propagate_on_container_swap = std::true_type;
+        //! Use default allocator for copy (same as construct+copy)
+        template<typename U = void>
+        HostAllocationPolicy select_on_container_copy_construction() const
+        {
+            return {};
+        }
     private:
-        /*! \brief Set the current pinning policy.
-         *
-         * Does not pin any current buffer. Use changePinningPolicy to
-         * orchestrate the necessary unpin, allocate, copy, pin for
-         * effectively changing the pinning policy of a HostVector.
-         *
-         * Does not throw.
-         */
-        void setPinningPolicy(PinningPolicy pinningPolicy);
-        /*! \brief Declare as a friend function the only supported way
-         * to change the pinning policy.
-         *
-         * When the pinning policy changes, we want the state of the
-         * allocation to match the new policy. However, that requires
-         * a copy and swap of the buffers, which can only take place
-         * at the level of the container. So we wrap the required
-         * operations in a helper friend function.
-         *
-         * Of course, if there is no allocation because the vector is
-         * empty, then nothing will change. */
-        template <class T> friend
-        void changePinningPolicy(HostVector<T> *v, PinningPolicy pinningPolicy);
-        //! Private implementation class.
-        class Impl;
-        /*! \brief State of the allocator.
-         *
-         * This could change through move- or copy-assignment of one
-         * policy to another, so isn't const. */
-        std::shared_ptr<Impl> impl_;
+        //! Pinning policy
+        PinningPolicy pinningPolicy_;
 };
 
 /*! \brief Helper function for changing the pinning policy of a HostVector.
@@ -242,26 +221,13 @@ class HostAllocationPolicy
 template <class T>
 void changePinningPolicy(HostVector<T> *v, PinningPolicy pinningPolicy)
 {
-    // Do we have anything to do?
-    HostAllocationPolicy vAllocationPolicy = v->get_allocator().getPolicy();
-    if (pinningPolicy == vAllocationPolicy.pinningPolicy())
+    if (v->get_allocator().pinningPolicy() == pinningPolicy)
     {
         return;
     }
-    // Make sure we never have two allocated buffers that are both pinned.
-    vAllocationPolicy.unpin();
-
-    // Construct a new vector that has the requested
-    // allocation+pinning policy, to swap into *v. If *v is empty,
-    // then no real work is done.
-    HostAllocator<T> newAllocator;
-    newAllocator.getPolicy().setPinningPolicy(pinningPolicy);
-    HostVector<T>    newV(v->begin(), v->end(), newAllocator);
-    // Replace the contents of *v, including the stateful allocator.
-    v->swap(newV);
-    // The destructor of newV cleans up the memory formerly managed by *v.
+    //Force reallocation by creating copy
+    *v = HostVector<T>(*v, {pinningPolicy});
 }
-
 }      // namespace gmx
 
 #endif
index 72b1e21bcbde0b4854e84fe983eaec571561e7a4..274e6f730c2ad838313b25295608f0c9b83c3fa4 100644 (file)
@@ -166,19 +166,20 @@ typedef ::testing::Types<int, real, RVec> TestTypes;
 
 //! Typed test fixture
 template <typename T>
-class HostAllocatorTest : public HostMemoryTest<T>
+struct HostAllocatorTest : HostMemoryTest<T>
 {
-    public:
-        //! Convenience type
-        using ValueType = T;
-        //! Convenience type
-        using AllocatorType = HostAllocator<T>;
-        //! Convenience type
-        using VectorType = std::vector<ValueType, AllocatorType>;
+    using VectorType = HostVector<T>; //!< HostVector of type tested
 };
-
 TYPED_TEST_CASE(HostAllocatorTest, TestTypes);
 
+//! Typed test fixture (no mem/gpu initializtion - much faster)
+template <typename T>
+struct HostAllocatorTestNoMem : ::testing::Test
+{
+    using VectorType = HostVector<T>; //!< HostVector of type tested
+};
+TYPED_TEST_CASE(HostAllocatorTestNoMem, TestTypes);
+
 // Note that in GoogleTest typed tests, the use of TestFixture:: and
 // this-> is sometimes required to get access to things in the fixture
 // class (or its base classes).
@@ -249,6 +250,71 @@ bool isPinned(const VectorType &v)
     return isHostMemoryPinned(data);
 }
 
+TYPED_TEST(HostAllocatorTestNoMem, CreateVector)
+{
+    typename TestFixture::VectorType input1;
+    EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+    typename TestFixture::VectorType input2({PinningPolicy::CanBePinned});
+    EXPECT_TRUE (input2.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+}
+
+TYPED_TEST(HostAllocatorTestNoMem, MoveAssignment)
+{
+    typename TestFixture::VectorType input1({PinningPolicy::CanBePinned});
+    input1 = typename TestFixture::VectorType();
+    EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+
+    typename TestFixture::VectorType input2;
+    input2 = typename TestFixture::VectorType({PinningPolicy::CanBePinned});
+    EXPECT_TRUE (input2.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+}
+
+TYPED_TEST(HostAllocatorTestNoMem, MoveConstruction)
+{
+    typename TestFixture::VectorType input1;
+    typename TestFixture::VectorType input2(std::move(input1));
+    EXPECT_FALSE(input2.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+
+    typename TestFixture::VectorType input3({PinningPolicy::CanBePinned});
+    typename TestFixture::VectorType input4(std::move(input3));
+    EXPECT_TRUE(input4.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+}
+
+TYPED_TEST(HostAllocatorTestNoMem, CopyAssignment)
+{
+    typename TestFixture::VectorType input1;
+    typename TestFixture::VectorType input2({PinningPolicy::CanBePinned});
+    input1 = input2;
+    EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+    EXPECT_TRUE (input2.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+    input2 = input1;
+    EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+    EXPECT_TRUE (input2.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+}
+
+TYPED_TEST(HostAllocatorTestNoMem, CopyConstruction)
+{
+    typename TestFixture::VectorType input1;
+    typename TestFixture::VectorType input2(input1);
+    EXPECT_FALSE(input2.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+
+    typename TestFixture::VectorType input3({PinningPolicy::CanBePinned});
+    typename TestFixture::VectorType input4(input3);
+    EXPECT_FALSE(input4.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+}
+
+TYPED_TEST(HostAllocatorTestNoMem, Swap)
+{
+    typename TestFixture::VectorType input1;
+    typename TestFixture::VectorType input2({PinningPolicy::CanBePinned});
+    swap(input1, input2);
+    EXPECT_TRUE (input1.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+    EXPECT_FALSE(input2.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+    swap(input2, input1);
+    EXPECT_FALSE(input1.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+    EXPECT_TRUE (input2.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
+}
+
 TYPED_TEST(HostAllocatorTest, ManualPinningOperationsWorkWithCuda)
 {
     if (!this->haveValidGpus())
@@ -258,14 +324,15 @@ TYPED_TEST(HostAllocatorTest, ManualPinningOperationsWorkWithCuda)
 
     typename TestFixture::VectorType input;
     changePinningPolicy(&input, PinningPolicy::CanBePinned);
+    EXPECT_TRUE(input.get_allocator().pinningPolicy() == PinningPolicy::CanBePinned);
     EXPECT_FALSE(isPinned(input));
 
     // Unpin before allocation is fine, but does nothing.
-    input.get_allocator().getPolicy().unpin();
+    input.get_allocator().getPolicy().unpin(input.data());
     EXPECT_FALSE(isPinned(input));
 
     // Pin with no contents is fine, but does nothing.
-    input.get_allocator().getPolicy().pin();
+    input.get_allocator().getPolicy().pin(input.data(), 0);
     EXPECT_FALSE(isPinned(input));
 
     // Fill some contents, which will be pinned because of the policy.
@@ -273,18 +340,19 @@ TYPED_TEST(HostAllocatorTest, ManualPinningOperationsWorkWithCuda)
     EXPECT_TRUE(isPinned(input));
 
     // Unpin after pin is fine.
-    input.get_allocator().getPolicy().unpin();
+    input.get_allocator().getPolicy().unpin(input.data());
     EXPECT_FALSE(isPinned(input));
 
     // Repeated unpin should be a no-op.
-    input.get_allocator().getPolicy().unpin();
+    input.get_allocator().getPolicy().unpin(input.data());
 
     // Pin after unpin is fine.
-    input.get_allocator().getPolicy().pin();
+    const size_t size = input.size() * sizeof(typename TestFixture::VectorType::value_type);
+    input.get_allocator().getPolicy().pin(input.data(), size);
     EXPECT_TRUE(isPinned(input));
 
     // Repeated pin should be a no-op, and still pinned.
-    input.get_allocator().getPolicy().pin();
+    input.get_allocator().getPolicy().pin(input.data(), size);
     EXPECT_TRUE(isPinned(input));
 
     // Switching policy to CannotBePinned must unpin the buffer (via
@@ -321,8 +389,9 @@ TYPED_TEST(HostAllocatorTest, ManualPinningOperationsWorkEvenWithoutCuda)
 
     // Since the buffer can't be pinned and isn't pinned, and the
     // calling code can't be unhappy about this, these are OK.
-    input.get_allocator().getPolicy().pin();
-    input.get_allocator().getPolicy().unpin();
+    input.get_allocator().getPolicy().pin(input.data(),
+                                          input.size()*sizeof(typename TestFixture::VectorType::value_type));
+    input.get_allocator().getPolicy().unpin(input.data());
 }
 
 #endif
index 5afb4bbf92edac78bca5949e43aa13e051403327..f67a0e9e213248f0c4d822a3efd878d97ede0bd1 100644 (file)
@@ -126,15 +126,6 @@ class Allocator : public AllocationPolicy
             typedef Allocator<U, AllocationPolicy> other; //!< Align class U with our alignment
         };
 
-        /*! \brief Templated copy constructor
-         *
-         * This template constructor cannot be auto-generated, and is
-         * normally unused, except e.g. MSVC2015 standard library uses
-         * it in debug mode, presumably to implement some checks.
-         */
-        template <class U>
-        explicit Allocator(const Allocator<U, AllocationPolicy> & /*unused*/) {}
-
         /*! \brief Constructor
          *
          * No constructor can be auto-generated in the presence of any
@@ -249,6 +240,12 @@ class Allocator : public AllocationPolicy
          */
         bool
         operator!=(const Allocator &rhs) const { return !operator==(rhs); }
+
+        //! Obtain allocator for copy construction
+        Allocator select_on_container_copy_construction() const
+        {
+            return Allocator(AllocationPolicy::select_on_container_copy_construction());
+        }
 };
 
 }      // namespace gmx