2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 * Note: this file was generated by the GROMACS avx_256_single kernel generator.
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
47 #include "gromacs/simd/math_x86_avx_256_single.h"
48 #include "kernelutil_x86_avx_256_single.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_VF_avx_256_single
52 * Electrostatics interaction: ReactionField
53 * VdW interaction: LennardJones
54 * Geometry: Particle-Particle
55 * Calculate force/pot: PotentialAndForce
58 nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_VF_avx_256_single
59 (t_nblist * gmx_restrict nlist,
60 rvec * gmx_restrict xx,
61 rvec * gmx_restrict ff,
62 t_forcerec * gmx_restrict fr,
63 t_mdatoms * gmx_restrict mdatoms,
64 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
65 t_nrnb * gmx_restrict nrnb)
67 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
68 * just 0 for non-waters.
69 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
70 * jnr indices corresponding to data put in the four positions in the SIMD register.
72 int i_shift_offset,i_coord_offset,outeriter,inneriter;
73 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
74 int jnrA,jnrB,jnrC,jnrD;
75 int jnrE,jnrF,jnrG,jnrH;
76 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
77 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
78 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
79 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
80 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
82 real *shiftvec,*fshift,*x,*f;
83 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
85 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
86 real * vdwioffsetptr0;
87 __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);
99 __m256 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
100 real rswitch_scalar,d_scalar;
101 __m256 dummy_mask,cutoff_mask;
102 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
103 __m256 one = _mm256_set1_ps(1.0);
104 __m256 two = _mm256_set1_ps(2.0);
110 jindex = nlist->jindex;
112 shiftidx = nlist->shift;
114 shiftvec = fr->shift_vec[0];
115 fshift = fr->fshift[0];
116 facel = _mm256_set1_ps(fr->epsfac);
117 charge = mdatoms->chargeA;
118 krf = _mm256_set1_ps(fr->ic->k_rf);
119 krf2 = _mm256_set1_ps(fr->ic->k_rf*2.0);
120 crf = _mm256_set1_ps(fr->ic->c_rf);
121 nvdwtype = fr->ntype;
123 vdwtype = mdatoms->typeA;
125 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
126 rcutoff_scalar = fr->rcoulomb;
127 rcutoff = _mm256_set1_ps(rcutoff_scalar);
128 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
130 rswitch_scalar = fr->rvdw_switch;
131 rswitch = _mm256_set1_ps(rswitch_scalar);
132 /* Setup switch parameters */
133 d_scalar = rcutoff_scalar-rswitch_scalar;
134 d = _mm256_set1_ps(d_scalar);
135 swV3 = _mm256_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
136 swV4 = _mm256_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
137 swV5 = _mm256_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
138 swF2 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
139 swF3 = _mm256_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
140 swF4 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
142 /* Avoid stupid compiler warnings */
143 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
156 for(iidx=0;iidx<4*DIM;iidx++)
161 /* Start outer loop over neighborlists */
162 for(iidx=0; iidx<nri; iidx++)
164 /* Load shift vector for this list */
165 i_shift_offset = DIM*shiftidx[iidx];
167 /* Load limits for loop over neighbors */
168 j_index_start = jindex[iidx];
169 j_index_end = jindex[iidx+1];
171 /* Get outer coordinate index */
173 i_coord_offset = DIM*inr;
175 /* Load i particle coords and add shift vector */
176 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
178 fix0 = _mm256_setzero_ps();
179 fiy0 = _mm256_setzero_ps();
180 fiz0 = _mm256_setzero_ps();
182 /* Load parameters for i particles */
183 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
184 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
186 /* Reset potential sums */
187 velecsum = _mm256_setzero_ps();
188 vvdwsum = _mm256_setzero_ps();
190 /* Start inner kernel loop */
191 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
194 /* Get j neighbor index, and coordinate index */
203 j_coord_offsetA = DIM*jnrA;
204 j_coord_offsetB = DIM*jnrB;
205 j_coord_offsetC = DIM*jnrC;
206 j_coord_offsetD = DIM*jnrD;
207 j_coord_offsetE = DIM*jnrE;
208 j_coord_offsetF = DIM*jnrF;
209 j_coord_offsetG = DIM*jnrG;
210 j_coord_offsetH = DIM*jnrH;
212 /* load j atom coordinates */
213 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
214 x+j_coord_offsetC,x+j_coord_offsetD,
215 x+j_coord_offsetE,x+j_coord_offsetF,
216 x+j_coord_offsetG,x+j_coord_offsetH,
219 /* Calculate displacement vector */
220 dx00 = _mm256_sub_ps(ix0,jx0);
221 dy00 = _mm256_sub_ps(iy0,jy0);
222 dz00 = _mm256_sub_ps(iz0,jz0);
224 /* Calculate squared distance and things based on it */
225 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
227 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
229 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
231 /* Load parameters for j particles */
232 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
233 charge+jnrC+0,charge+jnrD+0,
234 charge+jnrE+0,charge+jnrF+0,
235 charge+jnrG+0,charge+jnrH+0);
236 vdwjidx0A = 2*vdwtype[jnrA+0];
237 vdwjidx0B = 2*vdwtype[jnrB+0];
238 vdwjidx0C = 2*vdwtype[jnrC+0];
239 vdwjidx0D = 2*vdwtype[jnrD+0];
240 vdwjidx0E = 2*vdwtype[jnrE+0];
241 vdwjidx0F = 2*vdwtype[jnrF+0];
242 vdwjidx0G = 2*vdwtype[jnrG+0];
243 vdwjidx0H = 2*vdwtype[jnrH+0];
245 /**************************
246 * CALCULATE INTERACTIONS *
247 **************************/
249 if (gmx_mm256_any_lt(rsq00,rcutoff2))
252 r00 = _mm256_mul_ps(rsq00,rinv00);
254 /* Compute parameters for interactions between i and j atoms */
255 qq00 = _mm256_mul_ps(iq0,jq0);
256 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
257 vdwioffsetptr0+vdwjidx0B,
258 vdwioffsetptr0+vdwjidx0C,
259 vdwioffsetptr0+vdwjidx0D,
260 vdwioffsetptr0+vdwjidx0E,
261 vdwioffsetptr0+vdwjidx0F,
262 vdwioffsetptr0+vdwjidx0G,
263 vdwioffsetptr0+vdwjidx0H,
266 /* REACTION-FIELD ELECTROSTATICS */
267 velec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
268 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
270 /* LENNARD-JONES DISPERSION/REPULSION */
272 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
273 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
274 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
275 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
276 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
278 d = _mm256_sub_ps(r00,rswitch);
279 d = _mm256_max_ps(d,_mm256_setzero_ps());
280 d2 = _mm256_mul_ps(d,d);
281 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)))))));
283 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
285 /* Evaluate switch function */
286 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
287 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
288 vvdw = _mm256_mul_ps(vvdw,sw);
289 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
291 /* Update potential sum for this i atom from the interaction with this j atom. */
292 velec = _mm256_and_ps(velec,cutoff_mask);
293 velecsum = _mm256_add_ps(velecsum,velec);
294 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
295 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
297 fscal = _mm256_add_ps(felec,fvdw);
299 fscal = _mm256_and_ps(fscal,cutoff_mask);
301 /* Calculate temporary vectorial force */
302 tx = _mm256_mul_ps(fscal,dx00);
303 ty = _mm256_mul_ps(fscal,dy00);
304 tz = _mm256_mul_ps(fscal,dz00);
306 /* Update vectorial force */
307 fix0 = _mm256_add_ps(fix0,tx);
308 fiy0 = _mm256_add_ps(fiy0,ty);
309 fiz0 = _mm256_add_ps(fiz0,tz);
311 fjptrA = f+j_coord_offsetA;
312 fjptrB = f+j_coord_offsetB;
313 fjptrC = f+j_coord_offsetC;
314 fjptrD = f+j_coord_offsetD;
315 fjptrE = f+j_coord_offsetE;
316 fjptrF = f+j_coord_offsetF;
317 fjptrG = f+j_coord_offsetG;
318 fjptrH = f+j_coord_offsetH;
319 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
323 /* Inner loop uses 70 flops */
329 /* Get j neighbor index, and coordinate index */
330 jnrlistA = jjnr[jidx];
331 jnrlistB = jjnr[jidx+1];
332 jnrlistC = jjnr[jidx+2];
333 jnrlistD = jjnr[jidx+3];
334 jnrlistE = jjnr[jidx+4];
335 jnrlistF = jjnr[jidx+5];
336 jnrlistG = jjnr[jidx+6];
337 jnrlistH = jjnr[jidx+7];
338 /* Sign of each element will be negative for non-real atoms.
339 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
340 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
342 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
343 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
345 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
346 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
347 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
348 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
349 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
350 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
351 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
352 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
353 j_coord_offsetA = DIM*jnrA;
354 j_coord_offsetB = DIM*jnrB;
355 j_coord_offsetC = DIM*jnrC;
356 j_coord_offsetD = DIM*jnrD;
357 j_coord_offsetE = DIM*jnrE;
358 j_coord_offsetF = DIM*jnrF;
359 j_coord_offsetG = DIM*jnrG;
360 j_coord_offsetH = DIM*jnrH;
362 /* load j atom coordinates */
363 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
364 x+j_coord_offsetC,x+j_coord_offsetD,
365 x+j_coord_offsetE,x+j_coord_offsetF,
366 x+j_coord_offsetG,x+j_coord_offsetH,
369 /* Calculate displacement vector */
370 dx00 = _mm256_sub_ps(ix0,jx0);
371 dy00 = _mm256_sub_ps(iy0,jy0);
372 dz00 = _mm256_sub_ps(iz0,jz0);
374 /* Calculate squared distance and things based on it */
375 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
377 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
379 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
381 /* Load parameters for j particles */
382 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
383 charge+jnrC+0,charge+jnrD+0,
384 charge+jnrE+0,charge+jnrF+0,
385 charge+jnrG+0,charge+jnrH+0);
386 vdwjidx0A = 2*vdwtype[jnrA+0];
387 vdwjidx0B = 2*vdwtype[jnrB+0];
388 vdwjidx0C = 2*vdwtype[jnrC+0];
389 vdwjidx0D = 2*vdwtype[jnrD+0];
390 vdwjidx0E = 2*vdwtype[jnrE+0];
391 vdwjidx0F = 2*vdwtype[jnrF+0];
392 vdwjidx0G = 2*vdwtype[jnrG+0];
393 vdwjidx0H = 2*vdwtype[jnrH+0];
395 /**************************
396 * CALCULATE INTERACTIONS *
397 **************************/
399 if (gmx_mm256_any_lt(rsq00,rcutoff2))
402 r00 = _mm256_mul_ps(rsq00,rinv00);
403 r00 = _mm256_andnot_ps(dummy_mask,r00);
405 /* Compute parameters for interactions between i and j atoms */
406 qq00 = _mm256_mul_ps(iq0,jq0);
407 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
408 vdwioffsetptr0+vdwjidx0B,
409 vdwioffsetptr0+vdwjidx0C,
410 vdwioffsetptr0+vdwjidx0D,
411 vdwioffsetptr0+vdwjidx0E,
412 vdwioffsetptr0+vdwjidx0F,
413 vdwioffsetptr0+vdwjidx0G,
414 vdwioffsetptr0+vdwjidx0H,
417 /* REACTION-FIELD ELECTROSTATICS */
418 velec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
419 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
421 /* LENNARD-JONES DISPERSION/REPULSION */
423 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
424 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
425 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
426 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
427 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
429 d = _mm256_sub_ps(r00,rswitch);
430 d = _mm256_max_ps(d,_mm256_setzero_ps());
431 d2 = _mm256_mul_ps(d,d);
432 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)))))));
434 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
436 /* Evaluate switch function */
437 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
438 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
439 vvdw = _mm256_mul_ps(vvdw,sw);
440 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
442 /* Update potential sum for this i atom from the interaction with this j atom. */
443 velec = _mm256_and_ps(velec,cutoff_mask);
444 velec = _mm256_andnot_ps(dummy_mask,velec);
445 velecsum = _mm256_add_ps(velecsum,velec);
446 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
447 vvdw = _mm256_andnot_ps(dummy_mask,vvdw);
448 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
450 fscal = _mm256_add_ps(felec,fvdw);
452 fscal = _mm256_and_ps(fscal,cutoff_mask);
454 fscal = _mm256_andnot_ps(dummy_mask,fscal);
456 /* Calculate temporary vectorial force */
457 tx = _mm256_mul_ps(fscal,dx00);
458 ty = _mm256_mul_ps(fscal,dy00);
459 tz = _mm256_mul_ps(fscal,dz00);
461 /* Update vectorial force */
462 fix0 = _mm256_add_ps(fix0,tx);
463 fiy0 = _mm256_add_ps(fiy0,ty);
464 fiz0 = _mm256_add_ps(fiz0,tz);
466 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
467 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
468 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
469 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
470 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
471 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
472 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
473 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
474 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
478 /* Inner loop uses 71 flops */
481 /* End of innermost loop */
483 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
484 f+i_coord_offset,fshift+i_shift_offset);
487 /* Update potential energies */
488 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
489 gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
491 /* Increment number of inner iterations */
492 inneriter += j_index_end - j_index_start;
494 /* Outer loop uses 9 flops */
497 /* Increment number of outer iterations */
500 /* Update outer/inner flops */
502 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*71);
505 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_F_avx_256_single
506 * Electrostatics interaction: ReactionField
507 * VdW interaction: LennardJones
508 * Geometry: Particle-Particle
509 * Calculate force/pot: Force
512 nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_F_avx_256_single
513 (t_nblist * gmx_restrict nlist,
514 rvec * gmx_restrict xx,
515 rvec * gmx_restrict ff,
516 t_forcerec * gmx_restrict fr,
517 t_mdatoms * gmx_restrict mdatoms,
518 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
519 t_nrnb * gmx_restrict nrnb)
521 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
522 * just 0 for non-waters.
523 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
524 * jnr indices corresponding to data put in the four positions in the SIMD register.
526 int i_shift_offset,i_coord_offset,outeriter,inneriter;
527 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
528 int jnrA,jnrB,jnrC,jnrD;
529 int jnrE,jnrF,jnrG,jnrH;
530 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
531 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
532 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
533 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
534 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
536 real *shiftvec,*fshift,*x,*f;
537 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
539 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
540 real * vdwioffsetptr0;
541 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
542 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
543 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
544 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
545 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
548 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
551 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
552 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
553 __m256 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
554 real rswitch_scalar,d_scalar;
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->epsfac);
571 charge = mdatoms->chargeA;
572 krf = _mm256_set1_ps(fr->ic->k_rf);
573 krf2 = _mm256_set1_ps(fr->ic->k_rf*2.0);
574 crf = _mm256_set1_ps(fr->ic->c_rf);
575 nvdwtype = fr->ntype;
577 vdwtype = mdatoms->typeA;
579 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
580 rcutoff_scalar = fr->rcoulomb;
581 rcutoff = _mm256_set1_ps(rcutoff_scalar);
582 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
584 rswitch_scalar = fr->rvdw_switch;
585 rswitch = _mm256_set1_ps(rswitch_scalar);
586 /* Setup switch parameters */
587 d_scalar = rcutoff_scalar-rswitch_scalar;
588 d = _mm256_set1_ps(d_scalar);
589 swV3 = _mm256_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
590 swV4 = _mm256_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
591 swV5 = _mm256_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
592 swF2 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
593 swF3 = _mm256_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
594 swF4 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
596 /* Avoid stupid compiler warnings */
597 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
610 for(iidx=0;iidx<4*DIM;iidx++)
615 /* Start outer loop over neighborlists */
616 for(iidx=0; iidx<nri; iidx++)
618 /* Load shift vector for this list */
619 i_shift_offset = DIM*shiftidx[iidx];
621 /* Load limits for loop over neighbors */
622 j_index_start = jindex[iidx];
623 j_index_end = jindex[iidx+1];
625 /* Get outer coordinate index */
627 i_coord_offset = DIM*inr;
629 /* Load i particle coords and add shift vector */
630 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
632 fix0 = _mm256_setzero_ps();
633 fiy0 = _mm256_setzero_ps();
634 fiz0 = _mm256_setzero_ps();
636 /* Load parameters for i particles */
637 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
638 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
640 /* Start inner kernel loop */
641 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
644 /* Get j neighbor index, and coordinate index */
653 j_coord_offsetA = DIM*jnrA;
654 j_coord_offsetB = DIM*jnrB;
655 j_coord_offsetC = DIM*jnrC;
656 j_coord_offsetD = DIM*jnrD;
657 j_coord_offsetE = DIM*jnrE;
658 j_coord_offsetF = DIM*jnrF;
659 j_coord_offsetG = DIM*jnrG;
660 j_coord_offsetH = DIM*jnrH;
662 /* load j atom coordinates */
663 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
664 x+j_coord_offsetC,x+j_coord_offsetD,
665 x+j_coord_offsetE,x+j_coord_offsetF,
666 x+j_coord_offsetG,x+j_coord_offsetH,
669 /* Calculate displacement vector */
670 dx00 = _mm256_sub_ps(ix0,jx0);
671 dy00 = _mm256_sub_ps(iy0,jy0);
672 dz00 = _mm256_sub_ps(iz0,jz0);
674 /* Calculate squared distance and things based on it */
675 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
677 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
679 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
681 /* Load parameters for j particles */
682 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
683 charge+jnrC+0,charge+jnrD+0,
684 charge+jnrE+0,charge+jnrF+0,
685 charge+jnrG+0,charge+jnrH+0);
686 vdwjidx0A = 2*vdwtype[jnrA+0];
687 vdwjidx0B = 2*vdwtype[jnrB+0];
688 vdwjidx0C = 2*vdwtype[jnrC+0];
689 vdwjidx0D = 2*vdwtype[jnrD+0];
690 vdwjidx0E = 2*vdwtype[jnrE+0];
691 vdwjidx0F = 2*vdwtype[jnrF+0];
692 vdwjidx0G = 2*vdwtype[jnrG+0];
693 vdwjidx0H = 2*vdwtype[jnrH+0];
695 /**************************
696 * CALCULATE INTERACTIONS *
697 **************************/
699 if (gmx_mm256_any_lt(rsq00,rcutoff2))
702 r00 = _mm256_mul_ps(rsq00,rinv00);
704 /* Compute parameters for interactions between i and j atoms */
705 qq00 = _mm256_mul_ps(iq0,jq0);
706 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
707 vdwioffsetptr0+vdwjidx0B,
708 vdwioffsetptr0+vdwjidx0C,
709 vdwioffsetptr0+vdwjidx0D,
710 vdwioffsetptr0+vdwjidx0E,
711 vdwioffsetptr0+vdwjidx0F,
712 vdwioffsetptr0+vdwjidx0G,
713 vdwioffsetptr0+vdwjidx0H,
716 /* REACTION-FIELD ELECTROSTATICS */
717 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
719 /* LENNARD-JONES DISPERSION/REPULSION */
721 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
722 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
723 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
724 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
725 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
727 d = _mm256_sub_ps(r00,rswitch);
728 d = _mm256_max_ps(d,_mm256_setzero_ps());
729 d2 = _mm256_mul_ps(d,d);
730 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)))))));
732 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
734 /* Evaluate switch function */
735 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
736 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
737 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
739 fscal = _mm256_add_ps(felec,fvdw);
741 fscal = _mm256_and_ps(fscal,cutoff_mask);
743 /* Calculate temporary vectorial force */
744 tx = _mm256_mul_ps(fscal,dx00);
745 ty = _mm256_mul_ps(fscal,dy00);
746 tz = _mm256_mul_ps(fscal,dz00);
748 /* Update vectorial force */
749 fix0 = _mm256_add_ps(fix0,tx);
750 fiy0 = _mm256_add_ps(fiy0,ty);
751 fiz0 = _mm256_add_ps(fiz0,tz);
753 fjptrA = f+j_coord_offsetA;
754 fjptrB = f+j_coord_offsetB;
755 fjptrC = f+j_coord_offsetC;
756 fjptrD = f+j_coord_offsetD;
757 fjptrE = f+j_coord_offsetE;
758 fjptrF = f+j_coord_offsetF;
759 fjptrG = f+j_coord_offsetG;
760 fjptrH = f+j_coord_offsetH;
761 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
765 /* Inner loop uses 61 flops */
771 /* Get j neighbor index, and coordinate index */
772 jnrlistA = jjnr[jidx];
773 jnrlistB = jjnr[jidx+1];
774 jnrlistC = jjnr[jidx+2];
775 jnrlistD = jjnr[jidx+3];
776 jnrlistE = jjnr[jidx+4];
777 jnrlistF = jjnr[jidx+5];
778 jnrlistG = jjnr[jidx+6];
779 jnrlistH = jjnr[jidx+7];
780 /* Sign of each element will be negative for non-real atoms.
781 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
782 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
784 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
785 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
787 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
788 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
789 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
790 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
791 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
792 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
793 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
794 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
795 j_coord_offsetA = DIM*jnrA;
796 j_coord_offsetB = DIM*jnrB;
797 j_coord_offsetC = DIM*jnrC;
798 j_coord_offsetD = DIM*jnrD;
799 j_coord_offsetE = DIM*jnrE;
800 j_coord_offsetF = DIM*jnrF;
801 j_coord_offsetG = DIM*jnrG;
802 j_coord_offsetH = DIM*jnrH;
804 /* load j atom coordinates */
805 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
806 x+j_coord_offsetC,x+j_coord_offsetD,
807 x+j_coord_offsetE,x+j_coord_offsetF,
808 x+j_coord_offsetG,x+j_coord_offsetH,
811 /* Calculate displacement vector */
812 dx00 = _mm256_sub_ps(ix0,jx0);
813 dy00 = _mm256_sub_ps(iy0,jy0);
814 dz00 = _mm256_sub_ps(iz0,jz0);
816 /* Calculate squared distance and things based on it */
817 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
819 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
821 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
823 /* Load parameters for j particles */
824 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
825 charge+jnrC+0,charge+jnrD+0,
826 charge+jnrE+0,charge+jnrF+0,
827 charge+jnrG+0,charge+jnrH+0);
828 vdwjidx0A = 2*vdwtype[jnrA+0];
829 vdwjidx0B = 2*vdwtype[jnrB+0];
830 vdwjidx0C = 2*vdwtype[jnrC+0];
831 vdwjidx0D = 2*vdwtype[jnrD+0];
832 vdwjidx0E = 2*vdwtype[jnrE+0];
833 vdwjidx0F = 2*vdwtype[jnrF+0];
834 vdwjidx0G = 2*vdwtype[jnrG+0];
835 vdwjidx0H = 2*vdwtype[jnrH+0];
837 /**************************
838 * CALCULATE INTERACTIONS *
839 **************************/
841 if (gmx_mm256_any_lt(rsq00,rcutoff2))
844 r00 = _mm256_mul_ps(rsq00,rinv00);
845 r00 = _mm256_andnot_ps(dummy_mask,r00);
847 /* Compute parameters for interactions between i and j atoms */
848 qq00 = _mm256_mul_ps(iq0,jq0);
849 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
850 vdwioffsetptr0+vdwjidx0B,
851 vdwioffsetptr0+vdwjidx0C,
852 vdwioffsetptr0+vdwjidx0D,
853 vdwioffsetptr0+vdwjidx0E,
854 vdwioffsetptr0+vdwjidx0F,
855 vdwioffsetptr0+vdwjidx0G,
856 vdwioffsetptr0+vdwjidx0H,
859 /* REACTION-FIELD ELECTROSTATICS */
860 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
862 /* LENNARD-JONES DISPERSION/REPULSION */
864 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
865 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
866 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
867 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
868 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
870 d = _mm256_sub_ps(r00,rswitch);
871 d = _mm256_max_ps(d,_mm256_setzero_ps());
872 d2 = _mm256_mul_ps(d,d);
873 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)))))));
875 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
877 /* Evaluate switch function */
878 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
879 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
880 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
882 fscal = _mm256_add_ps(felec,fvdw);
884 fscal = _mm256_and_ps(fscal,cutoff_mask);
886 fscal = _mm256_andnot_ps(dummy_mask,fscal);
888 /* Calculate temporary vectorial force */
889 tx = _mm256_mul_ps(fscal,dx00);
890 ty = _mm256_mul_ps(fscal,dy00);
891 tz = _mm256_mul_ps(fscal,dz00);
893 /* Update vectorial force */
894 fix0 = _mm256_add_ps(fix0,tx);
895 fiy0 = _mm256_add_ps(fiy0,ty);
896 fiz0 = _mm256_add_ps(fiz0,tz);
898 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
899 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
900 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
901 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
902 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
903 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
904 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
905 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
906 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
910 /* Inner loop uses 62 flops */
913 /* End of innermost loop */
915 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
916 f+i_coord_offset,fshift+i_shift_offset);
918 /* Increment number of inner iterations */
919 inneriter += j_index_end - j_index_start;
921 /* Outer loop uses 7 flops */
924 /* Increment number of outer iterations */
927 /* Update outer/inner flops */
929 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*62);