Narrow-down Intel-specific #ifdef's
authorAndrey Alekseenko <al42and@gmail.com>
Thu, 29 Oct 2020 13:25:44 +0000 (13:25 +0000)
committerAndrey Alekseenko <al42and@gmail.com>
Thu, 29 Oct 2020 13:25:44 +0000 (13:25 +0000)
Some of the things hidden by `#ifdef __INTEL_COMPILER` does not seem to be relevant as of 2021.1-beta09/10 (ICC 2021.1 Beta 20200827).

Also, was getting a false-positive "warning: result of comparison of 0-bit unsigned value != 0 is always false" for line `if (grid_base[g_baseNR - 1] % 4 != 0)`. Replaced it with `static_assert`, because why not.

src/gromacs/fft/calcgrid.cpp
src/gromacs/listed_forces/tests/bonded.cpp
src/gromacs/pulling/pullutil.cpp
src/gromacs/simd/tests/simd_floatingpoint_util.cpp

index 8d293464b7b069a005aecb843a6689d8fbcfcc4f..de0bbb1b9a10eeb98b1d857a37e79e2238b5ea44 100644 (file)
  */
 
 /* Small grid size array */
-#define g_initNR 15
-const int grid_init[g_initNR] = { 6, 8, 10, 12, 14, 16, 20, 24, 25, 28, 32, 36, 40, 42, 44 };
+constexpr int g_initNR            = 15;
+constexpr int grid_init[g_initNR] = { 6, 8, 10, 12, 14, 16, 20, 24, 25, 28, 32, 36, 40, 42, 44 };
 
 /* For larger grid sizes, a prefactor with any power of 2 can be added.
  * Only sizes divisible by 4 should be used, 90 is allowed, 140 not.
  */
-#define g_baseNR 14
-const int grid_base[g_baseNR] = { 45, 48, 50, 52, 54, 56, 60, 64, 70, 72, 75, 80, 81, 84 };
+constexpr int g_baseNR            = 14;
+constexpr int grid_base[g_baseNR] = { 45, 48, 50, 52, 54, 56, 60, 64, 70, 72, 75, 80, 81, 84 };
 
 real calcFftGrid(FILE* fp, const matrix box, real gridSpacing, int minGridPointsPerDim, int* nx, int* ny, int* nz)
 {
@@ -74,10 +74,8 @@ real calcFftGrid(FILE* fp, const matrix box, real gridSpacing, int minGridPoints
         gmx_fatal(FARGS, "invalid fourier grid spacing: %g", gridSpacing);
     }
 
-    if (grid_base[g_baseNR - 1] % 4 != 0)
-    {
-        gmx_incons("the last entry in grid_base is not a multiple of 4");
-    }
+    static_assert(grid_base[g_baseNR - 1] % 4 == 0,
+                  "the last entry in grid_base is not a multiple of 4");
 
     /* New grid calculation setup:
      *
index 41c74f9294a23e40c02ce3fa6b5821afc4f4c1df..9cd1273cb9bc04a74ce27358273203233ca114b7 100644 (file)
@@ -770,8 +770,8 @@ std::vector<PaddedVector<RVec>> c_coordinatesForTestsZeroAngle = {
 //! PBC values for testing
 std::vector<PbcType> c_pbcForTests = { PbcType::No, PbcType::XY, PbcType::Xyz };
 
-// Those tests give errors with the intel compiler and nothing else, so we disable them only there.
-#ifndef __INTEL_COMPILER
+// Those tests give errors with the Intel compiler (as of October 2019) and nothing else, so we disable them only there.
+#if !defined(__INTEL_COMPILER) || (__INTEL_COMPILER >= 2021)
 INSTANTIATE_TEST_CASE_P(Bond,
                         ListedForcesTest,
                         ::testing::Combine(::testing::ValuesIn(c_InputBonds),
index a1c646c826ead2e847862b71783b96de6ff48b2a..875b535560fa85ffa334e41a22441c9b52d5f58f 100644 (file)
@@ -519,8 +519,9 @@ static void sum_com_part_cosweight(const pull_group_work_t* pgrp,
 
 /* calculates center of mass of selection index from all coordinates x */
 // Compiler segfault with 2019_update_5 and 2020_initial
-#if defined(__INTEL_COMPILER) \
-        && ((__INTEL_COMPILER == 1900 && __INTEL_COMPILER_UPDATE >= 5) || __INTEL_COMPILER >= 1910)
+#if defined(__INTEL_COMPILER)                                          \
+        && ((__INTEL_COMPILER == 1900 && __INTEL_COMPILER_UPDATE >= 5) \
+            || (__INTEL_COMPILER >= 1910 && __INTEL_COMPILER < 2021))
 #    pragma intel optimization_level 2
 #endif
 void pull_calc_coms(const t_commrec* cr, pull_t* pull, const real* masses, t_pbc* pbc, double t, const rvec x[], rvec* xp)
index ee43c3cada45242f2fe5c960048e3b7ae7c7dac1..502f5f0961f5ffcc3f192048c3e38970d25737fa 100644 (file)
@@ -466,7 +466,7 @@ TEST_F(SimdFloatingpointUtilTest, transposeScatterDecrU3Overlapping)
         mem0_[j] = refmem[j] = (1000.0 + j) * (1.0 + 100 * GMX_REAL_EPS);
     }
 
-#    ifdef __INTEL_COMPILER // Bug in (at least) 19u1 and 18u5 (03424712)
+#    if defined(__INTEL_COMPILER) && (__INTEL_COMPILER < 2021) // Bug in (at least) 19u1 and 18u5 (03424712)
 #        pragma novector
 #    endif
     for (std::size_t j = 0; j < GMX_SIMD_REAL_WIDTH; j++)