Fix several errors when compiling with MSVC
authorDmitry Morozov <aracsmd@gmail.com>
Wed, 7 Apr 2021 14:18:12 +0000 (17:18 +0300)
committerAndrey Alekseenko <al42and@gmail.com>
Sat, 10 Apr 2021 22:39:24 +0000 (22:39 +0000)
src/gromacs/mdlib/update.cpp
src/gromacs/simd/simd_math.h

index 466c6275ea033ca1042521e058bdca3e1c9ec79c..1ad97b2c5e9106d2d701f1c4dac71de68c37c17b 100644 (file)
@@ -357,7 +357,7 @@ enum class ApplyParrinelloRahmanVScaling
  * aligned (and padded) memory, possibly with some hints for the compilers.
  */
 template<StoreUpdatedVelocities storeUpdatedVelocities, NumTempScaleValues numTempScaleValues, ApplyParrinelloRahmanVScaling applyPRVScaling, typename VelocityType>
-static std::enable_if_t<std::is_same<VelocityType, rvec*>::value || std::is_same<VelocityType, const rvec*>::value, void>
+static std::enable_if_t<std::is_same<VelocityType, rvec>::value || std::is_same<VelocityType, const rvec>::value, void>
 updateMDLeapfrogSimple(int                                 start,
                        int                                 nrend,
                        real                                dt,
@@ -368,7 +368,7 @@ updateMDLeapfrogSimple(int                                 start,
                        const rvec                          pRVScaleMatrixDiagonal,
                        const rvec* gmx_restrict x,
                        rvec* gmx_restrict xprime,
-                       VelocityType gmx_restrict v,
+                       VelocityType* gmx_restrict v,
                        const rvec* gmx_restrict f)
 {
     real lambdaGroup;
@@ -479,7 +479,7 @@ static inline void simdStoreRvecs(rvec* r, int index, SimdReal r0, SimdReal r1,
  * \param[in]    f                      Forces
  */
 template<StoreUpdatedVelocities storeUpdatedVelocities, typename VelocityType>
-static std::enable_if_t<std::is_same<VelocityType, rvec*>::value || std::is_same<VelocityType, const rvec*>::value, void>
+static std::enable_if_t<std::is_same<VelocityType, rvec>::value || std::is_same<VelocityType, const rvec>::value, void>
 updateMDLeapfrogSimpleSimd(int                               start,
                            int                               nrend,
                            real                              dt,
@@ -487,7 +487,7 @@ updateMDLeapfrogSimpleSimd(int                               start,
                            gmx::ArrayRef<const t_grp_tcstat> tcstat,
                            const rvec* gmx_restrict x,
                            rvec* gmx_restrict xprime,
-                           VelocityType gmx_restrict v,
+                           VelocityType* gmx_restrict v,
                            const rvec* gmx_restrict f)
 {
     SimdReal timestep(dt);
index abd001c5f55878197167250cc3d429097b9352f4..5642c45e32ac8caa588919c61151bcdadfb43c08 100644 (file)
@@ -2,7 +2,7 @@
  * This file is part of the GROMACS molecular simulation package.
  *
  * Copyright (c) 2012,2013,2014,2015,2017 by the GROMACS development team.
- * Copyright (c) 2018,2019,2020, by the GROMACS development team, led by
+ * Copyright (c) 2018,2019,2020,2021, 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.
@@ -62,6 +62,7 @@
 
 #include <limits>
 
+#include "gromacs/math/units.h"
 #include "gromacs/math/utilities.h"
 #include "gromacs/simd/simd.h"
 #include "gromacs/utility/basedefinitions.h"