Fix multiple bugs detected on Windows
authorRoland Schulz <roland.schulz@intel.com>
Sun, 30 Sep 2018 19:09:03 +0000 (12:09 -0700)
committerSzilárd Páll <pall.szilard@gmail.com>
Fri, 5 Oct 2018 12:51:30 +0000 (14:51 +0200)
Change-Id: Ibcc1bcedcb292a0ec8f05069733adb2d127ddc5e

cmake/gmxManageSimd.cmake
cmake/gmxTestCXX11.cmake
docs/install-guide/index.rst
src/gromacs/CMakeLists.txt
src/gromacs/awh/grid.cpp
src/gromacs/awh/tests/bias.cpp
src/gromacs/ewald/pme-solve.cu
src/gromacs/gmxpreprocess/pdb2top.cpp
src/gromacs/gpu_utils/pinning.cu
src/gromacs/mdlib/nbnxn_cuda/nbnxn_cuda_kernel_utils.cuh
src/gromacs/utility/bitmask.h

index 5a6d91a0d872ba520ffd28e37d18e71a66fc6eea..4bfbe799c7d53da3672d1813afcb2ef2782d0ffc 100644 (file)
@@ -332,14 +332,16 @@ if(GMX_ENABLE_AVX512_TESTS AND
     # Since we might be overriding AVX2 architecture flags with the AVX512 flags for the
     # files where it is used, we also check for a flag not to warn about the first (unused) arch.
     # To avoid spamming the user with lots of gromacs tests we just call the CMake flag test directly.
     # Since we might be overriding AVX2 architecture flags with the AVX512 flags for the
     # files where it is used, we also check for a flag not to warn about the first (unused) arch.
     # To avoid spamming the user with lots of gromacs tests we just call the CMake flag test directly.
-    foreach(_testflag "-Wno-unused-command-line-argument" "-wd10121")
-        string(REGEX REPLACE "[^a-zA-Z0-9]+" "_" FLAG_ACCEPTED_VARIABLE "${_testflag}_FLAG_ACCEPTED")
-        check_cxx_compiler_flag("${_testflag}" ${FLAG_ACCEPTED_VARIABLE})
-        if(${FLAG_ACCEPTED_VARIABLE})
-            set(CXX_NO_UNUSED_OPTION_WARNING_FLAGS "${_testflag}")
-            break()
-        endif()
-    endforeach(_testflag)
+    if (NOT MSVC)
+        foreach(_testflag "-Wno-unused-command-line-argument" "-wd10121")
+            string(REGEX REPLACE "[^a-zA-Z0-9]+" "_" FLAG_ACCEPTED_VARIABLE "${_testflag}_FLAG_ACCEPTED")
+            check_cxx_compiler_flag("${_testflag}" ${FLAG_ACCEPTED_VARIABLE})
+            if(${FLAG_ACCEPTED_VARIABLE})
+                set(CXX_NO_UNUSED_OPTION_WARNING_FLAGS "${_testflag}")
+                break()
+            endif()
+        endforeach(_testflag)
+    endif()
 endif()
 
 # By default, 32-bit windows cannot pass SIMD (SSE/AVX) arguments in registers,
 endif()
 
 # By default, 32-bit windows cannot pass SIMD (SSE/AVX) arguments in registers,
index 2e02ffba742cc6901f7bc2e1d2aa5a3090638e58..9a2b296813491851cdc3e5183d634557baedea9b 100644 (file)
@@ -46,7 +46,7 @@ function(GMX_TEST_CXX11 CXX11_CXX_FLAG_NAME STDLIB_CXX_FLAG_NAME STDLIB_LIBRARIE
     # First check that the compiler is OK, and find the appropriate flag.
 
     if(WIN32 AND NOT MINGW)
     # First check that the compiler is OK, and find the appropriate flag.
 
     if(WIN32 AND NOT MINGW)
-        set(CXX11_CXX_FLAG "/Qstd=c++11")
+        set(CXX11_CXX_FLAG "/std:c++14 /Zc:__cplusplus")
     elseif(CYGWIN)
         set(CXX11_CXX_FLAG "-std=gnu++11") #required for strdup
     else()
     elseif(CYGWIN)
         set(CXX11_CXX_FLAG "-std=gnu++11") #required for strdup
     else()
index a0bee47c968eb6bf10e5f92dc476646d54bd2730..bdd61c3ef11b272d23b471c815b6958796ddea45 100644 (file)
@@ -117,7 +117,7 @@ compiler versions are
 * GNU (gcc) 4.8.1
 * Intel (icc) 17.0.1
 * LLVM (clang) 3.3
 * GNU (gcc) 4.8.1
 * Intel (icc) 17.0.1
 * LLVM (clang) 3.3
-* Microsoft (MSVC) 2017
+* Microsoft (MSVC) 2017 (C++14 is used)
 
 Other compilers may work (Cray, Pathscale, older clang) but do
 not offer competitive performance. We recommend against PGI because
 
 Other compilers may work (Cray, Pathscale, older clang) but do
 not offer competitive performance. We recommend against PGI because
index b2ee7eb1664bf396cbb1c2afbc2e236155a3ad23..b305b92326bf45e85b5baac028c4ba3d9d0080a9 100644 (file)
@@ -248,7 +248,9 @@ check_cxx_compiler_flag("-Wno-unused -Wno-unused-parameter" HAS_NO_UNUSED)
 check_cxx_compiler_flag(-Wno-missing-declarations HAS_NO_MISSING_DECL)
 check_cxx_compiler_flag(-Wno-missing-prototypes HAS_NO_MISSING_PROTO)
 check_cxx_compiler_flag(/wd4101 HAS_NO_MSVC_UNUSED)
 check_cxx_compiler_flag(-Wno-missing-declarations HAS_NO_MISSING_DECL)
 check_cxx_compiler_flag(-Wno-missing-prototypes HAS_NO_MISSING_PROTO)
 check_cxx_compiler_flag(/wd4101 HAS_NO_MSVC_UNUSED)
-check_cxx_compiler_flag(-wd1419 HAS_DECL_IN_SOURCE)
+if (NOT MSVC)
+    check_cxx_compiler_flag(-wd1419 HAS_DECL_IN_SOURCE)
+endif()
 if (HAS_NO_UNUSED)
     target_compile_options(libgromacs_generated PRIVATE "-Wno-unused;-Wno-unused-parameter")
 endif()
 if (HAS_NO_UNUSED)
     target_compile_options(libgromacs_generated PRIVATE "-Wno-unused;-Wno-unused-parameter")
 endif()
@@ -312,6 +314,7 @@ if (CMAKE_CXX_COMPILER_ID STREQUAL "MSVC")
      /wd28199 #uninitialized memory
      # For compile time constant (e.g. templates) the following warnings have flase postives
      /wd6239  #(<non-zero> && <expr>)
      /wd28199 #uninitialized memory
      # For compile time constant (e.g. templates) the following warnings have flase postives
      /wd6239  #(<non-zero> && <expr>)
+     /wd6240  #(<expr> && <non-zero>)
      /wd6294  #Ill-defined for-loop
      /wd6326  #comparison of constant with other constant
      /wd28020 #expression involving paramter is not true
      /wd6294  #Ill-defined for-loop
      /wd6326  #comparison of constant with other constant
      /wd28020 #expression involving paramter is not true
index 737b31f4f2268354a3f5c19a055989ffbfb56526..6cebf587febaaf5a17996ffc184263594b2fa042 100644 (file)
@@ -208,7 +208,6 @@ void linearArrayIndexToMultiDim(int indexLinear, int numDimensions, const awh_iv
     {
         int stride = 1;
 
     {
         int stride = 1;
 
-        /* Workaround for bug in clang */
         for (int k = d + 1; k < numDimensions; k++)
         {
             stride *= numPointsDim[k];
         for (int k = d + 1; k < numDimensions; k++)
         {
             stride *= numPointsDim[k];
index 84fa8bf22d6db2971ad17200dbcc3fc22818da2c..0a6ad40daadc407a07da92c904b36b91fdc490f7 100644 (file)
@@ -65,6 +65,15 @@ namespace test
  */
 struct AwhTestParameters
 {
  */
 struct AwhTestParameters
 {
+    AwhTestParameters() = default;
+    //!Move constructor
+    AwhTestParameters(AwhTestParameters &&o) noexcept : beta(o.beta), awhDimParams(o.awhDimParams),
+                                                        awhBiasParams(o.awhBiasParams), awhParams(o.awhParams),
+                                                        dimParams(std::move(o.dimParams))
+    {
+        awhBiasParams.dimParams = &awhDimParams;
+        awhParams.awhBiasParams = &awhBiasParams;
+    }
     double                 beta;          //!< 1/(kB*T)
 
     AwhDimParams           awhDimParams;  //!< Dimension parameters pointed to by \p awhBiasParams
     double                 beta;          //!< 1/(kB*T)
 
     AwhDimParams           awhDimParams;  //!< Dimension parameters pointed to by \p awhBiasParams
index 54cceab536caf938bb9949930645366bb146acc2..b163ddb3a84c9e069c4cb3f034f1cd1ee9dbd458 100644 (file)
@@ -43,6 +43,8 @@
 
 #include <cassert>
 
 
 #include <cassert>
 
+#include <math_constants.h>
+
 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
 
 #include "pme.cuh"
 #include "gromacs/gpu_utils/cuda_arch_utils.cuh"
 
 #include "pme.cuh"
@@ -203,7 +205,7 @@ __global__ void pme_solve_kernel(const struct PmeGpuCudaKernelParams kernelParam
             const float m2k        = mhxk * mhxk + mhyk * mhyk + mhzk * mhzk;
             assert(m2k != 0.0f);
             //TODO: use LDG/textures for gm_splineValue
             const float m2k        = mhxk * mhxk + mhyk * mhyk + mhzk * mhzk;
             assert(m2k != 0.0f);
             //TODO: use LDG/textures for gm_splineValue
-            float       denom = m2k * float(M_PI) * kernelParams.current.boxVolume * gm_splineValueMajor[kMajor] * gm_splineValueMiddle[kMiddle] * gm_splineValueMinor[kMinor];
+            float       denom = m2k * float(CUDART_PI_F) * kernelParams.current.boxVolume * gm_splineValueMajor[kMajor] * gm_splineValueMiddle[kMiddle] * gm_splineValueMinor[kMinor];
             assert(isfinite(denom));
             assert(denom != 0.0f);
             const float   tmp1   = expf(-kernelParams.grid.ewaldFactor * m2k);
             assert(isfinite(denom));
             assert(denom != 0.0f);
             const float   tmp1   = expf(-kernelParams.grid.ewaldFactor * m2k);
index 9eb8dcba84e1826263077449f88baa0ccf65d1e0..a505465eeb0e62a53a41c972dab092b773c77995 100644 (file)
@@ -642,25 +642,13 @@ void print_top_mols(FILE *out,
                     const char *title, const char *ffdir, const char *water,
                     int nincl, char **incls, int nmol, t_mols *mols)
 {
                     const char *title, const char *ffdir, const char *water,
                     int nincl, char **incls, int nmol, t_mols *mols)
 {
-    int   i;
-    char *incl;
 
     if (nincl > 0)
     {
         fprintf(out, "; Include chain topologies\n");
 
     if (nincl > 0)
     {
         fprintf(out, "; Include chain topologies\n");
-        for (i = 0; (i < nincl); i++)
+        for (int i = 0; i < nincl; i++)
         {
         {
-            incl = strrchr(incls[i], DIR_SEPARATOR);
-            if (incl == nullptr)
-            {
-                incl = incls[i];
-            }
-            else
-            {
-                /* Remove the path from the include name */
-                incl = incl + 1;
-            }
-            fprintf(out, "#include \"%s\"\n", incl);
+            fprintf(out, "#include \"%s\"\n", gmx::Path::getFilename(incls[i]).c_str());
         }
         fprintf(out, "\n");
     }
         }
         fprintf(out, "\n");
     }
@@ -675,7 +663,7 @@ void print_top_mols(FILE *out,
     {
         fprintf(out, "[ %s ]\n", dir2str(d_molecules));
         fprintf(out, "; %-15s %5s\n", "Compound", "#mols");
     {
         fprintf(out, "[ %s ]\n", dir2str(d_molecules));
         fprintf(out, "; %-15s %5s\n", "Compound", "#mols");
-        for (i = 0; (i < nmol); i++)
+        for (int i = 0; i < nmol; i++)
         {
             fprintf(out, "%-15s %5d\n", mols[i].name, mols[i].nr);
         }
         {
             fprintf(out, "%-15s %5d\n", mols[i].name, mols[i].nr);
         }
index bf140b064752dcceb66fe4367a0c0906621741d3..9aa6caf4e0a40b8d047f7b38b1e6be174b382053 100644 (file)
@@ -54,7 +54,7 @@
 
 namespace gmx
 {
 
 namespace gmx
 {
-
+#pragma diag_suppress 177 // No unused attribute which works with both MSVC and NVCC
 //! Is \c ptr aligned on a boundary that is a multiple of \c bytes.
 gmx_unused static inline bool isAligned(const void *ptr, size_t bytes)
 {
 //! Is \c ptr aligned on a boundary that is a multiple of \c bytes.
 gmx_unused static inline bool isAligned(const void *ptr, size_t bytes)
 {
index 082bab0ac76b69df964eb645808186a2070aa32f..3d0465d71caf222bb25db0087e324f28794aeb0c 100644 (file)
 
 /*! \brief Log of the i and j cluster size.
  *  change this together with c_clSize !*/
 
 /*! \brief Log of the i and j cluster size.
  *  change this together with c_clSize !*/
-static const int          c_clSizeLog2  = 3;
+static const int __device__          c_clSizeLog2  = 3;
 /*! \brief Square of cluster size. */
 /*! \brief Square of cluster size. */
-static const int          c_clSizeSq    = c_clSize*c_clSize;
+static const int __device__          c_clSizeSq    = c_clSize*c_clSize;
 /*! \brief j-cluster size after split (4 in the current implementation). */
 /*! \brief j-cluster size after split (4 in the current implementation). */
-static const int          c_splitClSize = c_clSize/c_nbnxnGpuClusterpairSplit;
+static const int __device__          c_splitClSize = c_clSize/c_nbnxnGpuClusterpairSplit;
 /*! \brief Stride in the force accumualation buffer */
 /*! \brief Stride in the force accumualation buffer */
-static const int          c_fbufStride  = c_clSizeSq;
+static const int __device__          c_fbufStride  = c_clSizeSq;
 /*! \brief i-cluster interaction mask for a super-cluster with all c_numClPerSupercl=8 bits set */
 /*! \brief i-cluster interaction mask for a super-cluster with all c_numClPerSupercl=8 bits set */
-static const unsigned     superClInteractionMask = ((1U << c_numClPerSupercl) - 1U);
+static const unsigned __device__     superClInteractionMask = ((1U << c_numClPerSupercl) - 1U);
 
 
-static const float        c_oneSixth    = 0.16666667f;
-static const float        c_oneTwelveth = 0.08333333f;
+static const float __device__        c_oneSixth    = 0.16666667f;
+static const float __device__        c_oneTwelveth = 0.08333333f;
 
 
 /*! Convert LJ sigma,epsilon parameters to C6,C12. */
 
 
 /*! Convert LJ sigma,epsilon parameters to C6,C12. */
index 1b3fe93f82f068b81f8a87809c7ef1b5544deb38..cfae126c940f70b567b510e52a549398fcd5d832 100644 (file)
@@ -147,7 +147,7 @@ inline static void bitmask_init_low_bits(gmx_bitmask_t* m, int b)
 {
     memset(m->data(), 255, b/64*8);
     (*m)[b/64] = (static_cast<uint64_t>(1) << (b%64)) - 1;
 {
     memset(m->data(), 255, b/64*8);
     (*m)[b/64] = (static_cast<uint64_t>(1) << (b%64)) - 1;
-    memset(&(*m)[b/64+1], 0, (BITMASK_ALEN-b/64-1)*8);
+    memset(m->data()+(b/64+1), 0, (BITMASK_ALEN-b/64-1)*8);
 }
 
 inline static bool bitmask_is_set(gmx_bitmask_t m, int b)
 }
 
 inline static bool bitmask_is_set(gmx_bitmask_t m, int b)