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_ElecRFCut_VdwLJSw_GeomP1P1_VF_avx_256_single
54 * Electrostatics interaction: ReactionField
55 * VdW interaction: LennardJones
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecRFCut_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);
101 __m256 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
102 real rswitch_scalar,d_scalar;
103 __m256 dummy_mask,cutoff_mask;
104 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
105 __m256 one = _mm256_set1_ps(1.0);
106 __m256 two = _mm256_set1_ps(2.0);
112 jindex = nlist->jindex;
114 shiftidx = nlist->shift;
116 shiftvec = fr->shift_vec[0];
117 fshift = fr->fshift[0];
118 facel = _mm256_set1_ps(fr->epsfac);
119 charge = mdatoms->chargeA;
120 krf = _mm256_set1_ps(fr->ic->k_rf);
121 krf2 = _mm256_set1_ps(fr->ic->k_rf*2.0);
122 crf = _mm256_set1_ps(fr->ic->c_rf);
123 nvdwtype = fr->ntype;
125 vdwtype = mdatoms->typeA;
127 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
128 rcutoff_scalar = fr->rcoulomb;
129 rcutoff = _mm256_set1_ps(rcutoff_scalar);
130 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
132 rswitch_scalar = fr->rvdw_switch;
133 rswitch = _mm256_set1_ps(rswitch_scalar);
134 /* Setup switch parameters */
135 d_scalar = rcutoff_scalar-rswitch_scalar;
136 d = _mm256_set1_ps(d_scalar);
137 swV3 = _mm256_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
138 swV4 = _mm256_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
139 swV5 = _mm256_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
140 swF2 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
141 swF3 = _mm256_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
142 swF4 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
144 /* Avoid stupid compiler warnings */
145 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
158 for(iidx=0;iidx<4*DIM;iidx++)
163 /* Start outer loop over neighborlists */
164 for(iidx=0; iidx<nri; iidx++)
166 /* Load shift vector for this list */
167 i_shift_offset = DIM*shiftidx[iidx];
169 /* Load limits for loop over neighbors */
170 j_index_start = jindex[iidx];
171 j_index_end = jindex[iidx+1];
173 /* Get outer coordinate index */
175 i_coord_offset = DIM*inr;
177 /* Load i particle coords and add shift vector */
178 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
180 fix0 = _mm256_setzero_ps();
181 fiy0 = _mm256_setzero_ps();
182 fiz0 = _mm256_setzero_ps();
184 /* Load parameters for i particles */
185 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
186 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
188 /* Reset potential sums */
189 velecsum = _mm256_setzero_ps();
190 vvdwsum = _mm256_setzero_ps();
192 /* Start inner kernel loop */
193 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
196 /* Get j neighbor index, and coordinate index */
205 j_coord_offsetA = DIM*jnrA;
206 j_coord_offsetB = DIM*jnrB;
207 j_coord_offsetC = DIM*jnrC;
208 j_coord_offsetD = DIM*jnrD;
209 j_coord_offsetE = DIM*jnrE;
210 j_coord_offsetF = DIM*jnrF;
211 j_coord_offsetG = DIM*jnrG;
212 j_coord_offsetH = DIM*jnrH;
214 /* load j atom coordinates */
215 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
216 x+j_coord_offsetC,x+j_coord_offsetD,
217 x+j_coord_offsetE,x+j_coord_offsetF,
218 x+j_coord_offsetG,x+j_coord_offsetH,
221 /* Calculate displacement vector */
222 dx00 = _mm256_sub_ps(ix0,jx0);
223 dy00 = _mm256_sub_ps(iy0,jy0);
224 dz00 = _mm256_sub_ps(iz0,jz0);
226 /* Calculate squared distance and things based on it */
227 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
229 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
231 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
233 /* Load parameters for j particles */
234 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
235 charge+jnrC+0,charge+jnrD+0,
236 charge+jnrE+0,charge+jnrF+0,
237 charge+jnrG+0,charge+jnrH+0);
238 vdwjidx0A = 2*vdwtype[jnrA+0];
239 vdwjidx0B = 2*vdwtype[jnrB+0];
240 vdwjidx0C = 2*vdwtype[jnrC+0];
241 vdwjidx0D = 2*vdwtype[jnrD+0];
242 vdwjidx0E = 2*vdwtype[jnrE+0];
243 vdwjidx0F = 2*vdwtype[jnrF+0];
244 vdwjidx0G = 2*vdwtype[jnrG+0];
245 vdwjidx0H = 2*vdwtype[jnrH+0];
247 /**************************
248 * CALCULATE INTERACTIONS *
249 **************************/
251 if (gmx_mm256_any_lt(rsq00,rcutoff2))
254 r00 = _mm256_mul_ps(rsq00,rinv00);
256 /* Compute parameters for interactions between i and j atoms */
257 qq00 = _mm256_mul_ps(iq0,jq0);
258 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
259 vdwioffsetptr0+vdwjidx0B,
260 vdwioffsetptr0+vdwjidx0C,
261 vdwioffsetptr0+vdwjidx0D,
262 vdwioffsetptr0+vdwjidx0E,
263 vdwioffsetptr0+vdwjidx0F,
264 vdwioffsetptr0+vdwjidx0G,
265 vdwioffsetptr0+vdwjidx0H,
268 /* REACTION-FIELD ELECTROSTATICS */
269 velec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
270 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
272 /* LENNARD-JONES DISPERSION/REPULSION */
274 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
275 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
276 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
277 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
278 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
280 d = _mm256_sub_ps(r00,rswitch);
281 d = _mm256_max_ps(d,_mm256_setzero_ps());
282 d2 = _mm256_mul_ps(d,d);
283 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)))))));
285 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
287 /* Evaluate switch function */
288 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
289 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
290 vvdw = _mm256_mul_ps(vvdw,sw);
291 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
293 /* Update potential sum for this i atom from the interaction with this j atom. */
294 velec = _mm256_and_ps(velec,cutoff_mask);
295 velecsum = _mm256_add_ps(velecsum,velec);
296 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
297 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
299 fscal = _mm256_add_ps(felec,fvdw);
301 fscal = _mm256_and_ps(fscal,cutoff_mask);
303 /* Calculate temporary vectorial force */
304 tx = _mm256_mul_ps(fscal,dx00);
305 ty = _mm256_mul_ps(fscal,dy00);
306 tz = _mm256_mul_ps(fscal,dz00);
308 /* Update vectorial force */
309 fix0 = _mm256_add_ps(fix0,tx);
310 fiy0 = _mm256_add_ps(fiy0,ty);
311 fiz0 = _mm256_add_ps(fiz0,tz);
313 fjptrA = f+j_coord_offsetA;
314 fjptrB = f+j_coord_offsetB;
315 fjptrC = f+j_coord_offsetC;
316 fjptrD = f+j_coord_offsetD;
317 fjptrE = f+j_coord_offsetE;
318 fjptrF = f+j_coord_offsetF;
319 fjptrG = f+j_coord_offsetG;
320 fjptrH = f+j_coord_offsetH;
321 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
325 /* Inner loop uses 70 flops */
331 /* Get j neighbor index, and coordinate index */
332 jnrlistA = jjnr[jidx];
333 jnrlistB = jjnr[jidx+1];
334 jnrlistC = jjnr[jidx+2];
335 jnrlistD = jjnr[jidx+3];
336 jnrlistE = jjnr[jidx+4];
337 jnrlistF = jjnr[jidx+5];
338 jnrlistG = jjnr[jidx+6];
339 jnrlistH = jjnr[jidx+7];
340 /* Sign of each element will be negative for non-real atoms.
341 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
342 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
344 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
345 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
347 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
348 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
349 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
350 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
351 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
352 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
353 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
354 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
355 j_coord_offsetA = DIM*jnrA;
356 j_coord_offsetB = DIM*jnrB;
357 j_coord_offsetC = DIM*jnrC;
358 j_coord_offsetD = DIM*jnrD;
359 j_coord_offsetE = DIM*jnrE;
360 j_coord_offsetF = DIM*jnrF;
361 j_coord_offsetG = DIM*jnrG;
362 j_coord_offsetH = DIM*jnrH;
364 /* load j atom coordinates */
365 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
366 x+j_coord_offsetC,x+j_coord_offsetD,
367 x+j_coord_offsetE,x+j_coord_offsetF,
368 x+j_coord_offsetG,x+j_coord_offsetH,
371 /* Calculate displacement vector */
372 dx00 = _mm256_sub_ps(ix0,jx0);
373 dy00 = _mm256_sub_ps(iy0,jy0);
374 dz00 = _mm256_sub_ps(iz0,jz0);
376 /* Calculate squared distance and things based on it */
377 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
379 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
381 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
383 /* Load parameters for j particles */
384 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
385 charge+jnrC+0,charge+jnrD+0,
386 charge+jnrE+0,charge+jnrF+0,
387 charge+jnrG+0,charge+jnrH+0);
388 vdwjidx0A = 2*vdwtype[jnrA+0];
389 vdwjidx0B = 2*vdwtype[jnrB+0];
390 vdwjidx0C = 2*vdwtype[jnrC+0];
391 vdwjidx0D = 2*vdwtype[jnrD+0];
392 vdwjidx0E = 2*vdwtype[jnrE+0];
393 vdwjidx0F = 2*vdwtype[jnrF+0];
394 vdwjidx0G = 2*vdwtype[jnrG+0];
395 vdwjidx0H = 2*vdwtype[jnrH+0];
397 /**************************
398 * CALCULATE INTERACTIONS *
399 **************************/
401 if (gmx_mm256_any_lt(rsq00,rcutoff2))
404 r00 = _mm256_mul_ps(rsq00,rinv00);
405 r00 = _mm256_andnot_ps(dummy_mask,r00);
407 /* Compute parameters for interactions between i and j atoms */
408 qq00 = _mm256_mul_ps(iq0,jq0);
409 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
410 vdwioffsetptr0+vdwjidx0B,
411 vdwioffsetptr0+vdwjidx0C,
412 vdwioffsetptr0+vdwjidx0D,
413 vdwioffsetptr0+vdwjidx0E,
414 vdwioffsetptr0+vdwjidx0F,
415 vdwioffsetptr0+vdwjidx0G,
416 vdwioffsetptr0+vdwjidx0H,
419 /* REACTION-FIELD ELECTROSTATICS */
420 velec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
421 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
423 /* LENNARD-JONES DISPERSION/REPULSION */
425 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
426 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
427 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
428 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
429 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
431 d = _mm256_sub_ps(r00,rswitch);
432 d = _mm256_max_ps(d,_mm256_setzero_ps());
433 d2 = _mm256_mul_ps(d,d);
434 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)))))));
436 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
438 /* Evaluate switch function */
439 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
440 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
441 vvdw = _mm256_mul_ps(vvdw,sw);
442 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
444 /* Update potential sum for this i atom from the interaction with this j atom. */
445 velec = _mm256_and_ps(velec,cutoff_mask);
446 velec = _mm256_andnot_ps(dummy_mask,velec);
447 velecsum = _mm256_add_ps(velecsum,velec);
448 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
449 vvdw = _mm256_andnot_ps(dummy_mask,vvdw);
450 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
452 fscal = _mm256_add_ps(felec,fvdw);
454 fscal = _mm256_and_ps(fscal,cutoff_mask);
456 fscal = _mm256_andnot_ps(dummy_mask,fscal);
458 /* Calculate temporary vectorial force */
459 tx = _mm256_mul_ps(fscal,dx00);
460 ty = _mm256_mul_ps(fscal,dy00);
461 tz = _mm256_mul_ps(fscal,dz00);
463 /* Update vectorial force */
464 fix0 = _mm256_add_ps(fix0,tx);
465 fiy0 = _mm256_add_ps(fiy0,ty);
466 fiz0 = _mm256_add_ps(fiz0,tz);
468 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
469 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
470 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
471 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
472 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
473 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
474 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
475 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
476 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
480 /* Inner loop uses 71 flops */
483 /* End of innermost loop */
485 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
486 f+i_coord_offset,fshift+i_shift_offset);
489 /* Update potential energies */
490 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
491 gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
493 /* Increment number of inner iterations */
494 inneriter += j_index_end - j_index_start;
496 /* Outer loop uses 9 flops */
499 /* Increment number of outer iterations */
502 /* Update outer/inner flops */
504 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*71);
507 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_F_avx_256_single
508 * Electrostatics interaction: ReactionField
509 * VdW interaction: LennardJones
510 * Geometry: Particle-Particle
511 * Calculate force/pot: Force
514 nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_F_avx_256_single
515 (t_nblist * gmx_restrict nlist,
516 rvec * gmx_restrict xx,
517 rvec * gmx_restrict ff,
518 t_forcerec * gmx_restrict fr,
519 t_mdatoms * gmx_restrict mdatoms,
520 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
521 t_nrnb * gmx_restrict nrnb)
523 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
524 * just 0 for non-waters.
525 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
526 * jnr indices corresponding to data put in the four positions in the SIMD register.
528 int i_shift_offset,i_coord_offset,outeriter,inneriter;
529 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
530 int jnrA,jnrB,jnrC,jnrD;
531 int jnrE,jnrF,jnrG,jnrH;
532 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
533 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
534 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
535 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
536 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
538 real *shiftvec,*fshift,*x,*f;
539 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
541 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
542 real * vdwioffsetptr0;
543 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
544 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
545 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
546 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
547 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
550 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
553 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
554 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
555 __m256 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
556 real rswitch_scalar,d_scalar;
557 __m256 dummy_mask,cutoff_mask;
558 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
559 __m256 one = _mm256_set1_ps(1.0);
560 __m256 two = _mm256_set1_ps(2.0);
566 jindex = nlist->jindex;
568 shiftidx = nlist->shift;
570 shiftvec = fr->shift_vec[0];
571 fshift = fr->fshift[0];
572 facel = _mm256_set1_ps(fr->epsfac);
573 charge = mdatoms->chargeA;
574 krf = _mm256_set1_ps(fr->ic->k_rf);
575 krf2 = _mm256_set1_ps(fr->ic->k_rf*2.0);
576 crf = _mm256_set1_ps(fr->ic->c_rf);
577 nvdwtype = fr->ntype;
579 vdwtype = mdatoms->typeA;
581 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
582 rcutoff_scalar = fr->rcoulomb;
583 rcutoff = _mm256_set1_ps(rcutoff_scalar);
584 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
586 rswitch_scalar = fr->rvdw_switch;
587 rswitch = _mm256_set1_ps(rswitch_scalar);
588 /* Setup switch parameters */
589 d_scalar = rcutoff_scalar-rswitch_scalar;
590 d = _mm256_set1_ps(d_scalar);
591 swV3 = _mm256_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
592 swV4 = _mm256_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
593 swV5 = _mm256_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
594 swF2 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
595 swF3 = _mm256_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
596 swF4 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
598 /* Avoid stupid compiler warnings */
599 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
612 for(iidx=0;iidx<4*DIM;iidx++)
617 /* Start outer loop over neighborlists */
618 for(iidx=0; iidx<nri; iidx++)
620 /* Load shift vector for this list */
621 i_shift_offset = DIM*shiftidx[iidx];
623 /* Load limits for loop over neighbors */
624 j_index_start = jindex[iidx];
625 j_index_end = jindex[iidx+1];
627 /* Get outer coordinate index */
629 i_coord_offset = DIM*inr;
631 /* Load i particle coords and add shift vector */
632 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
634 fix0 = _mm256_setzero_ps();
635 fiy0 = _mm256_setzero_ps();
636 fiz0 = _mm256_setzero_ps();
638 /* Load parameters for i particles */
639 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
640 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
642 /* Start inner kernel loop */
643 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
646 /* Get j neighbor index, and coordinate index */
655 j_coord_offsetA = DIM*jnrA;
656 j_coord_offsetB = DIM*jnrB;
657 j_coord_offsetC = DIM*jnrC;
658 j_coord_offsetD = DIM*jnrD;
659 j_coord_offsetE = DIM*jnrE;
660 j_coord_offsetF = DIM*jnrF;
661 j_coord_offsetG = DIM*jnrG;
662 j_coord_offsetH = DIM*jnrH;
664 /* load j atom coordinates */
665 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
666 x+j_coord_offsetC,x+j_coord_offsetD,
667 x+j_coord_offsetE,x+j_coord_offsetF,
668 x+j_coord_offsetG,x+j_coord_offsetH,
671 /* Calculate displacement vector */
672 dx00 = _mm256_sub_ps(ix0,jx0);
673 dy00 = _mm256_sub_ps(iy0,jy0);
674 dz00 = _mm256_sub_ps(iz0,jz0);
676 /* Calculate squared distance and things based on it */
677 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
679 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
681 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
683 /* Load parameters for j particles */
684 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
685 charge+jnrC+0,charge+jnrD+0,
686 charge+jnrE+0,charge+jnrF+0,
687 charge+jnrG+0,charge+jnrH+0);
688 vdwjidx0A = 2*vdwtype[jnrA+0];
689 vdwjidx0B = 2*vdwtype[jnrB+0];
690 vdwjidx0C = 2*vdwtype[jnrC+0];
691 vdwjidx0D = 2*vdwtype[jnrD+0];
692 vdwjidx0E = 2*vdwtype[jnrE+0];
693 vdwjidx0F = 2*vdwtype[jnrF+0];
694 vdwjidx0G = 2*vdwtype[jnrG+0];
695 vdwjidx0H = 2*vdwtype[jnrH+0];
697 /**************************
698 * CALCULATE INTERACTIONS *
699 **************************/
701 if (gmx_mm256_any_lt(rsq00,rcutoff2))
704 r00 = _mm256_mul_ps(rsq00,rinv00);
706 /* Compute parameters for interactions between i and j atoms */
707 qq00 = _mm256_mul_ps(iq0,jq0);
708 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
709 vdwioffsetptr0+vdwjidx0B,
710 vdwioffsetptr0+vdwjidx0C,
711 vdwioffsetptr0+vdwjidx0D,
712 vdwioffsetptr0+vdwjidx0E,
713 vdwioffsetptr0+vdwjidx0F,
714 vdwioffsetptr0+vdwjidx0G,
715 vdwioffsetptr0+vdwjidx0H,
718 /* REACTION-FIELD ELECTROSTATICS */
719 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
721 /* LENNARD-JONES DISPERSION/REPULSION */
723 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
724 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
725 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
726 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
727 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
729 d = _mm256_sub_ps(r00,rswitch);
730 d = _mm256_max_ps(d,_mm256_setzero_ps());
731 d2 = _mm256_mul_ps(d,d);
732 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)))))));
734 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
736 /* Evaluate switch function */
737 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
738 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
739 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
741 fscal = _mm256_add_ps(felec,fvdw);
743 fscal = _mm256_and_ps(fscal,cutoff_mask);
745 /* Calculate temporary vectorial force */
746 tx = _mm256_mul_ps(fscal,dx00);
747 ty = _mm256_mul_ps(fscal,dy00);
748 tz = _mm256_mul_ps(fscal,dz00);
750 /* Update vectorial force */
751 fix0 = _mm256_add_ps(fix0,tx);
752 fiy0 = _mm256_add_ps(fiy0,ty);
753 fiz0 = _mm256_add_ps(fiz0,tz);
755 fjptrA = f+j_coord_offsetA;
756 fjptrB = f+j_coord_offsetB;
757 fjptrC = f+j_coord_offsetC;
758 fjptrD = f+j_coord_offsetD;
759 fjptrE = f+j_coord_offsetE;
760 fjptrF = f+j_coord_offsetF;
761 fjptrG = f+j_coord_offsetG;
762 fjptrH = f+j_coord_offsetH;
763 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
767 /* Inner loop uses 61 flops */
773 /* Get j neighbor index, and coordinate index */
774 jnrlistA = jjnr[jidx];
775 jnrlistB = jjnr[jidx+1];
776 jnrlistC = jjnr[jidx+2];
777 jnrlistD = jjnr[jidx+3];
778 jnrlistE = jjnr[jidx+4];
779 jnrlistF = jjnr[jidx+5];
780 jnrlistG = jjnr[jidx+6];
781 jnrlistH = jjnr[jidx+7];
782 /* Sign of each element will be negative for non-real atoms.
783 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
784 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
786 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
787 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
789 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
790 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
791 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
792 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
793 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
794 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
795 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
796 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
797 j_coord_offsetA = DIM*jnrA;
798 j_coord_offsetB = DIM*jnrB;
799 j_coord_offsetC = DIM*jnrC;
800 j_coord_offsetD = DIM*jnrD;
801 j_coord_offsetE = DIM*jnrE;
802 j_coord_offsetF = DIM*jnrF;
803 j_coord_offsetG = DIM*jnrG;
804 j_coord_offsetH = DIM*jnrH;
806 /* load j atom coordinates */
807 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
808 x+j_coord_offsetC,x+j_coord_offsetD,
809 x+j_coord_offsetE,x+j_coord_offsetF,
810 x+j_coord_offsetG,x+j_coord_offsetH,
813 /* Calculate displacement vector */
814 dx00 = _mm256_sub_ps(ix0,jx0);
815 dy00 = _mm256_sub_ps(iy0,jy0);
816 dz00 = _mm256_sub_ps(iz0,jz0);
818 /* Calculate squared distance and things based on it */
819 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
821 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
823 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
825 /* Load parameters for j particles */
826 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
827 charge+jnrC+0,charge+jnrD+0,
828 charge+jnrE+0,charge+jnrF+0,
829 charge+jnrG+0,charge+jnrH+0);
830 vdwjidx0A = 2*vdwtype[jnrA+0];
831 vdwjidx0B = 2*vdwtype[jnrB+0];
832 vdwjidx0C = 2*vdwtype[jnrC+0];
833 vdwjidx0D = 2*vdwtype[jnrD+0];
834 vdwjidx0E = 2*vdwtype[jnrE+0];
835 vdwjidx0F = 2*vdwtype[jnrF+0];
836 vdwjidx0G = 2*vdwtype[jnrG+0];
837 vdwjidx0H = 2*vdwtype[jnrH+0];
839 /**************************
840 * CALCULATE INTERACTIONS *
841 **************************/
843 if (gmx_mm256_any_lt(rsq00,rcutoff2))
846 r00 = _mm256_mul_ps(rsq00,rinv00);
847 r00 = _mm256_andnot_ps(dummy_mask,r00);
849 /* Compute parameters for interactions between i and j atoms */
850 qq00 = _mm256_mul_ps(iq0,jq0);
851 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
852 vdwioffsetptr0+vdwjidx0B,
853 vdwioffsetptr0+vdwjidx0C,
854 vdwioffsetptr0+vdwjidx0D,
855 vdwioffsetptr0+vdwjidx0E,
856 vdwioffsetptr0+vdwjidx0F,
857 vdwioffsetptr0+vdwjidx0G,
858 vdwioffsetptr0+vdwjidx0H,
861 /* REACTION-FIELD ELECTROSTATICS */
862 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
864 /* LENNARD-JONES DISPERSION/REPULSION */
866 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
867 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
868 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
869 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
870 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
872 d = _mm256_sub_ps(r00,rswitch);
873 d = _mm256_max_ps(d,_mm256_setzero_ps());
874 d2 = _mm256_mul_ps(d,d);
875 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)))))));
877 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
879 /* Evaluate switch function */
880 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
881 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
882 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
884 fscal = _mm256_add_ps(felec,fvdw);
886 fscal = _mm256_and_ps(fscal,cutoff_mask);
888 fscal = _mm256_andnot_ps(dummy_mask,fscal);
890 /* Calculate temporary vectorial force */
891 tx = _mm256_mul_ps(fscal,dx00);
892 ty = _mm256_mul_ps(fscal,dy00);
893 tz = _mm256_mul_ps(fscal,dz00);
895 /* Update vectorial force */
896 fix0 = _mm256_add_ps(fix0,tx);
897 fiy0 = _mm256_add_ps(fiy0,ty);
898 fiz0 = _mm256_add_ps(fiz0,tz);
900 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
901 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
902 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
903 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
904 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
905 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
906 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
907 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
908 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
912 /* Inner loop uses 62 flops */
915 /* End of innermost loop */
917 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
918 f+i_coord_offset,fshift+i_shift_offset);
920 /* Increment number of inner iterations */
921 inneriter += j_index_end - j_index_start;
923 /* Outer loop uses 7 flops */
926 /* Increment number of outer iterations */
929 /* Update outer/inner flops */
931 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*62);