renamed gmx_x86_simd_macros.h to gmx_simd_macros.h
authorBerk Hess <hess@kth.se>
Wed, 12 Dec 2012 15:25:08 +0000 (16:25 +0100)
committerBerk Hess <hess@kth.se>
Tue, 18 Dec 2012 12:25:35 +0000 (13:25 +0100)
This is in preparation for non-x86 SIMD acceleration.

Change-Id: Idc652236f4b2e0f48e759d579c798c1b0dc6944f

include/gmx_simd_macros.h [moved from include/gmx_x86_simd_macros.h with 95% similarity]
src/mdlib/nbnxn_atomdata.c
src/mdlib/nbnxn_internal.h
src/mdlib/nbnxn_kernels/nbnxn_kernel_simd_2xnn.c
src/mdlib/nbnxn_kernels/nbnxn_kernel_simd_2xnn_outer.h
src/mdlib/nbnxn_kernels/nbnxn_kernel_simd_4xn.c
src/mdlib/nbnxn_kernels/nbnxn_kernel_simd_4xn_outer.h
src/mdlib/nbnxn_kernels/nbnxn_kernel_simd_utils.h
src/mdlib/nbnxn_search_simd_2xnn.h
src/mdlib/nbnxn_search_simd_4xn.h

similarity index 95%
rename from include/gmx_x86_simd_macros.h
rename to include/gmx_simd_macros.h
index e1b1403c24c16e702d7dd8a46039e2de76daf03c..2cd5ed91cb987e38ea40b554b22f042ace3c7c7c 100644 (file)
  * the research papers on the package. Check out http://www.gromacs.org.
  */
 
+/* The macros in this file are intended to be used for writing
+ * architecture independent SIMD intrinsics code.
+ * To support a new architecture, adding macros here should be (nearly)
+ * all that is needed.
+ */
+
 /* Undefine all defines used below so we can include this file multiple times
  * with different settings from the same source file.
  */
 
 /* NOTE: floor and blend are NOT available with SSE2 only acceleration */
 
-#undef GMX_X86_SIMD_WIDTH_HERE
+#undef GMX_SIMD_WIDTH_HERE
 
 #undef gmx_epi32
 
 
 #include "gmx_x86_simd_single.h"
 
-#define GMX_X86_SIMD_WIDTH_HERE  4
+#define GMX_SIMD_WIDTH_HERE  4
 
 #define gmx_mm_pr  __m128
 
 
 #include "gmx_x86_simd_double.h"
 
-#define GMX_X86_SIMD_WIDTH_HERE  2
+#define GMX_SIMD_WIDTH_HERE  2
 
 #define gmx_mm_pr  __m128d
 
 
 #include "gmx_x86_simd_single.h"
 
-#define GMX_X86_SIMD_WIDTH_HERE  8
+#define GMX_SIMD_WIDTH_HERE  8
 
 #define gmx_mm_pr  __m256
 
 
 #include "gmx_x86_simd_double.h"
 
-#define GMX_X86_SIMD_WIDTH_HERE  4
+#define GMX_SIMD_WIDTH_HERE  4
 
 #define gmx_mm_pr  __m256d
 
