Resolve C++14 TODOs
authorRoland Schulz <roland.schulz@intel.com>
Fri, 25 Jan 2019 04:23:15 +0000 (20:23 -0800)
committerMark Abraham <mark.j.abraham@gmail.com>
Mon, 28 Jan 2019 16:58:04 +0000 (17:58 +0100)
Also simplify constructor from fixed sized container.

Change-Id: I4983020939b1569949c39f7bf715fe2dbe69ab0d

src/gromacs/math/multidimarray.h
src/gromacs/mdspan/extents.h
src/gromacs/mdspan/layouts.h

index cc7f030446772eeed515d104d0fc4683ce5ccb19..fe79e08453388df916c791b7cc43c285631baf6f 100644 (file)
@@ -52,18 +52,21 @@ namespace gmx
 
 namespace detail
 {
-/*! \libinternal \brief
- * Determine static array size at compile time.
- * Statically evaluates the product of static array extents.
- * \tparam Extents Extents of an multidimensional array as in module_mdspan
- * \returns the product of the static extents
- */
-template <typename Extents>
-constexpr typename Extents::index_type staticExtentsProduct(size_t const i = 0)
-{
-    return (i < Extents::rank()) ? Extents::static_extent(i) * staticExtentsProduct<Extents>(i + 1) : 1;
-}
-}   // namespace detail
+//! Same as std::void_t from C++17
+template< class ... >
+using void_t = void;
+
+template<typename T, typename = void>
+struct is_resizable : std::false_type {};
+
+template<typename T>
+struct is_resizable<T, void_t<decltype(std::declval<T>().resize(size_t()))> > : std::true_type {};
+
+//! Type has a resize member function callable with size_t argument
+template<typename T>
+constexpr bool is_resizable_v = is_resizable<T>::value;
+}   //namespace detail
+
 /*! \libinternal \brief
  * Multidimensional array that manages its own memory.
  *
@@ -112,6 +115,11 @@ class MultiDimArray
          *  used, e.g., in free begin and end functions
          */
         using const_iterator = const typename ArrayRef<const value_type>::const_iterator;
+
+        static_assert(detail::is_resizable_v<TContainer> == (Extents::rank_dynamic() > 0),
+                      "Resizable container (e.g. std::vector) requires at least one dynamic rank. "
+                      "Non-resizable container (e.g. std::array) requires zero dynamic ranks.");
+
         /*! \brief
          * Allocate dynamic array data and set view with the dynamic extents.
          *
@@ -121,53 +129,27 @@ class MultiDimArray
          * \tparam IndexType        Parameter pack type holding the dynamic
          *                          extents of the multidimensional array
          */
-        template <class ... IndexType>
+        template <class ... IndexType, typename T = TContainer,
+                  typename = typename std::enable_if_t<detail::is_resizable_v<T> > >
         MultiDimArray(IndexType... dynamicExtent)
         {
-            static_assert(Extents::rank_dynamic() != 0,
-                          "Resizable container type not allowed in a static MultiDimArray.");
             resize(dynamicExtent ...);
         }
-        /*! \brief
-         * Construction with fixed sized array for storage if MultiDimArray size
-         * is static and layout policy allows compile time determination of the
-         * container size.
-         *
-         * Exampe use is MultiDimArray<std::array<float, 9>, extents<3,3>>
-         * \tparam StaticSizeAndLayout Template parameter for activation via SFINAE.
-         */
-        template <bool StaticSizeAndLayoutRight = Extents::rank_dynamic() == 0 &&
-                      std::is_same<LayoutPolicy, layout_right>::value,
-                  typename std::enable_if<StaticSizeAndLayoutRight, int>::type = 0>
-        constexpr MultiDimArray() noexcept : view_(data_.data())
-        {
-            // \todo replace staticExtentsProduct by the required_span_size
-            // of the LayoutPolicy, once required_span_size is a constexpr
-            // with C++14
-            // using tuple_size allows constexpr for this constructor (data_.size() is not available then)
-            static_assert(std::tuple_size<TContainer>() ==
-                          detail::staticExtentsProduct<Extents>(),
-                          "Non-resizable container type size must match static MultiDimArray size.");
-        }
         /*! \brief
          * Construction from fixed sized arrays if the array size is static and
          * layout policy allows compile time determination of the container size.
          *
          * Enables the expected initialization
          * MultiDimArray<std::array<float, 9>, extents<3,3>> arr = {{1,2...}}
-         * \tparam StaticSizeAndLayout Template parameter for activation via SFINAE.
+         * \tparam T Template parameter for activation via SFINAE.
          */
