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
44 #include "gromacs/mdlib/nbnxn_simd.h"
46 #ifdef GMX_NBNXN_SIMD_2XNN
48 /* Include the full-width SIMD macros */
49 #include "gromacs/simd/vector_operations.h"
51 #if !(GMX_SIMD_REAL_WIDTH == 8 || GMX_SIMD_REAL_WIDTH == 16)
52 #error "unsupported SIMD width"
55 #define GMX_SIMD_J_UNROLL_SIZE 2
56 #include "nbnxn_kernel_simd_2xnn.h"
57 #include "../nbnxn_kernel_common.h"
58 #include "gmx_omp_nthreads.h"
59 #include "types/force_flags.h"
60 #include "gromacs/utility/fatalerror.h"
62 /*! \brief Kinds of electrostatic treatments in SIMD Verlet kernels
65 coulktRF, coulktTAB, coulktTAB_TWIN, coulktEWALD, coulktEWALD_TWIN, coulktNR
68 /*! \brief Kinds of Van der Waals treatments in SIMD Verlet kernels
71 vdwktLJCUT_COMBGEOM, vdwktLJCUT_COMBLB, vdwktLJCUT_COMBNONE, vdwktLJFORCESWITCH, vdwktLJPOTSWITCH, vdwktLJEWALDCOMBGEOM, vdwktNR
74 /* Declare and define the kernel function pointer lookup tables.
75 * The minor index of the array goes over both the LJ combination rules,
76 * which is only supported by plain cut-off, and the LJ switch/PME functions.
78 static p_nbk_func_noener p_nbk_noener[coulktNR][vdwktNR] =
81 nbnxn_kernel_ElecRF_VdwLJCombGeom_F_2xnn,
82 nbnxn_kernel_ElecRF_VdwLJCombLB_F_2xnn,
83 nbnxn_kernel_ElecRF_VdwLJ_F_2xnn,
84 nbnxn_kernel_ElecRF_VdwLJFSw_F_2xnn,
85 nbnxn_kernel_ElecRF_VdwLJPSw_F_2xnn,
86 nbnxn_kernel_ElecRF_VdwLJEwCombGeom_F_2xnn,
89 nbnxn_kernel_ElecQSTab_VdwLJCombGeom_F_2xnn,
90 nbnxn_kernel_ElecQSTab_VdwLJCombLB_F_2xnn,
91 nbnxn_kernel_ElecQSTab_VdwLJ_F_2xnn,
92 nbnxn_kernel_ElecQSTab_VdwLJFSw_F_2xnn,
93 nbnxn_kernel_ElecQSTab_VdwLJPSw_F_2xnn,
94 nbnxn_kernel_ElecQSTab_VdwLJEwCombGeom_F_2xnn,
97 nbnxn_kernel_ElecQSTabTwinCut_VdwLJCombGeom_F_2xnn,
98 nbnxn_kernel_ElecQSTabTwinCut_VdwLJCombLB_F_2xnn,
99 nbnxn_kernel_ElecQSTabTwinCut_VdwLJ_F_2xnn,
100 nbnxn_kernel_ElecQSTabTwinCut_VdwLJFSw_F_2xnn,
101 nbnxn_kernel_ElecQSTabTwinCut_VdwLJPSw_F_2xnn,
102 nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_F_2xnn,
105 nbnxn_kernel_ElecEw_VdwLJCombGeom_F_2xnn,
106 nbnxn_kernel_ElecEw_VdwLJCombLB_F_2xnn,
107 nbnxn_kernel_ElecEw_VdwLJ_F_2xnn,
108 nbnxn_kernel_ElecEw_VdwLJFSw_F_2xnn,
109 nbnxn_kernel_ElecEw_VdwLJPSw_F_2xnn,
110 nbnxn_kernel_ElecEw_VdwLJEwCombGeom_F_2xnn,
113 nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_F_2xnn,
114 nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_F_2xnn,
115 nbnxn_kernel_ElecEwTwinCut_VdwLJ_F_2xnn,
116 nbnxn_kernel_ElecEwTwinCut_VdwLJFSw_F_2xnn,
117 nbnxn_kernel_ElecEwTwinCut_VdwLJPSw_F_2xnn,
118 nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_F_2xnn,
122 static p_nbk_func_ener p_nbk_ener[coulktNR][vdwktNR] =
125 nbnxn_kernel_ElecRF_VdwLJCombGeom_VF_2xnn,
126 nbnxn_kernel_ElecRF_VdwLJCombLB_VF_2xnn,
127 nbnxn_kernel_ElecRF_VdwLJ_VF_2xnn,
128 nbnxn_kernel_ElecRF_VdwLJFSw_VF_2xnn,
129 nbnxn_kernel_ElecRF_VdwLJPSw_VF_2xnn,
130 nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VF_2xnn,
133 nbnxn_kernel_ElecQSTab_VdwLJCombGeom_VF_2xnn,
134 nbnxn_kernel_ElecQSTab_VdwLJCombLB_VF_2xnn,
135 nbnxn_kernel_ElecQSTab_VdwLJ_VF_2xnn,
136 nbnxn_kernel_ElecQSTab_VdwLJFSw_VF_2xnn,
137 nbnxn_kernel_ElecQSTab_VdwLJPSw_VF_2xnn,
138 nbnxn_kernel_ElecQSTab_VdwLJEwCombGeom_VF_2xnn,
141 nbnxn_kernel_ElecQSTabTwinCut_VdwLJCombGeom_VF_2xnn,
142 nbnxn_kernel_ElecQSTabTwinCut_VdwLJCombLB_VF_2xnn,
143 nbnxn_kernel_ElecQSTabTwinCut_VdwLJ_VF_2xnn,
144 nbnxn_kernel_ElecQSTabTwinCut_VdwLJFSw_VF_2xnn,
145 nbnxn_kernel_ElecQSTabTwinCut_VdwLJPSw_VF_2xnn,
146 nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_VF_2xnn,
149 nbnxn_kernel_ElecEw_VdwLJCombGeom_VF_2xnn,
150 nbnxn_kernel_ElecEw_VdwLJCombLB_VF_2xnn,
151 nbnxn_kernel_ElecEw_VdwLJ_VF_2xnn,
152 nbnxn_kernel_ElecEw_VdwLJFSw_VF_2xnn,
153 nbnxn_kernel_ElecEw_VdwLJPSw_VF_2xnn,
154 nbnxn_kernel_ElecEw_VdwLJEwCombGeom_VF_2xnn,
157 nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_VF_2xnn,
158 nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_VF_2xnn,
159 nbnxn_kernel_ElecEwTwinCut_VdwLJ_VF_2xnn,
160 nbnxn_kernel_ElecEwTwinCut_VdwLJFSw_VF_2xnn,
161 nbnxn_kernel_ElecEwTwinCut_VdwLJPSw_VF_2xnn,
162 nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_VF_2xnn,
166 static p_nbk_func_ener p_nbk_energrp[coulktNR][vdwktNR] =
169 nbnxn_kernel_ElecRF_VdwLJCombGeom_VgrpF_2xnn,
170 nbnxn_kernel_ElecRF_VdwLJCombLB_VgrpF_2xnn,
171 nbnxn_kernel_ElecRF_VdwLJ_VgrpF_2xnn,
172 nbnxn_kernel_ElecRF_VdwLJFSw_VgrpF_2xnn,
173 nbnxn_kernel_ElecRF_VdwLJPSw_VgrpF_2xnn,
174 nbnxn_kernel_ElecRF_VdwLJEwCombGeom_VgrpF_2xnn,
177 nbnxn_kernel_ElecQSTab_VdwLJCombGeom_VgrpF_2xnn,
178 nbnxn_kernel_ElecQSTab_VdwLJCombLB_VgrpF_2xnn,
179 nbnxn_kernel_ElecQSTab_VdwLJ_VgrpF_2xnn,
180 nbnxn_kernel_ElecQSTab_VdwLJFSw_VgrpF_2xnn,
181 nbnxn_kernel_ElecQSTab_VdwLJPSw_VgrpF_2xnn,
182 nbnxn_kernel_ElecQSTab_VdwLJEwCombGeom_VgrpF_2xnn,
185 nbnxn_kernel_ElecQSTabTwinCut_VdwLJCombGeom_VgrpF_2xnn,
186 nbnxn_kernel_ElecQSTabTwinCut_VdwLJCombLB_VgrpF_2xnn,
187 nbnxn_kernel_ElecQSTabTwinCut_VdwLJ_VgrpF_2xnn,
188 nbnxn_kernel_ElecQSTabTwinCut_VdwLJFSw_VgrpF_2xnn,
189 nbnxn_kernel_ElecQSTabTwinCut_VdwLJPSw_VgrpF_2xnn,
190 nbnxn_kernel_ElecQSTabTwinCut_VdwLJEwCombGeom_VgrpF_2xnn,
193 nbnxn_kernel_ElecEw_VdwLJCombGeom_VgrpF_2xnn,
194 nbnxn_kernel_ElecEw_VdwLJCombLB_VgrpF_2xnn,
195 nbnxn_kernel_ElecEw_VdwLJ_VgrpF_2xnn,
196 nbnxn_kernel_ElecEw_VdwLJFSw_VgrpF_2xnn,
197 nbnxn_kernel_ElecEw_VdwLJPSw_VgrpF_2xnn,
198 nbnxn_kernel_ElecEw_VdwLJEwCombGeom_VgrpF_2xnn,
201 nbnxn_kernel_ElecEwTwinCut_VdwLJCombGeom_VgrpF_2xnn,
202 nbnxn_kernel_ElecEwTwinCut_VdwLJCombLB_VgrpF_2xnn,
203 nbnxn_kernel_ElecEwTwinCut_VdwLJ_VgrpF_2xnn,
204 nbnxn_kernel_ElecEwTwinCut_VdwLJFSw_VgrpF_2xnn,
205 nbnxn_kernel_ElecEwTwinCut_VdwLJPSw_VgrpF_2xnn,
206 nbnxn_kernel_ElecEwTwinCut_VdwLJEwCombGeom_VgrpF_2xnn,
212 reduce_group_energies(int ng, int ng_2log,
213 const real *VSvdw, const real *VSc,
214 real *Vvdw, real *Vc)
216 const int unrollj = GMX_SIMD_REAL_WIDTH/GMX_SIMD_J_UNROLL_SIZE;
217 const int unrollj_half = unrollj/2;
218 int ng_p2, i, j, j0, j1, c, s;
220 ng_p2 = (1<<ng_2log);
222 /* The size of the x86 SIMD energy group buffer array is:
223 * ng*ng*ng_p2*unrollj_half*simd_width
225 for (i = 0; i < ng; i++)
227 for (j = 0; j < ng; j++)
233 for (j1 = 0; j1 < ng; j1++)
235 for (j0 = 0; j0 < ng; j0++)
237 c = ((i*ng + j1)*ng_p2 + j0)*unrollj_half*unrollj;
238 for (s = 0; s < unrollj_half; s++)
240 Vvdw[i*ng+j0] += VSvdw[c+0];
241 Vvdw[i*ng+j1] += VSvdw[c+1];
242 Vc [i*ng+j0] += VSc [c+0];
243 Vc [i*ng+j1] += VSc [c+1];
251 #else /* GMX_NBNXN_SIMD_2XNN */
253 #include "gromacs/utility/fatalerror.h"
255 #endif /* GMX_NBNXN_SIMD_2XNN */
258 nbnxn_kernel_simd_2xnn(nbnxn_pairlist_set_t gmx_unused *nbl_list,
259 const nbnxn_atomdata_t gmx_unused *nbat,
260 const interaction_const_t gmx_unused *ic,
261 int gmx_unused ewald_excl,
262 rvec gmx_unused *shift_vec,
263 int gmx_unused force_flags,
264 int gmx_unused clearF,
265 real gmx_unused *fshift,
267 real gmx_unused *Vvdw)
268 #ifdef GMX_NBNXN_SIMD_2XNN
271 nbnxn_pairlist_t **nbl;
272 int coulkt, vdwkt = 0;
275 nnbl = nbl_list->nnbl;
278 if (EEL_RF(ic->eeltype) || ic->eeltype == eelCUT)
284 if (ewald_excl == ewaldexclTable)
286 if (ic->rcoulomb == ic->rvdw)
292 coulkt = coulktTAB_TWIN;
297 if (ic->rcoulomb == ic->rvdw)
299 coulkt = coulktEWALD;
303 coulkt = coulktEWALD_TWIN;
308 if (ic->vdwtype == evdwCUT)
310 switch (ic->vdw_modifier)
313 case eintmodPOTSHIFT:
314 switch (nbat->comb_rule)
316 case ljcrGEOM: vdwkt = vdwktLJCUT_COMBGEOM; break;
317 case ljcrLB: vdwkt = vdwktLJCUT_COMBLB; break;
318 case ljcrNONE: vdwkt = vdwktLJCUT_COMBNONE; break;
319 default: gmx_incons("Unknown combination rule");
322 case eintmodFORCESWITCH:
323 vdwkt = vdwktLJFORCESWITCH;
325 case eintmodPOTSWITCH:
326 vdwkt = vdwktLJPOTSWITCH;
329 gmx_incons("Unsupported VdW interaction modifier");
332 else if (ic->vdwtype == evdwPME)
334 if (ic->ljpme_comb_rule == eljpmeLB)
336 gmx_incons("The nbnxn SIMD kernels don't suport LJ-PME with LB");
338 vdwkt = vdwktLJEWALDCOMBGEOM;
342 gmx_incons("Unsupported VdW interaction type");
345 #pragma omp parallel for schedule(static) num_threads(gmx_omp_nthreads_get(emntNonbonded))
346 for (nb = 0; nb < nnbl; nb++)
348 nbnxn_atomdata_output_t *out;
351 out = &nbat->out[nb];
353 if (clearF == enbvClearFYes)
355 clear_f(nbat, nb, out->f);
358 if ((force_flags & GMX_FORCE_VIRIAL) && nnbl == 1)
364 fshift_p = out->fshift;
366 if (clearF == enbvClearFYes)
368 clear_fshift(fshift_p);
372 if (!(force_flags & GMX_FORCE_ENERGY))
374 /* Don't calculate energies */
375 p_nbk_noener[coulkt][vdwkt](nbl[nb], nbat,
381 else if (out->nV == 1)
383 /* No energy groups */
387 p_nbk_ener[coulkt][vdwkt](nbl[nb], nbat,
397 /* Calculate energy group contributions */
400 for (i = 0; i < out->nVS; i++)
404 for (i = 0; i < out->nVS; i++)
409 p_nbk_energrp[coulkt][vdwkt](nbl[nb], nbat,
417 reduce_group_energies(nbat->nenergrp, nbat->neg_2log,
418 out->VSvdw, out->VSc,
423 if (force_flags & GMX_FORCE_ENERGY)
425 reduce_energies_over_lists(nbat, nnbl, Vvdw, Vc);
430 gmx_incons("nbnxn_kernel_simd_2xnn called when such kernels "
431 " are not enabled.");
434 #undef GMX_SIMD_J_UNROLL_SIZE