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_ElecEwSh_VdwLJSh_GeomP1P1_VF_avx_256_single
51 * Electrostatics interaction: Ewald
52 * VdW interaction: LennardJones
53 * Geometry: Particle-Particle
54 * Calculate force/pot: PotentialAndForce
57 nb_kernel_ElecEwSh_VdwLJSh_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 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
87 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
88 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
89 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
90 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
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);
99 __m128i ewitab_lo,ewitab_hi;
100 __m256 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
101 __m256 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
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 facel = _mm256_set1_ps(fr->ic->epsfac);
119 charge = mdatoms->chargeA;
120 nvdwtype = fr->ntype;
122 vdwtype = mdatoms->typeA;
124 sh_ewald = _mm256_set1_ps(fr->ic->sh_ewald);
125 beta = _mm256_set1_ps(fr->ic->ewaldcoeff_q);
126 beta2 = _mm256_mul_ps(beta,beta);
127 beta3 = _mm256_mul_ps(beta,beta2);
129 ewtab = fr->ic->tabq_coul_FDV0;
130 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
131 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
133 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
134 rcutoff_scalar = fr->ic->rcoulomb;
135 rcutoff = _mm256_set1_ps(rcutoff_scalar);
136 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
138 sh_vdw_invrcut6 = _mm256_set1_ps(fr->ic->sh_invrc6);
139 rvdw = _mm256_set1_ps(fr->ic->rvdw);
141 /* Avoid stupid compiler warnings */
142 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
155 for(iidx=0;iidx<4*DIM;iidx++)
160 /* Start outer loop over neighborlists */
161 for(iidx=0; iidx<nri; iidx++)
163 /* Load shift vector for this list */
164 i_shift_offset = DIM*shiftidx[iidx];
166 /* Load limits for loop over neighbors */
167 j_index_start = jindex[iidx];
168 j_index_end = jindex[iidx+1];
170 /* Get outer coordinate index */
172 i_coord_offset = DIM*inr;
174 /* Load i particle coords and add shift vector */
175 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
177 fix0 = _mm256_setzero_ps();
178 fiy0 = _mm256_setzero_ps();
179 fiz0 = _mm256_setzero_ps();
181 /* Load parameters for i particles */
182 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
183 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
185 /* Reset potential sums */
186 velecsum = _mm256_setzero_ps();
187 vvdwsum = _mm256_setzero_ps();
189 /* Start inner kernel loop */
190 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
193 /* Get j neighbor index, and coordinate index */
202 j_coord_offsetA = DIM*jnrA;
203 j_coord_offsetB = DIM*jnrB;
204 j_coord_offsetC = DIM*jnrC;
205 j_coord_offsetD = DIM*jnrD;
206 j_coord_offsetE = DIM*jnrE;
207 j_coord_offsetF = DIM*jnrF;
208 j_coord_offsetG = DIM*jnrG;
209 j_coord_offsetH = DIM*jnrH;
211 /* load j atom coordinates */
212 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
213 x+j_coord_offsetC,x+j_coord_offsetD,
214 x+j_coord_offsetE,x+j_coord_offsetF,
215 x+j_coord_offsetG,x+j_coord_offsetH,
218 /* Calculate displacement vector */
219 dx00 = _mm256_sub_ps(ix0,jx0);
220 dy00 = _mm256_sub_ps(iy0,jy0);
221 dz00 = _mm256_sub_ps(iz0,jz0);
223 /* Calculate squared distance and things based on it */
224 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
226 rinv00 = avx256_invsqrt_f(rsq00);
228 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
230 /* Load parameters for j particles */
231 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
232 charge+jnrC+0,charge+jnrD+0,
233 charge+jnrE+0,charge+jnrF+0,
234 charge+jnrG+0,charge+jnrH+0);
235 vdwjidx0A = 2*vdwtype[jnrA+0];
236 vdwjidx0B = 2*vdwtype[jnrB+0];
237 vdwjidx0C = 2*vdwtype[jnrC+0];
238 vdwjidx0D = 2*vdwtype[jnrD+0];
239 vdwjidx0E = 2*vdwtype[jnrE+0];
240 vdwjidx0F = 2*vdwtype[jnrF+0];
241 vdwjidx0G = 2*vdwtype[jnrG+0];
242 vdwjidx0H = 2*vdwtype[jnrH+0];
244 /**************************
245 * CALCULATE INTERACTIONS *
246 **************************/
248 if (gmx_mm256_any_lt(rsq00,rcutoff2))
251 r00 = _mm256_mul_ps(rsq00,rinv00);
253 /* Compute parameters for interactions between i and j atoms */
254 qq00 = _mm256_mul_ps(iq0,jq0);
255 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
256 vdwioffsetptr0+vdwjidx0B,
257 vdwioffsetptr0+vdwjidx0C,
258 vdwioffsetptr0+vdwjidx0D,
259 vdwioffsetptr0+vdwjidx0E,
260 vdwioffsetptr0+vdwjidx0F,
261 vdwioffsetptr0+vdwjidx0G,
262 vdwioffsetptr0+vdwjidx0H,
265 /* EWALD ELECTROSTATICS */
267 /* Analytical PME correction */
268 zeta2 = _mm256_mul_ps(beta2,rsq00);
269 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
270 pmecorrF = avx256_pmecorrF_f(zeta2);
271 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
272 felec = _mm256_mul_ps(qq00,felec);
273 pmecorrV = avx256_pmecorrV_f(zeta2);
274 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
275 velec = _mm256_sub_ps(_mm256_sub_ps(rinv00,sh_ewald),pmecorrV);
276 velec = _mm256_mul_ps(qq00,velec);
278 /* LENNARD-JONES DISPERSION/REPULSION */
280 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
281 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
282 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
283 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) ,
284 _mm256_mul_ps( _mm256_sub_ps(vvdw6,_mm256_mul_ps(c6_00,sh_vdw_invrcut6)),one_sixth));
285 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
287 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
289 /* Update potential sum for this i atom from the interaction with this j atom. */
290 velec = _mm256_and_ps(velec,cutoff_mask);
291 velecsum = _mm256_add_ps(velecsum,velec);
292 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
293 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
295 fscal = _mm256_add_ps(felec,fvdw);
297 fscal = _mm256_and_ps(fscal,cutoff_mask);
299 /* Calculate temporary vectorial force */
300 tx = _mm256_mul_ps(fscal,dx00);
301 ty = _mm256_mul_ps(fscal,dy00);
302 tz = _mm256_mul_ps(fscal,dz00);
304 /* Update vectorial force */
305 fix0 = _mm256_add_ps(fix0,tx);
306 fiy0 = _mm256_add_ps(fiy0,ty);
307 fiz0 = _mm256_add_ps(fiz0,tz);
309 fjptrA = f+j_coord_offsetA;
310 fjptrB = f+j_coord_offsetB;
311 fjptrC = f+j_coord_offsetC;
312 fjptrD = f+j_coord_offsetD;
313 fjptrE = f+j_coord_offsetE;
314 fjptrF = f+j_coord_offsetF;
315 fjptrG = f+j_coord_offsetG;
316 fjptrH = f+j_coord_offsetH;
317 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
321 /* Inner loop uses 127 flops */
327 /* Get j neighbor index, and coordinate index */
328 jnrlistA = jjnr[jidx];
329 jnrlistB = jjnr[jidx+1];
330 jnrlistC = jjnr[jidx+2];
331 jnrlistD = jjnr[jidx+3];
332 jnrlistE = jjnr[jidx+4];
333 jnrlistF = jjnr[jidx+5];
334 jnrlistG = jjnr[jidx+6];
335 jnrlistH = jjnr[jidx+7];
336 /* Sign of each element will be negative for non-real atoms.
337 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
338 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
340 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
341 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
343 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
344 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
345 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
346 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
347 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
348 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
349 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
350 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
351 j_coord_offsetA = DIM*jnrA;
352 j_coord_offsetB = DIM*jnrB;
353 j_coord_offsetC = DIM*jnrC;
354 j_coord_offsetD = DIM*jnrD;
355 j_coord_offsetE = DIM*jnrE;
356 j_coord_offsetF = DIM*jnrF;
357 j_coord_offsetG = DIM*jnrG;
358 j_coord_offsetH = DIM*jnrH;
360 /* load j atom coordinates */
361 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
362 x+j_coord_offsetC,x+j_coord_offsetD,
363 x+j_coord_offsetE,x+j_coord_offsetF,
364 x+j_coord_offsetG,x+j_coord_offsetH,
367 /* Calculate displacement vector */
368 dx00 = _mm256_sub_ps(ix0,jx0);
369 dy00 = _mm256_sub_ps(iy0,jy0);
370 dz00 = _mm256_sub_ps(iz0,jz0);
372 /* Calculate squared distance and things based on it */
373 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
375 rinv00 = avx256_invsqrt_f(rsq00);
377 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
379 /* Load parameters for j particles */
380 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
381 charge+jnrC+0,charge+jnrD+0,
382 charge+jnrE+0,charge+jnrF+0,
383 charge+jnrG+0,charge+jnrH+0);
384 vdwjidx0A = 2*vdwtype[jnrA+0];
385 vdwjidx0B = 2*vdwtype[jnrB+0];
386 vdwjidx0C = 2*vdwtype[jnrC+0];
387 vdwjidx0D = 2*vdwtype[jnrD+0];
388 vdwjidx0E = 2*vdwtype[jnrE+0];
389 vdwjidx0F = 2*vdwtype[jnrF+0];
390 vdwjidx0G = 2*vdwtype[jnrG+0];
391 vdwjidx0H = 2*vdwtype[jnrH+0];
393 /**************************
394 * CALCULATE INTERACTIONS *
395 **************************/
397 if (gmx_mm256_any_lt(rsq00,rcutoff2))
400 r00 = _mm256_mul_ps(rsq00,rinv00);
401 r00 = _mm256_andnot_ps(dummy_mask,r00);
403 /* Compute parameters for interactions between i and j atoms */
404 qq00 = _mm256_mul_ps(iq0,jq0);
405 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
406 vdwioffsetptr0+vdwjidx0B,
407 vdwioffsetptr0+vdwjidx0C,
408 vdwioffsetptr0+vdwjidx0D,
409 vdwioffsetptr0+vdwjidx0E,
410 vdwioffsetptr0+vdwjidx0F,
411 vdwioffsetptr0+vdwjidx0G,
412 vdwioffsetptr0+vdwjidx0H,
415 /* EWALD ELECTROSTATICS */
417 /* Analytical PME correction */
418 zeta2 = _mm256_mul_ps(beta2,rsq00);
419 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
420 pmecorrF = avx256_pmecorrF_f(zeta2);
421 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
422 felec = _mm256_mul_ps(qq00,felec);
423 pmecorrV = avx256_pmecorrV_f(zeta2);
424 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
425 velec = _mm256_sub_ps(_mm256_sub_ps(rinv00,sh_ewald),pmecorrV);
426 velec = _mm256_mul_ps(qq00,velec);
428 /* LENNARD-JONES DISPERSION/REPULSION */
430 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
431 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
432 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
433 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) ,
434 _mm256_mul_ps( _mm256_sub_ps(vvdw6,_mm256_mul_ps(c6_00,sh_vdw_invrcut6)),one_sixth));
435 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
437 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
439 /* Update potential sum for this i atom from the interaction with this j atom. */
440 velec = _mm256_and_ps(velec,cutoff_mask);
441 velec = _mm256_andnot_ps(dummy_mask,velec);
442 velecsum = _mm256_add_ps(velecsum,velec);
443 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
444 vvdw = _mm256_andnot_ps(dummy_mask,vvdw);
445 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
447 fscal = _mm256_add_ps(felec,fvdw);
449 fscal = _mm256_and_ps(fscal,cutoff_mask);
451 fscal = _mm256_andnot_ps(dummy_mask,fscal);
453 /* Calculate temporary vectorial force */
454 tx = _mm256_mul_ps(fscal,dx00);
455 ty = _mm256_mul_ps(fscal,dy00);
456 tz = _mm256_mul_ps(fscal,dz00);
458 /* Update vectorial force */
459 fix0 = _mm256_add_ps(fix0,tx);
460 fiy0 = _mm256_add_ps(fiy0,ty);
461 fiz0 = _mm256_add_ps(fiz0,tz);
463 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
464 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
465 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
466 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
467 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
468 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
469 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
470 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
471 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
475 /* Inner loop uses 128 flops */
478 /* End of innermost loop */
480 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
481 f+i_coord_offset,fshift+i_shift_offset);
484 /* Update potential energies */
485 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
486 gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
488 /* Increment number of inner iterations */
489 inneriter += j_index_end - j_index_start;
491 /* Outer loop uses 9 flops */
494 /* Increment number of outer iterations */
497 /* Update outer/inner flops */
499 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*128);
502 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJSh_GeomP1P1_F_avx_256_single
503 * Electrostatics interaction: Ewald
504 * VdW interaction: LennardJones
505 * Geometry: Particle-Particle
506 * Calculate force/pot: Force
509 nb_kernel_ElecEwSh_VdwLJSh_GeomP1P1_F_avx_256_single
510 (t_nblist * gmx_restrict nlist,
511 rvec * gmx_restrict xx,
512 rvec * gmx_restrict ff,
513 struct t_forcerec * gmx_restrict fr,
514 t_mdatoms * gmx_restrict mdatoms,
515 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
516 t_nrnb * gmx_restrict nrnb)
518 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
519 * just 0 for non-waters.
520 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
521 * jnr indices corresponding to data put in the four positions in the SIMD register.
523 int i_shift_offset,i_coord_offset,outeriter,inneriter;
524 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
525 int jnrA,jnrB,jnrC,jnrD;
526 int jnrE,jnrF,jnrG,jnrH;
527 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
528 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
529 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
530 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
531 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
533 real *shiftvec,*fshift,*x,*f;
534 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
536 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
537 real * vdwioffsetptr0;
538 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
539 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
540 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
541 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
542 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
545 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
548 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
549 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
551 __m128i ewitab_lo,ewitab_hi;
552 __m256 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
553 __m256 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
555 __m256 dummy_mask,cutoff_mask;
556 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
557 __m256 one = _mm256_set1_ps(1.0);
558 __m256 two = _mm256_set1_ps(2.0);
564 jindex = nlist->jindex;
566 shiftidx = nlist->shift;
568 shiftvec = fr->shift_vec[0];
569 fshift = fr->fshift[0];
570 facel = _mm256_set1_ps(fr->ic->epsfac);
571 charge = mdatoms->chargeA;
572 nvdwtype = fr->ntype;
574 vdwtype = mdatoms->typeA;
576 sh_ewald = _mm256_set1_ps(fr->ic->sh_ewald);
577 beta = _mm256_set1_ps(fr->ic->ewaldcoeff_q);
578 beta2 = _mm256_mul_ps(beta,beta);
579 beta3 = _mm256_mul_ps(beta,beta2);
581 ewtab = fr->ic->tabq_coul_F;
582 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
583 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
585 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
586 rcutoff_scalar = fr->ic->rcoulomb;
587 rcutoff = _mm256_set1_ps(rcutoff_scalar);
588 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
590 sh_vdw_invrcut6 = _mm256_set1_ps(fr->ic->sh_invrc6);
591 rvdw = _mm256_set1_ps(fr->ic->rvdw);
593 /* Avoid stupid compiler warnings */
594 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
607 for(iidx=0;iidx<4*DIM;iidx++)
612 /* Start outer loop over neighborlists */
613 for(iidx=0; iidx<nri; iidx++)
615 /* Load shift vector for this list */
616 i_shift_offset = DIM*shiftidx[iidx];
618 /* Load limits for loop over neighbors */
619 j_index_start = jindex[iidx];
620 j_index_end = jindex[iidx+1];
622 /* Get outer coordinate index */
624 i_coord_offset = DIM*inr;
626 /* Load i particle coords and add shift vector */
627 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
629 fix0 = _mm256_setzero_ps();
630 fiy0 = _mm256_setzero_ps();
631 fiz0 = _mm256_setzero_ps();
633 /* Load parameters for i particles */
634 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
635 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
637 /* Start inner kernel loop */
638 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
641 /* Get j neighbor index, and coordinate index */
650 j_coord_offsetA = DIM*jnrA;
651 j_coord_offsetB = DIM*jnrB;
652 j_coord_offsetC = DIM*jnrC;
653 j_coord_offsetD = DIM*jnrD;
654 j_coord_offsetE = DIM*jnrE;
655 j_coord_offsetF = DIM*jnrF;
656 j_coord_offsetG = DIM*jnrG;
657 j_coord_offsetH = DIM*jnrH;
659 /* load j atom coordinates */
660 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
661 x+j_coord_offsetC,x+j_coord_offsetD,
662 x+j_coord_offsetE,x+j_coord_offsetF,
663 x+j_coord_offsetG,x+j_coord_offsetH,
666 /* Calculate displacement vector */
667 dx00 = _mm256_sub_ps(ix0,jx0);
668 dy00 = _mm256_sub_ps(iy0,jy0);
669 dz00 = _mm256_sub_ps(iz0,jz0);
671 /* Calculate squared distance and things based on it */
672 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
674 rinv00 = avx256_invsqrt_f(rsq00);
676 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
678 /* Load parameters for j particles */
679 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
680 charge+jnrC+0,charge+jnrD+0,
681 charge+jnrE+0,charge+jnrF+0,
682 charge+jnrG+0,charge+jnrH+0);
683 vdwjidx0A = 2*vdwtype[jnrA+0];
684 vdwjidx0B = 2*vdwtype[jnrB+0];
685 vdwjidx0C = 2*vdwtype[jnrC+0];
686 vdwjidx0D = 2*vdwtype[jnrD+0];
687 vdwjidx0E = 2*vdwtype[jnrE+0];
688 vdwjidx0F = 2*vdwtype[jnrF+0];
689 vdwjidx0G = 2*vdwtype[jnrG+0];
690 vdwjidx0H = 2*vdwtype[jnrH+0];
692 /**************************
693 * CALCULATE INTERACTIONS *
694 **************************/
696 if (gmx_mm256_any_lt(rsq00,rcutoff2))
699 r00 = _mm256_mul_ps(rsq00,rinv00);
701 /* Compute parameters for interactions between i and j atoms */
702 qq00 = _mm256_mul_ps(iq0,jq0);
703 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
704 vdwioffsetptr0+vdwjidx0B,
705 vdwioffsetptr0+vdwjidx0C,
706 vdwioffsetptr0+vdwjidx0D,
707 vdwioffsetptr0+vdwjidx0E,
708 vdwioffsetptr0+vdwjidx0F,
709 vdwioffsetptr0+vdwjidx0G,
710 vdwioffsetptr0+vdwjidx0H,
713 /* EWALD ELECTROSTATICS */
715 /* Analytical PME correction */
716 zeta2 = _mm256_mul_ps(beta2,rsq00);
717 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
718 pmecorrF = avx256_pmecorrF_f(zeta2);
719 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
720 felec = _mm256_mul_ps(qq00,felec);
722 /* LENNARD-JONES DISPERSION/REPULSION */
724 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
725 fvdw = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_00,rinvsix),c6_00),_mm256_mul_ps(rinvsix,rinvsq00));
727 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
729 fscal = _mm256_add_ps(felec,fvdw);
731 fscal = _mm256_and_ps(fscal,cutoff_mask);
733 /* Calculate temporary vectorial force */
734 tx = _mm256_mul_ps(fscal,dx00);
735 ty = _mm256_mul_ps(fscal,dy00);
736 tz = _mm256_mul_ps(fscal,dz00);
738 /* Update vectorial force */
739 fix0 = _mm256_add_ps(fix0,tx);
740 fiy0 = _mm256_add_ps(fiy0,ty);
741 fiz0 = _mm256_add_ps(fiz0,tz);
743 fjptrA = f+j_coord_offsetA;
744 fjptrB = f+j_coord_offsetB;
745 fjptrC = f+j_coord_offsetC;
746 fjptrD = f+j_coord_offsetD;
747 fjptrE = f+j_coord_offsetE;
748 fjptrF = f+j_coord_offsetF;
749 fjptrG = f+j_coord_offsetG;
750 fjptrH = f+j_coord_offsetH;
751 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
755 /* Inner loop uses 66 flops */
761 /* Get j neighbor index, and coordinate index */
762 jnrlistA = jjnr[jidx];
763 jnrlistB = jjnr[jidx+1];
764 jnrlistC = jjnr[jidx+2];
765 jnrlistD = jjnr[jidx+3];
766 jnrlistE = jjnr[jidx+4];
767 jnrlistF = jjnr[jidx+5];
768 jnrlistG = jjnr[jidx+6];
769 jnrlistH = jjnr[jidx+7];
770 /* Sign of each element will be negative for non-real atoms.
771 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
772 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
774 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
775 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
777 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
778 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
779 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
780 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
781 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
782 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
783 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
784 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
785 j_coord_offsetA = DIM*jnrA;
786 j_coord_offsetB = DIM*jnrB;
787 j_coord_offsetC = DIM*jnrC;
788 j_coord_offsetD = DIM*jnrD;
789 j_coord_offsetE = DIM*jnrE;
790 j_coord_offsetF = DIM*jnrF;
791 j_coord_offsetG = DIM*jnrG;
792 j_coord_offsetH = DIM*jnrH;
794 /* load j atom coordinates */
795 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
796 x+j_coord_offsetC,x+j_coord_offsetD,
797 x+j_coord_offsetE,x+j_coord_offsetF,
798 x+j_coord_offsetG,x+j_coord_offsetH,
801 /* Calculate displacement vector */
802 dx00 = _mm256_sub_ps(ix0,jx0);
803 dy00 = _mm256_sub_ps(iy0,jy0);
804 dz00 = _mm256_sub_ps(iz0,jz0);
806 /* Calculate squared distance and things based on it */
807 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
809 rinv00 = avx256_invsqrt_f(rsq00);
811 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
813 /* Load parameters for j particles */
814 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
815 charge+jnrC+0,charge+jnrD+0,
816 charge+jnrE+0,charge+jnrF+0,
817 charge+jnrG+0,charge+jnrH+0);
818 vdwjidx0A = 2*vdwtype[jnrA+0];
819 vdwjidx0B = 2*vdwtype[jnrB+0];
820 vdwjidx0C = 2*vdwtype[jnrC+0];
821 vdwjidx0D = 2*vdwtype[jnrD+0];
822 vdwjidx0E = 2*vdwtype[jnrE+0];
823 vdwjidx0F = 2*vdwtype[jnrF+0];
824 vdwjidx0G = 2*vdwtype[jnrG+0];
825 vdwjidx0H = 2*vdwtype[jnrH+0];
827 /**************************
828 * CALCULATE INTERACTIONS *
829 **************************/
831 if (gmx_mm256_any_lt(rsq00,rcutoff2))
834 r00 = _mm256_mul_ps(rsq00,rinv00);
835 r00 = _mm256_andnot_ps(dummy_mask,r00);
837 /* Compute parameters for interactions between i and j atoms */
838 qq00 = _mm256_mul_ps(iq0,jq0);
839 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
840 vdwioffsetptr0+vdwjidx0B,
841 vdwioffsetptr0+vdwjidx0C,
842 vdwioffsetptr0+vdwjidx0D,
843 vdwioffsetptr0+vdwjidx0E,
844 vdwioffsetptr0+vdwjidx0F,
845 vdwioffsetptr0+vdwjidx0G,
846 vdwioffsetptr0+vdwjidx0H,
849 /* EWALD ELECTROSTATICS */
851 /* Analytical PME correction */
852 zeta2 = _mm256_mul_ps(beta2,rsq00);
853 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
854 pmecorrF = avx256_pmecorrF_f(zeta2);
855 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
856 felec = _mm256_mul_ps(qq00,felec);
858 /* LENNARD-JONES DISPERSION/REPULSION */
860 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
861 fvdw = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_00,rinvsix),c6_00),_mm256_mul_ps(rinvsix,rinvsq00));
863 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
865 fscal = _mm256_add_ps(felec,fvdw);
867 fscal = _mm256_and_ps(fscal,cutoff_mask);
869 fscal = _mm256_andnot_ps(dummy_mask,fscal);
871 /* Calculate temporary vectorial force */
872 tx = _mm256_mul_ps(fscal,dx00);
873 ty = _mm256_mul_ps(fscal,dy00);
874 tz = _mm256_mul_ps(fscal,dz00);
876 /* Update vectorial force */
877 fix0 = _mm256_add_ps(fix0,tx);
878 fiy0 = _mm256_add_ps(fiy0,ty);
879 fiz0 = _mm256_add_ps(fiz0,tz);
881 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
882 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
883 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
884 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
885 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
886 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
887 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
888 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
889 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
893 /* Inner loop uses 67 flops */
896 /* End of innermost loop */
898 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
899 f+i_coord_offset,fshift+i_shift_offset);
901 /* Increment number of inner iterations */
902 inneriter += j_index_end - j_index_start;
904 /* Outer loop uses 7 flops */
907 /* Increment number of outer iterations */
910 /* Update outer/inner flops */
912 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*67);