2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, 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.
36 * Note: this file was generated by the Verlet kernel generator for
46 #include "gromacs/mdlib/nbnxn_simd.h"
48 #ifdef GMX_NBNXN_SIMD_4XN
50 #include "gromacs/simd/vector_operations.h"
52 #if !(GMX_SIMD_REAL_WIDTH == 2 || GMX_SIMD_REAL_WIDTH == 4 || GMX_SIMD_REAL_WIDTH == 8)
53 #error "unsupported SIMD width"
56 #define GMX_SIMD_J_UNROLL_SIZE 1
57 #include "nbnxn_kernel_simd_4xn.h"
58 #include "../nbnxn_kernel_common.h"
59 #include "gmx_omp_nthreads.h"
60 #include "types/force_flags.h"
61 #include "gmx_fatal.h"
63 /*! \brief Kinds of electrostatic treatments in SIMD Verlet kernels
66 coulktRF, coulktTAB, coulktTAB_TWIN, coulktEWALD, coulktEWALD_TWIN, coulktNR
69 /*! \brief Kinds of Van der Waals treatments in SIMD Verlet kernels
72 vdwktLJCUT_COMBGEOM, vdwktLJCUT_COMBLB, vdwktLJCUT_COMBNONE, vdwktLJFORCESWITCH, vdwktLJPOTSWITCH, vdwktLJEWALDCOMBGEOM, vdwktNR
75 /* Declare and define the kernel function pointer lookup tables.
76 * The minor index of the array goes over both the LJ combination rules,
77 * which is only supported by plain cut-off, and the LJ switch/PME functions.
79 static p_nbk_func_noener p_nbk_noener[coulktNR][vdwktNR] =
82 nbnxn_kernel_ElecRF_VdwLJCombGeom_F_4xn,
83 nbnxn_kernel_ElecRF_VdwLJCombLB_F_4xn,
84 nbnxn_kernel_ElecRF_VdwLJ_F_4xn,
85 nbnxn_kernel_ElecRF_VdwLJFSw_F_4xn,
86 nbnxn_kernel_ElecRF_VdwLJPSw_F_4xn,
87 nbnxn_kernel_ElecRF_VdwLJEwCombGeom_F_4xn,
90 nbnxn_kernel_ElecQSTab_VdwLJCombGeom_F_4xn,
91 nbnxn_kernel_ElecQSTab_VdwLJCombLB_F_4xn,
92 nbnxn_kernel_ElecQSTab_VdwLJ_F_4xn,
93 nbnxn_kernel_ElecQSTab_VdwLJFSw_F_4xn,
94 nbnxn_kernel_ElecQSTab_VdwLJPSw_F_4xn,
95 nbnxn_kernel_ElecQSTab_VdwLJEwCombGeom_F_4xn,
98 nbnxn_kernel_ElecQSTabTwinCut_VdwLJCombGeom_F_4xn,
99 nbnxn_kernel_ElecQSTabTwinCut_VdwLJCombLB_F_4xn,
100 nbnxn_kernel_ElecQSTabTwinCut_VdwLJ_F_4xn,
101 nbnxn_kernel_ElecQSTabTwinCut_VdwLJFSw_F_4xn,
102 nbnxn_kernel_ElecQSTabTwinCut_VdwLJPSw_F_4xn,
103 nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_F_4xn,
106 nbnxn_kernel_ElecEw_VdwLJCombGeom_F_4xn,
107 nbnxn_kernel_ElecEw_VdwLJCombLB_F_4xn,
108 nbnxn_kernel_ElecEw_VdwLJ_F_4xn,
109 nbnxn_kernel_ElecEw_VdwLJFSw_F_4xn,
110 nbnxn_kernel_ElecEw_VdwLJPSw_F_4xn,
111 nbnxn_kernel_ElecEw_VdwLJEwCombGeom_F_4xn,
114 nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_F_4xn,
115 nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_F_4xn,
116 nbnxn_kernel_ElecEwTwinCut_VdwLJ_F_4xn,
117 nbnxn_kernel_ElecEwTwinCut_VdwLJFSw_F_4xn,
118 nbnxn_kernel_ElecEwTwinCut_VdwLJPSw_F_4xn,
119 nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_F_4xn,
123 static p_nbk_func_ener p_nbk_ener[coulktNR][vdwktNR] =
126 nbnxn_kernel_ElecRF_VdwLJCombGeom_VF_4xn,
127 nbnxn_kernel_ElecRF_VdwLJCombLB_VF_4xn,
128 nbnxn_kernel_ElecRF_VdwLJ_VF_4xn,
129 nbnxn_kernel_ElecRF_VdwLJFSw_VF_4xn,
130 nbnxn_kernel_ElecRF_VdwLJPSw_VF_4xn,
131 nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VF_4xn,
134 nbnxn_kernel_ElecQSTab_VdwLJCombGeom_VF_4xn,
135 nbnxn_kernel_ElecQSTab_VdwLJCombLB_VF_4xn,
136 nbnxn_kernel_ElecQSTab_VdwLJ_VF_4xn,
137 nbnxn_kernel_ElecQSTab_VdwLJFSw_VF_4xn,
138 nbnxn_kernel_ElecQSTab_VdwLJPSw_VF_4xn,
139 nbnxn_kernel_ElecQSTab_VdwLJEwCombGeom_VF_4xn,
142 nbnxn_kernel_ElecQSTabTwinCut_VdwLJCombGeom_VF_4xn,
143 nbnxn_kernel_ElecQSTabTwinCut_VdwLJCombLB_VF_4xn,
144 nbnxn_kernel_ElecQSTabTwinCut_VdwLJ_VF_4xn,
145 nbnxn_kernel_ElecQSTabTwinCut_VdwLJFSw_VF_4xn,
146 nbnxn_kernel_ElecQSTabTwinCut_VdwLJPSw_VF_4xn,
147 nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_VF_4xn,
150 nbnxn_kernel_ElecEw_VdwLJCombGeom_VF_4xn,
151 nbnxn_kernel_ElecEw_VdwLJCombLB_VF_4xn,
152 nbnxn_kernel_ElecEw_VdwLJ_VF_4xn,
153 nbnxn_kernel_ElecEw_VdwLJFSw_VF_4xn,
154 nbnxn_kernel_ElecEw_VdwLJPSw_VF_4xn,
155 nbnxn_kernel_ElecEw_VdwLJEwCombGeom_VF_4xn,
158 nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_VF_4xn,
159 nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_VF_4xn,
160 nbnxn_kernel_ElecEwTwinCut_VdwLJ_VF_4xn,
161 nbnxn_kernel_ElecEwTwinCut_VdwLJFSw_VF_4xn,
162 nbnxn_kernel_ElecEwTwinCut_VdwLJPSw_VF_4xn,
163 nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_VF_4xn,
167 static p_nbk_func_ener p_nbk_energrp[coulktNR][vdwktNR] =
170 nbnxn_kernel_ElecRF_VdwLJCombGeom_VgrpF_4xn,
171 nbnxn_kernel_ElecRF_VdwLJCombLB_VgrpF_4xn,
172 nbnxn_kernel_ElecRF_VdwLJ_VgrpF_4xn,
173 nbnxn_kernel_ElecRF_VdwLJFSw_VgrpF_4xn,
174 nbnxn_kernel_ElecRF_VdwLJPSw_VgrpF_4xn,
175 nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VgrpF_4xn,
178 nbnxn_kernel_ElecQSTab_VdwLJCombGeom_VgrpF_4xn,
179 nbnxn_kernel_ElecQSTab_VdwLJCombLB_VgrpF_4xn,
180 nbnxn_kernel_ElecQSTab_VdwLJ_VgrpF_4xn,
181 nbnxn_kernel_ElecQSTab_VdwLJFSw_VgrpF_4xn,
182 nbnxn_kernel_ElecQSTab_VdwLJPSw_VgrpF_4xn,
183 nbnxn_kernel_ElecQSTab_VdwLJEwCombGeom_VgrpF_4xn,
186 nbnxn_kernel_ElecQSTabTwinCut_VdwLJCombGeom_VgrpF_4xn,
187 nbnxn_kernel_ElecQSTabTwinCut_VdwLJCombLB_VgrpF_4xn,
188 nbnxn_kernel_ElecQSTabTwinCut_VdwLJ_VgrpF_4xn,
189 nbnxn_kernel_ElecQSTabTwinCut_VdwLJFSw_VgrpF_4xn,
190 nbnxn_kernel_ElecQSTabTwinCut_VdwLJPSw_VgrpF_4xn,
191 nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_VgrpF_4xn,
194 nbnxn_kernel_ElecEw_VdwLJCombGeom_VgrpF_4xn,
195 nbnxn_kernel_ElecEw_VdwLJCombLB_VgrpF_4xn,
196 nbnxn_kernel_ElecEw_VdwLJ_VgrpF_4xn,
197 nbnxn_kernel_ElecEw_VdwLJFSw_VgrpF_4xn,
198 nbnxn_kernel_ElecEw_VdwLJPSw_VgrpF_4xn,
199 nbnxn_kernel_ElecEw_VdwLJEwCombGeom_VgrpF_4xn,
202 nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_VgrpF_4xn,
203 nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_VgrpF_4xn,
204 nbnxn_kernel_ElecEwTwinCut_VdwLJ_VgrpF_4xn,
205 nbnxn_kernel_ElecEwTwinCut_VdwLJFSw_VgrpF_4xn,
206 nbnxn_kernel_ElecEwTwinCut_VdwLJPSw_VgrpF_4xn,
207 nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_VgrpF_4xn,
213 reduce_group_energies(int ng, int ng_2log,
214 const real *VSvdw, const real *VSc,
215 real *Vvdw, real *Vc)
217 const int unrollj = GMX_SIMD_REAL_WIDTH/GMX_SIMD_J_UNROLL_SIZE;
218 const int unrollj_half = unrollj/2;
219 int ng_p2, i, j, j0, j1, c, s;
221 ng_p2 = (1<<ng_2log);
223 /* The size of the x86 SIMD energy group buffer array is:
224 * ng*ng*ng_p2*unrollj_half*simd_width
226 for (i = 0; i < ng; i++)
228 for (j = 0; j < ng; j++)
234 for (j1 = 0; j1 < ng; j1++)
236 for (j0 = 0; j0 < ng; j0++)
238 c = ((i*ng + j1)*ng_p2 + j0)*unrollj_half*unrollj;
239 for (s = 0; s < unrollj_half; s++)
241 Vvdw[i*ng+j0] += VSvdw[c+0];
242 Vvdw[i*ng+j1] += VSvdw[c+1];
243 Vc [i*ng+j0] += VSc [c+0];
244 Vc [i*ng+j1] += VSc [c+1];
252 #else /* GMX_NBNXN_SIMD_4XN */
254 #include "gmx_fatal.h"
256 #endif /* GMX_NBNXN_SIMD_4XN */
259 nbnxn_kernel_simd_4xn(nbnxn_pairlist_set_t gmx_unused *nbl_list,
260 const nbnxn_atomdata_t gmx_unused *nbat,
261 const interaction_const_t gmx_unused *ic,
262 int gmx_unused ewald_excl,
263 rvec gmx_unused *shift_vec,
264 int gmx_unused force_flags,
265 int gmx_unused clearF,
266 real gmx_unused *fshift,
268 real gmx_unused *Vvdw)
269 #ifdef GMX_NBNXN_SIMD_4XN
272 nbnxn_pairlist_t **nbl;
273 int coulkt, vdwkt = 0;
275 int nthreads gmx_unused;
277 nnbl = nbl_list->nnbl;
280 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
286 if (ewald_excl == ewaldexclTable)
288 if (ic->rcoulomb == ic->rvdw)
294 coulkt = coulktTAB_TWIN;
299 if (ic->rcoulomb == ic->rvdw)
301 coulkt = coulktEWALD;
305 coulkt = coulktEWALD_TWIN;
310 if (ic->vdwtype == evdwCUT)
312 switch (ic->vdw_modifier)
315 case eintmodPOTSHIFT:
316 switch (nbat->comb_rule)
318 case ljcrGEOM: vdwkt = vdwktLJCUT_COMBGEOM; break;
319 case ljcrLB: vdwkt = vdwktLJCUT_COMBLB; break;
320 case ljcrNONE: vdwkt = vdwktLJCUT_COMBNONE; break;
321 default: gmx_incons("Unknown combination rule");
324 case eintmodFORCESWITCH:
325 vdwkt = vdwktLJFORCESWITCH;
327 case eintmodPOTSWITCH:
328 vdwkt = vdwktLJPOTSWITCH;
331 gmx_incons("Unsupported VdW interaction modifier");
334 else if (ic->vdwtype == evdwPME)
336 if (ic->ljpme_comb_rule == eljpmeLB)
338 gmx_incons("The nbnxn SIMD kernels don't suport LJ-PME with LB");
340 vdwkt = vdwktLJEWALDCOMBGEOM;
344 gmx_incons("Unsupported VdW interaction type");
347 nthreads = gmx_omp_nthreads_get(emntNonbonded);
348 #pragma omp parallel for schedule(static) num_threads(nthreads)
349 for (nb = 0; nb < nnbl; nb++)
351 nbnxn_atomdata_output_t *out;
354 out = &nbat->out[nb];
356 if (clearF == enbvClearFYes)
358 clear_f(nbat, nb, out->f);
361 if ((force_flags & GMX_FORCE_VIRIAL) && nnbl == 1)
367 fshift_p = out->fshift;
369 if (clearF == enbvClearFYes)
371 clear_fshift(fshift_p);
375 if (!(force_flags & GMX_FORCE_ENERGY))
377 /* Don't calculate energies */
378 p_nbk_noener[coulkt][vdwkt](nbl[nb], nbat,
384 else if (out->nV == 1)
386 /* No energy groups */
390 p_nbk_ener[coulkt][vdwkt](nbl[nb], nbat,
400 /* Calculate energy group contributions */
403 for (i = 0; i < out->nVS; i++)
407 for (i = 0; i < out->nVS; i++)
412 p_nbk_energrp[coulkt][vdwkt](nbl[nb], nbat,
420 reduce_group_energies(nbat->nenergrp, nbat->neg_2log,
421 out->VSvdw, out->VSc,
426 if (force_flags & GMX_FORCE_ENERGY)
428 reduce_energies_over_lists(nbat, nnbl, Vvdw, Vc);
433 gmx_incons("nbnxn_kernel_simd_4xn called when such kernels "
434 " are not enabled.");
437 #undef GMX_SIMD_J_UNROLL_SIZE