2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2017, 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_VdwLJEwSh_GeomP1P1_VF_avx_256_single
51 * Electrostatics interaction: Ewald
52 * VdW interaction: LJEwald
53 * Geometry: Particle-Particle
54 * Calculate force/pot: PotentialAndForce
57 nb_kernel_ElecEwSh_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;
91 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
94 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
97 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
98 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
101 __m256 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
102 __m256 one_half = _mm256_set1_ps(0.5);
103 __m256 minus_one = _mm256_set1_ps(-1.0);
105 __m128i ewitab_lo,ewitab_hi;
106 __m256 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
107 __m256 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
109 __m256 dummy_mask,cutoff_mask;
110 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
111 __m256 one = _mm256_set1_ps(1.0);
112 __m256 two = _mm256_set1_ps(2.0);
118 jindex = nlist->jindex;
120 shiftidx = nlist->shift;
122 shiftvec = fr->shift_vec[0];
123 fshift = fr->fshift[0];
124 facel = _mm256_set1_ps(fr->ic->epsfac);
125 charge = mdatoms->chargeA;
126 nvdwtype = fr->ntype;
128 vdwtype = mdatoms->typeA;
129 vdwgridparam = fr->ljpme_c6grid;
130 sh_lj_ewald = _mm256_set1_ps(fr->ic->sh_lj_ewald);
131 ewclj = _mm256_set1_ps(fr->ic->ewaldcoeff_lj);
132 ewclj2 = _mm256_mul_ps(minus_one,_mm256_mul_ps(ewclj,ewclj));
134 sh_ewald = _mm256_set1_ps(fr->ic->sh_ewald);
135 beta = _mm256_set1_ps(fr->ic->ewaldcoeff_q);
136 beta2 = _mm256_mul_ps(beta,beta);
137 beta3 = _mm256_mul_ps(beta,beta2);
139 ewtab = fr->ic->tabq_coul_FDV0;
140 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
141 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
143 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
144 rcutoff_scalar = fr->ic->rcoulomb;
145 rcutoff = _mm256_set1_ps(rcutoff_scalar);
146 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
148 sh_vdw_invrcut6 = _mm256_set1_ps(fr->ic->sh_invrc6);
149 rvdw = _mm256_set1_ps(fr->ic->rvdw);
151 /* Avoid stupid compiler warnings */
152 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
165 for(iidx=0;iidx<4*DIM;iidx++)
170 /* Start outer loop over neighborlists */
171 for(iidx=0; iidx<nri; iidx++)
173 /* Load shift vector for this list */
174 i_shift_offset = DIM*shiftidx[iidx];
176 /* Load limits for loop over neighbors */
177 j_index_start = jindex[iidx];
178 j_index_end = jindex[iidx+1];
180 /* Get outer coordinate index */
182 i_coord_offset = DIM*inr;
184 /* Load i particle coords and add shift vector */
185 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
187 fix0 = _mm256_setzero_ps();
188 fiy0 = _mm256_setzero_ps();
189 fiz0 = _mm256_setzero_ps();
191 /* Load parameters for i particles */
192 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
193 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
194 vdwgridioffsetptr0 = vdwgridparam+2*nvdwtype*vdwtype[inr+0];
196 /* Reset potential sums */
197 velecsum = _mm256_setzero_ps();
198 vvdwsum = _mm256_setzero_ps();
200 /* Start inner kernel loop */
201 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
204 /* Get j neighbor index, and coordinate index */
213 j_coord_offsetA = DIM*jnrA;
214 j_coord_offsetB = DIM*jnrB;
215 j_coord_offsetC = DIM*jnrC;
216 j_coord_offsetD = DIM*jnrD;
217 j_coord_offsetE = DIM*jnrE;
218 j_coord_offsetF = DIM*jnrF;
219 j_coord_offsetG = DIM*jnrG;
220 j_coord_offsetH = DIM*jnrH;
222 /* load j atom coordinates */
223 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
224 x+j_coord_offsetC,x+j_coord_offsetD,
225 x+j_coord_offsetE,x+j_coord_offsetF,
226 x+j_coord_offsetG,x+j_coord_offsetH,
229 /* Calculate displacement vector */
230 dx00 = _mm256_sub_ps(ix0,jx0);
231 dy00 = _mm256_sub_ps(iy0,jy0);
232 dz00 = _mm256_sub_ps(iz0,jz0);
234 /* Calculate squared distance and things based on it */
235 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
237 rinv00 = avx256_invsqrt_f(rsq00);
239 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
241 /* Load parameters for j particles */
242 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
243 charge+jnrC+0,charge+jnrD+0,
244 charge+jnrE+0,charge+jnrF+0,
245 charge+jnrG+0,charge+jnrH+0);
246 vdwjidx0A = 2*vdwtype[jnrA+0];
247 vdwjidx0B = 2*vdwtype[jnrB+0];
248 vdwjidx0C = 2*vdwtype[jnrC+0];
249 vdwjidx0D = 2*vdwtype[jnrD+0];
250 vdwjidx0E = 2*vdwtype[jnrE+0];
251 vdwjidx0F = 2*vdwtype[jnrF+0];
252 vdwjidx0G = 2*vdwtype[jnrG+0];
253 vdwjidx0H = 2*vdwtype[jnrH+0];
255 /**************************
256 * CALCULATE INTERACTIONS *
257 **************************/
259 if (gmx_mm256_any_lt(rsq00,rcutoff2))
262 r00 = _mm256_mul_ps(rsq00,rinv00);
264 /* Compute parameters for interactions between i and j atoms */
265 qq00 = _mm256_mul_ps(iq0,jq0);
266 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
267 vdwioffsetptr0+vdwjidx0B,
268 vdwioffsetptr0+vdwjidx0C,
269 vdwioffsetptr0+vdwjidx0D,
270 vdwioffsetptr0+vdwjidx0E,
271 vdwioffsetptr0+vdwjidx0F,
272 vdwioffsetptr0+vdwjidx0G,
273 vdwioffsetptr0+vdwjidx0H,
276 c6grid_00 = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
277 vdwgridioffsetptr0+vdwjidx0B,
278 vdwgridioffsetptr0+vdwjidx0C,
279 vdwgridioffsetptr0+vdwjidx0D,
280 vdwgridioffsetptr0+vdwjidx0E,
281 vdwgridioffsetptr0+vdwjidx0F,
282 vdwgridioffsetptr0+vdwjidx0G,
283 vdwgridioffsetptr0+vdwjidx0H);
285 /* EWALD ELECTROSTATICS */
287 /* Analytical PME correction */
288 zeta2 = _mm256_mul_ps(beta2,rsq00);
289 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
290 pmecorrF = avx256_pmecorrF_f(zeta2);
291 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
292 felec = _mm256_mul_ps(qq00,felec);
293 pmecorrV = avx256_pmecorrV_f(zeta2);
294 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
295 velec = _mm256_sub_ps(_mm256_sub_ps(rinv00,sh_ewald),pmecorrV);
296 velec = _mm256_mul_ps(qq00,velec);
298 /* Analytical LJ-PME */
299 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
300 ewcljrsq = _mm256_mul_ps(ewclj2,rsq00);
301 ewclj6 = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
302 exponent = avx256_exp_f(ewcljrsq);
303 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
304 poly = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
305 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
306 vvdw6 = _mm256_mul_ps(_mm256_sub_ps(c6_00,_mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly))),rinvsix);
307 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
308 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) ,
309 _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));
310 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
311 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);
313 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
315 /* Update potential sum for this i atom from the interaction with this j atom. */
316 velec = _mm256_and_ps(velec,cutoff_mask);
317 velecsum = _mm256_add_ps(velecsum,velec);
318 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
319 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
321 fscal = _mm256_add_ps(felec,fvdw);
323 fscal = _mm256_and_ps(fscal,cutoff_mask);
325 /* Calculate temporary vectorial force */
326 tx = _mm256_mul_ps(fscal,dx00);
327 ty = _mm256_mul_ps(fscal,dy00);
328 tz = _mm256_mul_ps(fscal,dz00);
330 /* Update vectorial force */
331 fix0 = _mm256_add_ps(fix0,tx);
332 fiy0 = _mm256_add_ps(fiy0,ty);
333 fiz0 = _mm256_add_ps(fiz0,tz);
335 fjptrA = f+j_coord_offsetA;
336 fjptrB = f+j_coord_offsetB;
337 fjptrC = f+j_coord_offsetC;
338 fjptrD = f+j_coord_offsetD;
339 fjptrE = f+j_coord_offsetE;
340 fjptrF = f+j_coord_offsetF;
341 fjptrG = f+j_coord_offsetG;
342 fjptrH = f+j_coord_offsetH;
343 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
347 /* Inner loop uses 145 flops */
353 /* Get j neighbor index, and coordinate index */
354 jnrlistA = jjnr[jidx];
355 jnrlistB = jjnr[jidx+1];
356 jnrlistC = jjnr[jidx+2];
357 jnrlistD = jjnr[jidx+3];
358 jnrlistE = jjnr[jidx+4];
359 jnrlistF = jjnr[jidx+5];
360 jnrlistG = jjnr[jidx+6];
361 jnrlistH = jjnr[jidx+7];
362 /* Sign of each element will be negative for non-real atoms.
363 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
364 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
366 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
367 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
369 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
370 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
371 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
372 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
373 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
374 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
375 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
376 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
377 j_coord_offsetA = DIM*jnrA;
378 j_coord_offsetB = DIM*jnrB;
379 j_coord_offsetC = DIM*jnrC;
380 j_coord_offsetD = DIM*jnrD;
381 j_coord_offsetE = DIM*jnrE;
382 j_coord_offsetF = DIM*jnrF;
383 j_coord_offsetG = DIM*jnrG;
384 j_coord_offsetH = DIM*jnrH;
386 /* load j atom coordinates */
387 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
388 x+j_coord_offsetC,x+j_coord_offsetD,
389 x+j_coord_offsetE,x+j_coord_offsetF,
390 x+j_coord_offsetG,x+j_coord_offsetH,
393 /* Calculate displacement vector */
394 dx00 = _mm256_sub_ps(ix0,jx0);
395 dy00 = _mm256_sub_ps(iy0,jy0);
396 dz00 = _mm256_sub_ps(iz0,jz0);
398 /* Calculate squared distance and things based on it */
399 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
401 rinv00 = avx256_invsqrt_f(rsq00);
403 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
405 /* Load parameters for j particles */
406 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
407 charge+jnrC+0,charge+jnrD+0,
408 charge+jnrE+0,charge+jnrF+0,
409 charge+jnrG+0,charge+jnrH+0);
410 vdwjidx0A = 2*vdwtype[jnrA+0];
411 vdwjidx0B = 2*vdwtype[jnrB+0];
412 vdwjidx0C = 2*vdwtype[jnrC+0];
413 vdwjidx0D = 2*vdwtype[jnrD+0];
414 vdwjidx0E = 2*vdwtype[jnrE+0];
415 vdwjidx0F = 2*vdwtype[jnrF+0];
416 vdwjidx0G = 2*vdwtype[jnrG+0];
417 vdwjidx0H = 2*vdwtype[jnrH+0];
419 /**************************
420 * CALCULATE INTERACTIONS *
421 **************************/
423 if (gmx_mm256_any_lt(rsq00,rcutoff2))
426 r00 = _mm256_mul_ps(rsq00,rinv00);
427 r00 = _mm256_andnot_ps(dummy_mask,r00);
429 /* Compute parameters for interactions between i and j atoms */
430 qq00 = _mm256_mul_ps(iq0,jq0);
431 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
432 vdwioffsetptr0+vdwjidx0B,
433 vdwioffsetptr0+vdwjidx0C,
434 vdwioffsetptr0+vdwjidx0D,
435 vdwioffsetptr0+vdwjidx0E,
436 vdwioffsetptr0+vdwjidx0F,
437 vdwioffsetptr0+vdwjidx0G,
438 vdwioffsetptr0+vdwjidx0H,
441 c6grid_00 = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
442 vdwgridioffsetptr0+vdwjidx0B,
443 vdwgridioffsetptr0+vdwjidx0C,
444 vdwgridioffsetptr0+vdwjidx0D,
445 vdwgridioffsetptr0+vdwjidx0E,
446 vdwgridioffsetptr0+vdwjidx0F,
447 vdwgridioffsetptr0+vdwjidx0G,
448 vdwgridioffsetptr0+vdwjidx0H);
450 /* EWALD ELECTROSTATICS */
452 /* Analytical PME correction */
453 zeta2 = _mm256_mul_ps(beta2,rsq00);
454 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
455 pmecorrF = avx256_pmecorrF_f(zeta2);
456 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
457 felec = _mm256_mul_ps(qq00,felec);
458 pmecorrV = avx256_pmecorrV_f(zeta2);
459 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
460 velec = _mm256_sub_ps(_mm256_sub_ps(rinv00,sh_ewald),pmecorrV);
461 velec = _mm256_mul_ps(qq00,velec);
463 /* Analytical LJ-PME */
464 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
465 ewcljrsq = _mm256_mul_ps(ewclj2,rsq00);
466 ewclj6 = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
467 exponent = avx256_exp_f(ewcljrsq);
468 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
469 poly = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
470 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
471 vvdw6 = _mm256_mul_ps(_mm256_sub_ps(c6_00,_mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly))),rinvsix);
472 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
473 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) ,
474 _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));
475 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
476 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);
478 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
480 /* Update potential sum for this i atom from the interaction with this j atom. */
481 velec = _mm256_and_ps(velec,cutoff_mask);
482 velec = _mm256_andnot_ps(dummy_mask,velec);
483 velecsum = _mm256_add_ps(velecsum,velec);
484 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
485 vvdw = _mm256_andnot_ps(dummy_mask,vvdw);
486 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
488 fscal = _mm256_add_ps(felec,fvdw);
490 fscal = _mm256_and_ps(fscal,cutoff_mask);
492 fscal = _mm256_andnot_ps(dummy_mask,fscal);
494 /* Calculate temporary vectorial force */
495 tx = _mm256_mul_ps(fscal,dx00);
496 ty = _mm256_mul_ps(fscal,dy00);
497 tz = _mm256_mul_ps(fscal,dz00);
499 /* Update vectorial force */
500 fix0 = _mm256_add_ps(fix0,tx);
501 fiy0 = _mm256_add_ps(fiy0,ty);
502 fiz0 = _mm256_add_ps(fiz0,tz);
504 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
505 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
506 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
507 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
508 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
509 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
510 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
511 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
512 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
516 /* Inner loop uses 146 flops */
519 /* End of innermost loop */
521 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
522 f+i_coord_offset,fshift+i_shift_offset);
525 /* Update potential energies */
526 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
527 gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
529 /* Increment number of inner iterations */
530 inneriter += j_index_end - j_index_start;
532 /* Outer loop uses 9 flops */
535 /* Increment number of outer iterations */
538 /* Update outer/inner flops */
540 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*146);
543 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_F_avx_256_single
544 * Electrostatics interaction: Ewald
545 * VdW interaction: LJEwald
546 * Geometry: Particle-Particle
547 * Calculate force/pot: Force
550 nb_kernel_ElecEwSh_VdwLJEwSh_GeomP1P1_F_avx_256_single
551 (t_nblist * gmx_restrict nlist,
552 rvec * gmx_restrict xx,
553 rvec * gmx_restrict ff,
554 struct t_forcerec * gmx_restrict fr,
555 t_mdatoms * gmx_restrict mdatoms,
556 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
557 t_nrnb * gmx_restrict nrnb)
559 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
560 * just 0 for non-waters.
561 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
562 * jnr indices corresponding to data put in the four positions in the SIMD register.
564 int i_shift_offset,i_coord_offset,outeriter,inneriter;
565 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
566 int jnrA,jnrB,jnrC,jnrD;
567 int jnrE,jnrF,jnrG,jnrH;
568 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
569 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
570 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
571 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
572 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
574 real *shiftvec,*fshift,*x,*f;
575 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
577 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
578 real * vdwioffsetptr0;
579 real * vdwgridioffsetptr0;
580 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
581 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
582 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
583 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
584 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
587 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
590 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
591 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
594 __m256 ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
595 __m256 one_half = _mm256_set1_ps(0.5);
596 __m256 minus_one = _mm256_set1_ps(-1.0);
598 __m128i ewitab_lo,ewitab_hi;
599 __m256 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
600 __m256 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
602 __m256 dummy_mask,cutoff_mask;
603 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
604 __m256 one = _mm256_set1_ps(1.0);
605 __m256 two = _mm256_set1_ps(2.0);
611 jindex = nlist->jindex;
613 shiftidx = nlist->shift;
615 shiftvec = fr->shift_vec[0];
616 fshift = fr->fshift[0];
617 facel = _mm256_set1_ps(fr->ic->epsfac);
618 charge = mdatoms->chargeA;
619 nvdwtype = fr->ntype;
621 vdwtype = mdatoms->typeA;
622 vdwgridparam = fr->ljpme_c6grid;
623 sh_lj_ewald = _mm256_set1_ps(fr->ic->sh_lj_ewald);
624 ewclj = _mm256_set1_ps(fr->ic->ewaldcoeff_lj);
625 ewclj2 = _mm256_mul_ps(minus_one,_mm256_mul_ps(ewclj,ewclj));
627 sh_ewald = _mm256_set1_ps(fr->ic->sh_ewald);
628 beta = _mm256_set1_ps(fr->ic->ewaldcoeff_q);
629 beta2 = _mm256_mul_ps(beta,beta);
630 beta3 = _mm256_mul_ps(beta,beta2);
632 ewtab = fr->ic->tabq_coul_F;
633 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
634 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
636 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
637 rcutoff_scalar = fr->ic->rcoulomb;
638 rcutoff = _mm256_set1_ps(rcutoff_scalar);
639 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
641 sh_vdw_invrcut6 = _mm256_set1_ps(fr->ic->sh_invrc6);
642 rvdw = _mm256_set1_ps(fr->ic->rvdw);
644 /* Avoid stupid compiler warnings */
645 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
658 for(iidx=0;iidx<4*DIM;iidx++)
663 /* Start outer loop over neighborlists */
664 for(iidx=0; iidx<nri; iidx++)
666 /* Load shift vector for this list */
667 i_shift_offset = DIM*shiftidx[iidx];
669 /* Load limits for loop over neighbors */
670 j_index_start = jindex[iidx];
671 j_index_end = jindex[iidx+1];
673 /* Get outer coordinate index */
675 i_coord_offset = DIM*inr;
677 /* Load i particle coords and add shift vector */
678 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
680 fix0 = _mm256_setzero_ps();
681 fiy0 = _mm256_setzero_ps();
682 fiz0 = _mm256_setzero_ps();
684 /* Load parameters for i particles */
685 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
686 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
687 vdwgridioffsetptr0 = vdwgridparam+2*nvdwtype*vdwtype[inr+0];
689 /* Start inner kernel loop */
690 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
693 /* Get j neighbor index, and coordinate index */
702 j_coord_offsetA = DIM*jnrA;
703 j_coord_offsetB = DIM*jnrB;
704 j_coord_offsetC = DIM*jnrC;
705 j_coord_offsetD = DIM*jnrD;
706 j_coord_offsetE = DIM*jnrE;
707 j_coord_offsetF = DIM*jnrF;
708 j_coord_offsetG = DIM*jnrG;
709 j_coord_offsetH = DIM*jnrH;
711 /* load j atom coordinates */
712 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
713 x+j_coord_offsetC,x+j_coord_offsetD,
714 x+j_coord_offsetE,x+j_coord_offsetF,
715 x+j_coord_offsetG,x+j_coord_offsetH,
718 /* Calculate displacement vector */
719 dx00 = _mm256_sub_ps(ix0,jx0);
720 dy00 = _mm256_sub_ps(iy0,jy0);
721 dz00 = _mm256_sub_ps(iz0,jz0);
723 /* Calculate squared distance and things based on it */
724 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
726 rinv00 = avx256_invsqrt_f(rsq00);
728 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
730 /* Load parameters for j particles */
731 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
732 charge+jnrC+0,charge+jnrD+0,
733 charge+jnrE+0,charge+jnrF+0,
734 charge+jnrG+0,charge+jnrH+0);
735 vdwjidx0A = 2*vdwtype[jnrA+0];
736 vdwjidx0B = 2*vdwtype[jnrB+0];
737 vdwjidx0C = 2*vdwtype[jnrC+0];
738 vdwjidx0D = 2*vdwtype[jnrD+0];
739 vdwjidx0E = 2*vdwtype[jnrE+0];
740 vdwjidx0F = 2*vdwtype[jnrF+0];
741 vdwjidx0G = 2*vdwtype[jnrG+0];
742 vdwjidx0H = 2*vdwtype[jnrH+0];
744 /**************************
745 * CALCULATE INTERACTIONS *
746 **************************/
748 if (gmx_mm256_any_lt(rsq00,rcutoff2))
751 r00 = _mm256_mul_ps(rsq00,rinv00);
753 /* Compute parameters for interactions between i and j atoms */
754 qq00 = _mm256_mul_ps(iq0,jq0);
755 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
756 vdwioffsetptr0+vdwjidx0B,
757 vdwioffsetptr0+vdwjidx0C,
758 vdwioffsetptr0+vdwjidx0D,
759 vdwioffsetptr0+vdwjidx0E,
760 vdwioffsetptr0+vdwjidx0F,
761 vdwioffsetptr0+vdwjidx0G,
762 vdwioffsetptr0+vdwjidx0H,
765 c6grid_00 = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
766 vdwgridioffsetptr0+vdwjidx0B,
767 vdwgridioffsetptr0+vdwjidx0C,
768 vdwgridioffsetptr0+vdwjidx0D,
769 vdwgridioffsetptr0+vdwjidx0E,
770 vdwgridioffsetptr0+vdwjidx0F,
771 vdwgridioffsetptr0+vdwjidx0G,
772 vdwgridioffsetptr0+vdwjidx0H);
774 /* EWALD ELECTROSTATICS */
776 /* Analytical PME correction */
777 zeta2 = _mm256_mul_ps(beta2,rsq00);
778 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
779 pmecorrF = avx256_pmecorrF_f(zeta2);
780 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
781 felec = _mm256_mul_ps(qq00,felec);
783 /* Analytical LJ-PME */
784 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
785 ewcljrsq = _mm256_mul_ps(ewclj2,rsq00);
786 ewclj6 = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
787 exponent = avx256_exp_f(ewcljrsq);
788 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
789 poly = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
790 /* f6A = 6 * C6grid * (1 - poly) */
791 f6A = _mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly));
792 /* f6B = C6grid * exponent * beta^6 */
793 f6B = _mm256_mul_ps(_mm256_mul_ps(c6grid_00,one_sixth),_mm256_mul_ps(exponent,ewclj6));
794 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
795 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);
797 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
799 fscal = _mm256_add_ps(felec,fvdw);
801 fscal = _mm256_and_ps(fscal,cutoff_mask);
803 /* Calculate temporary vectorial force */
804 tx = _mm256_mul_ps(fscal,dx00);
805 ty = _mm256_mul_ps(fscal,dy00);
806 tz = _mm256_mul_ps(fscal,dz00);
808 /* Update vectorial force */
809 fix0 = _mm256_add_ps(fix0,tx);
810 fiy0 = _mm256_add_ps(fiy0,ty);
811 fiz0 = _mm256_add_ps(fiz0,tz);
813 fjptrA = f+j_coord_offsetA;
814 fjptrB = f+j_coord_offsetB;
815 fjptrC = f+j_coord_offsetC;
816 fjptrD = f+j_coord_offsetD;
817 fjptrE = f+j_coord_offsetE;
818 fjptrF = f+j_coord_offsetF;
819 fjptrG = f+j_coord_offsetG;
820 fjptrH = f+j_coord_offsetH;
821 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
825 /* Inner loop uses 82 flops */
831 /* Get j neighbor index, and coordinate index */
832 jnrlistA = jjnr[jidx];
833 jnrlistB = jjnr[jidx+1];
834 jnrlistC = jjnr[jidx+2];
835 jnrlistD = jjnr[jidx+3];
836 jnrlistE = jjnr[jidx+4];
837 jnrlistF = jjnr[jidx+5];
838 jnrlistG = jjnr[jidx+6];
839 jnrlistH = jjnr[jidx+7];
840 /* Sign of each element will be negative for non-real atoms.
841 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
842 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
844 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
845 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
847 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
848 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
849 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
850 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
851 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
852 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
853 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
854 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
855 j_coord_offsetA = DIM*jnrA;
856 j_coord_offsetB = DIM*jnrB;
857 j_coord_offsetC = DIM*jnrC;
858 j_coord_offsetD = DIM*jnrD;
859 j_coord_offsetE = DIM*jnrE;
860 j_coord_offsetF = DIM*jnrF;
861 j_coord_offsetG = DIM*jnrG;
862 j_coord_offsetH = DIM*jnrH;
864 /* load j atom coordinates */
865 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
866 x+j_coord_offsetC,x+j_coord_offsetD,
867 x+j_coord_offsetE,x+j_coord_offsetF,
868 x+j_coord_offsetG,x+j_coord_offsetH,
871 /* Calculate displacement vector */
872 dx00 = _mm256_sub_ps(ix0,jx0);
873 dy00 = _mm256_sub_ps(iy0,jy0);
874 dz00 = _mm256_sub_ps(iz0,jz0);
876 /* Calculate squared distance and things based on it */
877 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
879 rinv00 = avx256_invsqrt_f(rsq00);
881 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
883 /* Load parameters for j particles */
884 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
885 charge+jnrC+0,charge+jnrD+0,
886 charge+jnrE+0,charge+jnrF+0,
887 charge+jnrG+0,charge+jnrH+0);
888 vdwjidx0A = 2*vdwtype[jnrA+0];
889 vdwjidx0B = 2*vdwtype[jnrB+0];
890 vdwjidx0C = 2*vdwtype[jnrC+0];
891 vdwjidx0D = 2*vdwtype[jnrD+0];
892 vdwjidx0E = 2*vdwtype[jnrE+0];
893 vdwjidx0F = 2*vdwtype[jnrF+0];
894 vdwjidx0G = 2*vdwtype[jnrG+0];
895 vdwjidx0H = 2*vdwtype[jnrH+0];
897 /**************************
898 * CALCULATE INTERACTIONS *
899 **************************/
901 if (gmx_mm256_any_lt(rsq00,rcutoff2))
904 r00 = _mm256_mul_ps(rsq00,rinv00);
905 r00 = _mm256_andnot_ps(dummy_mask,r00);
907 /* Compute parameters for interactions between i and j atoms */
908 qq00 = _mm256_mul_ps(iq0,jq0);
909 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
910 vdwioffsetptr0+vdwjidx0B,
911 vdwioffsetptr0+vdwjidx0C,
912 vdwioffsetptr0+vdwjidx0D,
913 vdwioffsetptr0+vdwjidx0E,
914 vdwioffsetptr0+vdwjidx0F,
915 vdwioffsetptr0+vdwjidx0G,
916 vdwioffsetptr0+vdwjidx0H,
919 c6grid_00 = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
920 vdwgridioffsetptr0+vdwjidx0B,
921 vdwgridioffsetptr0+vdwjidx0C,
922 vdwgridioffsetptr0+vdwjidx0D,
923 vdwgridioffsetptr0+vdwjidx0E,
924 vdwgridioffsetptr0+vdwjidx0F,
925 vdwgridioffsetptr0+vdwjidx0G,
926 vdwgridioffsetptr0+vdwjidx0H);
928 /* EWALD ELECTROSTATICS */
930 /* Analytical PME correction */
931 zeta2 = _mm256_mul_ps(beta2,rsq00);
932 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
933 pmecorrF = avx256_pmecorrF_f(zeta2);
934 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
935 felec = _mm256_mul_ps(qq00,felec);
937 /* Analytical LJ-PME */
938 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
939 ewcljrsq = _mm256_mul_ps(ewclj2,rsq00);
940 ewclj6 = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
941 exponent = avx256_exp_f(ewcljrsq);
942 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
943 poly = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
944 /* f6A = 6 * C6grid * (1 - poly) */
945 f6A = _mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly));
946 /* f6B = C6grid * exponent * beta^6 */
947 f6B = _mm256_mul_ps(_mm256_mul_ps(c6grid_00,one_sixth),_mm256_mul_ps(exponent,ewclj6));
948 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
949 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);
951 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
953 fscal = _mm256_add_ps(felec,fvdw);
955 fscal = _mm256_and_ps(fscal,cutoff_mask);
957 fscal = _mm256_andnot_ps(dummy_mask,fscal);
959 /* Calculate temporary vectorial force */
960 tx = _mm256_mul_ps(fscal,dx00);
961 ty = _mm256_mul_ps(fscal,dy00);
962 tz = _mm256_mul_ps(fscal,dz00);
964 /* Update vectorial force */
965 fix0 = _mm256_add_ps(fix0,tx);
966 fiy0 = _mm256_add_ps(fiy0,ty);
967 fiz0 = _mm256_add_ps(fiz0,tz);
969 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
970 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
971 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
972 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
973 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
974 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
975 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
976 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
977 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
981 /* Inner loop uses 83 flops */
984 /* End of innermost loop */
986 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
987 f+i_coord_offset,fshift+i_shift_offset);
989 /* Increment number of inner iterations */
990 inneriter += j_index_end - j_index_start;
992 /* Outer loop uses 7 flops */
995 /* Increment number of outer iterations */
998 /* Update outer/inner flops */
1000 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*83);