2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 * Note: this file was generated by the GROMACS avx_256_single kernel generator.
44 #include "../nb_kernel.h"
45 #include "gromacs/legacyheaders/types/simple.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/legacyheaders/nrnb.h"
49 #include "gromacs/simd/math_x86_avx_256_single.h"
50 #include "kernelutil_x86_avx_256_single.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecEwSw_VdwLJSw_GeomP1P1_VF_avx_256_single
54 * Electrostatics interaction: Ewald
55 * VdW interaction: LennardJones
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecEwSw_VdwLJSw_GeomP1P1_VF_avx_256_single
61 (t_nblist * gmx_restrict nlist,
62 rvec * gmx_restrict xx,
63 rvec * gmx_restrict ff,
64 t_forcerec * gmx_restrict fr,
65 t_mdatoms * gmx_restrict mdatoms,
66 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67 t_nrnb * gmx_restrict nrnb)
69 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
70 * just 0 for non-waters.
71 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
72 * jnr indices corresponding to data put in the four positions in the SIMD register.
74 int i_shift_offset,i_coord_offset,outeriter,inneriter;
75 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
76 int jnrA,jnrB,jnrC,jnrD;
77 int jnrE,jnrF,jnrG,jnrH;
78 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
79 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
80 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
81 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
82 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
84 real *shiftvec,*fshift,*x,*f;
85 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
87 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
88 real * vdwioffsetptr0;
89 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
90 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
91 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
92 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
93 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
96 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
99 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
100 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
102 __m128i ewitab_lo,ewitab_hi;
103 __m256 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
104 __m256 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
106 __m256 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
107 real rswitch_scalar,d_scalar;
108 __m256 dummy_mask,cutoff_mask;
109 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
110 __m256 one = _mm256_set1_ps(1.0);
111 __m256 two = _mm256_set1_ps(2.0);
117 jindex = nlist->jindex;
119 shiftidx = nlist->shift;
121 shiftvec = fr->shift_vec[0];
122 fshift = fr->fshift[0];
123 facel = _mm256_set1_ps(fr->epsfac);
124 charge = mdatoms->chargeA;
125 nvdwtype = fr->ntype;
127 vdwtype = mdatoms->typeA;
129 sh_ewald = _mm256_set1_ps(fr->ic->sh_ewald);
130 beta = _mm256_set1_ps(fr->ic->ewaldcoeff_q);
131 beta2 = _mm256_mul_ps(beta,beta);
132 beta3 = _mm256_mul_ps(beta,beta2);
134 ewtab = fr->ic->tabq_coul_FDV0;
135 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
136 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
138 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
139 rcutoff_scalar = fr->rcoulomb;
140 rcutoff = _mm256_set1_ps(rcutoff_scalar);
141 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
143 rswitch_scalar = fr->rcoulomb_switch;
144 rswitch = _mm256_set1_ps(rswitch_scalar);
145 /* Setup switch parameters */
146 d_scalar = rcutoff_scalar-rswitch_scalar;
147 d = _mm256_set1_ps(d_scalar);
148 swV3 = _mm256_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
149 swV4 = _mm256_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
150 swV5 = _mm256_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
151 swF2 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
152 swF3 = _mm256_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
153 swF4 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
155 /* Avoid stupid compiler warnings */
156 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
169 for(iidx=0;iidx<4*DIM;iidx++)
174 /* Start outer loop over neighborlists */
175 for(iidx=0; iidx<nri; iidx++)
177 /* Load shift vector for this list */
178 i_shift_offset = DIM*shiftidx[iidx];
180 /* Load limits for loop over neighbors */
181 j_index_start = jindex[iidx];
182 j_index_end = jindex[iidx+1];
184 /* Get outer coordinate index */
186 i_coord_offset = DIM*inr;
188 /* Load i particle coords and add shift vector */
189 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
191 fix0 = _mm256_setzero_ps();
192 fiy0 = _mm256_setzero_ps();
193 fiz0 = _mm256_setzero_ps();
195 /* Load parameters for i particles */
196 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
197 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
199 /* Reset potential sums */
200 velecsum = _mm256_setzero_ps();
201 vvdwsum = _mm256_setzero_ps();
203 /* Start inner kernel loop */
204 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
207 /* Get j neighbor index, and coordinate index */
216 j_coord_offsetA = DIM*jnrA;
217 j_coord_offsetB = DIM*jnrB;
218 j_coord_offsetC = DIM*jnrC;
219 j_coord_offsetD = DIM*jnrD;
220 j_coord_offsetE = DIM*jnrE;
221 j_coord_offsetF = DIM*jnrF;
222 j_coord_offsetG = DIM*jnrG;
223 j_coord_offsetH = DIM*jnrH;
225 /* load j atom coordinates */
226 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
227 x+j_coord_offsetC,x+j_coord_offsetD,
228 x+j_coord_offsetE,x+j_coord_offsetF,
229 x+j_coord_offsetG,x+j_coord_offsetH,
232 /* Calculate displacement vector */
233 dx00 = _mm256_sub_ps(ix0,jx0);
234 dy00 = _mm256_sub_ps(iy0,jy0);
235 dz00 = _mm256_sub_ps(iz0,jz0);
237 /* Calculate squared distance and things based on it */
238 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
240 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
242 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
244 /* Load parameters for j particles */
245 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
246 charge+jnrC+0,charge+jnrD+0,
247 charge+jnrE+0,charge+jnrF+0,
248 charge+jnrG+0,charge+jnrH+0);
249 vdwjidx0A = 2*vdwtype[jnrA+0];
250 vdwjidx0B = 2*vdwtype[jnrB+0];
251 vdwjidx0C = 2*vdwtype[jnrC+0];
252 vdwjidx0D = 2*vdwtype[jnrD+0];
253 vdwjidx0E = 2*vdwtype[jnrE+0];
254 vdwjidx0F = 2*vdwtype[jnrF+0];
255 vdwjidx0G = 2*vdwtype[jnrG+0];
256 vdwjidx0H = 2*vdwtype[jnrH+0];
258 /**************************
259 * CALCULATE INTERACTIONS *
260 **************************/
262 if (gmx_mm256_any_lt(rsq00,rcutoff2))
265 r00 = _mm256_mul_ps(rsq00,rinv00);
267 /* Compute parameters for interactions between i and j atoms */
268 qq00 = _mm256_mul_ps(iq0,jq0);
269 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
270 vdwioffsetptr0+vdwjidx0B,
271 vdwioffsetptr0+vdwjidx0C,
272 vdwioffsetptr0+vdwjidx0D,
273 vdwioffsetptr0+vdwjidx0E,
274 vdwioffsetptr0+vdwjidx0F,
275 vdwioffsetptr0+vdwjidx0G,
276 vdwioffsetptr0+vdwjidx0H,
279 /* EWALD ELECTROSTATICS */
281 /* Analytical PME correction */
282 zeta2 = _mm256_mul_ps(beta2,rsq00);
283 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
284 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
285 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
286 felec = _mm256_mul_ps(qq00,felec);
287 pmecorrV = gmx_mm256_pmecorrV_ps(zeta2);
288 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
289 velec = _mm256_sub_ps(rinv00,pmecorrV);
290 velec = _mm256_mul_ps(qq00,velec);
292 /* LENNARD-JONES DISPERSION/REPULSION */
294 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
295 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
296 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
297 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
298 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
300 d = _mm256_sub_ps(r00,rswitch);
301 d = _mm256_max_ps(d,_mm256_setzero_ps());
302 d2 = _mm256_mul_ps(d,d);
303 sw = _mm256_add_ps(one,_mm256_mul_ps(d2,_mm256_mul_ps(d,_mm256_add_ps(swV3,_mm256_mul_ps(d,_mm256_add_ps(swV4,_mm256_mul_ps(d,swV5)))))));
305 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
307 /* Evaluate switch function */
308 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
309 felec = _mm256_sub_ps( _mm256_mul_ps(felec,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(velec,dsw)) );
310 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
311 velec = _mm256_mul_ps(velec,sw);
312 vvdw = _mm256_mul_ps(vvdw,sw);
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 126 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 = gmx_mm256_invsqrt_ps(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 /* EWALD ELECTROSTATICS */
443 /* Analytical PME correction */
444 zeta2 = _mm256_mul_ps(beta2,rsq00);
445 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
446 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
447 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
448 felec = _mm256_mul_ps(qq00,felec);
449 pmecorrV = gmx_mm256_pmecorrV_ps(zeta2);
450 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
451 velec = _mm256_sub_ps(rinv00,pmecorrV);
452 velec = _mm256_mul_ps(qq00,velec);
454 /* LENNARD-JONES DISPERSION/REPULSION */
456 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
457 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
458 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
459 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
460 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
462 d = _mm256_sub_ps(r00,rswitch);
463 d = _mm256_max_ps(d,_mm256_setzero_ps());
464 d2 = _mm256_mul_ps(d,d);
465 sw = _mm256_add_ps(one,_mm256_mul_ps(d2,_mm256_mul_ps(d,_mm256_add_ps(swV3,_mm256_mul_ps(d,_mm256_add_ps(swV4,_mm256_mul_ps(d,swV5)))))));
467 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
469 /* Evaluate switch function */
470 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
471 felec = _mm256_sub_ps( _mm256_mul_ps(felec,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(velec,dsw)) );
472 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
473 velec = _mm256_mul_ps(velec,sw);
474 vvdw = _mm256_mul_ps(vvdw,sw);
475 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
477 /* Update potential sum for this i atom from the interaction with this j atom. */
478 velec = _mm256_and_ps(velec,cutoff_mask);
479 velec = _mm256_andnot_ps(dummy_mask,velec);
480 velecsum = _mm256_add_ps(velecsum,velec);
481 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
482 vvdw = _mm256_andnot_ps(dummy_mask,vvdw);
483 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
485 fscal = _mm256_add_ps(felec,fvdw);
487 fscal = _mm256_and_ps(fscal,cutoff_mask);
489 fscal = _mm256_andnot_ps(dummy_mask,fscal);
491 /* Calculate temporary vectorial force */
492 tx = _mm256_mul_ps(fscal,dx00);
493 ty = _mm256_mul_ps(fscal,dy00);
494 tz = _mm256_mul_ps(fscal,dz00);
496 /* Update vectorial force */
497 fix0 = _mm256_add_ps(fix0,tx);
498 fiy0 = _mm256_add_ps(fiy0,ty);
499 fiz0 = _mm256_add_ps(fiz0,tz);
501 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
502 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
503 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
504 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
505 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
506 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
507 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
508 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
509 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
513 /* Inner loop uses 127 flops */
516 /* End of innermost loop */
518 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
519 f+i_coord_offset,fshift+i_shift_offset);
522 /* Update potential energies */
523 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
524 gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
526 /* Increment number of inner iterations */
527 inneriter += j_index_end - j_index_start;
529 /* Outer loop uses 9 flops */
532 /* Increment number of outer iterations */
535 /* Update outer/inner flops */
537 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*127);
540 * Gromacs nonbonded kernel: nb_kernel_ElecEwSw_VdwLJSw_GeomP1P1_F_avx_256_single
541 * Electrostatics interaction: Ewald
542 * VdW interaction: LennardJones
543 * Geometry: Particle-Particle
544 * Calculate force/pot: Force
547 nb_kernel_ElecEwSw_VdwLJSw_GeomP1P1_F_avx_256_single
548 (t_nblist * gmx_restrict nlist,
549 rvec * gmx_restrict xx,
550 rvec * gmx_restrict ff,
551 t_forcerec * gmx_restrict fr,
552 t_mdatoms * gmx_restrict mdatoms,
553 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
554 t_nrnb * gmx_restrict nrnb)
556 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
557 * just 0 for non-waters.
558 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
559 * jnr indices corresponding to data put in the four positions in the SIMD register.
561 int i_shift_offset,i_coord_offset,outeriter,inneriter;
562 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
563 int jnrA,jnrB,jnrC,jnrD;
564 int jnrE,jnrF,jnrG,jnrH;
565 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
566 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
567 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
568 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
569 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
571 real *shiftvec,*fshift,*x,*f;
572 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
574 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
575 real * vdwioffsetptr0;
576 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
577 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
578 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
579 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
580 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
583 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
586 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
587 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
589 __m128i ewitab_lo,ewitab_hi;
590 __m256 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
591 __m256 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
593 __m256 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
594 real rswitch_scalar,d_scalar;
595 __m256 dummy_mask,cutoff_mask;
596 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
597 __m256 one = _mm256_set1_ps(1.0);
598 __m256 two = _mm256_set1_ps(2.0);
604 jindex = nlist->jindex;
606 shiftidx = nlist->shift;
608 shiftvec = fr->shift_vec[0];
609 fshift = fr->fshift[0];
610 facel = _mm256_set1_ps(fr->epsfac);
611 charge = mdatoms->chargeA;
612 nvdwtype = fr->ntype;
614 vdwtype = mdatoms->typeA;
616 sh_ewald = _mm256_set1_ps(fr->ic->sh_ewald);
617 beta = _mm256_set1_ps(fr->ic->ewaldcoeff_q);
618 beta2 = _mm256_mul_ps(beta,beta);
619 beta3 = _mm256_mul_ps(beta,beta2);
621 ewtab = fr->ic->tabq_coul_FDV0;
622 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
623 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
625 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
626 rcutoff_scalar = fr->rcoulomb;
627 rcutoff = _mm256_set1_ps(rcutoff_scalar);
628 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
630 rswitch_scalar = fr->rcoulomb_switch;
631 rswitch = _mm256_set1_ps(rswitch_scalar);
632 /* Setup switch parameters */
633 d_scalar = rcutoff_scalar-rswitch_scalar;
634 d = _mm256_set1_ps(d_scalar);
635 swV3 = _mm256_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
636 swV4 = _mm256_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
637 swV5 = _mm256_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
638 swF2 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
639 swF3 = _mm256_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
640 swF4 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
642 /* Avoid stupid compiler warnings */
643 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
656 for(iidx=0;iidx<4*DIM;iidx++)
661 /* Start outer loop over neighborlists */
662 for(iidx=0; iidx<nri; iidx++)
664 /* Load shift vector for this list */
665 i_shift_offset = DIM*shiftidx[iidx];
667 /* Load limits for loop over neighbors */
668 j_index_start = jindex[iidx];
669 j_index_end = jindex[iidx+1];
671 /* Get outer coordinate index */
673 i_coord_offset = DIM*inr;
675 /* Load i particle coords and add shift vector */
676 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
678 fix0 = _mm256_setzero_ps();
679 fiy0 = _mm256_setzero_ps();
680 fiz0 = _mm256_setzero_ps();
682 /* Load parameters for i particles */
683 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
684 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
686 /* Start inner kernel loop */
687 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
690 /* Get j neighbor index, and coordinate index */
699 j_coord_offsetA = DIM*jnrA;
700 j_coord_offsetB = DIM*jnrB;
701 j_coord_offsetC = DIM*jnrC;
702 j_coord_offsetD = DIM*jnrD;
703 j_coord_offsetE = DIM*jnrE;
704 j_coord_offsetF = DIM*jnrF;
705 j_coord_offsetG = DIM*jnrG;
706 j_coord_offsetH = DIM*jnrH;
708 /* load j atom coordinates */
709 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
710 x+j_coord_offsetC,x+j_coord_offsetD,
711 x+j_coord_offsetE,x+j_coord_offsetF,
712 x+j_coord_offsetG,x+j_coord_offsetH,
715 /* Calculate displacement vector */
716 dx00 = _mm256_sub_ps(ix0,jx0);
717 dy00 = _mm256_sub_ps(iy0,jy0);
718 dz00 = _mm256_sub_ps(iz0,jz0);
720 /* Calculate squared distance and things based on it */
721 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
723 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
725 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
727 /* Load parameters for j particles */
728 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
729 charge+jnrC+0,charge+jnrD+0,
730 charge+jnrE+0,charge+jnrF+0,
731 charge+jnrG+0,charge+jnrH+0);
732 vdwjidx0A = 2*vdwtype[jnrA+0];
733 vdwjidx0B = 2*vdwtype[jnrB+0];
734 vdwjidx0C = 2*vdwtype[jnrC+0];
735 vdwjidx0D = 2*vdwtype[jnrD+0];
736 vdwjidx0E = 2*vdwtype[jnrE+0];
737 vdwjidx0F = 2*vdwtype[jnrF+0];
738 vdwjidx0G = 2*vdwtype[jnrG+0];
739 vdwjidx0H = 2*vdwtype[jnrH+0];
741 /**************************
742 * CALCULATE INTERACTIONS *
743 **************************/
745 if (gmx_mm256_any_lt(rsq00,rcutoff2))
748 r00 = _mm256_mul_ps(rsq00,rinv00);
750 /* Compute parameters for interactions between i and j atoms */
751 qq00 = _mm256_mul_ps(iq0,jq0);
752 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
753 vdwioffsetptr0+vdwjidx0B,
754 vdwioffsetptr0+vdwjidx0C,
755 vdwioffsetptr0+vdwjidx0D,
756 vdwioffsetptr0+vdwjidx0E,
757 vdwioffsetptr0+vdwjidx0F,
758 vdwioffsetptr0+vdwjidx0G,
759 vdwioffsetptr0+vdwjidx0H,
762 /* EWALD ELECTROSTATICS */
764 /* Analytical PME correction */
765 zeta2 = _mm256_mul_ps(beta2,rsq00);
766 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
767 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
768 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
769 felec = _mm256_mul_ps(qq00,felec);
770 pmecorrV = gmx_mm256_pmecorrV_ps(zeta2);
771 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
772 velec = _mm256_sub_ps(rinv00,pmecorrV);
773 velec = _mm256_mul_ps(qq00,velec);
775 /* LENNARD-JONES DISPERSION/REPULSION */
777 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
778 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
779 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
780 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
781 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
783 d = _mm256_sub_ps(r00,rswitch);
784 d = _mm256_max_ps(d,_mm256_setzero_ps());
785 d2 = _mm256_mul_ps(d,d);
786 sw = _mm256_add_ps(one,_mm256_mul_ps(d2,_mm256_mul_ps(d,_mm256_add_ps(swV3,_mm256_mul_ps(d,_mm256_add_ps(swV4,_mm256_mul_ps(d,swV5)))))));
788 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
790 /* Evaluate switch function */
791 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
792 felec = _mm256_sub_ps( _mm256_mul_ps(felec,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(velec,dsw)) );
793 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
794 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
796 fscal = _mm256_add_ps(felec,fvdw);
798 fscal = _mm256_and_ps(fscal,cutoff_mask);
800 /* Calculate temporary vectorial force */
801 tx = _mm256_mul_ps(fscal,dx00);
802 ty = _mm256_mul_ps(fscal,dy00);
803 tz = _mm256_mul_ps(fscal,dz00);
805 /* Update vectorial force */
806 fix0 = _mm256_add_ps(fix0,tx);
807 fiy0 = _mm256_add_ps(fiy0,ty);
808 fiz0 = _mm256_add_ps(fiz0,tz);
810 fjptrA = f+j_coord_offsetA;
811 fjptrB = f+j_coord_offsetB;
812 fjptrC = f+j_coord_offsetC;
813 fjptrD = f+j_coord_offsetD;
814 fjptrE = f+j_coord_offsetE;
815 fjptrF = f+j_coord_offsetF;
816 fjptrG = f+j_coord_offsetG;
817 fjptrH = f+j_coord_offsetH;
818 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
822 /* Inner loop uses 120 flops */
828 /* Get j neighbor index, and coordinate index */
829 jnrlistA = jjnr[jidx];
830 jnrlistB = jjnr[jidx+1];
831 jnrlistC = jjnr[jidx+2];
832 jnrlistD = jjnr[jidx+3];
833 jnrlistE = jjnr[jidx+4];
834 jnrlistF = jjnr[jidx+5];
835 jnrlistG = jjnr[jidx+6];
836 jnrlistH = jjnr[jidx+7];
837 /* Sign of each element will be negative for non-real atoms.
838 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
839 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
841 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
842 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
844 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
845 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
846 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
847 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
848 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
849 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
850 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
851 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
852 j_coord_offsetA = DIM*jnrA;
853 j_coord_offsetB = DIM*jnrB;
854 j_coord_offsetC = DIM*jnrC;
855 j_coord_offsetD = DIM*jnrD;
856 j_coord_offsetE = DIM*jnrE;
857 j_coord_offsetF = DIM*jnrF;
858 j_coord_offsetG = DIM*jnrG;
859 j_coord_offsetH = DIM*jnrH;
861 /* load j atom coordinates */
862 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
863 x+j_coord_offsetC,x+j_coord_offsetD,
864 x+j_coord_offsetE,x+j_coord_offsetF,
865 x+j_coord_offsetG,x+j_coord_offsetH,
868 /* Calculate displacement vector */
869 dx00 = _mm256_sub_ps(ix0,jx0);
870 dy00 = _mm256_sub_ps(iy0,jy0);
871 dz00 = _mm256_sub_ps(iz0,jz0);
873 /* Calculate squared distance and things based on it */
874 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
876 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
878 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
880 /* Load parameters for j particles */
881 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
882 charge+jnrC+0,charge+jnrD+0,
883 charge+jnrE+0,charge+jnrF+0,
884 charge+jnrG+0,charge+jnrH+0);
885 vdwjidx0A = 2*vdwtype[jnrA+0];
886 vdwjidx0B = 2*vdwtype[jnrB+0];
887 vdwjidx0C = 2*vdwtype[jnrC+0];
888 vdwjidx0D = 2*vdwtype[jnrD+0];
889 vdwjidx0E = 2*vdwtype[jnrE+0];
890 vdwjidx0F = 2*vdwtype[jnrF+0];
891 vdwjidx0G = 2*vdwtype[jnrG+0];
892 vdwjidx0H = 2*vdwtype[jnrH+0];
894 /**************************
895 * CALCULATE INTERACTIONS *
896 **************************/
898 if (gmx_mm256_any_lt(rsq00,rcutoff2))
901 r00 = _mm256_mul_ps(rsq00,rinv00);
902 r00 = _mm256_andnot_ps(dummy_mask,r00);
904 /* Compute parameters for interactions between i and j atoms */
905 qq00 = _mm256_mul_ps(iq0,jq0);
906 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
907 vdwioffsetptr0+vdwjidx0B,
908 vdwioffsetptr0+vdwjidx0C,
909 vdwioffsetptr0+vdwjidx0D,
910 vdwioffsetptr0+vdwjidx0E,
911 vdwioffsetptr0+vdwjidx0F,
912 vdwioffsetptr0+vdwjidx0G,
913 vdwioffsetptr0+vdwjidx0H,
916 /* EWALD ELECTROSTATICS */
918 /* Analytical PME correction */
919 zeta2 = _mm256_mul_ps(beta2,rsq00);
920 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
921 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
922 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
923 felec = _mm256_mul_ps(qq00,felec);
924 pmecorrV = gmx_mm256_pmecorrV_ps(zeta2);
925 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
926 velec = _mm256_sub_ps(rinv00,pmecorrV);
927 velec = _mm256_mul_ps(qq00,velec);
929 /* LENNARD-JONES DISPERSION/REPULSION */
931 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
932 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
933 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
934 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
935 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
937 d = _mm256_sub_ps(r00,rswitch);
938 d = _mm256_max_ps(d,_mm256_setzero_ps());
939 d2 = _mm256_mul_ps(d,d);
940 sw = _mm256_add_ps(one,_mm256_mul_ps(d2,_mm256_mul_ps(d,_mm256_add_ps(swV3,_mm256_mul_ps(d,_mm256_add_ps(swV4,_mm256_mul_ps(d,swV5)))))));
942 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
944 /* Evaluate switch function */
945 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
946 felec = _mm256_sub_ps( _mm256_mul_ps(felec,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(velec,dsw)) );
947 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
948 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
950 fscal = _mm256_add_ps(felec,fvdw);
952 fscal = _mm256_and_ps(fscal,cutoff_mask);
954 fscal = _mm256_andnot_ps(dummy_mask,fscal);
956 /* Calculate temporary vectorial force */
957 tx = _mm256_mul_ps(fscal,dx00);
958 ty = _mm256_mul_ps(fscal,dy00);
959 tz = _mm256_mul_ps(fscal,dz00);
961 /* Update vectorial force */
962 fix0 = _mm256_add_ps(fix0,tx);
963 fiy0 = _mm256_add_ps(fiy0,ty);
964 fiz0 = _mm256_add_ps(fiz0,tz);
966 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
967 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
968 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
969 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
970 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
971 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
972 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
973 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
974 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
978 /* Inner loop uses 121 flops */
981 /* End of innermost loop */
983 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
984 f+i_coord_offset,fshift+i_shift_offset);
986 /* Increment number of inner iterations */
987 inneriter += j_index_end - j_index_start;
989 /* Outer loop uses 7 flops */
992 /* Increment number of outer iterations */
995 /* Update outer/inner flops */
997 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*121);