Fix nbnxn SIMD 2xNN PME tabulated kernel
authorBerk Hess <hess@kth.se>
Thu, 14 Sep 2017 14:34:35 +0000 (16:34 +0200)
committerBerk Hess <hess@kth.se>
Thu, 14 Sep 2017 14:34:35 +0000 (16:34 +0200)
This kernel produced completety incorrect forces and slightly
incorrect energies.
This kernel was not (yet) selected by default in a 2016 release.

Fixes #2247

Change-Id: I297ad257932eaabe6cf84b17f34ad555921f48b0

src/gromacs/mdlib/nbnxn_kernels/simd_2xnn/nbnxn_kernel_simd_2xnn_inner.h

index ac30b0b0ef516c909f1ca107cdaada3820154334..a680aff199e26c6d0b045f20ed5d3444dd6868be 100644 (file)
 #else
     gatherLoadUBySimdIntTranspose<1>(tab_coul_F, ti_S0, &ctab0_S0, &ctab1_S0);
     gatherLoadUBySimdIntTranspose<1>(tab_coul_F, ti_S2, &ctab0_S2, &ctab1_S2);
+    ctab1_S0 = ctab1_S0 - ctab0_S0;
+    ctab1_S2 = ctab1_S2 - ctab0_S2;
 #endif
 #else
 #ifdef TAB_FDV0
 #endif
     fsub_S0     = fma(frac_S0, ctab1_S0, ctab0_S0);
     fsub_S2     = fma(frac_S2, ctab1_S2, ctab0_S2);
-    frcoul_S0   = qq_S0, fnma(fsub_S0, r_S0, rinv_ex_S0);
-    frcoul_S2   = qq_S2, fnma(fsub_S2, r_S2, rinv_ex_S2);
+    frcoul_S0   = qq_S0 * fnma(fsub_S0, r_S0, rinv_ex_S0);
+    frcoul_S2   = qq_S2 * fnma(fsub_S2, r_S2, rinv_ex_S2);
 
 #ifdef CALC_ENERGIES
-    vc_sub_S0   = ctabv_S0 + (mhalfsp_S * frac_S0 * (ctab0_S0 + fsub_S0));
-    vc_sub_S2   = ctabv_S2 + (mhalfsp_S * frac_S2 * (ctab0_S2 + fsub_S2));
+    vc_sub_S0   = fma((mhalfsp_S * frac_S0), (ctab0_S0 + fsub_S0), ctabv_S0);
+    vc_sub_S2   = fma((mhalfsp_S * frac_S2), (ctab0_S2 + fsub_S2), ctabv_S2);
 #endif
 #endif /* CALC_COUL_TAB */
 
 
     vcoul_S0    = qq_S0 * (rinv_ex_S0 - vc_sub_S0);
     vcoul_S2    = qq_S2 * (rinv_ex_S2 - vc_sub_S2);
+
 #endif
 
 #ifdef CALC_ENERGIES