2 * Note: this file was generated by the Gromacs avx_256_single kernel generator.
4 * This source code is part of
8 * Copyright (c) 2001-2012, The GROMACS Development Team
10 * Gromacs is a library for molecular simulation and trajectory analysis,
11 * written by Erik Lindahl, David van der Spoel, Berk Hess, and others - for
12 * a full list of developers and information, check out http://www.gromacs.org
14 * This program is free software; you can redistribute it and/or modify it under
15 * the terms of the GNU Lesser General Public License as published by the Free
16 * Software Foundation; either version 2 of the License, or (at your option) any
19 * To help fund GROMACS development, we humbly ask that you cite
20 * the papers people have written on it - you can find them on the website.
28 #include "../nb_kernel.h"
29 #include "types/simple.h"
33 #include "gmx_math_x86_avx_256_single.h"
34 #include "kernelutil_x86_avx_256_single.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_VF_avx_256_single
38 * Electrostatics interaction: ReactionField
39 * VdW interaction: LennardJones
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_VF_avx_256_single
45 (t_nblist * gmx_restrict nlist,
46 rvec * gmx_restrict xx,
47 rvec * gmx_restrict ff,
48 t_forcerec * gmx_restrict fr,
49 t_mdatoms * gmx_restrict mdatoms,
50 nb_kernel_data_t * gmx_restrict kernel_data,
51 t_nrnb * gmx_restrict nrnb)
53 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
54 * just 0 for non-waters.
55 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
56 * jnr indices corresponding to data put in the four positions in the SIMD register.
58 int i_shift_offset,i_coord_offset,outeriter,inneriter;
59 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
60 int jnrA,jnrB,jnrC,jnrD;
61 int jnrE,jnrF,jnrG,jnrH;
62 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
63 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
64 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
65 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
66 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
68 real *shiftvec,*fshift,*x,*f;
69 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
71 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
72 real * vdwioffsetptr0;
73 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
74 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
75 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
76 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
77 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
80 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
83 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
84 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
85 __m256 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
86 real rswitch_scalar,d_scalar;
87 __m256 dummy_mask,cutoff_mask;
88 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
89 __m256 one = _mm256_set1_ps(1.0);
90 __m256 two = _mm256_set1_ps(2.0);
96 jindex = nlist->jindex;
98 shiftidx = nlist->shift;
100 shiftvec = fr->shift_vec[0];
101 fshift = fr->fshift[0];
102 facel = _mm256_set1_ps(fr->epsfac);
103 charge = mdatoms->chargeA;
104 krf = _mm256_set1_ps(fr->ic->k_rf);
105 krf2 = _mm256_set1_ps(fr->ic->k_rf*2.0);
106 crf = _mm256_set1_ps(fr->ic->c_rf);
107 nvdwtype = fr->ntype;
109 vdwtype = mdatoms->typeA;
111 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
112 rcutoff_scalar = fr->rcoulomb;
113 rcutoff = _mm256_set1_ps(rcutoff_scalar);
114 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
116 rswitch_scalar = fr->rvdw_switch;
117 rswitch = _mm256_set1_ps(rswitch_scalar);
118 /* Setup switch parameters */
119 d_scalar = rcutoff_scalar-rswitch_scalar;
120 d = _mm256_set1_ps(d_scalar);
121 swV3 = _mm256_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
122 swV4 = _mm256_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
123 swV5 = _mm256_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
124 swF2 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
125 swF3 = _mm256_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
126 swF4 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
128 /* Avoid stupid compiler warnings */
129 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
142 for(iidx=0;iidx<4*DIM;iidx++)
147 /* Start outer loop over neighborlists */
148 for(iidx=0; iidx<nri; iidx++)
150 /* Load shift vector for this list */
151 i_shift_offset = DIM*shiftidx[iidx];
153 /* Load limits for loop over neighbors */
154 j_index_start = jindex[iidx];
155 j_index_end = jindex[iidx+1];
157 /* Get outer coordinate index */
159 i_coord_offset = DIM*inr;
161 /* Load i particle coords and add shift vector */
162 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
164 fix0 = _mm256_setzero_ps();
165 fiy0 = _mm256_setzero_ps();
166 fiz0 = _mm256_setzero_ps();
168 /* Load parameters for i particles */
169 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
170 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
172 /* Reset potential sums */
173 velecsum = _mm256_setzero_ps();
174 vvdwsum = _mm256_setzero_ps();
176 /* Start inner kernel loop */
177 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
180 /* Get j neighbor index, and coordinate index */
189 j_coord_offsetA = DIM*jnrA;
190 j_coord_offsetB = DIM*jnrB;
191 j_coord_offsetC = DIM*jnrC;
192 j_coord_offsetD = DIM*jnrD;
193 j_coord_offsetE = DIM*jnrE;
194 j_coord_offsetF = DIM*jnrF;
195 j_coord_offsetG = DIM*jnrG;
196 j_coord_offsetH = DIM*jnrH;
198 /* load j atom coordinates */
199 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
200 x+j_coord_offsetC,x+j_coord_offsetD,
201 x+j_coord_offsetE,x+j_coord_offsetF,
202 x+j_coord_offsetG,x+j_coord_offsetH,
205 /* Calculate displacement vector */
206 dx00 = _mm256_sub_ps(ix0,jx0);
207 dy00 = _mm256_sub_ps(iy0,jy0);
208 dz00 = _mm256_sub_ps(iz0,jz0);
210 /* Calculate squared distance and things based on it */
211 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
213 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
215 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
217 /* Load parameters for j particles */
218 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
219 charge+jnrC+0,charge+jnrD+0,
220 charge+jnrE+0,charge+jnrF+0,
221 charge+jnrG+0,charge+jnrH+0);
222 vdwjidx0A = 2*vdwtype[jnrA+0];
223 vdwjidx0B = 2*vdwtype[jnrB+0];
224 vdwjidx0C = 2*vdwtype[jnrC+0];
225 vdwjidx0D = 2*vdwtype[jnrD+0];
226 vdwjidx0E = 2*vdwtype[jnrE+0];
227 vdwjidx0F = 2*vdwtype[jnrF+0];
228 vdwjidx0G = 2*vdwtype[jnrG+0];
229 vdwjidx0H = 2*vdwtype[jnrH+0];
231 /**************************
232 * CALCULATE INTERACTIONS *
233 **************************/
235 if (gmx_mm256_any_lt(rsq00,rcutoff2))
238 r00 = _mm256_mul_ps(rsq00,rinv00);
240 /* Compute parameters for interactions between i and j atoms */
241 qq00 = _mm256_mul_ps(iq0,jq0);
242 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
243 vdwioffsetptr0+vdwjidx0B,
244 vdwioffsetptr0+vdwjidx0C,
245 vdwioffsetptr0+vdwjidx0D,
246 vdwioffsetptr0+vdwjidx0E,
247 vdwioffsetptr0+vdwjidx0F,
248 vdwioffsetptr0+vdwjidx0G,
249 vdwioffsetptr0+vdwjidx0H,
252 /* REACTION-FIELD ELECTROSTATICS */
253 velec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
254 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
256 /* LENNARD-JONES DISPERSION/REPULSION */
258 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
259 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
260 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
261 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
262 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
264 d = _mm256_sub_ps(r00,rswitch);
265 d = _mm256_max_ps(d,_mm256_setzero_ps());
266 d2 = _mm256_mul_ps(d,d);
267 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)))))));
269 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
271 /* Evaluate switch function */
272 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
273 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
274 vvdw = _mm256_mul_ps(vvdw,sw);
275 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
277 /* Update potential sum for this i atom from the interaction with this j atom. */
278 velec = _mm256_and_ps(velec,cutoff_mask);
279 velecsum = _mm256_add_ps(velecsum,velec);
280 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
281 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
283 fscal = _mm256_add_ps(felec,fvdw);
285 fscal = _mm256_and_ps(fscal,cutoff_mask);
287 /* Calculate temporary vectorial force */
288 tx = _mm256_mul_ps(fscal,dx00);
289 ty = _mm256_mul_ps(fscal,dy00);
290 tz = _mm256_mul_ps(fscal,dz00);
292 /* Update vectorial force */
293 fix0 = _mm256_add_ps(fix0,tx);
294 fiy0 = _mm256_add_ps(fiy0,ty);
295 fiz0 = _mm256_add_ps(fiz0,tz);
297 fjptrA = f+j_coord_offsetA;
298 fjptrB = f+j_coord_offsetB;
299 fjptrC = f+j_coord_offsetC;
300 fjptrD = f+j_coord_offsetD;
301 fjptrE = f+j_coord_offsetE;
302 fjptrF = f+j_coord_offsetF;
303 fjptrG = f+j_coord_offsetG;
304 fjptrH = f+j_coord_offsetH;
305 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
309 /* Inner loop uses 70 flops */
315 /* Get j neighbor index, and coordinate index */
316 jnrlistA = jjnr[jidx];
317 jnrlistB = jjnr[jidx+1];
318 jnrlistC = jjnr[jidx+2];
319 jnrlistD = jjnr[jidx+3];
320 jnrlistE = jjnr[jidx+4];
321 jnrlistF = jjnr[jidx+5];
322 jnrlistG = jjnr[jidx+6];
323 jnrlistH = jjnr[jidx+7];
324 /* Sign of each element will be negative for non-real atoms.
325 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
326 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
328 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
329 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
331 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
332 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
333 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
334 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
335 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
336 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
337 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
338 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
339 j_coord_offsetA = DIM*jnrA;
340 j_coord_offsetB = DIM*jnrB;
341 j_coord_offsetC = DIM*jnrC;
342 j_coord_offsetD = DIM*jnrD;
343 j_coord_offsetE = DIM*jnrE;
344 j_coord_offsetF = DIM*jnrF;
345 j_coord_offsetG = DIM*jnrG;
346 j_coord_offsetH = DIM*jnrH;
348 /* load j atom coordinates */
349 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
350 x+j_coord_offsetC,x+j_coord_offsetD,
351 x+j_coord_offsetE,x+j_coord_offsetF,
352 x+j_coord_offsetG,x+j_coord_offsetH,
355 /* Calculate displacement vector */
356 dx00 = _mm256_sub_ps(ix0,jx0);
357 dy00 = _mm256_sub_ps(iy0,jy0);
358 dz00 = _mm256_sub_ps(iz0,jz0);
360 /* Calculate squared distance and things based on it */
361 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
363 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
365 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
367 /* Load parameters for j particles */
368 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
369 charge+jnrC+0,charge+jnrD+0,
370 charge+jnrE+0,charge+jnrF+0,
371 charge+jnrG+0,charge+jnrH+0);
372 vdwjidx0A = 2*vdwtype[jnrA+0];
373 vdwjidx0B = 2*vdwtype[jnrB+0];
374 vdwjidx0C = 2*vdwtype[jnrC+0];
375 vdwjidx0D = 2*vdwtype[jnrD+0];
376 vdwjidx0E = 2*vdwtype[jnrE+0];
377 vdwjidx0F = 2*vdwtype[jnrF+0];
378 vdwjidx0G = 2*vdwtype[jnrG+0];
379 vdwjidx0H = 2*vdwtype[jnrH+0];
381 /**************************
382 * CALCULATE INTERACTIONS *
383 **************************/
385 if (gmx_mm256_any_lt(rsq00,rcutoff2))
388 r00 = _mm256_mul_ps(rsq00,rinv00);
389 r00 = _mm256_andnot_ps(dummy_mask,r00);
391 /* Compute parameters for interactions between i and j atoms */
392 qq00 = _mm256_mul_ps(iq0,jq0);
393 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
394 vdwioffsetptr0+vdwjidx0B,
395 vdwioffsetptr0+vdwjidx0C,
396 vdwioffsetptr0+vdwjidx0D,
397 vdwioffsetptr0+vdwjidx0E,
398 vdwioffsetptr0+vdwjidx0F,
399 vdwioffsetptr0+vdwjidx0G,
400 vdwioffsetptr0+vdwjidx0H,
403 /* REACTION-FIELD ELECTROSTATICS */
404 velec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
405 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
407 /* LENNARD-JONES DISPERSION/REPULSION */
409 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
410 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
411 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
412 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
413 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
415 d = _mm256_sub_ps(r00,rswitch);
416 d = _mm256_max_ps(d,_mm256_setzero_ps());
417 d2 = _mm256_mul_ps(d,d);
418 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)))))));
420 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
422 /* Evaluate switch function */
423 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
424 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
425 vvdw = _mm256_mul_ps(vvdw,sw);
426 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
428 /* Update potential sum for this i atom from the interaction with this j atom. */
429 velec = _mm256_and_ps(velec,cutoff_mask);
430 velec = _mm256_andnot_ps(dummy_mask,velec);
431 velecsum = _mm256_add_ps(velecsum,velec);
432 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
433 vvdw = _mm256_andnot_ps(dummy_mask,vvdw);
434 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
436 fscal = _mm256_add_ps(felec,fvdw);
438 fscal = _mm256_and_ps(fscal,cutoff_mask);
440 fscal = _mm256_andnot_ps(dummy_mask,fscal);
442 /* Calculate temporary vectorial force */
443 tx = _mm256_mul_ps(fscal,dx00);
444 ty = _mm256_mul_ps(fscal,dy00);
445 tz = _mm256_mul_ps(fscal,dz00);
447 /* Update vectorial force */
448 fix0 = _mm256_add_ps(fix0,tx);
449 fiy0 = _mm256_add_ps(fiy0,ty);
450 fiz0 = _mm256_add_ps(fiz0,tz);
452 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
453 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
454 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
455 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
456 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
457 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
458 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
459 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
460 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
464 /* Inner loop uses 71 flops */
467 /* End of innermost loop */
469 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
470 f+i_coord_offset,fshift+i_shift_offset);
473 /* Update potential energies */
474 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
475 gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
477 /* Increment number of inner iterations */
478 inneriter += j_index_end - j_index_start;
480 /* Outer loop uses 9 flops */
483 /* Increment number of outer iterations */
486 /* Update outer/inner flops */
488 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*71);
491 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_F_avx_256_single
492 * Electrostatics interaction: ReactionField
493 * VdW interaction: LennardJones
494 * Geometry: Particle-Particle
495 * Calculate force/pot: Force
498 nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_F_avx_256_single
499 (t_nblist * gmx_restrict nlist,
500 rvec * gmx_restrict xx,
501 rvec * gmx_restrict ff,
502 t_forcerec * gmx_restrict fr,
503 t_mdatoms * gmx_restrict mdatoms,
504 nb_kernel_data_t * gmx_restrict kernel_data,
505 t_nrnb * gmx_restrict nrnb)
507 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
508 * just 0 for non-waters.
509 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
510 * jnr indices corresponding to data put in the four positions in the SIMD register.
512 int i_shift_offset,i_coord_offset,outeriter,inneriter;
513 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
514 int jnrA,jnrB,jnrC,jnrD;
515 int jnrE,jnrF,jnrG,jnrH;
516 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
517 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
518 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
519 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
520 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
522 real *shiftvec,*fshift,*x,*f;
523 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
525 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
526 real * vdwioffsetptr0;
527 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
528 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
529 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
530 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
531 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
534 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
537 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
538 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
539 __m256 rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
540 real rswitch_scalar,d_scalar;
541 __m256 dummy_mask,cutoff_mask;
542 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
543 __m256 one = _mm256_set1_ps(1.0);
544 __m256 two = _mm256_set1_ps(2.0);
550 jindex = nlist->jindex;
552 shiftidx = nlist->shift;
554 shiftvec = fr->shift_vec[0];
555 fshift = fr->fshift[0];
556 facel = _mm256_set1_ps(fr->epsfac);
557 charge = mdatoms->chargeA;
558 krf = _mm256_set1_ps(fr->ic->k_rf);
559 krf2 = _mm256_set1_ps(fr->ic->k_rf*2.0);
560 crf = _mm256_set1_ps(fr->ic->c_rf);
561 nvdwtype = fr->ntype;
563 vdwtype = mdatoms->typeA;
565 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
566 rcutoff_scalar = fr->rcoulomb;
567 rcutoff = _mm256_set1_ps(rcutoff_scalar);
568 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
570 rswitch_scalar = fr->rvdw_switch;
571 rswitch = _mm256_set1_ps(rswitch_scalar);
572 /* Setup switch parameters */
573 d_scalar = rcutoff_scalar-rswitch_scalar;
574 d = _mm256_set1_ps(d_scalar);
575 swV3 = _mm256_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
576 swV4 = _mm256_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
577 swV5 = _mm256_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
578 swF2 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
579 swF3 = _mm256_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
580 swF4 = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
582 /* Avoid stupid compiler warnings */
583 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
596 for(iidx=0;iidx<4*DIM;iidx++)
601 /* Start outer loop over neighborlists */
602 for(iidx=0; iidx<nri; iidx++)
604 /* Load shift vector for this list */
605 i_shift_offset = DIM*shiftidx[iidx];
607 /* Load limits for loop over neighbors */
608 j_index_start = jindex[iidx];
609 j_index_end = jindex[iidx+1];
611 /* Get outer coordinate index */
613 i_coord_offset = DIM*inr;
615 /* Load i particle coords and add shift vector */
616 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
618 fix0 = _mm256_setzero_ps();
619 fiy0 = _mm256_setzero_ps();
620 fiz0 = _mm256_setzero_ps();
622 /* Load parameters for i particles */
623 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
624 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
626 /* Start inner kernel loop */
627 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
630 /* Get j neighbor index, and coordinate index */
639 j_coord_offsetA = DIM*jnrA;
640 j_coord_offsetB = DIM*jnrB;
641 j_coord_offsetC = DIM*jnrC;
642 j_coord_offsetD = DIM*jnrD;
643 j_coord_offsetE = DIM*jnrE;
644 j_coord_offsetF = DIM*jnrF;
645 j_coord_offsetG = DIM*jnrG;
646 j_coord_offsetH = DIM*jnrH;
648 /* load j atom coordinates */
649 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
650 x+j_coord_offsetC,x+j_coord_offsetD,
651 x+j_coord_offsetE,x+j_coord_offsetF,
652 x+j_coord_offsetG,x+j_coord_offsetH,
655 /* Calculate displacement vector */
656 dx00 = _mm256_sub_ps(ix0,jx0);
657 dy00 = _mm256_sub_ps(iy0,jy0);
658 dz00 = _mm256_sub_ps(iz0,jz0);
660 /* Calculate squared distance and things based on it */
661 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
663 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
665 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
667 /* Load parameters for j particles */
668 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
669 charge+jnrC+0,charge+jnrD+0,
670 charge+jnrE+0,charge+jnrF+0,
671 charge+jnrG+0,charge+jnrH+0);
672 vdwjidx0A = 2*vdwtype[jnrA+0];
673 vdwjidx0B = 2*vdwtype[jnrB+0];
674 vdwjidx0C = 2*vdwtype[jnrC+0];
675 vdwjidx0D = 2*vdwtype[jnrD+0];
676 vdwjidx0E = 2*vdwtype[jnrE+0];
677 vdwjidx0F = 2*vdwtype[jnrF+0];
678 vdwjidx0G = 2*vdwtype[jnrG+0];
679 vdwjidx0H = 2*vdwtype[jnrH+0];
681 /**************************
682 * CALCULATE INTERACTIONS *
683 **************************/
685 if (gmx_mm256_any_lt(rsq00,rcutoff2))
688 r00 = _mm256_mul_ps(rsq00,rinv00);
690 /* Compute parameters for interactions between i and j atoms */
691 qq00 = _mm256_mul_ps(iq0,jq0);
692 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
693 vdwioffsetptr0+vdwjidx0B,
694 vdwioffsetptr0+vdwjidx0C,
695 vdwioffsetptr0+vdwjidx0D,
696 vdwioffsetptr0+vdwjidx0E,
697 vdwioffsetptr0+vdwjidx0F,
698 vdwioffsetptr0+vdwjidx0G,
699 vdwioffsetptr0+vdwjidx0H,
702 /* REACTION-FIELD ELECTROSTATICS */
703 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
705 /* LENNARD-JONES DISPERSION/REPULSION */
707 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
708 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
709 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
710 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
711 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
713 d = _mm256_sub_ps(r00,rswitch);
714 d = _mm256_max_ps(d,_mm256_setzero_ps());
715 d2 = _mm256_mul_ps(d,d);
716 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)))))));
718 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
720 /* Evaluate switch function */
721 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
722 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
723 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
725 fscal = _mm256_add_ps(felec,fvdw);
727 fscal = _mm256_and_ps(fscal,cutoff_mask);
729 /* Calculate temporary vectorial force */
730 tx = _mm256_mul_ps(fscal,dx00);
731 ty = _mm256_mul_ps(fscal,dy00);
732 tz = _mm256_mul_ps(fscal,dz00);
734 /* Update vectorial force */
735 fix0 = _mm256_add_ps(fix0,tx);
736 fiy0 = _mm256_add_ps(fiy0,ty);
737 fiz0 = _mm256_add_ps(fiz0,tz);
739 fjptrA = f+j_coord_offsetA;
740 fjptrB = f+j_coord_offsetB;
741 fjptrC = f+j_coord_offsetC;
742 fjptrD = f+j_coord_offsetD;
743 fjptrE = f+j_coord_offsetE;
744 fjptrF = f+j_coord_offsetF;
745 fjptrG = f+j_coord_offsetG;
746 fjptrH = f+j_coord_offsetH;
747 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
751 /* Inner loop uses 61 flops */
757 /* Get j neighbor index, and coordinate index */
758 jnrlistA = jjnr[jidx];
759 jnrlistB = jjnr[jidx+1];
760 jnrlistC = jjnr[jidx+2];
761 jnrlistD = jjnr[jidx+3];
762 jnrlistE = jjnr[jidx+4];
763 jnrlistF = jjnr[jidx+5];
764 jnrlistG = jjnr[jidx+6];
765 jnrlistH = jjnr[jidx+7];
766 /* Sign of each element will be negative for non-real atoms.
767 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
768 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
770 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
771 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
773 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
774 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
775 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
776 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
777 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
778 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
779 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
780 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
781 j_coord_offsetA = DIM*jnrA;
782 j_coord_offsetB = DIM*jnrB;
783 j_coord_offsetC = DIM*jnrC;
784 j_coord_offsetD = DIM*jnrD;
785 j_coord_offsetE = DIM*jnrE;
786 j_coord_offsetF = DIM*jnrF;
787 j_coord_offsetG = DIM*jnrG;
788 j_coord_offsetH = DIM*jnrH;
790 /* load j atom coordinates */
791 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
792 x+j_coord_offsetC,x+j_coord_offsetD,
793 x+j_coord_offsetE,x+j_coord_offsetF,
794 x+j_coord_offsetG,x+j_coord_offsetH,
797 /* Calculate displacement vector */
798 dx00 = _mm256_sub_ps(ix0,jx0);
799 dy00 = _mm256_sub_ps(iy0,jy0);
800 dz00 = _mm256_sub_ps(iz0,jz0);
802 /* Calculate squared distance and things based on it */
803 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
805 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
807 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
809 /* Load parameters for j particles */
810 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
811 charge+jnrC+0,charge+jnrD+0,
812 charge+jnrE+0,charge+jnrF+0,
813 charge+jnrG+0,charge+jnrH+0);
814 vdwjidx0A = 2*vdwtype[jnrA+0];
815 vdwjidx0B = 2*vdwtype[jnrB+0];
816 vdwjidx0C = 2*vdwtype[jnrC+0];
817 vdwjidx0D = 2*vdwtype[jnrD+0];
818 vdwjidx0E = 2*vdwtype[jnrE+0];
819 vdwjidx0F = 2*vdwtype[jnrF+0];
820 vdwjidx0G = 2*vdwtype[jnrG+0];
821 vdwjidx0H = 2*vdwtype[jnrH+0];
823 /**************************
824 * CALCULATE INTERACTIONS *
825 **************************/
827 if (gmx_mm256_any_lt(rsq00,rcutoff2))
830 r00 = _mm256_mul_ps(rsq00,rinv00);
831 r00 = _mm256_andnot_ps(dummy_mask,r00);
833 /* Compute parameters for interactions between i and j atoms */
834 qq00 = _mm256_mul_ps(iq0,jq0);
835 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
836 vdwioffsetptr0+vdwjidx0B,
837 vdwioffsetptr0+vdwjidx0C,
838 vdwioffsetptr0+vdwjidx0D,
839 vdwioffsetptr0+vdwjidx0E,
840 vdwioffsetptr0+vdwjidx0F,
841 vdwioffsetptr0+vdwjidx0G,
842 vdwioffsetptr0+vdwjidx0H,
845 /* REACTION-FIELD ELECTROSTATICS */
846 felec = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
848 /* LENNARD-JONES DISPERSION/REPULSION */
850 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
851 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
852 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
853 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
854 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
856 d = _mm256_sub_ps(r00,rswitch);
857 d = _mm256_max_ps(d,_mm256_setzero_ps());
858 d2 = _mm256_mul_ps(d,d);
859 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)))))));
861 dsw = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
863 /* Evaluate switch function */
864 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
865 fvdw = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv00,_mm256_mul_ps(vvdw,dsw)) );
866 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
868 fscal = _mm256_add_ps(felec,fvdw);
870 fscal = _mm256_and_ps(fscal,cutoff_mask);
872 fscal = _mm256_andnot_ps(dummy_mask,fscal);
874 /* Calculate temporary vectorial force */
875 tx = _mm256_mul_ps(fscal,dx00);
876 ty = _mm256_mul_ps(fscal,dy00);
877 tz = _mm256_mul_ps(fscal,dz00);
879 /* Update vectorial force */
880 fix0 = _mm256_add_ps(fix0,tx);
881 fiy0 = _mm256_add_ps(fiy0,ty);
882 fiz0 = _mm256_add_ps(fiz0,tz);
884 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
885 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
886 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
887 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
888 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
889 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
890 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
891 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
892 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
896 /* Inner loop uses 62 flops */
899 /* End of innermost loop */
901 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
902 f+i_coord_offset,fshift+i_shift_offset);
904 /* Increment number of inner iterations */
905 inneriter += j_index_end - j_index_start;
907 /* Outer loop uses 7 flops */
910 /* Increment number of outer iterations */
913 /* Update outer/inner flops */
915 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*62);