Add isIntegralConstant type trait
authorRoland Schulz <roland.schulz@intel.com>
Mon, 31 Jul 2017 23:56:00 +0000 (16:56 -0700)
committerMark Abraham <mark.j.abraham@gmail.com>
Thu, 7 Sep 2017 15:00:28 +0000 (17:00 +0200)
Add static assert to pme-gather to make overload safer

Change-Id: I71c67f3ee40185e31796752d090d2aa9f0918ec8

src/gromacs/ewald/pme-gather.cpp
src/gromacs/utility/tests/CMakeLists.txt
src/gromacs/utility/tests/typetraits.cpp [new file with mode: 0644]
src/gromacs/utility/typetraits.h [new file with mode: 0644]

index e1d6059c21a9ed526f865067dc3ef498cff80516..a1112f99bc8cf8df0f2e31545f452c43e51a1433 100644 (file)
@@ -44,6 +44,7 @@
 #include "gromacs/utility/basedefinitions.h"
 #include "gromacs/utility/gmxassert.h"
 #include "gromacs/utility/smalloc.h"
+#include "gromacs/utility/typetraits.h"
 
 #include "pme-internal.h"
 #include "pme-simd.h"
@@ -86,6 +87,9 @@ struct do_fspline
     template <typename Int>
     RVec operator()(Int order) const
     {
+        static_assert(isIntegralConstant<Int, int>::value || std::is_same<Int, int>::value,
+                      "'order' needs to be either of type integral_constant<int,N> or int.");
+
         const int  norder = nn*order;
 
         /* Pointer arithmetic alert, next six statements */
index 5fa46959d7a97c242c5156d29dac4df26c206bf7..57cf2a19f199326019a87d5197a1f361855f8fa5 100644 (file)
@@ -45,4 +45,5 @@ gmx_add_unit_test(UtilityUnitTests utility-test
                   stringutil.cpp
                   textreader.cpp
                   textwriter.cpp
+                  typetraits.cpp
                   )
diff --git a/src/gromacs/utility/tests/typetraits.cpp b/src/gromacs/utility/tests/typetraits.cpp
new file mode 100644 (file)
index 0000000..d6c8abd
--- /dev/null
@@ -0,0 +1,53 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2017, 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.
+ */
+#include "gmxpre.h"
+
+#include "gromacs/utility/typetraits.h"
+
+#include <gtest/gtest.h>
+
+namespace gmx
+{
+namespace test
+{
+TEST(TypeTraitsTest, IsIntegralConstant)
+{
+    EXPECT_FALSE(isIntegralConstant<int>::value);
+    EXPECT_TRUE((isIntegralConstant<std::integral_constant<int, 5> >::value));
+    EXPECT_FALSE((isIntegralConstant<std::integral_constant<short, 5>, long>::value));
+    EXPECT_TRUE((isIntegralConstant<std::integral_constant<long, 5>, long>::value));
+}
+}
+}
diff --git a/src/gromacs/utility/typetraits.h b/src/gromacs/utility/typetraits.h
new file mode 100644 (file)
index 0000000..2ba69ba
--- /dev/null
@@ -0,0 +1,71 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2017, 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.
+ */
+/*! \libinternal \file
+ * \brief
+ * Declares type traits
+ *
+ * \inlibraryapi
+ * \ingroup module_utility
+ */
+
+#ifndef GMX_UTILITY_TYPETRAITS_H
+#define GMX_UTILITY_TYPETRAITS_H
+
+#include <type_traits>
+
+namespace gmx
+{
+
+/*! \libinternal \brief
+ * Is true if type is a std::integral_constant
+ *
+ * If the optional integral type is given, than it is only true for a
+ * std::integral_constant of that integral type.
+ *
+ * \tparam T type to check
+ * \tparam Int optional integral type
+ */
+template<typename T, typename Int = void>
+struct isIntegralConstant : public std::false_type {};
+
+template<typename Int, Int N>
+struct isIntegralConstant<std::integral_constant<Int, N>, void> : public std::true_type {};
+
+template<typename Int, Int N>
+struct isIntegralConstant<std::integral_constant<Int, N>, Int> : public std::true_type {};
+
+} // namespace gmx
+
+#endif