index ad231287cf06f34a0cfbe19c9107ffcb29f12b03..a11f764914917d3c6019c117815c1567a7f14f95 100644 (file)
@@ -1042,14 +1042,14 @@ nbnxn_atomdata_reduce_reals_x86_simd(real * gmx_restrict dest,
 #else
 #define GMX_MM128_HERE
 #endif
-#include "gmx_x86_simd_macros.h"
+#include "gmx_simd_macros.h"
 
     int       i,s;
     gmx_mm_pr dest_SSE,src_SSE;
 
     if (bDestSet)
     {
-        for(i=i0; i<i1; i+=GMX_X86_SIMD_WIDTH_HERE)
+        for(i=i0; i<i1; i+=GMX_SIMD_WIDTH_HERE)
         {
             dest_SSE = gmx_load_pr(dest+i);
             for(s=0; s<nsrc; s++)
@@ -1062,7 +1062,7 @@ nbnxn_atomdata_reduce_reals_x86_simd(real * gmx_restrict dest,
     }
     else
     {
-        for(i=i0; i<i1; i+=GMX_X86_SIMD_WIDTH_HERE)
+        for(i=i0; i<i1; i+=GMX_SIMD_WIDTH_HERE)
         {
             dest_SSE = gmx_load_pr(src[0]+i);
             for(s=1; s<nsrc; s++)
index 8e9f005b7c68c87395f609208f029c42b8004234..1a7f5bfcea634b6722f5851c9198f2b9f0a692d7 100644 (file)
@@ -105,7 +105,7 @@ typedef struct {
 #error "unsupported GMX_NBNXN_SIMD_BITWIDTH"
 #endif
 #endif
-#include "gmx_x86_simd_macros.h"
+#include "gmx_simd_macros.h"
 
 typedef struct nbnxn_x_ci_simd_4xn {
     /* The i-cluster coordinates for simple search */
index a294bad9c05ba54f23eaf19bcc9ea13c2dee15b3..0106d4d7801b1d716ae1199a0cea87c61b143a48 100644 (file)
@@ -149,8 +149,8 @@ static void reduce_group_energies(int ng,int ng_2log,
                                   const real *VSvdw,const real *VSc,
                                   real *Vvdw,real *Vc)
 {
-    const int simd_width   = GMX_X86_SIMD_WIDTH_HERE;
-    const int unrollj_half = GMX_X86_SIMD_WIDTH_HERE/4;
+    const int simd_width   = GMX_SIMD_WIDTH_HERE;
+    const int unrollj_half = GMX_SIMD_WIDTH_HERE/4;
     int ng_p2,i,j,j0,j1,c,s;
 
     ng_p2 = (1<<ng_2log);
index d36a99f15ca94facd7ffa2210f640c92553e0572..faa445efbfb1fb2ef9d67becdbc88ed1a9b46de2 100644 (file)
  */
 
 /* GMX_MM128_HERE or GMX_MM256_HERE should be set before including this file */
-#include "gmx_x86_simd_macros.h"
+#include "gmx_simd_macros.h"
 
 #define SUM_SIMD4(x) (x[0]+x[1]+x[2]+x[3])
 
 #define UNROLLI    NBNXN_CPU_CLUSTER_I_SIZE
-#define UNROLLJ    (GMX_X86_SIMD_WIDTH_HERE/2)
+#define UNROLLJ    (GMX_SIMD_WIDTH_HERE/2)
 
 #if defined GMX_MM128_HERE || defined GMX_DOUBLE
 #define STRIDE     4
index 640509d75c33a06f80bbe9ad0f6435f10377b3ee..470d27dcf7f96b6412979f1ecbbd3fd7f7dc4319 100644 (file)
@@ -149,8 +149,8 @@ static void reduce_group_energies(int ng,int ng_2log,
                                   const real *VSvdw,const real *VSc,
                                   real *Vvdw,real *Vc)
 {
-    const int simd_width   = GMX_X86_SIMD_WIDTH_HERE;
-    const int unrollj_half = GMX_X86_SIMD_WIDTH_HERE/2;
+    const int simd_width   = GMX_SIMD_WIDTH_HERE;
+    const int unrollj_half = GMX_SIMD_WIDTH_HERE/2;
     int ng_p2,i,j,j0,j1,c,s;
 
     ng_p2 = (1<<ng_2log);
index b347b3785f1f91fb42cd95a36267050a34455444..1ab915deaecc3e5670ac0de2fffbd6b7da69bf78 100644 (file)
  */
 
 /* GMX_MM128_HERE or GMX_MM256_HERE should be set before including this file */
-#include "gmx_x86_simd_macros.h"
+#include "gmx_simd_macros.h"
 
 #define SUM_SIMD4(x) (x[0]+x[1]+x[2]+x[3])
 
 #define UNROLLI    NBNXN_CPU_CLUSTER_I_SIZE
-#define UNROLLJ    GMX_X86_SIMD_WIDTH_HERE
+#define UNROLLJ    GMX_SIMD_WIDTH_HERE
 
 #if defined GMX_MM128_HERE || defined GMX_DOUBLE
 #define STRIDE     4
index d4f964ea645f49c7d01c1f2a94a6201b0bd26d10..45ab2aedcc9205f3d7d86d41cb3fbef186d09ef2 100644 (file)
@@ -511,12 +511,12 @@ static inline void add_ener_grp(gmx_mm_pr e_SSE,real *v,const int *offset_jj)
     {
         gmx_mm_pr v_SSE;
 
-        v_SSE = gmx_load_pr(v+offset_jj[jj]+jj*GMX_X86_SIMD_WIDTH_HERE);
-        gmx_store_pr(v+offset_jj[jj]+jj*GMX_X86_SIMD_WIDTH_HERE,gmx_add_pr(v_SSE,e_SSE));
+        v_SSE = gmx_load_pr(v+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE);
+        gmx_store_pr(v+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE,gmx_add_pr(v_SSE,e_SSE));
     }
 }
 
-#if defined GMX_X86_AVX_256 && GMX_X86_SIMD_WIDTH_HERE == 8
+#if defined GMX_X86_AVX_256 && GMX_SIMD_WIDTH_HERE == 8
 /* As add_ener_grp above, but for two groups of UNROLLJ/2 stored in
  * a single SIMD register.
  */
@@ -533,15 +533,15 @@ static inline void add_ener_grp_halves(gmx_mm_pr e_SSE,
     {
         gmx_mm_hpr v_SSE;
 
-        v_SSE = gmx_load_hpr(v0+offset_jj[jj]+jj*GMX_X86_SIMD_WIDTH_HERE/2);
-        gmx_store_hpr(v0+offset_jj[jj]+jj*GMX_X86_SIMD_WIDTH_HERE/2,gmx_add_hpr(v_SSE,e_SSE0));
+        v_SSE = gmx_load_hpr(v0+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2);
+        gmx_store_hpr(v0+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2,gmx_add_hpr(v_SSE,e_SSE0));
     }
     for(jj=0; jj<(UNROLLJ/2); jj++)
     {
         gmx_mm_hpr v_SSE;
 
-        v_SSE = gmx_load_hpr(v1+offset_jj[jj]+jj*GMX_X86_SIMD_WIDTH_HERE/2);
-        gmx_store_hpr(v1+offset_jj[jj]+jj*GMX_X86_SIMD_WIDTH_HERE/2,gmx_add_hpr(v_SSE,e_SSE1));
+        v_SSE = gmx_load_hpr(v1+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2);
+        gmx_store_hpr(v1+offset_jj[jj]+jj*GMX_SIMD_WIDTH_HERE/2,gmx_add_hpr(v_SSE,e_SSE1));
     }
 }
 #endif
index baf6f7bbb94920169c1910956d2a2cb6010968e8..04dd50105d212827f568d11f7b38cac2b5056c2a 100644 (file)
 #error "unsupported GMX_NBNXN_SIMD_BITWIDTH"
 #endif
 #endif
-#include "gmx_x86_simd_macros.h"
+#include "gmx_simd_macros.h"
 
-#if GMX_X86_SIMD_WIDTH_HERE >= 2*NBNXN_CPU_CLUSTER_I_SIZE
-#define STRIDE_S  (GMX_X86_SIMD_WIDTH_HERE/2)
+#if GMX_SIMD_WIDTH_HERE >= 2*NBNXN_CPU_CLUSTER_I_SIZE
+#define STRIDE_S  (GMX_SIMD_WIDTH_HERE/2)
 #else
 #define STRIDE_S  NBNXN_CPU_CLUSTER_I_SIZE
 #endif
@@ -181,7 +181,7 @@ make_cluster_list_simd_2xnn(const nbnxn_grid_t *gridj,
 
             InRange            = gmx_movemask_pr(wco_any_SSE);
 
-            *ndistc += 2*GMX_X86_SIMD_WIDTH_HERE;
+            *ndistc += 2*GMX_SIMD_WIDTH_HERE;
         }
         if (!InRange)
         {
@@ -235,7 +235,7 @@ make_cluster_list_simd_2xnn(const nbnxn_grid_t *gridj,
 
             InRange            = gmx_movemask_pr(wco_any_SSE);
 
-            *ndistc += 2*GMX_X86_SIMD_WIDTH_HERE;
+            *ndistc += 2*GMX_SIMD_WIDTH_HERE;
         }
         if (!InRange)
         {
index d754d07bc09325f9fc310a97281d3c775fb435c9..60742fb7196a16f90df3601808f4827995dbd96c 100644 (file)
 #error "unsupported GMX_NBNXN_SIMD_BITWIDTH"
 #endif
 #endif
-#include "gmx_x86_simd_macros.h"
+#include "gmx_simd_macros.h"
 
-#if GMX_X86_SIMD_WIDTH_HERE >= NBNXN_CPU_CLUSTER_I_SIZE
-#define STRIDE_S  (GMX_X86_SIMD_WIDTH_HERE)
+#if GMX_SIMD_WIDTH_HERE >= NBNXN_CPU_CLUSTER_I_SIZE
+#define STRIDE_S  (GMX_SIMD_WIDTH_HERE)
 #else
 #define STRIDE_S  NBNXN_CPU_CLUSTER_I_SIZE
 #endif
@@ -187,7 +187,7 @@ make_cluster_list_simd_4xn(const nbnxn_grid_t *gridj,
             
             InRange            = gmx_movemask_pr(wco_any_SSE);
 
-            *ndistc += 4*GMX_X86_SIMD_WIDTH_HERE;
+            *ndistc += 4*GMX_SIMD_WIDTH_HERE;
         }
         if (!InRange)
         {
@@ -253,7 +253,7 @@ make_cluster_list_simd_4xn(const nbnxn_grid_t *gridj,
             
             InRange            = gmx_movemask_pr(wco_any_SSE);
 
-            *ndistc += 4*GMX_X86_SIMD_WIDTH_HERE;
+            *ndistc += 4*GMX_SIMD_WIDTH_HERE;
         }
         if (!InRange)
         {