From 2755334da0b19be6d478eae8d6cceb81b0c0916c Mon Sep 17 00:00:00 2001 From: Berk Hess Date: Fri, 16 Jan 2015 10:05:14 +0100 Subject: [PATCH] Re-enabled SIMD for angles and dihedrals The SIMD acceleration for angles and dihedrals was accidentally disabled recently during code reorganization, since the SIMD header file was not includes and thus the SIMD macro was not set. Fixes #1673 Change-Id: I04ddde71c9dd1a846ecf088a5c973e9a23846a52 --- src/gromacs/listed-forces/bonded.cpp | 6 ++-- src/gromacs/listed-forces/bonded.h | 32 ++++++++++++++++++++- src/gromacs/listed-forces/listed-forces.cpp | 1 + 3 files changed, 35 insertions(+), 4 deletions(-) diff --git a/src/gromacs/listed-forces/bonded.cpp b/src/gromacs/listed-forces/bonded.cpp index 7e4619e82b..63604f2032 100644 --- a/src/gromacs/listed-forces/bonded.cpp +++ b/src/gromacs/listed-forces/bonded.cpp @@ -1037,10 +1037,10 @@ real angles(int nbonds, #ifdef GMX_SIMD_HAVE_REAL -/* As angles, but using SIMD to calculate many dihedrals at once. +/* As angles, but using SIMD to calculate many angles at once. * This routines does not calculate energies and shift forces. */ -static gmx_inline void +void angles_noener_simd(int nbonds, const t_iatom forceatoms[], const t_iparams forceparams[], const rvec x[], rvec f[], @@ -1966,7 +1966,7 @@ pdihs_noener(int nbonds, #ifdef GMX_SIMD_HAVE_REAL /* As pdihs_noner above, but using SIMD to calculate many dihedrals at once */ -static void +void pdihs_noener_simd(int nbonds, const t_iatom forceatoms[], const t_iparams forceparams[], const rvec x[], rvec f[], diff --git a/src/gromacs/listed-forces/bonded.h b/src/gromacs/listed-forces/bonded.h index b8abfdd0c6..93e8ed343f 100644 --- a/src/gromacs/listed-forces/bonded.h +++ b/src/gromacs/listed-forces/bonded.h @@ -56,6 +56,7 @@ #include "gromacs/legacyheaders/types/mdatom.h" #include "gromacs/legacyheaders/types/nrnb.h" #include "gromacs/math/vec.h" +#include "gromacs/simd/simd.h" #include "gromacs/topology/idef.h" #include "gromacs/utility/basedefinitions.h" @@ -115,7 +116,7 @@ t_ifunc restrdihs, cbtdihs; t_ifunc tab_bonds, tab_angles, tab_dihs; t_ifunc polarize, anharm_polarize, water_pol, thole_pol, angres, angresz, dihres, unimplemented; -/* As pdihs above, but without calculating energies and shift forces */ +/* As pdihs(), but without calculating energies and shift forces */ void pdihs_noener(int nbonds, const t_iatom forceatoms[], const t_iparams forceparams[], @@ -125,6 +126,35 @@ void real lambda, const t_mdatoms gmx_unused *md, t_fcdata gmx_unused *fcd, int gmx_unused *global_atom_index); + +#ifdef GMX_SIMD_HAVE_REAL + +/* As angles(), but using SIMD to calculate many angles at once. + * This routines does not calculate energies and shift forces. + */ +void + angles_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); + +/* As pdihs_noener(), but using SIMD to calculate many dihedrals at once. */ +void + pdihs_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 #ifdef __cplusplus diff --git a/src/gromacs/listed-forces/listed-forces.cpp b/src/gromacs/listed-forces/listed-forces.cpp index 1864494bf4..bc5dc6bbd1 100644 --- a/src/gromacs/listed-forces/listed-forces.cpp +++ b/src/gromacs/listed-forces/listed-forces.cpp @@ -63,6 +63,7 @@ #include "gromacs/mdlib/forcerec-threading.h" #include "gromacs/pbcutil/ishift.h" #include "gromacs/pbcutil/pbc.h" +#include "gromacs/simd/simd.h" #include "gromacs/timing/wallcycle.h" #include "gromacs/utility/smalloc.h" -- 2.22.0