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_ElecCSTab_VdwLJ_GeomP1P1_VF_avx_256_single
38 * Electrostatics interaction: CubicSplineTable
39 * VdW interaction: LennardJones
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecCSTab_VdwLJ_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);
86 __m128i vfitab_lo,vfitab_hi;
87 __m128i ifour = _mm_set1_epi32(4);
88 __m256 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
90 __m256 dummy_mask,cutoff_mask;
91 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
92 __m256 one = _mm256_set1_ps(1.0);
93 __m256 two = _mm256_set1_ps(2.0);
99 jindex = nlist->jindex;
101 shiftidx = nlist->shift;
103 shiftvec = fr->shift_vec[0];
104 fshift = fr->fshift[0];
105 facel = _mm256_set1_ps(fr->epsfac);
106 charge = mdatoms->chargeA;
107 nvdwtype = fr->ntype;
109 vdwtype = mdatoms->typeA;
111 vftab = kernel_data->table_elec->data;
112 vftabscale = _mm256_set1_ps(kernel_data->table_elec->scale);
114 /* Avoid stupid compiler warnings */
115 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
128 for(iidx=0;iidx<4*DIM;iidx++)
133 /* Start outer loop over neighborlists */
134 for(iidx=0; iidx<nri; iidx++)
136 /* Load shift vector for this list */
137 i_shift_offset = DIM*shiftidx[iidx];
139 /* Load limits for loop over neighbors */
140 j_index_start = jindex[iidx];
141 j_index_end = jindex[iidx+1];
143 /* Get outer coordinate index */
145 i_coord_offset = DIM*inr;
147 /* Load i particle coords and add shift vector */
148 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
150 fix0 = _mm256_setzero_ps();
151 fiy0 = _mm256_setzero_ps();
152 fiz0 = _mm256_setzero_ps();
154 /* Load parameters for i particles */
155 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
156 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
158 /* Reset potential sums */
159 velecsum = _mm256_setzero_ps();
160 vvdwsum = _mm256_setzero_ps();
162 /* Start inner kernel loop */
163 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
166 /* Get j neighbor index, and coordinate index */
175 j_coord_offsetA = DIM*jnrA;
176 j_coord_offsetB = DIM*jnrB;
177 j_coord_offsetC = DIM*jnrC;
178 j_coord_offsetD = DIM*jnrD;
179 j_coord_offsetE = DIM*jnrE;
180 j_coord_offsetF = DIM*jnrF;
181 j_coord_offsetG = DIM*jnrG;
182 j_coord_offsetH = DIM*jnrH;
184 /* load j atom coordinates */
185 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
186 x+j_coord_offsetC,x+j_coord_offsetD,
187 x+j_coord_offsetE,x+j_coord_offsetF,
188 x+j_coord_offsetG,x+j_coord_offsetH,
191 /* Calculate displacement vector */
192 dx00 = _mm256_sub_ps(ix0,jx0);
193 dy00 = _mm256_sub_ps(iy0,jy0);
194 dz00 = _mm256_sub_ps(iz0,jz0);
196 /* Calculate squared distance and things based on it */
197 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
199 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
201 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
203 /* Load parameters for j particles */
204 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
205 charge+jnrC+0,charge+jnrD+0,
206 charge+jnrE+0,charge+jnrF+0,
207 charge+jnrG+0,charge+jnrH+0);
208 vdwjidx0A = 2*vdwtype[jnrA+0];
209 vdwjidx0B = 2*vdwtype[jnrB+0];
210 vdwjidx0C = 2*vdwtype[jnrC+0];
211 vdwjidx0D = 2*vdwtype[jnrD+0];
212 vdwjidx0E = 2*vdwtype[jnrE+0];
213 vdwjidx0F = 2*vdwtype[jnrF+0];
214 vdwjidx0G = 2*vdwtype[jnrG+0];
215 vdwjidx0H = 2*vdwtype[jnrH+0];
217 /**************************
218 * CALCULATE INTERACTIONS *
219 **************************/
221 r00 = _mm256_mul_ps(rsq00,rinv00);
223 /* Compute parameters for interactions between i and j atoms */
224 qq00 = _mm256_mul_ps(iq0,jq0);
225 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
226 vdwioffsetptr0+vdwjidx0B,
227 vdwioffsetptr0+vdwjidx0C,
228 vdwioffsetptr0+vdwjidx0D,
229 vdwioffsetptr0+vdwjidx0E,
230 vdwioffsetptr0+vdwjidx0F,
231 vdwioffsetptr0+vdwjidx0G,
232 vdwioffsetptr0+vdwjidx0H,
235 /* Calculate table index by multiplying r with table scale and truncate to integer */
236 rt = _mm256_mul_ps(r00,vftabscale);
237 vfitab = _mm256_cvttps_epi32(rt);
238 vfeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
239 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
240 vfitab_lo = _mm256_extractf128_si256(vfitab,0x0);
241 vfitab_hi = _mm256_extractf128_si256(vfitab,0x1);
242 vfitab_lo = _mm_slli_epi32(vfitab_lo,2);
243 vfitab_hi = _mm_slli_epi32(vfitab_hi,2);
245 /* CUBIC SPLINE TABLE ELECTROSTATICS */
246 Y = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
247 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
248 F = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
249 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
250 G = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
251 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
252 H = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
253 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
254 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
255 Heps = _mm256_mul_ps(vfeps,H);
256 Fp = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
257 VV = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
258 velec = _mm256_mul_ps(qq00,VV);
259 FF = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
260 felec = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_mul_ps(qq00,FF),_mm256_mul_ps(vftabscale,rinv00)));
262 /* LENNARD-JONES DISPERSION/REPULSION */
264 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
265 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
266 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
267 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
268 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
270 /* Update potential sum for this i atom from the interaction with this j atom. */
271 velecsum = _mm256_add_ps(velecsum,velec);
272 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
274 fscal = _mm256_add_ps(felec,fvdw);
276 /* Calculate temporary vectorial force */
277 tx = _mm256_mul_ps(fscal,dx00);
278 ty = _mm256_mul_ps(fscal,dy00);
279 tz = _mm256_mul_ps(fscal,dz00);
281 /* Update vectorial force */
282 fix0 = _mm256_add_ps(fix0,tx);
283 fiy0 = _mm256_add_ps(fiy0,ty);
284 fiz0 = _mm256_add_ps(fiz0,tz);
286 fjptrA = f+j_coord_offsetA;
287 fjptrB = f+j_coord_offsetB;
288 fjptrC = f+j_coord_offsetC;
289 fjptrD = f+j_coord_offsetD;
290 fjptrE = f+j_coord_offsetE;
291 fjptrF = f+j_coord_offsetF;
292 fjptrG = f+j_coord_offsetG;
293 fjptrH = f+j_coord_offsetH;
294 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
296 /* Inner loop uses 56 flops */
302 /* Get j neighbor index, and coordinate index */
303 jnrlistA = jjnr[jidx];
304 jnrlistB = jjnr[jidx+1];
305 jnrlistC = jjnr[jidx+2];
306 jnrlistD = jjnr[jidx+3];
307 jnrlistE = jjnr[jidx+4];
308 jnrlistF = jjnr[jidx+5];
309 jnrlistG = jjnr[jidx+6];
310 jnrlistH = jjnr[jidx+7];
311 /* Sign of each element will be negative for non-real atoms.
312 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
313 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
315 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
316 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
318 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
319 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
320 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
321 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
322 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
323 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
324 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
325 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
326 j_coord_offsetA = DIM*jnrA;
327 j_coord_offsetB = DIM*jnrB;
328 j_coord_offsetC = DIM*jnrC;
329 j_coord_offsetD = DIM*jnrD;
330 j_coord_offsetE = DIM*jnrE;
331 j_coord_offsetF = DIM*jnrF;
332 j_coord_offsetG = DIM*jnrG;
333 j_coord_offsetH = DIM*jnrH;
335 /* load j atom coordinates */
336 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
337 x+j_coord_offsetC,x+j_coord_offsetD,
338 x+j_coord_offsetE,x+j_coord_offsetF,
339 x+j_coord_offsetG,x+j_coord_offsetH,
342 /* Calculate displacement vector */
343 dx00 = _mm256_sub_ps(ix0,jx0);
344 dy00 = _mm256_sub_ps(iy0,jy0);
345 dz00 = _mm256_sub_ps(iz0,jz0);
347 /* Calculate squared distance and things based on it */
348 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
350 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
352 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
354 /* Load parameters for j particles */
355 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
356 charge+jnrC+0,charge+jnrD+0,
357 charge+jnrE+0,charge+jnrF+0,
358 charge+jnrG+0,charge+jnrH+0);
359 vdwjidx0A = 2*vdwtype[jnrA+0];
360 vdwjidx0B = 2*vdwtype[jnrB+0];
361 vdwjidx0C = 2*vdwtype[jnrC+0];
362 vdwjidx0D = 2*vdwtype[jnrD+0];
363 vdwjidx0E = 2*vdwtype[jnrE+0];
364 vdwjidx0F = 2*vdwtype[jnrF+0];
365 vdwjidx0G = 2*vdwtype[jnrG+0];
366 vdwjidx0H = 2*vdwtype[jnrH+0];
368 /**************************
369 * CALCULATE INTERACTIONS *
370 **************************/
372 r00 = _mm256_mul_ps(rsq00,rinv00);
373 r00 = _mm256_andnot_ps(dummy_mask,r00);
375 /* Compute parameters for interactions between i and j atoms */
376 qq00 = _mm256_mul_ps(iq0,jq0);
377 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
378 vdwioffsetptr0+vdwjidx0B,
379 vdwioffsetptr0+vdwjidx0C,
380 vdwioffsetptr0+vdwjidx0D,
381 vdwioffsetptr0+vdwjidx0E,
382 vdwioffsetptr0+vdwjidx0F,
383 vdwioffsetptr0+vdwjidx0G,
384 vdwioffsetptr0+vdwjidx0H,
387 /* Calculate table index by multiplying r with table scale and truncate to integer */
388 rt = _mm256_mul_ps(r00,vftabscale);
389 vfitab = _mm256_cvttps_epi32(rt);
390 vfeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
391 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
392 vfitab_lo = _mm256_extractf128_si256(vfitab,0x0);
393 vfitab_hi = _mm256_extractf128_si256(vfitab,0x1);
394 vfitab_lo = _mm_slli_epi32(vfitab_lo,2);
395 vfitab_hi = _mm_slli_epi32(vfitab_hi,2);
397 /* CUBIC SPLINE TABLE ELECTROSTATICS */
398 Y = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
399 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
400 F = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
401 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
402 G = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
403 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
404 H = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
405 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
406 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
407 Heps = _mm256_mul_ps(vfeps,H);
408 Fp = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
409 VV = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
410 velec = _mm256_mul_ps(qq00,VV);
411 FF = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
412 felec = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_mul_ps(qq00,FF),_mm256_mul_ps(vftabscale,rinv00)));
414 /* LENNARD-JONES DISPERSION/REPULSION */
416 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
417 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
418 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
419 vvdw = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
420 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
422 /* Update potential sum for this i atom from the interaction with this j atom. */
423 velec = _mm256_andnot_ps(dummy_mask,velec);
424 velecsum = _mm256_add_ps(velecsum,velec);
425 vvdw = _mm256_andnot_ps(dummy_mask,vvdw);
426 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
428 fscal = _mm256_add_ps(felec,fvdw);
430 fscal = _mm256_andnot_ps(dummy_mask,fscal);
432 /* Calculate temporary vectorial force */
433 tx = _mm256_mul_ps(fscal,dx00);
434 ty = _mm256_mul_ps(fscal,dy00);
435 tz = _mm256_mul_ps(fscal,dz00);
437 /* Update vectorial force */
438 fix0 = _mm256_add_ps(fix0,tx);
439 fiy0 = _mm256_add_ps(fiy0,ty);
440 fiz0 = _mm256_add_ps(fiz0,tz);
442 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
443 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
444 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
445 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
446 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
447 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
448 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
449 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
450 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
452 /* Inner loop uses 57 flops */
455 /* End of innermost loop */
457 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
458 f+i_coord_offset,fshift+i_shift_offset);
461 /* Update potential energies */
462 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
463 gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
465 /* Increment number of inner iterations */
466 inneriter += j_index_end - j_index_start;
468 /* Outer loop uses 9 flops */
471 /* Increment number of outer iterations */
474 /* Update outer/inner flops */
476 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*57);
479 * Gromacs nonbonded kernel: nb_kernel_ElecCSTab_VdwLJ_GeomP1P1_F_avx_256_single
480 * Electrostatics interaction: CubicSplineTable
481 * VdW interaction: LennardJones
482 * Geometry: Particle-Particle
483 * Calculate force/pot: Force
486 nb_kernel_ElecCSTab_VdwLJ_GeomP1P1_F_avx_256_single
487 (t_nblist * gmx_restrict nlist,
488 rvec * gmx_restrict xx,
489 rvec * gmx_restrict ff,
490 t_forcerec * gmx_restrict fr,
491 t_mdatoms * gmx_restrict mdatoms,
492 nb_kernel_data_t * gmx_restrict kernel_data,
493 t_nrnb * gmx_restrict nrnb)
495 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
496 * just 0 for non-waters.
497 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
498 * jnr indices corresponding to data put in the four positions in the SIMD register.
500 int i_shift_offset,i_coord_offset,outeriter,inneriter;
501 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
502 int jnrA,jnrB,jnrC,jnrD;
503 int jnrE,jnrF,jnrG,jnrH;
504 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
505 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
506 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
507 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
508 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
510 real *shiftvec,*fshift,*x,*f;
511 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
513 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
514 real * vdwioffsetptr0;
515 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
516 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
517 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
518 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
519 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
522 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
525 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
526 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
528 __m128i vfitab_lo,vfitab_hi;
529 __m128i ifour = _mm_set1_epi32(4);
530 __m256 rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
532 __m256 dummy_mask,cutoff_mask;
533 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
534 __m256 one = _mm256_set1_ps(1.0);
535 __m256 two = _mm256_set1_ps(2.0);
541 jindex = nlist->jindex;
543 shiftidx = nlist->shift;
545 shiftvec = fr->shift_vec[0];
546 fshift = fr->fshift[0];
547 facel = _mm256_set1_ps(fr->epsfac);
548 charge = mdatoms->chargeA;
549 nvdwtype = fr->ntype;
551 vdwtype = mdatoms->typeA;
553 vftab = kernel_data->table_elec->data;
554 vftabscale = _mm256_set1_ps(kernel_data->table_elec->scale);
556 /* Avoid stupid compiler warnings */
557 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
570 for(iidx=0;iidx<4*DIM;iidx++)
575 /* Start outer loop over neighborlists */
576 for(iidx=0; iidx<nri; iidx++)
578 /* Load shift vector for this list */
579 i_shift_offset = DIM*shiftidx[iidx];
581 /* Load limits for loop over neighbors */
582 j_index_start = jindex[iidx];
583 j_index_end = jindex[iidx+1];
585 /* Get outer coordinate index */
587 i_coord_offset = DIM*inr;
589 /* Load i particle coords and add shift vector */
590 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
592 fix0 = _mm256_setzero_ps();
593 fiy0 = _mm256_setzero_ps();
594 fiz0 = _mm256_setzero_ps();
596 /* Load parameters for i particles */
597 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
598 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
600 /* Start inner kernel loop */
601 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
604 /* Get j neighbor index, and coordinate index */
613 j_coord_offsetA = DIM*jnrA;
614 j_coord_offsetB = DIM*jnrB;
615 j_coord_offsetC = DIM*jnrC;
616 j_coord_offsetD = DIM*jnrD;
617 j_coord_offsetE = DIM*jnrE;
618 j_coord_offsetF = DIM*jnrF;
619 j_coord_offsetG = DIM*jnrG;
620 j_coord_offsetH = DIM*jnrH;
622 /* load j atom coordinates */
623 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
624 x+j_coord_offsetC,x+j_coord_offsetD,
625 x+j_coord_offsetE,x+j_coord_offsetF,
626 x+j_coord_offsetG,x+j_coord_offsetH,
629 /* Calculate displacement vector */
630 dx00 = _mm256_sub_ps(ix0,jx0);
631 dy00 = _mm256_sub_ps(iy0,jy0);
632 dz00 = _mm256_sub_ps(iz0,jz0);
634 /* Calculate squared distance and things based on it */
635 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
637 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
639 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
641 /* Load parameters for j particles */
642 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
643 charge+jnrC+0,charge+jnrD+0,
644 charge+jnrE+0,charge+jnrF+0,
645 charge+jnrG+0,charge+jnrH+0);
646 vdwjidx0A = 2*vdwtype[jnrA+0];
647 vdwjidx0B = 2*vdwtype[jnrB+0];
648 vdwjidx0C = 2*vdwtype[jnrC+0];
649 vdwjidx0D = 2*vdwtype[jnrD+0];
650 vdwjidx0E = 2*vdwtype[jnrE+0];
651 vdwjidx0F = 2*vdwtype[jnrF+0];
652 vdwjidx0G = 2*vdwtype[jnrG+0];
653 vdwjidx0H = 2*vdwtype[jnrH+0];
655 /**************************
656 * CALCULATE INTERACTIONS *
657 **************************/
659 r00 = _mm256_mul_ps(rsq00,rinv00);
661 /* Compute parameters for interactions between i and j atoms */
662 qq00 = _mm256_mul_ps(iq0,jq0);
663 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
664 vdwioffsetptr0+vdwjidx0B,
665 vdwioffsetptr0+vdwjidx0C,
666 vdwioffsetptr0+vdwjidx0D,
667 vdwioffsetptr0+vdwjidx0E,
668 vdwioffsetptr0+vdwjidx0F,
669 vdwioffsetptr0+vdwjidx0G,
670 vdwioffsetptr0+vdwjidx0H,
673 /* Calculate table index by multiplying r with table scale and truncate to integer */
674 rt = _mm256_mul_ps(r00,vftabscale);
675 vfitab = _mm256_cvttps_epi32(rt);
676 vfeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
677 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
678 vfitab_lo = _mm256_extractf128_si256(vfitab,0x0);
679 vfitab_hi = _mm256_extractf128_si256(vfitab,0x1);
680 vfitab_lo = _mm_slli_epi32(vfitab_lo,2);
681 vfitab_hi = _mm_slli_epi32(vfitab_hi,2);
683 /* CUBIC SPLINE TABLE ELECTROSTATICS */
684 Y = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
685 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
686 F = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
687 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
688 G = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
689 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
690 H = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
691 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
692 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
693 Heps = _mm256_mul_ps(vfeps,H);
694 Fp = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
695 FF = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
696 felec = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_mul_ps(qq00,FF),_mm256_mul_ps(vftabscale,rinv00)));
698 /* LENNARD-JONES DISPERSION/REPULSION */
700 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
701 fvdw = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_00,rinvsix),c6_00),_mm256_mul_ps(rinvsix,rinvsq00));
703 fscal = _mm256_add_ps(felec,fvdw);
705 /* Calculate temporary vectorial force */
706 tx = _mm256_mul_ps(fscal,dx00);
707 ty = _mm256_mul_ps(fscal,dy00);
708 tz = _mm256_mul_ps(fscal,dz00);
710 /* Update vectorial force */
711 fix0 = _mm256_add_ps(fix0,tx);
712 fiy0 = _mm256_add_ps(fiy0,ty);
713 fiz0 = _mm256_add_ps(fiz0,tz);
715 fjptrA = f+j_coord_offsetA;
716 fjptrB = f+j_coord_offsetB;
717 fjptrC = f+j_coord_offsetC;
718 fjptrD = f+j_coord_offsetD;
719 fjptrE = f+j_coord_offsetE;
720 fjptrF = f+j_coord_offsetF;
721 fjptrG = f+j_coord_offsetG;
722 fjptrH = f+j_coord_offsetH;
723 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
725 /* Inner loop uses 47 flops */
731 /* Get j neighbor index, and coordinate index */
732 jnrlistA = jjnr[jidx];
733 jnrlistB = jjnr[jidx+1];
734 jnrlistC = jjnr[jidx+2];
735 jnrlistD = jjnr[jidx+3];
736 jnrlistE = jjnr[jidx+4];
737 jnrlistF = jjnr[jidx+5];
738 jnrlistG = jjnr[jidx+6];
739 jnrlistH = jjnr[jidx+7];
740 /* Sign of each element will be negative for non-real atoms.
741 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
742 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
744 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
745 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
747 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
748 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
749 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
750 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
751 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
752 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
753 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
754 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
755 j_coord_offsetA = DIM*jnrA;
756 j_coord_offsetB = DIM*jnrB;
757 j_coord_offsetC = DIM*jnrC;
758 j_coord_offsetD = DIM*jnrD;
759 j_coord_offsetE = DIM*jnrE;
760 j_coord_offsetF = DIM*jnrF;
761 j_coord_offsetG = DIM*jnrG;
762 j_coord_offsetH = DIM*jnrH;
764 /* load j atom coordinates */
765 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
766 x+j_coord_offsetC,x+j_coord_offsetD,
767 x+j_coord_offsetE,x+j_coord_offsetF,
768 x+j_coord_offsetG,x+j_coord_offsetH,
771 /* Calculate displacement vector */
772 dx00 = _mm256_sub_ps(ix0,jx0);
773 dy00 = _mm256_sub_ps(iy0,jy0);
774 dz00 = _mm256_sub_ps(iz0,jz0);
776 /* Calculate squared distance and things based on it */
777 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
779 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
781 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
783 /* Load parameters for j particles */
784 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
785 charge+jnrC+0,charge+jnrD+0,
786 charge+jnrE+0,charge+jnrF+0,
787 charge+jnrG+0,charge+jnrH+0);
788 vdwjidx0A = 2*vdwtype[jnrA+0];
789 vdwjidx0B = 2*vdwtype[jnrB+0];
790 vdwjidx0C = 2*vdwtype[jnrC+0];
791 vdwjidx0D = 2*vdwtype[jnrD+0];
792 vdwjidx0E = 2*vdwtype[jnrE+0];
793 vdwjidx0F = 2*vdwtype[jnrF+0];
794 vdwjidx0G = 2*vdwtype[jnrG+0];
795 vdwjidx0H = 2*vdwtype[jnrH+0];
797 /**************************
798 * CALCULATE INTERACTIONS *
799 **************************/
801 r00 = _mm256_mul_ps(rsq00,rinv00);
802 r00 = _mm256_andnot_ps(dummy_mask,r00);
804 /* Compute parameters for interactions between i and j atoms */
805 qq00 = _mm256_mul_ps(iq0,jq0);
806 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
807 vdwioffsetptr0+vdwjidx0B,
808 vdwioffsetptr0+vdwjidx0C,
809 vdwioffsetptr0+vdwjidx0D,
810 vdwioffsetptr0+vdwjidx0E,
811 vdwioffsetptr0+vdwjidx0F,
812 vdwioffsetptr0+vdwjidx0G,
813 vdwioffsetptr0+vdwjidx0H,
816 /* Calculate table index by multiplying r with table scale and truncate to integer */
817 rt = _mm256_mul_ps(r00,vftabscale);
818 vfitab = _mm256_cvttps_epi32(rt);
819 vfeps = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
820 /* AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
821 vfitab_lo = _mm256_extractf128_si256(vfitab,0x0);
822 vfitab_hi = _mm256_extractf128_si256(vfitab,0x1);
823 vfitab_lo = _mm_slli_epi32(vfitab_lo,2);
824 vfitab_hi = _mm_slli_epi32(vfitab_hi,2);
826 /* CUBIC SPLINE TABLE ELECTROSTATICS */
827 Y = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
828 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
829 F = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
830 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
831 G = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
832 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
833 H = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
834 _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
835 GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
836 Heps = _mm256_mul_ps(vfeps,H);
837 Fp = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
838 FF = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
839 felec = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_mul_ps(qq00,FF),_mm256_mul_ps(vftabscale,rinv00)));
841 /* LENNARD-JONES DISPERSION/REPULSION */
843 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
844 fvdw = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_00,rinvsix),c6_00),_mm256_mul_ps(rinvsix,rinvsq00));
846 fscal = _mm256_add_ps(felec,fvdw);
848 fscal = _mm256_andnot_ps(dummy_mask,fscal);
850 /* Calculate temporary vectorial force */
851 tx = _mm256_mul_ps(fscal,dx00);
852 ty = _mm256_mul_ps(fscal,dy00);
853 tz = _mm256_mul_ps(fscal,dz00);
855 /* Update vectorial force */
856 fix0 = _mm256_add_ps(fix0,tx);
857 fiy0 = _mm256_add_ps(fiy0,ty);
858 fiz0 = _mm256_add_ps(fiz0,tz);
860 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
861 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
862 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
863 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
864 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
865 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
866 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
867 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
868 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
870 /* Inner loop uses 48 flops */
873 /* End of innermost loop */
875 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
876 f+i_coord_offset,fshift+i_shift_offset);
878 /* Increment number of inner iterations */
879 inneriter += j_index_end - j_index_start;
881 /* Outer loop uses 7 flops */
884 /* Increment number of outer iterations */
887 /* Update outer/inner flops */
889 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*48);