Reenable SIMD for R-B dihedrals
authorMark Abraham <mark.j.abraham@gmail.com>
Sun, 15 Feb 2015 18:24:50 +0000 (19:24 +0100)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Wed, 18 Feb 2015 14:26:30 +0000 (15:26 +0100)
Fixes #1686

Change-Id: If612a5afde4bcfd55319d2b6e38b521aa43a9973

src/gromacs/listed-forces/bonded.cpp
src/gromacs/listed-forces/bonded.h
src/gromacs/listed-forces/listed-forces.cpp

index 63604f20321f0f7af10fdd56611f31ee5de48be9..d3a8c24a9eac44175a049246701b7e17ffc9d34e 100644 (file)
@@ -2093,7 +2093,7 @@ pdihs_noener_simd(int nbonds,
  * the RB potential instead of a harmonic potential.
  * This function can replace rbdihs() when no energy and virial are needed.
  */
-static void
+void
 rbdihs_noener_simd(int nbonds,
                    const t_iatom forceatoms[], const t_iparams forceparams[],
                    const rvec x[], rvec f[],
index 93e8ed343f17160125cf6ca986baf427b627bf6c..4c9e185b2ec1ee3ae3cb36297301a40ffaafab8a 100644 (file)
@@ -153,6 +153,17 @@ void
                       const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
                       int gmx_unused *global_atom_index);
 
+/* As rbdihs(), when not needing energy or shift force, using SIMD to calculate many dihedrals at once. */
+void
+    rbdihs_noener_simd(int nbonds,
+                       const t_iatom forceatoms[], const t_iparams forceparams[],
+                       const rvec x[], rvec f[],
+                       const struct t_pbc *pbc,
+                       const struct t_graph gmx_unused *g,
+                       real gmx_unused lambda,
+                       const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd,
+                       int gmx_unused *global_atom_index);
+
 #endif
 
 //! \endcond
index bc5dc6bbd1e1ab8fe4a238e9d8dd2acc856e7bc1..417d224e7375b3bf2615f7330482ae7e69ff9e92 100644 (file)
@@ -319,6 +319,19 @@ calc_one_bond(int thread,
                 global_atom_index);
             v = 0;
         }
+#ifdef GMX_SIMD_HAVE_REAL
+        else if (ftype == F_RBDIHS &&
+                 !bCalcEnerVir && fr->efep == efepNO)
+        {
+            /* No energies, shift forces, dvdl */
+            rbdihs_noener_simd(nbn, idef->il[ftype].iatoms+nb0,
+                               idef->iparams,
+                               (const rvec*)x, f,
+                               pbc, g, lambda[efptFTYPE], md, fcd,
+                               global_atom_index);
+            v = 0;
+        }
+#endif
         else
         {
             v = interaction_function[ftype].ifunc(nbn, iatoms+nb0,