Improved check of if only forces are computed in do_pairs()
authorMagnus Lundborg <lundborg.magnus@gmail.com>
Thu, 1 Mar 2018 07:39:18 +0000 (08:39 +0100)
committerMark Abraham <mark.j.abraham@gmail.com>
Thu, 29 Mar 2018 12:10:37 +0000 (14:10 +0200)
Use the check introduced in I67bdb0cba9c6ec20b to see if the
SIMD version can be used.

Change-Id: I3dfd401aa1dc421757bf243809ac85e17a9062b6

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

index 32a8bc97ecb4f044e94b0b3bd1788dd24b1a1a1e..d0fc37e7bb50a083067603331ef69527c833c78a 100644 (file)
@@ -386,7 +386,7 @@ calc_one_bond(int thread,
            extended to support calling from multiple threads. */
         do_pairs(ftype, nbn, iatoms+nb0, idef->iparams, x, f, fshift,
                  pbc, g, lambda, dvdl, md, fr,
-                 bCalcEnerVir, grpp, global_atom_index);
+                 computeForcesOnly, grpp, global_atom_index);
         v = 0;
     }
 
index fdc2279c20983c9a7eeae6f7195316ae3c190ac5..01bf6130a45b3e001e9a86b45825ba42db7a4e2e 100644 (file)
@@ -646,12 +646,13 @@ do_pairs(int ftype, int nbonds,
          const real *lambda, real *dvdl,
          const t_mdatoms *md,
          const t_forcerec *fr,
-         gmx_bool bCalcEnergyAndVirial, gmx_grppairener_t *grppener,
+         const gmx_bool computeForcesOnly,
+         gmx_grppairener_t *grppener,
          int *global_atom_index)
 {
     if (ftype == F_LJ14 &&
         fr->ic->vdwtype != evdwUSER && !EEL_USER(fr->ic->eeltype) &&
-        !bCalcEnergyAndVirial && fr->efep == efepNO)
+        computeForcesOnly)
     {
         /* We use a fast code-path for plain LJ 1-4 without FEP.
          *
index c6dbcc256b035cdbc1cdde507b05b50f77745ddf..d3c4bc9b0897c0708f6a866712cdfd66a65791d2 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2014,2015,2016,2017, by the GROMACS development team, led by
+ * Copyright (c) 2014,2015,2016,2017,2018, 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.
@@ -64,7 +64,7 @@ do_pairs(int ftype, int nbonds, const t_iatom iatoms[], const t_iparams iparams[
          const rvec x[], rvec4 f[], rvec fshift[],
          const struct t_pbc *pbc, const struct t_graph *g,
          const real *lambda, real *dvdl, const t_mdatoms *md, const t_forcerec *fr,
-         gmx_bool bCalcEnergyAndVirial, gmx_grppairener_t *grppener,
+         const gmx_bool computeForcesOnly, gmx_grppairener_t *grppener,
          int *global_atom_index);
 
 #endif