-        template <bool StaticSizeAndLayoutRight = Extents::rank_dynamic() == 0 &&
-                      std::is_same<LayoutPolicy, layout_right>::value,
-                  typename std::enable_if<StaticSizeAndLayoutRight, int>::type = 0>
-        constexpr MultiDimArray(const TContainer &data) noexcept : data_(data), view_(data_.data())
+        // SFINAE required because std::vector::size isn't constexpr and is_constexpr doesn't exist.
+        template <typename T = TContainer,
+                  typename   = typename std::enable_if_t<!detail::is_resizable_v<T> > >
+        constexpr MultiDimArray(const TContainer &data = {}
+                                ) noexcept : data_(data), view_(data_.data())
         {
-            // \todo replace staticExtentsProduct by the required_span_size
-            // of the LayoutPolicy, once required_span_size is a constexpr
-            // with C++14
-            // using tuple_size allows constexpr for this constructor (data_.size() is not available then)
-            static_assert(std::tuple_size<TContainer>() ==
-                          detail::staticExtentsProduct<Extents>(),
+            static_assert(TContainer().size() == typename view_type::mapping_type().required_span_size(),
                           "Non-resizable container type size must match static MultiDimArray size.");
         }
         //! Copy constructor
index a7a99b31ac1cb85080cafcc6a9fcdd22426095bd..beff43f94b523b5b164fff1b54162423a36ee670 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2018, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019, 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.
@@ -105,8 +105,8 @@ template< std::ptrdiff_t ... StaticExtents >
 class extents;
 
 template<std::ptrdiff_t... LHS, std::ptrdiff_t... RHS>
-bool operator==(const extents<LHS...> &lhs,
-                const extents<RHS...> &rhs) noexcept;
+constexpr bool operator==(const extents<LHS...> &lhs,
+                          const extents<RHS...> &rhs) noexcept;
 
 template<std::ptrdiff_t... LHS, std::ptrdiff_t... RHS>
 constexpr bool operator!=(const extents<LHS...> &lhs,
@@ -142,14 +142,14 @@ struct extents_analyse<R, E0, StaticExtents...> {
     next_extents_analyse next;
 
     //! Trivial constructor.
-    extents_analyse() : next() {}
+    constexpr extents_analyse() : next() {}
 
     /*! \brief Construction from dynamic extents hands the extents down
      * to the next extents analysis of lower rank.
      * \param[in] de dynamic extents
      */
     template<class ... DynamicExtents>
-    extents_analyse(DynamicExtents... de) : next(de ...) {}
+    constexpr extents_analyse(DynamicExtents... de) : next(de ...) {}
 
     /*! \brief Construct from an array of dynamic extentes and rank.
      * Hand down the dynamic rank parameters to the next extents analysis rank
@@ -157,7 +157,7 @@ struct extents_analyse<R, E0, StaticExtents...> {
      * \param[in] r rank to read from the dynamic extent
      */
     template<std::size_t Rank>
-    extents_analyse(const std::array<std::ptrdiff_t, Rank> &de, const std::size_t r) : next(de, r) {}
+    constexpr extents_analyse(const std::array<std::ptrdiff_t, Rank> &de, const std::size_t r) : next(de, r) {}
 
     //! Copy constructor.
     template<std::ptrdiff_t... OtherStaticExtents>
@@ -177,8 +177,6 @@ struct extents_analyse<R, E0, StaticExtents...> {
      */
     constexpr std::ptrdiff_t extent(const std::size_t r) const noexcept
     {
-        // TODO use early return and if instead of ternary operator
-        // after bumping requirements to C++14
         return (r == R) ? E0 : next.extent(r);
     }
     /*! \brief Report the static extent of dimension r.
@@ -187,8 +185,6 @@ struct extents_analyse<R, E0, StaticExtents...> {
      */
     static constexpr std::ptrdiff_t static_extent(const std::size_t r) noexcept
     {
-        // TODO use early return and if instead of ternary operator
-        // after bumping requirements to C++14
         return (r == R) ? E0 : next_extents_analyse::static_extent(r);
     }
 };
@@ -252,7 +248,7 @@ struct extents_analyse<R, dynamic_extent, StaticExtents...> {
      */
     constexpr std::ptrdiff_t extent(const std::size_t r) const noexcept
     {
-        return (r == R) ? this_extent : next.extent(r); // TODO use early return and if instead of ternary operator after bumping requirements to C++14
+        return (r == R) ? this_extent : next.extent(r);
     }
     /*! \brief Report the static extent of dimension r.
      * \param[in] r the dimension to query
@@ -260,7 +256,7 @@ struct extents_analyse<R, dynamic_extent, StaticExtents...> {
      */
     static constexpr std::ptrdiff_t static_extent(const std::size_t r) noexcept
     {
-        return (r == R) ? dynamic_extent : next_extents_analyse::static_extent(r); // TODO use early return and if instead of ternary operator after bumping requirements to C++14
+        return (r == R) ? dynamic_extent : next_extents_analyse::static_extent(r);
     }
 };
 
@@ -279,7 +275,7 @@ struct extents_analyse<0> {
     static constexpr std::size_t rank_dynamic() noexcept { return 0; }
 
     //! Trivial constructor.
-    extents_analyse() {}
+    constexpr extents_analyse() {}
 
     //! Construct from array and rank, doing nothing.
     template<std::size_t Rank>
@@ -402,11 +398,10 @@ class extents
 
 /*! \brief Comparison operator.
  * \returns true if extents are equal
- * TODO add constexpr when C++14 is required
  */
 template<std::ptrdiff_t... LHS, std::ptrdiff_t... RHS>
-bool operator==(const extents<LHS...> &lhs,
-                const extents<RHS...> &rhs) noexcept
+constexpr bool operator==(const extents<LHS...> &lhs,
+                          const extents<RHS...> &rhs) noexcept
 {
     bool equal = lhs.rank() == rhs.rank();
     for (std::size_t r = 0; r < lhs.rank(); r++)
index 7addb27e7bab5da27dff15dbf4aefb0f4981d0a8..4e1cda739dde4949e6a259ac44787177a5d7cdaf 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2018, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019, 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.
@@ -166,11 +166,10 @@ class layout_right
             public:
                 /*! \brief Return the size of the underlying one-dimensional
                  * data structure, so that the mapping is always valid.
-                 * \todo add constexpr when C++14 is required
                  *
                  * \returns number of span elements
                  */
-                index_type required_span_size() const noexcept
+                constexpr index_type required_span_size() const noexcept
                 {
                     index_type size = 1;
                     for (size_t r = 0; r < m_extents.rank(); r++)
@@ -185,11 +184,10 @@ class layout_right
                  * \tparam Indices type of the indices to be mapped
                  * \param[in] indices the indices to be mapped
                  * \returns One-dimensional integer index.
-                 * \todo C++14 activation of constexpr
                  */
                 template<class ... Indices >
                 typename std::enable_if<sizeof ... (Indices) == Extents::rank(), index_type>::type
-                operator()( Indices ... indices ) const noexcept
+                constexpr operator()( Indices ... indices ) const noexcept
                 { return offset( 0, 0, indices ... ); }
 
                 //! Report that this mapping is always unique.
@@ -208,9 +206,8 @@ class layout_right
                 /*!\brief Return the stride of dimension r.
                  * \param[in] R rank of the stride to be queried.
                  * \returns the stride along dimension r.
-                 * \todo C++14 activation of constexpr
                  */
-                index_type stride(const size_t R) const noexcept
+                constexpr index_type stride(const size_t R) const noexcept
                 {
                     ptrdiff_t stride = 1;
                     for (size_t r = m_extents.rank()-1; r > R; r--)