2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
45 #include "gmx_omp_nthreads.h"
46 #include "nbnxn_kernel_ref.h"
47 #include "../nbnxn_consts.h"
48 #include "nbnxn_kernel_common.h"
50 /*! \brief Typedefs for declaring lookup tables of kernel functions.
52 typedef void (*p_nbk_func_ener)(const nbnxn_pairlist_t *nbl,
53 const nbnxn_atomdata_t *nbat,
54 const interaction_const_t *ic,
61 typedef void (*p_nbk_func_noener)(const nbnxn_pairlist_t *nbl,
62 const nbnxn_atomdata_t *nbat,
63 const interaction_const_t *ic,
68 /* Analytical reaction-field kernels */
71 /* Include the force+energy kernels */
73 #include "nbnxn_kernel_ref_outer.h"
76 /* Include the force+energygroups kernels */
79 #include "nbnxn_kernel_ref_outer.h"
83 /* Include the force only kernels */
84 #include "nbnxn_kernel_ref_outer.h"
89 /* Tabulated exclusion interaction electrostatics kernels */
92 /* Include the force+energy kernels */
94 #include "nbnxn_kernel_ref_outer.h"
97 /* Include the force+energygroups kernels */
100 #include "nbnxn_kernel_ref_outer.h"
104 /* Include the force only kernels */
105 #include "nbnxn_kernel_ref_outer.h"
107 /* Twin-range cut-off kernels */
108 #define VDW_CUTOFF_CHECK
110 /* Include the force+energy kernels */
111 #define CALC_ENERGIES
112 #include "nbnxn_kernel_ref_outer.h"
115 /* Include the force+energygroups kernels */
116 #define CALC_ENERGIES
117 #define ENERGY_GROUPS
118 #include "nbnxn_kernel_ref_outer.h"
122 /* Include the force only kernels */
123 #include "nbnxn_kernel_ref_outer.h"
125 #undef VDW_CUTOFF_CHECK
131 coultRF, coultTAB, coultTAB_TWIN, coultNR
134 p_nbk_func_ener p_nbk_c_ener[coultNR] =
136 nbnxn_kernel_ref_rf_ener,
137 nbnxn_kernel_ref_tab_ener,
138 nbnxn_kernel_ref_tab_twin_ener
141 p_nbk_func_ener p_nbk_c_energrp[coultNR] =
143 nbnxn_kernel_ref_rf_energrp,
144 nbnxn_kernel_ref_tab_energrp,
145 nbnxn_kernel_ref_tab_twin_energrp
148 p_nbk_func_noener p_nbk_c_noener[coultNR] =
150 nbnxn_kernel_ref_rf_noener,
151 nbnxn_kernel_ref_tab_noener,
152 nbnxn_kernel_ref_tab_twin_noener
156 nbnxn_kernel_ref(const nbnxn_pairlist_set_t *nbl_list,
157 const nbnxn_atomdata_t *nbat,
158 const interaction_const_t *ic,
167 nbnxn_pairlist_t **nbl;
171 nnbl = nbl_list->nnbl;
174 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
180 if (ic->rcoulomb == ic->rvdw)
186 coult = coultTAB_TWIN;
190 #pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
191 for (nb = 0; nb < nnbl; nb++)
193 nbnxn_atomdata_output_t *out;
196 out = &nbat->out[nb];
198 if (clearF == enbvClearFYes)
200 clear_f(nbat, nb, out->f);
203 if ((force_flags & GMX_FORCE_VIRIAL) && nnbl == 1)
209 fshift_p = out->fshift;
211 if (clearF == enbvClearFYes)
213 clear_fshift(fshift_p);
217 if (!(force_flags & GMX_FORCE_ENERGY))
219 /* Don't calculate energies */
220 p_nbk_c_noener[coult](nbl[nb], nbat,
226 else if (out->nV == 1)
228 /* No energy groups */
232 p_nbk_c_ener[coult](nbl[nb], nbat,
242 /* Calculate energy group contributions */
245 for (i = 0; i < out->nV; i++)
249 for (i = 0; i < out->nV; i++)
254 p_nbk_c_energrp[coult](nbl[nb], nbat,
264 if (force_flags & GMX_FORCE_ENERGY)
266 reduce_energies_over_lists(nbat, nnbl, Vvdw, Vc);