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 GROMACS avx_256_single kernel generator.
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
47 #include "gromacs/simd/math_x86_avx_256_single.h"
48 #include "kernelutil_x86_avx_256_single.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_avx_256_single
52 * Electrostatics interaction: None
53 * VdW interaction: LJEwald
54 * Geometry: Particle-Particle
55 * Calculate force/pot: PotentialAndForce
58 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_avx_256_single
59 (t_nblist
* gmx_restrict nlist
,
60 rvec
* gmx_restrict xx
,
61 rvec
* gmx_restrict ff
,
62 t_forcerec
* gmx_restrict fr
,
63 t_mdatoms
* gmx_restrict mdatoms
,
64 nb_kernel_data_t gmx_unused
* gmx_restrict kernel_data
,
65 t_nrnb
* gmx_restrict nrnb
)
67 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
68 * just 0 for non-waters.
69 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
70 * jnr indices corresponding to data put in the four positions in the SIMD register.
72 int i_shift_offset
,i_coord_offset
,outeriter
,inneriter
;
73 int j_index_start
,j_index_end
,jidx
,nri
,inr
,ggid
,iidx
;
74 int jnrA
,jnrB
,jnrC
,jnrD
;
75 int jnrE
,jnrF
,jnrG
,jnrH
;
76 int jnrlistA
,jnrlistB
,jnrlistC
,jnrlistD
;
77 int jnrlistE
,jnrlistF
,jnrlistG
,jnrlistH
;
78 int j_coord_offsetA
,j_coord_offsetB
,j_coord_offsetC
,j_coord_offsetD
;
79 int j_coord_offsetE
,j_coord_offsetF
,j_coord_offsetG
,j_coord_offsetH
;
80 int *iinr
,*jindex
,*jjnr
,*shiftidx
,*gid
;
82 real
*shiftvec
,*fshift
,*x
,*f
;
83 real
*fjptrA
,*fjptrB
,*fjptrC
,*fjptrD
,*fjptrE
,*fjptrF
,*fjptrG
,*fjptrH
;
85 __m256 tx
,ty
,tz
,fscal
,rcutoff
,rcutoff2
,jidxall
;
86 real
* vdwioffsetptr0
;
87 real
* vdwgridioffsetptr0
;
88 __m256 ix0
,iy0
,iz0
,fix0
,fiy0
,fiz0
,iq0
,isai0
;
89 int vdwjidx0A
,vdwjidx0B
,vdwjidx0C
,vdwjidx0D
,vdwjidx0E
,vdwjidx0F
,vdwjidx0G
,vdwjidx0H
;
90 __m256 jx0
,jy0
,jz0
,fjx0
,fjy0
,fjz0
,jq0
,isaj0
;
91 __m256 dx00
,dy00
,dz00
,rsq00
,rinv00
,rinvsq00
,r00
,qq00
,c6_00
,c12_00
;
93 __m256 rinvsix
,rvdw
,vvdw
,vvdw6
,vvdw12
,fvdw
,fvdw6
,fvdw12
,vvdwsum
,sh_vdw_invrcut6
;
96 __m256 one_sixth
= _mm256_set1_ps(1.0/6.0);
97 __m256 one_twelfth
= _mm256_set1_ps(1.0/12.0);
100 __m256 ewclj
,ewclj2
,ewclj6
,ewcljrsq
,poly
,exponent
,f6A
,f6B
,sh_lj_ewald
;
101 __m256 one_half
= _mm256_set1_ps(0.5);
102 __m256 minus_one
= _mm256_set1_ps(-1.0);
103 __m256 dummy_mask
,cutoff_mask
;
104 __m256 signbit
= _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
105 __m256 one
= _mm256_set1_ps(1.0);
106 __m256 two
= _mm256_set1_ps(2.0);
112 jindex
= nlist
->jindex
;
114 shiftidx
= nlist
->shift
;
116 shiftvec
= fr
->shift_vec
[0];
117 fshift
= fr
->fshift
[0];
118 nvdwtype
= fr
->ntype
;
120 vdwtype
= mdatoms
->typeA
;
121 vdwgridparam
= fr
->ljpme_c6grid
;
122 sh_lj_ewald
= _mm256_set1_ps(fr
->ic
->sh_lj_ewald
);
123 ewclj
= _mm256_set1_ps(fr
->ewaldcoeff_lj
);
124 ewclj2
= _mm256_mul_ps(minus_one
,_mm256_mul_ps(ewclj
,ewclj
));
126 rcutoff_scalar
= fr
->rvdw
;
127 rcutoff
= _mm256_set1_ps(rcutoff_scalar
);
128 rcutoff2
= _mm256_mul_ps(rcutoff
,rcutoff
);
130 sh_vdw_invrcut6
= _mm256_set1_ps(fr
->ic
->sh_invrc6
);
131 rvdw
= _mm256_set1_ps(fr
->rvdw
);
133 /* Avoid stupid compiler warnings */
134 jnrA
= jnrB
= jnrC
= jnrD
= jnrE
= jnrF
= jnrG
= jnrH
= 0;
147 for(iidx
=0;iidx
<4*DIM
;iidx
++)
152 /* Start outer loop over neighborlists */
153 for(iidx
=0; iidx
<nri
; iidx
++)
155 /* Load shift vector for this list */
156 i_shift_offset
= DIM
*shiftidx
[iidx
];
158 /* Load limits for loop over neighbors */
159 j_index_start
= jindex
[iidx
];
160 j_index_end
= jindex
[iidx
+1];
162 /* Get outer coordinate index */
164 i_coord_offset
= DIM
*inr
;
166 /* Load i particle coords and add shift vector */
167 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec
+i_shift_offset
,x
+i_coord_offset
,&ix0
,&iy0
,&iz0
);
169 fix0
= _mm256_setzero_ps();
170 fiy0
= _mm256_setzero_ps();
171 fiz0
= _mm256_setzero_ps();
173 /* Load parameters for i particles */
174 vdwioffsetptr0
= vdwparam
+2*nvdwtype
*vdwtype
[inr
+0];
175 vdwgridioffsetptr0
= vdwgridparam
+2*nvdwtype
*vdwtype
[inr
+0];
177 /* Reset potential sums */
178 vvdwsum
= _mm256_setzero_ps();
180 /* Start inner kernel loop */
181 for(jidx
=j_index_start
; jidx
<j_index_end
&& jjnr
[jidx
+7]>=0; jidx
+=8)
184 /* Get j neighbor index, and coordinate index */
193 j_coord_offsetA
= DIM
*jnrA
;
194 j_coord_offsetB
= DIM
*jnrB
;
195 j_coord_offsetC
= DIM
*jnrC
;
196 j_coord_offsetD
= DIM
*jnrD
;
197 j_coord_offsetE
= DIM
*jnrE
;
198 j_coord_offsetF
= DIM
*jnrF
;
199 j_coord_offsetG
= DIM
*jnrG
;
200 j_coord_offsetH
= DIM
*jnrH
;
202 /* load j atom coordinates */
203 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x
+j_coord_offsetA
,x
+j_coord_offsetB
,
204 x
+j_coord_offsetC
,x
+j_coord_offsetD
,
205 x
+j_coord_offsetE
,x
+j_coord_offsetF
,
206 x
+j_coord_offsetG
,x
+j_coord_offsetH
,
209 /* Calculate displacement vector */
210 dx00
= _mm256_sub_ps(ix0
,jx0
);
211 dy00
= _mm256_sub_ps(iy0
,jy0
);
212 dz00
= _mm256_sub_ps(iz0
,jz0
);
214 /* Calculate squared distance and things based on it */
215 rsq00
= gmx_mm256_calc_rsq_ps(dx00
,dy00
,dz00
);
217 rinv00
= gmx_mm256_invsqrt_ps(rsq00
);
219 rinvsq00
= _mm256_mul_ps(rinv00
,rinv00
);
221 /* Load parameters for j particles */
222 vdwjidx0A
= 2*vdwtype
[jnrA
+0];
223 vdwjidx0B
= 2*vdwtype
[jnrB
+0];
224 vdwjidx0C
= 2*vdwtype
[jnrC
+0];
225 vdwjidx0D
= 2*vdwtype
[jnrD
+0];
226 vdwjidx0E
= 2*vdwtype
[jnrE
+0];
227 vdwjidx0F
= 2*vdwtype
[jnrF
+0];
228 vdwjidx0G
= 2*vdwtype
[jnrG
+0];
229 vdwjidx0H
= 2*vdwtype
[jnrH
+0];
231 /**************************
232 * CALCULATE INTERACTIONS *
233 **************************/
235 if (gmx_mm256_any_lt(rsq00
,rcutoff2
))
238 r00
= _mm256_mul_ps(rsq00
,rinv00
);
240 /* Compute parameters for interactions between i and j atoms */
241 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0
+vdwjidx0A
,
242 vdwioffsetptr0
+vdwjidx0B
,
243 vdwioffsetptr0
+vdwjidx0C
,
244 vdwioffsetptr0
+vdwjidx0D
,
245 vdwioffsetptr0
+vdwjidx0E
,
246 vdwioffsetptr0
+vdwjidx0F
,
247 vdwioffsetptr0
+vdwjidx0G
,
248 vdwioffsetptr0
+vdwjidx0H
,
251 c6grid_00
= gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0
+vdwjidx0A
,
252 vdwgridioffsetptr0
+vdwjidx0B
,
253 vdwgridioffsetptr0
+vdwjidx0C
,
254 vdwgridioffsetptr0
+vdwjidx0D
,
255 vdwgridioffsetptr0
+vdwjidx0E
,
256 vdwgridioffsetptr0
+vdwjidx0F
,
257 vdwgridioffsetptr0
+vdwjidx0G
,
258 vdwgridioffsetptr0
+vdwjidx0H
);
260 /* Analytical LJ-PME */
261 rinvsix
= _mm256_mul_ps(_mm256_mul_ps(rinvsq00
,rinvsq00
),rinvsq00
);
262 ewcljrsq
= _mm256_mul_ps(ewclj2
,rsq00
);
263 ewclj6
= _mm256_mul_ps(ewclj2
,_mm256_mul_ps(ewclj2
,ewclj2
));
264 exponent
= gmx_simd_exp_r(ewcljrsq
);
265 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
266 poly
= _mm256_mul_ps(exponent
,_mm256_add_ps(_mm256_sub_ps(one
,ewcljrsq
),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq
,ewcljrsq
),one_half
)));
267 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
268 vvdw6
= _mm256_mul_ps(_mm256_sub_ps(c6_00
,_mm256_mul_ps(c6grid_00
,_mm256_sub_ps(one
,poly
))),rinvsix
);
269 vvdw12
= _mm256_mul_ps(c12_00
,_mm256_mul_ps(rinvsix
,rinvsix
));
270 vvdw
= _mm256_sub_ps(_mm256_mul_ps( _mm256_sub_ps(vvdw12
, _mm256_mul_ps(c12_00
,_mm256_mul_ps(sh_vdw_invrcut6
,sh_vdw_invrcut6
))), one_twelfth
) ,
271 _mm256_mul_ps( _mm256_sub_ps(vvdw6
,_mm256_add_ps(_mm256_mul_ps(c6_00
,sh_vdw_invrcut6
),_mm256_mul_ps(c6grid_00
,sh_lj_ewald
))),one_sixth
));
272 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
273 fvdw
= _mm256_mul_ps(_mm256_sub_ps(vvdw12
,_mm256_sub_ps(vvdw6
,_mm256_mul_ps(_mm256_mul_ps(c6grid_00
,one_sixth
),_mm256_mul_ps(exponent
,ewclj6
)))),rinvsq00
);
275 cutoff_mask
= _mm256_cmp_ps(rsq00
,rcutoff2
,_CMP_LT_OQ
);
277 /* Update potential sum for this i atom from the interaction with this j atom. */
278 vvdw
= _mm256_and_ps(vvdw
,cutoff_mask
);
279 vvdwsum
= _mm256_add_ps(vvdwsum
,vvdw
);
283 fscal
= _mm256_and_ps(fscal
,cutoff_mask
);
285 /* Calculate temporary vectorial force */
286 tx
= _mm256_mul_ps(fscal
,dx00
);
287 ty
= _mm256_mul_ps(fscal
,dy00
);
288 tz
= _mm256_mul_ps(fscal
,dz00
);
290 /* Update vectorial force */
291 fix0
= _mm256_add_ps(fix0
,tx
);
292 fiy0
= _mm256_add_ps(fiy0
,ty
);
293 fiz0
= _mm256_add_ps(fiz0
,tz
);
295 fjptrA
= f
+j_coord_offsetA
;
296 fjptrB
= f
+j_coord_offsetB
;
297 fjptrC
= f
+j_coord_offsetC
;
298 fjptrD
= f
+j_coord_offsetD
;
299 fjptrE
= f
+j_coord_offsetE
;
300 fjptrF
= f
+j_coord_offsetF
;
301 fjptrG
= f
+j_coord_offsetG
;
302 fjptrH
= f
+j_coord_offsetH
;
303 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA
,fjptrB
,fjptrC
,fjptrD
,fjptrE
,fjptrF
,fjptrG
,fjptrH
,tx
,ty
,tz
);
307 /* Inner loop uses 62 flops */
313 /* Get j neighbor index, and coordinate index */
314 jnrlistA
= jjnr
[jidx
];
315 jnrlistB
= jjnr
[jidx
+1];
316 jnrlistC
= jjnr
[jidx
+2];
317 jnrlistD
= jjnr
[jidx
+3];
318 jnrlistE
= jjnr
[jidx
+4];
319 jnrlistF
= jjnr
[jidx
+5];
320 jnrlistG
= jjnr
[jidx
+6];
321 jnrlistH
= jjnr
[jidx
+7];
322 /* Sign of each element will be negative for non-real atoms.
323 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
324 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
326 dummy_mask
= gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i
*)(jjnr
+jidx
+4)),_mm_setzero_si128())),
327 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i
*)(jjnr
+jidx
)),_mm_setzero_si128())));
329 jnrA
= (jnrlistA
>=0) ? jnrlistA
: 0;
330 jnrB
= (jnrlistB
>=0) ? jnrlistB
: 0;
331 jnrC
= (jnrlistC
>=0) ? jnrlistC
: 0;
332 jnrD
= (jnrlistD
>=0) ? jnrlistD
: 0;
333 jnrE
= (jnrlistE
>=0) ? jnrlistE
: 0;
334 jnrF
= (jnrlistF
>=0) ? jnrlistF
: 0;
335 jnrG
= (jnrlistG
>=0) ? jnrlistG
: 0;
336 jnrH
= (jnrlistH
>=0) ? jnrlistH
: 0;
337 j_coord_offsetA
= DIM
*jnrA
;
338 j_coord_offsetB
= DIM
*jnrB
;
339 j_coord_offsetC
= DIM
*jnrC
;
340 j_coord_offsetD
= DIM
*jnrD
;
341 j_coord_offsetE
= DIM
*jnrE
;
342 j_coord_offsetF
= DIM
*jnrF
;
343 j_coord_offsetG
= DIM
*jnrG
;
344 j_coord_offsetH
= DIM
*jnrH
;
346 /* load j atom coordinates */
347 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x
+j_coord_offsetA
,x
+j_coord_offsetB
,
348 x
+j_coord_offsetC
,x
+j_coord_offsetD
,
349 x
+j_coord_offsetE
,x
+j_coord_offsetF
,
350 x
+j_coord_offsetG
,x
+j_coord_offsetH
,
353 /* Calculate displacement vector */
354 dx00
= _mm256_sub_ps(ix0
,jx0
);
355 dy00
= _mm256_sub_ps(iy0
,jy0
);
356 dz00
= _mm256_sub_ps(iz0
,jz0
);
358 /* Calculate squared distance and things based on it */
359 rsq00
= gmx_mm256_calc_rsq_ps(dx00
,dy00
,dz00
);
361 rinv00
= gmx_mm256_invsqrt_ps(rsq00
);
363 rinvsq00
= _mm256_mul_ps(rinv00
,rinv00
);
365 /* Load parameters for j particles */
366 vdwjidx0A
= 2*vdwtype
[jnrA
+0];
367 vdwjidx0B
= 2*vdwtype
[jnrB
+0];
368 vdwjidx0C
= 2*vdwtype
[jnrC
+0];
369 vdwjidx0D
= 2*vdwtype
[jnrD
+0];
370 vdwjidx0E
= 2*vdwtype
[jnrE
+0];
371 vdwjidx0F
= 2*vdwtype
[jnrF
+0];
372 vdwjidx0G
= 2*vdwtype
[jnrG
+0];
373 vdwjidx0H
= 2*vdwtype
[jnrH
+0];
375 /**************************
376 * CALCULATE INTERACTIONS *
377 **************************/
379 if (gmx_mm256_any_lt(rsq00
,rcutoff2
))
382 r00
= _mm256_mul_ps(rsq00
,rinv00
);
383 r00
= _mm256_andnot_ps(dummy_mask
,r00
);
385 /* Compute parameters for interactions between i and j atoms */
386 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0
+vdwjidx0A
,
387 vdwioffsetptr0
+vdwjidx0B
,
388 vdwioffsetptr0
+vdwjidx0C
,
389 vdwioffsetptr0
+vdwjidx0D
,
390 vdwioffsetptr0
+vdwjidx0E
,
391 vdwioffsetptr0
+vdwjidx0F
,
392 vdwioffsetptr0
+vdwjidx0G
,
393 vdwioffsetptr0
+vdwjidx0H
,
396 c6grid_00
= gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0
+vdwjidx0A
,
397 vdwgridioffsetptr0
+vdwjidx0B
,
398 vdwgridioffsetptr0
+vdwjidx0C
,
399 vdwgridioffsetptr0
+vdwjidx0D
,
400 vdwgridioffsetptr0
+vdwjidx0E
,
401 vdwgridioffsetptr0
+vdwjidx0F
,
402 vdwgridioffsetptr0
+vdwjidx0G
,
403 vdwgridioffsetptr0
+vdwjidx0H
);
405 /* Analytical LJ-PME */
406 rinvsix
= _mm256_mul_ps(_mm256_mul_ps(rinvsq00
,rinvsq00
),rinvsq00
);
407 ewcljrsq
= _mm256_mul_ps(ewclj2
,rsq00
);
408 ewclj6
= _mm256_mul_ps(ewclj2
,_mm256_mul_ps(ewclj2
,ewclj2
));
409 exponent
= gmx_simd_exp_r(ewcljrsq
);
410 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
411 poly
= _mm256_mul_ps(exponent
,_mm256_add_ps(_mm256_sub_ps(one
,ewcljrsq
),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq
,ewcljrsq
),one_half
)));
412 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
413 vvdw6
= _mm256_mul_ps(_mm256_sub_ps(c6_00
,_mm256_mul_ps(c6grid_00
,_mm256_sub_ps(one
,poly
))),rinvsix
);
414 vvdw12
= _mm256_mul_ps(c12_00
,_mm256_mul_ps(rinvsix
,rinvsix
));
415 vvdw
= _mm256_sub_ps(_mm256_mul_ps( _mm256_sub_ps(vvdw12
, _mm256_mul_ps(c12_00
,_mm256_mul_ps(sh_vdw_invrcut6
,sh_vdw_invrcut6
))), one_twelfth
) ,
416 _mm256_mul_ps( _mm256_sub_ps(vvdw6
,_mm256_add_ps(_mm256_mul_ps(c6_00
,sh_vdw_invrcut6
),_mm256_mul_ps(c6grid_00
,sh_lj_ewald
))),one_sixth
));
417 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
418 fvdw
= _mm256_mul_ps(_mm256_sub_ps(vvdw12
,_mm256_sub_ps(vvdw6
,_mm256_mul_ps(_mm256_mul_ps(c6grid_00
,one_sixth
),_mm256_mul_ps(exponent
,ewclj6
)))),rinvsq00
);
420 cutoff_mask
= _mm256_cmp_ps(rsq00
,rcutoff2
,_CMP_LT_OQ
);
422 /* Update potential sum for this i atom from the interaction with this j atom. */
423 vvdw
= _mm256_and_ps(vvdw
,cutoff_mask
);
424 vvdw
= _mm256_andnot_ps(dummy_mask
,vvdw
);
425 vvdwsum
= _mm256_add_ps(vvdwsum
,vvdw
);
429 fscal
= _mm256_and_ps(fscal
,cutoff_mask
);
431 fscal
= _mm256_andnot_ps(dummy_mask
,fscal
);
433 /* Calculate temporary vectorial force */
434 tx
= _mm256_mul_ps(fscal
,dx00
);
435 ty
= _mm256_mul_ps(fscal
,dy00
);
436 tz
= _mm256_mul_ps(fscal
,dz00
);
438 /* Update vectorial force */
439 fix0
= _mm256_add_ps(fix0
,tx
);
440 fiy0
= _mm256_add_ps(fiy0
,ty
);
441 fiz0
= _mm256_add_ps(fiz0
,tz
);
443 fjptrA
= (jnrlistA
>=0) ? f
+j_coord_offsetA
: scratch
;
444 fjptrB
= (jnrlistB
>=0) ? f
+j_coord_offsetB
: scratch
;
445 fjptrC
= (jnrlistC
>=0) ? f
+j_coord_offsetC
: scratch
;
446 fjptrD
= (jnrlistD
>=0) ? f
+j_coord_offsetD
: scratch
;
447 fjptrE
= (jnrlistE
>=0) ? f
+j_coord_offsetE
: scratch
;
448 fjptrF
= (jnrlistF
>=0) ? f
+j_coord_offsetF
: scratch
;
449 fjptrG
= (jnrlistG
>=0) ? f
+j_coord_offsetG
: scratch
;
450 fjptrH
= (jnrlistH
>=0) ? f
+j_coord_offsetH
: scratch
;
451 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA
,fjptrB
,fjptrC
,fjptrD
,fjptrE
,fjptrF
,fjptrG
,fjptrH
,tx
,ty
,tz
);
455 /* Inner loop uses 63 flops */
458 /* End of innermost loop */
460 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0
,fiy0
,fiz0
,
461 f
+i_coord_offset
,fshift
+i_shift_offset
);
464 /* Update potential energies */
465 gmx_mm256_update_1pot_ps(vvdwsum
,kernel_data
->energygrp_vdw
+ggid
);
467 /* Increment number of inner iterations */
468 inneriter
+= j_index_end
- j_index_start
;
470 /* Outer loop uses 7 flops */
473 /* Increment number of outer iterations */
476 /* Update outer/inner flops */
478 inc_nrnb(nrnb
,eNR_NBKERNEL_VDW_VF
,outeriter
*7 + inneriter
*63);
481 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_avx_256_single
482 * Electrostatics interaction: None
483 * VdW interaction: LJEwald
484 * Geometry: Particle-Particle
485 * Calculate force/pot: Force
488 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_avx_256_single
489 (t_nblist
* gmx_restrict nlist
,
490 rvec
* gmx_restrict xx
,
491 rvec
* gmx_restrict ff
,
492 t_forcerec
* gmx_restrict fr
,
493 t_mdatoms
* gmx_restrict mdatoms
,
494 nb_kernel_data_t gmx_unused
* gmx_restrict kernel_data
,
495 t_nrnb
* gmx_restrict nrnb
)
497 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
498 * just 0 for non-waters.
499 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
500 * jnr indices corresponding to data put in the four positions in the SIMD register.
502 int i_shift_offset
,i_coord_offset
,outeriter
,inneriter
;
503 int j_index_start
,j_index_end
,jidx
,nri
,inr
,ggid
,iidx
;
504 int jnrA
,jnrB
,jnrC
,jnrD
;
505 int jnrE
,jnrF
,jnrG
,jnrH
;
506 int jnrlistA
,jnrlistB
,jnrlistC
,jnrlistD
;
507 int jnrlistE
,jnrlistF
,jnrlistG
,jnrlistH
;
508 int j_coord_offsetA
,j_coord_offsetB
,j_coord_offsetC
,j_coord_offsetD
;
509 int j_coord_offsetE
,j_coord_offsetF
,j_coord_offsetG
,j_coord_offsetH
;
510 int *iinr
,*jindex
,*jjnr
,*shiftidx
,*gid
;
512 real
*shiftvec
,*fshift
,*x
,*f
;
513 real
*fjptrA
,*fjptrB
,*fjptrC
,*fjptrD
,*fjptrE
,*fjptrF
,*fjptrG
,*fjptrH
;
515 __m256 tx
,ty
,tz
,fscal
,rcutoff
,rcutoff2
,jidxall
;
516 real
* vdwioffsetptr0
;
517 real
* vdwgridioffsetptr0
;
518 __m256 ix0
,iy0
,iz0
,fix0
,fiy0
,fiz0
,iq0
,isai0
;
519 int vdwjidx0A
,vdwjidx0B
,vdwjidx0C
,vdwjidx0D
,vdwjidx0E
,vdwjidx0F
,vdwjidx0G
,vdwjidx0H
;
520 __m256 jx0
,jy0
,jz0
,fjx0
,fjy0
,fjz0
,jq0
,isaj0
;
521 __m256 dx00
,dy00
,dz00
,rsq00
,rinv00
,rinvsq00
,r00
,qq00
,c6_00
,c12_00
;
523 __m256 rinvsix
,rvdw
,vvdw
,vvdw6
,vvdw12
,fvdw
,fvdw6
,fvdw12
,vvdwsum
,sh_vdw_invrcut6
;
526 __m256 one_sixth
= _mm256_set1_ps(1.0/6.0);
527 __m256 one_twelfth
= _mm256_set1_ps(1.0/12.0);
530 __m256 ewclj
,ewclj2
,ewclj6
,ewcljrsq
,poly
,exponent
,f6A
,f6B
,sh_lj_ewald
;
531 __m256 one_half
= _mm256_set1_ps(0.5);
532 __m256 minus_one
= _mm256_set1_ps(-1.0);
533 __m256 dummy_mask
,cutoff_mask
;
534 __m256 signbit
= _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
535 __m256 one
= _mm256_set1_ps(1.0);
536 __m256 two
= _mm256_set1_ps(2.0);
542 jindex
= nlist
->jindex
;
544 shiftidx
= nlist
->shift
;
546 shiftvec
= fr
->shift_vec
[0];
547 fshift
= fr
->fshift
[0];
548 nvdwtype
= fr
->ntype
;
550 vdwtype
= mdatoms
->typeA
;
551 vdwgridparam
= fr
->ljpme_c6grid
;
552 sh_lj_ewald
= _mm256_set1_ps(fr
->ic
->sh_lj_ewald
);
553 ewclj
= _mm256_set1_ps(fr
->ewaldcoeff_lj
);
554 ewclj2
= _mm256_mul_ps(minus_one
,_mm256_mul_ps(ewclj
,ewclj
));
556 rcutoff_scalar
= fr
->rvdw
;
557 rcutoff
= _mm256_set1_ps(rcutoff_scalar
);
558 rcutoff2
= _mm256_mul_ps(rcutoff
,rcutoff
);
560 sh_vdw_invrcut6
= _mm256_set1_ps(fr
->ic
->sh_invrc6
);
561 rvdw
= _mm256_set1_ps(fr
->rvdw
);
563 /* Avoid stupid compiler warnings */
564 jnrA
= jnrB
= jnrC
= jnrD
= jnrE
= jnrF
= jnrG
= jnrH
= 0;
577 for(iidx
=0;iidx
<4*DIM
;iidx
++)
582 /* Start outer loop over neighborlists */
583 for(iidx
=0; iidx
<nri
; iidx
++)
585 /* Load shift vector for this list */
586 i_shift_offset
= DIM
*shiftidx
[iidx
];
588 /* Load limits for loop over neighbors */
589 j_index_start
= jindex
[iidx
];
590 j_index_end
= jindex
[iidx
+1];
592 /* Get outer coordinate index */
594 i_coord_offset
= DIM
*inr
;
596 /* Load i particle coords and add shift vector */
597 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec
+i_shift_offset
,x
+i_coord_offset
,&ix0
,&iy0
,&iz0
);
599 fix0
= _mm256_setzero_ps();
600 fiy0
= _mm256_setzero_ps();
601 fiz0
= _mm256_setzero_ps();
603 /* Load parameters for i particles */
604 vdwioffsetptr0
= vdwparam
+2*nvdwtype
*vdwtype
[inr
+0];
605 vdwgridioffsetptr0
= vdwgridparam
+2*nvdwtype
*vdwtype
[inr
+0];
607 /* Start inner kernel loop */
608 for(jidx
=j_index_start
; jidx
<j_index_end
&& jjnr
[jidx
+7]>=0; jidx
+=8)
611 /* Get j neighbor index, and coordinate index */
620 j_coord_offsetA
= DIM
*jnrA
;
621 j_coord_offsetB
= DIM
*jnrB
;
622 j_coord_offsetC
= DIM
*jnrC
;
623 j_coord_offsetD
= DIM
*jnrD
;
624 j_coord_offsetE
= DIM
*jnrE
;
625 j_coord_offsetF
= DIM
*jnrF
;
626 j_coord_offsetG
= DIM
*jnrG
;
627 j_coord_offsetH
= DIM
*jnrH
;
629 /* load j atom coordinates */
630 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x
+j_coord_offsetA
,x
+j_coord_offsetB
,
631 x
+j_coord_offsetC
,x
+j_coord_offsetD
,
632 x
+j_coord_offsetE
,x
+j_coord_offsetF
,
633 x
+j_coord_offsetG
,x
+j_coord_offsetH
,
636 /* Calculate displacement vector */
637 dx00
= _mm256_sub_ps(ix0
,jx0
);
638 dy00
= _mm256_sub_ps(iy0
,jy0
);
639 dz00
= _mm256_sub_ps(iz0
,jz0
);
641 /* Calculate squared distance and things based on it */
642 rsq00
= gmx_mm256_calc_rsq_ps(dx00
,dy00
,dz00
);
644 rinv00
= gmx_mm256_invsqrt_ps(rsq00
);
646 rinvsq00
= _mm256_mul_ps(rinv00
,rinv00
);
648 /* Load parameters for j particles */
649 vdwjidx0A
= 2*vdwtype
[jnrA
+0];
650 vdwjidx0B
= 2*vdwtype
[jnrB
+0];
651 vdwjidx0C
= 2*vdwtype
[jnrC
+0];
652 vdwjidx0D
= 2*vdwtype
[jnrD
+0];
653 vdwjidx0E
= 2*vdwtype
[jnrE
+0];
654 vdwjidx0F
= 2*vdwtype
[jnrF
+0];
655 vdwjidx0G
= 2*vdwtype
[jnrG
+0];
656 vdwjidx0H
= 2*vdwtype
[jnrH
+0];
658 /**************************
659 * CALCULATE INTERACTIONS *
660 **************************/
662 if (gmx_mm256_any_lt(rsq00
,rcutoff2
))
665 r00
= _mm256_mul_ps(rsq00
,rinv00
);
667 /* Compute parameters for interactions between i and j atoms */
668 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0
+vdwjidx0A
,
669 vdwioffsetptr0
+vdwjidx0B
,
670 vdwioffsetptr0
+vdwjidx0C
,
671 vdwioffsetptr0
+vdwjidx0D
,
672 vdwioffsetptr0
+vdwjidx0E
,
673 vdwioffsetptr0
+vdwjidx0F
,
674 vdwioffsetptr0
+vdwjidx0G
,
675 vdwioffsetptr0
+vdwjidx0H
,
678 c6grid_00
= gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0
+vdwjidx0A
,
679 vdwgridioffsetptr0
+vdwjidx0B
,
680 vdwgridioffsetptr0
+vdwjidx0C
,
681 vdwgridioffsetptr0
+vdwjidx0D
,
682 vdwgridioffsetptr0
+vdwjidx0E
,
683 vdwgridioffsetptr0
+vdwjidx0F
,
684 vdwgridioffsetptr0
+vdwjidx0G
,
685 vdwgridioffsetptr0
+vdwjidx0H
);
687 /* Analytical LJ-PME */
688 rinvsix
= _mm256_mul_ps(_mm256_mul_ps(rinvsq00
,rinvsq00
),rinvsq00
);
689 ewcljrsq
= _mm256_mul_ps(ewclj2
,rsq00
);
690 ewclj6
= _mm256_mul_ps(ewclj2
,_mm256_mul_ps(ewclj2
,ewclj2
));
691 exponent
= gmx_simd_exp_r(ewcljrsq
);
692 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
693 poly
= _mm256_mul_ps(exponent
,_mm256_add_ps(_mm256_sub_ps(one
,ewcljrsq
),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq
,ewcljrsq
),one_half
)));
694 /* f6A = 6 * C6grid * (1 - poly) */
695 f6A
= _mm256_mul_ps(c6grid_00
,_mm256_sub_ps(one
,poly
));
696 /* f6B = C6grid * exponent * beta^6 */
697 f6B
= _mm256_mul_ps(_mm256_mul_ps(c6grid_00
,one_sixth
),_mm256_mul_ps(exponent
,ewclj6
));
698 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
699 fvdw
= _mm256_mul_ps(_mm256_add_ps(_mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_00
,rinvsix
),_mm256_sub_ps(c6_00
,f6A
)),rinvsix
),f6B
),rinvsq00
);
701 cutoff_mask
= _mm256_cmp_ps(rsq00
,rcutoff2
,_CMP_LT_OQ
);
705 fscal
= _mm256_and_ps(fscal
,cutoff_mask
);
707 /* Calculate temporary vectorial force */
708 tx
= _mm256_mul_ps(fscal
,dx00
);
709 ty
= _mm256_mul_ps(fscal
,dy00
);
710 tz
= _mm256_mul_ps(fscal
,dz00
);
712 /* Update vectorial force */
713 fix0
= _mm256_add_ps(fix0
,tx
);
714 fiy0
= _mm256_add_ps(fiy0
,ty
);
715 fiz0
= _mm256_add_ps(fiz0
,tz
);
717 fjptrA
= f
+j_coord_offsetA
;
718 fjptrB
= f
+j_coord_offsetB
;
719 fjptrC
= f
+j_coord_offsetC
;
720 fjptrD
= f
+j_coord_offsetD
;
721 fjptrE
= f
+j_coord_offsetE
;
722 fjptrF
= f
+j_coord_offsetF
;
723 fjptrG
= f
+j_coord_offsetG
;
724 fjptrH
= f
+j_coord_offsetH
;
725 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA
,fjptrB
,fjptrC
,fjptrD
,fjptrE
,fjptrF
,fjptrG
,fjptrH
,tx
,ty
,tz
);
729 /* Inner loop uses 49 flops */
735 /* Get j neighbor index, and coordinate index */
736 jnrlistA
= jjnr
[jidx
];
737 jnrlistB
= jjnr
[jidx
+1];
738 jnrlistC
= jjnr
[jidx
+2];
739 jnrlistD
= jjnr
[jidx
+3];
740 jnrlistE
= jjnr
[jidx
+4];
741 jnrlistF
= jjnr
[jidx
+5];
742 jnrlistG
= jjnr
[jidx
+6];
743 jnrlistH
= jjnr
[jidx
+7];
744 /* Sign of each element will be negative for non-real atoms.
745 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
746 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
748 dummy_mask
= gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i
*)(jjnr
+jidx
+4)),_mm_setzero_si128())),
749 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i
*)(jjnr
+jidx
)),_mm_setzero_si128())));
751 jnrA
= (jnrlistA
>=0) ? jnrlistA
: 0;
752 jnrB
= (jnrlistB
>=0) ? jnrlistB
: 0;
753 jnrC
= (jnrlistC
>=0) ? jnrlistC
: 0;
754 jnrD
= (jnrlistD
>=0) ? jnrlistD
: 0;
755 jnrE
= (jnrlistE
>=0) ? jnrlistE
: 0;
756 jnrF
= (jnrlistF
>=0) ? jnrlistF
: 0;
757 jnrG
= (jnrlistG
>=0) ? jnrlistG
: 0;
758 jnrH
= (jnrlistH
>=0) ? jnrlistH
: 0;
759 j_coord_offsetA
= DIM
*jnrA
;
760 j_coord_offsetB
= DIM
*jnrB
;
761 j_coord_offsetC
= DIM
*jnrC
;
762 j_coord_offsetD
= DIM
*jnrD
;
763 j_coord_offsetE
= DIM
*jnrE
;
764 j_coord_offsetF
= DIM
*jnrF
;
765 j_coord_offsetG
= DIM
*jnrG
;
766 j_coord_offsetH
= DIM
*jnrH
;
768 /* load j atom coordinates */
769 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x
+j_coord_offsetA
,x
+j_coord_offsetB
,
770 x
+j_coord_offsetC
,x
+j_coord_offsetD
,
771 x
+j_coord_offsetE
,x
+j_coord_offsetF
,
772 x
+j_coord_offsetG
,x
+j_coord_offsetH
,
775 /* Calculate displacement vector */
776 dx00
= _mm256_sub_ps(ix0
,jx0
);
777 dy00
= _mm256_sub_ps(iy0
,jy0
);
778 dz00
= _mm256_sub_ps(iz0
,jz0
);
780 /* Calculate squared distance and things based on it */
781 rsq00
= gmx_mm256_calc_rsq_ps(dx00
,dy00
,dz00
);
783 rinv00
= gmx_mm256_invsqrt_ps(rsq00
);
785 rinvsq00
= _mm256_mul_ps(rinv00
,rinv00
);
787 /* Load parameters for j particles */
788 vdwjidx0A
= 2*vdwtype
[jnrA
+0];
789 vdwjidx0B
= 2*vdwtype
[jnrB
+0];
790 vdwjidx0C
= 2*vdwtype
[jnrC
+0];
791 vdwjidx0D
= 2*vdwtype
[jnrD
+0];
792 vdwjidx0E
= 2*vdwtype
[jnrE
+0];
793 vdwjidx0F
= 2*vdwtype
[jnrF
+0];
794 vdwjidx0G
= 2*vdwtype
[jnrG
+0];
795 vdwjidx0H
= 2*vdwtype
[jnrH
+0];
797 /**************************
798 * CALCULATE INTERACTIONS *
799 **************************/
801 if (gmx_mm256_any_lt(rsq00
,rcutoff2
))
804 r00
= _mm256_mul_ps(rsq00
,rinv00
);
805 r00
= _mm256_andnot_ps(dummy_mask
,r00
);
807 /* Compute parameters for interactions between i and j atoms */
808 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0
+vdwjidx0A
,
809 vdwioffsetptr0
+vdwjidx0B
,
810 vdwioffsetptr0
+vdwjidx0C
,
811 vdwioffsetptr0
+vdwjidx0D
,
812 vdwioffsetptr0
+vdwjidx0E
,
813 vdwioffsetptr0
+vdwjidx0F
,
814 vdwioffsetptr0
+vdwjidx0G
,
815 vdwioffsetptr0
+vdwjidx0H
,
818 c6grid_00
= gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0
+vdwjidx0A
,
819 vdwgridioffsetptr0
+vdwjidx0B
,
820 vdwgridioffsetptr0
+vdwjidx0C
,
821 vdwgridioffsetptr0
+vdwjidx0D
,
822 vdwgridioffsetptr0
+vdwjidx0E
,
823 vdwgridioffsetptr0
+vdwjidx0F
,
824 vdwgridioffsetptr0
+vdwjidx0G
,
825 vdwgridioffsetptr0
+vdwjidx0H
);
827 /* Analytical LJ-PME */
828 rinvsix
= _mm256_mul_ps(_mm256_mul_ps(rinvsq00
,rinvsq00
),rinvsq00
);
829 ewcljrsq
= _mm256_mul_ps(ewclj2
,rsq00
);
830 ewclj6
= _mm256_mul_ps(ewclj2
,_mm256_mul_ps(ewclj2
,ewclj2
));
831 exponent
= gmx_simd_exp_r(ewcljrsq
);
832 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
833 poly
= _mm256_mul_ps(exponent
,_mm256_add_ps(_mm256_sub_ps(one
,ewcljrsq
),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq
,ewcljrsq
),one_half
)));
834 /* f6A = 6 * C6grid * (1 - poly) */
835 f6A
= _mm256_mul_ps(c6grid_00
,_mm256_sub_ps(one
,poly
));
836 /* f6B = C6grid * exponent * beta^6 */
837 f6B
= _mm256_mul_ps(_mm256_mul_ps(c6grid_00
,one_sixth
),_mm256_mul_ps(exponent
,ewclj6
));
838 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
839 fvdw
= _mm256_mul_ps(_mm256_add_ps(_mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_00
,rinvsix
),_mm256_sub_ps(c6_00
,f6A
)),rinvsix
),f6B
),rinvsq00
);
841 cutoff_mask
= _mm256_cmp_ps(rsq00
,rcutoff2
,_CMP_LT_OQ
);
845 fscal
= _mm256_and_ps(fscal
,cutoff_mask
);
847 fscal
= _mm256_andnot_ps(dummy_mask
,fscal
);
849 /* Calculate temporary vectorial force */
850 tx
= _mm256_mul_ps(fscal
,dx00
);
851 ty
= _mm256_mul_ps(fscal
,dy00
);
852 tz
= _mm256_mul_ps(fscal
,dz00
);
854 /* Update vectorial force */
855 fix0
= _mm256_add_ps(fix0
,tx
);
856 fiy0
= _mm256_add_ps(fiy0
,ty
);
857 fiz0
= _mm256_add_ps(fiz0
,tz
);
859 fjptrA
= (jnrlistA
>=0) ? f
+j_coord_offsetA
: scratch
;
860 fjptrB
= (jnrlistB
>=0) ? f
+j_coord_offsetB
: scratch
;
861 fjptrC
= (jnrlistC
>=0) ? f
+j_coord_offsetC
: scratch
;
862 fjptrD
= (jnrlistD
>=0) ? f
+j_coord_offsetD
: scratch
;
863 fjptrE
= (jnrlistE
>=0) ? f
+j_coord_offsetE
: scratch
;
864 fjptrF
= (jnrlistF
>=0) ? f
+j_coord_offsetF
: scratch
;
865 fjptrG
= (jnrlistG
>=0) ? f
+j_coord_offsetG
: scratch
;
866 fjptrH
= (jnrlistH
>=0) ? f
+j_coord_offsetH
: scratch
;
867 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA
,fjptrB
,fjptrC
,fjptrD
,fjptrE
,fjptrF
,fjptrG
,fjptrH
,tx
,ty
,tz
);
871 /* Inner loop uses 50 flops */
874 /* End of innermost loop */
876 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0
,fiy0
,fiz0
,
877 f
+i_coord_offset
,fshift
+i_shift_offset
);
879 /* Increment number of inner iterations */
880 inneriter
+= j_index_end
- j_index_start
;
882 /* Outer loop uses 6 flops */
885 /* Increment number of outer iterations */
888 /* Update outer/inner flops */
890 inc_nrnb(nrnb
,eNR_NBKERNEL_VDW_F
,outeriter
*6 + inneriter
*50);