2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2017,2018, 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.
44 #include "../nb_kernel.h"
45 #include "gromacs/gmxlib/nrnb.h"
47 #include "kernelutil_x86_avx_256_single.h"
50 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_avx_256_single
51 * Electrostatics interaction: None
52 * VdW interaction: LJEwald
53 * Geometry: Particle-Particle
54 * Calculate force/pot: PotentialAndForce
57 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_avx_256_single
58 (t_nblist * gmx_restrict nlist,
59 rvec * gmx_restrict xx,
60 rvec * gmx_restrict ff,
61 struct t_forcerec * gmx_restrict fr,
62 t_mdatoms * gmx_restrict mdatoms,
63 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
64 t_nrnb * gmx_restrict nrnb)
66 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
67 * just 0 for non-waters.
68 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
69 * jnr indices corresponding to data put in the four positions in the SIMD register.
71 int i_shift_offset,i_coord_offset,outeriter,inneriter;
72 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
73 int jnrA,jnrB,jnrC,jnrD;
74 int jnrE,jnrF,jnrG,jnrH;
75 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
76 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
77 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
78 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
79 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
81 real *shiftvec,*fshift,*x,*f;
82 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
84 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
85 real * vdwioffsetptr0;
86 real * vdwgridioffsetptr0;
87 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
88 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
89 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
90 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
92 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
95 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
96 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
99 __m256 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
100 __m256 one_half = _mm256_set1_ps(0.5);
101 __m256 minus_one = _mm256_set1_ps(-1.0);
102 __m256 dummy_mask,cutoff_mask;
103 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
104 __m256 one = _mm256_set1_ps(1.0);
105 __m256 two = _mm256_set1_ps(2.0);
111 jindex = nlist->jindex;
113 shiftidx = nlist->shift;
115 shiftvec = fr->shift_vec[0];
116 fshift = fr->fshift[0];
117 nvdwtype = fr->ntype;
119 vdwtype = mdatoms->typeA;
120 vdwgridparam = fr->ljpme_c6grid;
121 sh_lj_ewald = _mm256_set1_ps(fr->ic->sh_lj_ewald);
122 ewclj = _mm256_set1_ps(fr->ic->ewaldcoeff_lj);
123 ewclj2 = _mm256_mul_ps(minus_one,_mm256_mul_ps(ewclj,ewclj));
125 rcutoff_scalar = fr->ic->rvdw;
126 rcutoff = _mm256_set1_ps(rcutoff_scalar);
127 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
129 sh_vdw_invrcut6 = _mm256_set1_ps(fr->ic->sh_invrc6);
130 rvdw = _mm256_set1_ps(fr->ic->rvdw);
132 /* Avoid stupid compiler warnings */
133 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
146 for(iidx=0;iidx<4*DIM;iidx++)
151 /* Start outer loop over neighborlists */
152 for(iidx=0; iidx<nri; iidx++)
154 /* Load shift vector for this list */
155 i_shift_offset = DIM*shiftidx[iidx];
157 /* Load limits for loop over neighbors */
158 j_index_start = jindex[iidx];
159 j_index_end = jindex[iidx+1];
161 /* Get outer coordinate index */
163 i_coord_offset = DIM*inr;
165 /* Load i particle coords and add shift vector */
166 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
168 fix0 = _mm256_setzero_ps();
169 fiy0 = _mm256_setzero_ps();
170 fiz0 = _mm256_setzero_ps();
172 /* Load parameters for i particles */
173 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
174 vdwgridioffsetptr0 = vdwgridparam+2*nvdwtype*vdwtype[inr+0];
176 /* Reset potential sums */
177 vvdwsum = _mm256_setzero_ps();
179 /* Start inner kernel loop */
180 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
183 /* Get j neighbor index, and coordinate index */
192 j_coord_offsetA = DIM*jnrA;
193 j_coord_offsetB = DIM*jnrB;
194 j_coord_offsetC = DIM*jnrC;
195 j_coord_offsetD = DIM*jnrD;
196 j_coord_offsetE = DIM*jnrE;
197 j_coord_offsetF = DIM*jnrF;
198 j_coord_offsetG = DIM*jnrG;
199 j_coord_offsetH = DIM*jnrH;
201 /* load j atom coordinates */
202 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
203 x+j_coord_offsetC,x+j_coord_offsetD,
204 x+j_coord_offsetE,x+j_coord_offsetF,
205 x+j_coord_offsetG,x+j_coord_offsetH,
208 /* Calculate displacement vector */
209 dx00 = _mm256_sub_ps(ix0,jx0);
210 dy00 = _mm256_sub_ps(iy0,jy0);
211 dz00 = _mm256_sub_ps(iz0,jz0);
213 /* Calculate squared distance and things based on it */
214 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
216 rinv00 = avx256_invsqrt_f(rsq00);
218 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
220 /* Load parameters for j particles */
221 vdwjidx0A = 2*vdwtype[jnrA+0];
222 vdwjidx0B = 2*vdwtype[jnrB+0];
223 vdwjidx0C = 2*vdwtype[jnrC+0];
224 vdwjidx0D = 2*vdwtype[jnrD+0];
225 vdwjidx0E = 2*vdwtype[jnrE+0];
226 vdwjidx0F = 2*vdwtype[jnrF+0];
227 vdwjidx0G = 2*vdwtype[jnrG+0];
228 vdwjidx0H = 2*vdwtype[jnrH+0];
230 /**************************
231 * CALCULATE INTERACTIONS *
232 **************************/
234 if (gmx_mm256_any_lt(rsq00,rcutoff2))
237 r00 = _mm256_mul_ps(rsq00,rinv00);
239 /* Compute parameters for interactions between i and j atoms */
240 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
241 vdwioffsetptr0+vdwjidx0B,
242 vdwioffsetptr0+vdwjidx0C,
243 vdwioffsetptr0+vdwjidx0D,
244 vdwioffsetptr0+vdwjidx0E,
245 vdwioffsetptr0+vdwjidx0F,
246 vdwioffsetptr0+vdwjidx0G,
247 vdwioffsetptr0+vdwjidx0H,
250 c6grid_00 = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
251 vdwgridioffsetptr0+vdwjidx0B,
252 vdwgridioffsetptr0+vdwjidx0C,
253 vdwgridioffsetptr0+vdwjidx0D,
254 vdwgridioffsetptr0+vdwjidx0E,
255 vdwgridioffsetptr0+vdwjidx0F,
256 vdwgridioffsetptr0+vdwjidx0G,
257 vdwgridioffsetptr0+vdwjidx0H);
259 /* Analytical LJ-PME */
260 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
261 ewcljrsq = _mm256_mul_ps(ewclj2,rsq00);
262 ewclj6 = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
263 exponent = avx256_exp_f(ewcljrsq);
264 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
265 poly = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
266 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
267 vvdw6 = _mm256_mul_ps(_mm256_sub_ps(c6_00,_mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly))),rinvsix);
268 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
269 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) ,
270 _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));
271 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
272 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);
274 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
276 /* Update potential sum for this i atom from the interaction with this j atom. */
277 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
278 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
282 fscal = _mm256_and_ps(fscal,cutoff_mask);
284 /* Calculate temporary vectorial force */
285 tx = _mm256_mul_ps(fscal,dx00);
286 ty = _mm256_mul_ps(fscal,dy00);
287 tz = _mm256_mul_ps(fscal,dz00);
289 /* Update vectorial force */
290 fix0 = _mm256_add_ps(fix0,tx);
291 fiy0 = _mm256_add_ps(fiy0,ty);
292 fiz0 = _mm256_add_ps(fiz0,tz);
294 fjptrA = f+j_coord_offsetA;
295 fjptrB = f+j_coord_offsetB;
296 fjptrC = f+j_coord_offsetC;
297 fjptrD = f+j_coord_offsetD;
298 fjptrE = f+j_coord_offsetE;
299 fjptrF = f+j_coord_offsetF;
300 fjptrG = f+j_coord_offsetG;
301 fjptrH = f+j_coord_offsetH;
302 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
306 /* Inner loop uses 62 flops */
312 /* Get j neighbor index, and coordinate index */
313 jnrlistA = jjnr[jidx];
314 jnrlistB = jjnr[jidx+1];
315 jnrlistC = jjnr[jidx+2];
316 jnrlistD = jjnr[jidx+3];
317 jnrlistE = jjnr[jidx+4];
318 jnrlistF = jjnr[jidx+5];
319 jnrlistG = jjnr[jidx+6];
320 jnrlistH = jjnr[jidx+7];
321 /* Sign of each element will be negative for non-real atoms.
322 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
323 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
325 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
326 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
328 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
329 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
330 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
331 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
332 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
333 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
334 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
335 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
336 j_coord_offsetA = DIM*jnrA;
337 j_coord_offsetB = DIM*jnrB;
338 j_coord_offsetC = DIM*jnrC;
339 j_coord_offsetD = DIM*jnrD;
340 j_coord_offsetE = DIM*jnrE;
341 j_coord_offsetF = DIM*jnrF;
342 j_coord_offsetG = DIM*jnrG;
343 j_coord_offsetH = DIM*jnrH;
345 /* load j atom coordinates */
346 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
347 x+j_coord_offsetC,x+j_coord_offsetD,
348 x+j_coord_offsetE,x+j_coord_offsetF,
349 x+j_coord_offsetG,x+j_coord_offsetH,
352 /* Calculate displacement vector */
353 dx00 = _mm256_sub_ps(ix0,jx0);
354 dy00 = _mm256_sub_ps(iy0,jy0);
355 dz00 = _mm256_sub_ps(iz0,jz0);
357 /* Calculate squared distance and things based on it */
358 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
360 rinv00 = avx256_invsqrt_f(rsq00);
362 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
364 /* Load parameters for j particles */
365 vdwjidx0A = 2*vdwtype[jnrA+0];
366 vdwjidx0B = 2*vdwtype[jnrB+0];
367 vdwjidx0C = 2*vdwtype[jnrC+0];
368 vdwjidx0D = 2*vdwtype[jnrD+0];
369 vdwjidx0E = 2*vdwtype[jnrE+0];
370 vdwjidx0F = 2*vdwtype[jnrF+0];
371 vdwjidx0G = 2*vdwtype[jnrG+0];
372 vdwjidx0H = 2*vdwtype[jnrH+0];
374 /**************************
375 * CALCULATE INTERACTIONS *
376 **************************/
378 if (gmx_mm256_any_lt(rsq00,rcutoff2))
381 r00 = _mm256_mul_ps(rsq00,rinv00);
382 r00 = _mm256_andnot_ps(dummy_mask,r00);
384 /* Compute parameters for interactions between i and j atoms */
385 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
386 vdwioffsetptr0+vdwjidx0B,
387 vdwioffsetptr0+vdwjidx0C,
388 vdwioffsetptr0+vdwjidx0D,
389 vdwioffsetptr0+vdwjidx0E,
390 vdwioffsetptr0+vdwjidx0F,
391 vdwioffsetptr0+vdwjidx0G,
392 vdwioffsetptr0+vdwjidx0H,
395 c6grid_00 = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
396 vdwgridioffsetptr0+vdwjidx0B,
397 vdwgridioffsetptr0+vdwjidx0C,
398 vdwgridioffsetptr0+vdwjidx0D,
399 vdwgridioffsetptr0+vdwjidx0E,
400 vdwgridioffsetptr0+vdwjidx0F,
401 vdwgridioffsetptr0+vdwjidx0G,
402 vdwgridioffsetptr0+vdwjidx0H);
404 /* Analytical LJ-PME */
405 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
406 ewcljrsq = _mm256_mul_ps(ewclj2,rsq00);
407 ewclj6 = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
408 exponent = avx256_exp_f(ewcljrsq);
409 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
410 poly = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
411 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
412 vvdw6 = _mm256_mul_ps(_mm256_sub_ps(c6_00,_mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly))),rinvsix);
413 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
414 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) ,
415 _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));
416 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
417 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);
419 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
421 /* Update potential sum for this i atom from the interaction with this j atom. */
422 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
423 vvdw = _mm256_andnot_ps(dummy_mask,vvdw);
424 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
428 fscal = _mm256_and_ps(fscal,cutoff_mask);
430 fscal = _mm256_andnot_ps(dummy_mask,fscal);
432 /* Calculate temporary vectorial force */
433 tx = _mm256_mul_ps(fscal,dx00);
434 ty = _mm256_mul_ps(fscal,dy00);
435 tz = _mm256_mul_ps(fscal,dz00);
437 /* Update vectorial force */
438 fix0 = _mm256_add_ps(fix0,tx);
439 fiy0 = _mm256_add_ps(fiy0,ty);
440 fiz0 = _mm256_add_ps(fiz0,tz);
442 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
443 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
444 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
445 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
446 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
447 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
448 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
449 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
450 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
454 /* Inner loop uses 63 flops */
457 /* End of innermost loop */
459 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
460 f+i_coord_offset,fshift+i_shift_offset);
463 /* Update potential energies */
464 gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
466 /* Increment number of inner iterations */
467 inneriter += j_index_end - j_index_start;
469 /* Outer loop uses 7 flops */
472 /* Increment number of outer iterations */
475 /* Update outer/inner flops */
477 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_VF,outeriter*7 + inneriter*63);
480 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_avx_256_single
481 * Electrostatics interaction: None
482 * VdW interaction: LJEwald
483 * Geometry: Particle-Particle
484 * Calculate force/pot: Force
487 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_avx_256_single
488 (t_nblist * gmx_restrict nlist,
489 rvec * gmx_restrict xx,
490 rvec * gmx_restrict ff,
491 struct t_forcerec * gmx_restrict fr,
492 t_mdatoms * gmx_restrict mdatoms,
493 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
494 t_nrnb * gmx_restrict nrnb)
496 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
497 * just 0 for non-waters.
498 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
499 * jnr indices corresponding to data put in the four positions in the SIMD register.
501 int i_shift_offset,i_coord_offset,outeriter,inneriter;
502 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
503 int jnrA,jnrB,jnrC,jnrD;
504 int jnrE,jnrF,jnrG,jnrH;
505 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
506 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
507 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
508 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
509 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
511 real *shiftvec,*fshift,*x,*f;
512 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
514 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
515 real * vdwioffsetptr0;
516 real * vdwgridioffsetptr0;
517 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
518 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
519 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
520 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
522 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
525 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
526 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
529 __m256 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
530 __m256 one_half = _mm256_set1_ps(0.5);
531 __m256 minus_one = _mm256_set1_ps(-1.0);
532 __m256 dummy_mask,cutoff_mask;
533 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
534 __m256 one = _mm256_set1_ps(1.0);
535 __m256 two = _mm256_set1_ps(2.0);
541 jindex = nlist->jindex;
543 shiftidx = nlist->shift;
545 shiftvec = fr->shift_vec[0];
546 fshift = fr->fshift[0];
547 nvdwtype = fr->ntype;
549 vdwtype = mdatoms->typeA;
550 vdwgridparam = fr->ljpme_c6grid;
551 sh_lj_ewald = _mm256_set1_ps(fr->ic->sh_lj_ewald);
552 ewclj = _mm256_set1_ps(fr->ic->ewaldcoeff_lj);
553 ewclj2 = _mm256_mul_ps(minus_one,_mm256_mul_ps(ewclj,ewclj));
555 rcutoff_scalar = fr->ic->rvdw;
556 rcutoff = _mm256_set1_ps(rcutoff_scalar);
557 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
559 sh_vdw_invrcut6 = _mm256_set1_ps(fr->ic->sh_invrc6);
560 rvdw = _mm256_set1_ps(fr->ic->rvdw);
562 /* Avoid stupid compiler warnings */
563 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
576 for(iidx=0;iidx<4*DIM;iidx++)
581 /* Start outer loop over neighborlists */
582 for(iidx=0; iidx<nri; iidx++)
584 /* Load shift vector for this list */
585 i_shift_offset = DIM*shiftidx[iidx];
587 /* Load limits for loop over neighbors */
588 j_index_start = jindex[iidx];
589 j_index_end = jindex[iidx+1];
591 /* Get outer coordinate index */
593 i_coord_offset = DIM*inr;
595 /* Load i particle coords and add shift vector */
596 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
598 fix0 = _mm256_setzero_ps();
599 fiy0 = _mm256_setzero_ps();
600 fiz0 = _mm256_setzero_ps();
602 /* Load parameters for i particles */
603 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
604 vdwgridioffsetptr0 = vdwgridparam+2*nvdwtype*vdwtype[inr+0];
606 /* Start inner kernel loop */
607 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
610 /* Get j neighbor index, and coordinate index */
619 j_coord_offsetA = DIM*jnrA;
620 j_coord_offsetB = DIM*jnrB;
621 j_coord_offsetC = DIM*jnrC;
622 j_coord_offsetD = DIM*jnrD;
623 j_coord_offsetE = DIM*jnrE;
624 j_coord_offsetF = DIM*jnrF;
625 j_coord_offsetG = DIM*jnrG;
626 j_coord_offsetH = DIM*jnrH;
628 /* load j atom coordinates */
629 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
630 x+j_coord_offsetC,x+j_coord_offsetD,
631 x+j_coord_offsetE,x+j_coord_offsetF,
632 x+j_coord_offsetG,x+j_coord_offsetH,
635 /* Calculate displacement vector */
636 dx00 = _mm256_sub_ps(ix0,jx0);
637 dy00 = _mm256_sub_ps(iy0,jy0);
638 dz00 = _mm256_sub_ps(iz0,jz0);
640 /* Calculate squared distance and things based on it */
641 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
643 rinv00 = avx256_invsqrt_f(rsq00);
645 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
647 /* Load parameters for j particles */
648 vdwjidx0A = 2*vdwtype[jnrA+0];
649 vdwjidx0B = 2*vdwtype[jnrB+0];
650 vdwjidx0C = 2*vdwtype[jnrC+0];
651 vdwjidx0D = 2*vdwtype[jnrD+0];
652 vdwjidx0E = 2*vdwtype[jnrE+0];
653 vdwjidx0F = 2*vdwtype[jnrF+0];
654 vdwjidx0G = 2*vdwtype[jnrG+0];
655 vdwjidx0H = 2*vdwtype[jnrH+0];
657 /**************************
658 * CALCULATE INTERACTIONS *
659 **************************/
661 if (gmx_mm256_any_lt(rsq00,rcutoff2))
664 r00 = _mm256_mul_ps(rsq00,rinv00);
666 /* Compute parameters for interactions between i and j atoms */
667 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
668 vdwioffsetptr0+vdwjidx0B,
669 vdwioffsetptr0+vdwjidx0C,
670 vdwioffsetptr0+vdwjidx0D,
671 vdwioffsetptr0+vdwjidx0E,
672 vdwioffsetptr0+vdwjidx0F,
673 vdwioffsetptr0+vdwjidx0G,
674 vdwioffsetptr0+vdwjidx0H,
677 c6grid_00 = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
678 vdwgridioffsetptr0+vdwjidx0B,
679 vdwgridioffsetptr0+vdwjidx0C,
680 vdwgridioffsetptr0+vdwjidx0D,
681 vdwgridioffsetptr0+vdwjidx0E,
682 vdwgridioffsetptr0+vdwjidx0F,
683 vdwgridioffsetptr0+vdwjidx0G,
684 vdwgridioffsetptr0+vdwjidx0H);
686 /* Analytical LJ-PME */
687 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
688 ewcljrsq = _mm256_mul_ps(ewclj2,rsq00);
689 ewclj6 = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
690 exponent = avx256_exp_f(ewcljrsq);
691 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
692 poly = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
693 /* f6A = 6 * C6grid * (1 - poly) */
694 f6A = _mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly));
695 /* f6B = C6grid * exponent * beta^6 */
696 f6B = _mm256_mul_ps(_mm256_mul_ps(c6grid_00,one_sixth),_mm256_mul_ps(exponent,ewclj6));
697 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
698 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);
700 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
704 fscal = _mm256_and_ps(fscal,cutoff_mask);
706 /* Calculate temporary vectorial force */
707 tx = _mm256_mul_ps(fscal,dx00);
708 ty = _mm256_mul_ps(fscal,dy00);
709 tz = _mm256_mul_ps(fscal,dz00);
711 /* Update vectorial force */
712 fix0 = _mm256_add_ps(fix0,tx);
713 fiy0 = _mm256_add_ps(fiy0,ty);
714 fiz0 = _mm256_add_ps(fiz0,tz);
716 fjptrA = f+j_coord_offsetA;
717 fjptrB = f+j_coord_offsetB;
718 fjptrC = f+j_coord_offsetC;
719 fjptrD = f+j_coord_offsetD;
720 fjptrE = f+j_coord_offsetE;
721 fjptrF = f+j_coord_offsetF;
722 fjptrG = f+j_coord_offsetG;
723 fjptrH = f+j_coord_offsetH;
724 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
728 /* Inner loop uses 49 flops */
734 /* Get j neighbor index, and coordinate index */
735 jnrlistA = jjnr[jidx];
736 jnrlistB = jjnr[jidx+1];
737 jnrlistC = jjnr[jidx+2];
738 jnrlistD = jjnr[jidx+3];
739 jnrlistE = jjnr[jidx+4];
740 jnrlistF = jjnr[jidx+5];
741 jnrlistG = jjnr[jidx+6];
742 jnrlistH = jjnr[jidx+7];
743 /* Sign of each element will be negative for non-real atoms.
744 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
745 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
747 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
748 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
750 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
751 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
752 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
753 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
754 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
755 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
756 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
757 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
758 j_coord_offsetA = DIM*jnrA;
759 j_coord_offsetB = DIM*jnrB;
760 j_coord_offsetC = DIM*jnrC;
761 j_coord_offsetD = DIM*jnrD;
762 j_coord_offsetE = DIM*jnrE;
763 j_coord_offsetF = DIM*jnrF;
764 j_coord_offsetG = DIM*jnrG;
765 j_coord_offsetH = DIM*jnrH;
767 /* load j atom coordinates */
768 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
769 x+j_coord_offsetC,x+j_coord_offsetD,
770 x+j_coord_offsetE,x+j_coord_offsetF,
771 x+j_coord_offsetG,x+j_coord_offsetH,
774 /* Calculate displacement vector */
775 dx00 = _mm256_sub_ps(ix0,jx0);
776 dy00 = _mm256_sub_ps(iy0,jy0);
777 dz00 = _mm256_sub_ps(iz0,jz0);
779 /* Calculate squared distance and things based on it */
780 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
782 rinv00 = avx256_invsqrt_f(rsq00);
784 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
786 /* Load parameters for j particles */
787 vdwjidx0A = 2*vdwtype[jnrA+0];
788 vdwjidx0B = 2*vdwtype[jnrB+0];
789 vdwjidx0C = 2*vdwtype[jnrC+0];
790 vdwjidx0D = 2*vdwtype[jnrD+0];
791 vdwjidx0E = 2*vdwtype[jnrE+0];
792 vdwjidx0F = 2*vdwtype[jnrF+0];
793 vdwjidx0G = 2*vdwtype[jnrG+0];
794 vdwjidx0H = 2*vdwtype[jnrH+0];
796 /**************************
797 * CALCULATE INTERACTIONS *
798 **************************/
800 if (gmx_mm256_any_lt(rsq00,rcutoff2))
803 r00 = _mm256_mul_ps(rsq00,rinv00);
804 r00 = _mm256_andnot_ps(dummy_mask,r00);
806 /* Compute parameters for interactions between i and j atoms */
807 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
808 vdwioffsetptr0+vdwjidx0B,
809 vdwioffsetptr0+vdwjidx0C,
810 vdwioffsetptr0+vdwjidx0D,
811 vdwioffsetptr0+vdwjidx0E,
812 vdwioffsetptr0+vdwjidx0F,
813 vdwioffsetptr0+vdwjidx0G,
814 vdwioffsetptr0+vdwjidx0H,
817 c6grid_00 = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
818 vdwgridioffsetptr0+vdwjidx0B,
819 vdwgridioffsetptr0+vdwjidx0C,
820 vdwgridioffsetptr0+vdwjidx0D,
821 vdwgridioffsetptr0+vdwjidx0E,
822 vdwgridioffsetptr0+vdwjidx0F,
823 vdwgridioffsetptr0+vdwjidx0G,
824 vdwgridioffsetptr0+vdwjidx0H);
826 /* Analytical LJ-PME */
827 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
828 ewcljrsq = _mm256_mul_ps(ewclj2,rsq00);
829 ewclj6 = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
830 exponent = avx256_exp_f(ewcljrsq);
831 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
832 poly = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
833 /* f6A = 6 * C6grid * (1 - poly) */
834 f6A = _mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly));
835 /* f6B = C6grid * exponent * beta^6 */
836 f6B = _mm256_mul_ps(_mm256_mul_ps(c6grid_00,one_sixth),_mm256_mul_ps(exponent,ewclj6));
837 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
838 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);
840 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
844 fscal = _mm256_and_ps(fscal,cutoff_mask);
846 fscal = _mm256_andnot_ps(dummy_mask,fscal);
848 /* Calculate temporary vectorial force */
849 tx = _mm256_mul_ps(fscal,dx00);
850 ty = _mm256_mul_ps(fscal,dy00);
851 tz = _mm256_mul_ps(fscal,dz00);
853 /* Update vectorial force */
854 fix0 = _mm256_add_ps(fix0,tx);
855 fiy0 = _mm256_add_ps(fiy0,ty);
856 fiz0 = _mm256_add_ps(fiz0,tz);
858 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
859 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
860 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
861 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
862 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
863 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
864 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
865 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
866 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
870 /* Inner loop uses 50 flops */
873 /* End of innermost loop */
875 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
876 f+i_coord_offset,fshift+i_shift_offset);
878 /* Increment number of inner iterations */
879 inneriter += j_index_end - j_index_start;
881 /* Outer loop uses 6 flops */
884 /* Increment number of outer iterations */
887 /* Update outer/inner flops */
889 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_F,outeriter*6 + inneriter*50);