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_ElecEwSh_VdwLJSh_GeomP1P1_VF_avx_256_single
38 * Electrostatics interaction: Ewald
39 * VdW interaction: LennardJones
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecEwSh_VdwLJSh_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 ewitab_lo,ewitab_hi;
87 __m256 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
88 __m256 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
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 sh_ewald = _mm256_set1_ps(fr->ic->sh_ewald);
112 beta = _mm256_set1_ps(fr->ic->ewaldcoeff);
113 beta2 = _mm256_mul_ps(beta,beta);
114 beta3 = _mm256_mul_ps(beta,beta2);
116 ewtab = fr->ic->tabq_coul_FDV0;
117 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
118 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
120 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
121 rcutoff_scalar = fr->rcoulomb;
122 rcutoff = _mm256_set1_ps(rcutoff_scalar);
123 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
125 sh_vdw_invrcut6 = _mm256_set1_ps(fr->ic->sh_invrc6);
126 rvdw = _mm256_set1_ps(fr->rvdw);
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 /* EWALD ELECTROSTATICS */
254 /* Analytical PME correction */
255 zeta2 = _mm256_mul_ps(beta2,rsq00);
256 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
257 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
258 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
259 felec = _mm256_mul_ps(qq00,felec);
260 pmecorrV = gmx_mm256_pmecorrV_ps(zeta2);
261 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
262 velec = _mm256_sub_ps(_mm256_sub_ps(rinv00,sh_ewald),pmecorrV);
263 velec = _mm256_mul_ps(qq00,velec);
265 /* LENNARD-JONES DISPERSION/REPULSION */
267 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
268 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
269 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
270 vvdw = _mm256_sub_ps(_mm256_mul_ps( _mm256_sub_ps(vvdw12 , _mm256_mul_ps(c12_00,_mm256_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
271 _mm256_mul_ps( _mm256_sub_ps(vvdw6,_mm256_mul_ps(c6_00,sh_vdw_invrcut6)),one_sixth));
272 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
274 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
276 /* Update potential sum for this i atom from the interaction with this j atom. */
277 velec = _mm256_and_ps(velec,cutoff_mask);
278 velecsum = _mm256_add_ps(velecsum,velec);
279 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
280 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
282 fscal = _mm256_add_ps(felec,fvdw);
284 fscal = _mm256_and_ps(fscal,cutoff_mask);
286 /* Calculate temporary vectorial force */
287 tx = _mm256_mul_ps(fscal,dx00);
288 ty = _mm256_mul_ps(fscal,dy00);
289 tz = _mm256_mul_ps(fscal,dz00);
291 /* Update vectorial force */
292 fix0 = _mm256_add_ps(fix0,tx);
293 fiy0 = _mm256_add_ps(fiy0,ty);
294 fiz0 = _mm256_add_ps(fiz0,tz);
296 fjptrA = f+j_coord_offsetA;
297 fjptrB = f+j_coord_offsetB;
298 fjptrC = f+j_coord_offsetC;
299 fjptrD = f+j_coord_offsetD;
300 fjptrE = f+j_coord_offsetE;
301 fjptrF = f+j_coord_offsetF;
302 fjptrG = f+j_coord_offsetG;
303 fjptrH = f+j_coord_offsetH;
304 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
308 /* Inner loop uses 127 flops */
314 /* Get j neighbor index, and coordinate index */
315 jnrlistA = jjnr[jidx];
316 jnrlistB = jjnr[jidx+1];
317 jnrlistC = jjnr[jidx+2];
318 jnrlistD = jjnr[jidx+3];
319 jnrlistE = jjnr[jidx+4];
320 jnrlistF = jjnr[jidx+5];
321 jnrlistG = jjnr[jidx+6];
322 jnrlistH = jjnr[jidx+7];
323 /* Sign of each element will be negative for non-real atoms.
324 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
325 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
327 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
328 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
330 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
331 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
332 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
333 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
334 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
335 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
336 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
337 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
338 j_coord_offsetA = DIM*jnrA;
339 j_coord_offsetB = DIM*jnrB;
340 j_coord_offsetC = DIM*jnrC;
341 j_coord_offsetD = DIM*jnrD;
342 j_coord_offsetE = DIM*jnrE;
343 j_coord_offsetF = DIM*jnrF;
344 j_coord_offsetG = DIM*jnrG;
345 j_coord_offsetH = DIM*jnrH;
347 /* load j atom coordinates */
348 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
349 x+j_coord_offsetC,x+j_coord_offsetD,
350 x+j_coord_offsetE,x+j_coord_offsetF,
351 x+j_coord_offsetG,x+j_coord_offsetH,
354 /* Calculate displacement vector */
355 dx00 = _mm256_sub_ps(ix0,jx0);
356 dy00 = _mm256_sub_ps(iy0,jy0);
357 dz00 = _mm256_sub_ps(iz0,jz0);
359 /* Calculate squared distance and things based on it */
360 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
362 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
364 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
366 /* Load parameters for j particles */
367 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
368 charge+jnrC+0,charge+jnrD+0,
369 charge+jnrE+0,charge+jnrF+0,
370 charge+jnrG+0,charge+jnrH+0);
371 vdwjidx0A = 2*vdwtype[jnrA+0];
372 vdwjidx0B = 2*vdwtype[jnrB+0];
373 vdwjidx0C = 2*vdwtype[jnrC+0];
374 vdwjidx0D = 2*vdwtype[jnrD+0];
375 vdwjidx0E = 2*vdwtype[jnrE+0];
376 vdwjidx0F = 2*vdwtype[jnrF+0];
377 vdwjidx0G = 2*vdwtype[jnrG+0];
378 vdwjidx0H = 2*vdwtype[jnrH+0];
380 /**************************
381 * CALCULATE INTERACTIONS *
382 **************************/
384 if (gmx_mm256_any_lt(rsq00,rcutoff2))
387 r00 = _mm256_mul_ps(rsq00,rinv00);
388 r00 = _mm256_andnot_ps(dummy_mask,r00);
390 /* Compute parameters for interactions between i and j atoms */
391 qq00 = _mm256_mul_ps(iq0,jq0);
392 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
393 vdwioffsetptr0+vdwjidx0B,
394 vdwioffsetptr0+vdwjidx0C,
395 vdwioffsetptr0+vdwjidx0D,
396 vdwioffsetptr0+vdwjidx0E,
397 vdwioffsetptr0+vdwjidx0F,
398 vdwioffsetptr0+vdwjidx0G,
399 vdwioffsetptr0+vdwjidx0H,
402 /* EWALD ELECTROSTATICS */
404 /* Analytical PME correction */
405 zeta2 = _mm256_mul_ps(beta2,rsq00);
406 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
407 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
408 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
409 felec = _mm256_mul_ps(qq00,felec);
410 pmecorrV = gmx_mm256_pmecorrV_ps(zeta2);
411 pmecorrV = _mm256_mul_ps(pmecorrV,beta);
412 velec = _mm256_sub_ps(_mm256_sub_ps(rinv00,sh_ewald),pmecorrV);
413 velec = _mm256_mul_ps(qq00,velec);
415 /* LENNARD-JONES DISPERSION/REPULSION */
417 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
418 vvdw6 = _mm256_mul_ps(c6_00,rinvsix);
419 vvdw12 = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
420 vvdw = _mm256_sub_ps(_mm256_mul_ps( _mm256_sub_ps(vvdw12 , _mm256_mul_ps(c12_00,_mm256_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
421 _mm256_mul_ps( _mm256_sub_ps(vvdw6,_mm256_mul_ps(c6_00,sh_vdw_invrcut6)),one_sixth));
422 fvdw = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
424 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
426 /* Update potential sum for this i atom from the interaction with this j atom. */
427 velec = _mm256_and_ps(velec,cutoff_mask);
428 velec = _mm256_andnot_ps(dummy_mask,velec);
429 velecsum = _mm256_add_ps(velecsum,velec);
430 vvdw = _mm256_and_ps(vvdw,cutoff_mask);
431 vvdw = _mm256_andnot_ps(dummy_mask,vvdw);
432 vvdwsum = _mm256_add_ps(vvdwsum,vvdw);
434 fscal = _mm256_add_ps(felec,fvdw);
436 fscal = _mm256_and_ps(fscal,cutoff_mask);
438 fscal = _mm256_andnot_ps(dummy_mask,fscal);
440 /* Calculate temporary vectorial force */
441 tx = _mm256_mul_ps(fscal,dx00);
442 ty = _mm256_mul_ps(fscal,dy00);
443 tz = _mm256_mul_ps(fscal,dz00);
445 /* Update vectorial force */
446 fix0 = _mm256_add_ps(fix0,tx);
447 fiy0 = _mm256_add_ps(fiy0,ty);
448 fiz0 = _mm256_add_ps(fiz0,tz);
450 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
451 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
452 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
453 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
454 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
455 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
456 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
457 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
458 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
462 /* Inner loop uses 128 flops */
465 /* End of innermost loop */
467 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
468 f+i_coord_offset,fshift+i_shift_offset);
471 /* Update potential energies */
472 gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
473 gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
475 /* Increment number of inner iterations */
476 inneriter += j_index_end - j_index_start;
478 /* Outer loop uses 9 flops */
481 /* Increment number of outer iterations */
484 /* Update outer/inner flops */
486 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*128);
489 * Gromacs nonbonded kernel: nb_kernel_ElecEwSh_VdwLJSh_GeomP1P1_F_avx_256_single
490 * Electrostatics interaction: Ewald
491 * VdW interaction: LennardJones
492 * Geometry: Particle-Particle
493 * Calculate force/pot: Force
496 nb_kernel_ElecEwSh_VdwLJSh_GeomP1P1_F_avx_256_single
497 (t_nblist * gmx_restrict nlist,
498 rvec * gmx_restrict xx,
499 rvec * gmx_restrict ff,
500 t_forcerec * gmx_restrict fr,
501 t_mdatoms * gmx_restrict mdatoms,
502 nb_kernel_data_t * gmx_restrict kernel_data,
503 t_nrnb * gmx_restrict nrnb)
505 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
506 * just 0 for non-waters.
507 * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
508 * jnr indices corresponding to data put in the four positions in the SIMD register.
510 int i_shift_offset,i_coord_offset,outeriter,inneriter;
511 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
512 int jnrA,jnrB,jnrC,jnrD;
513 int jnrE,jnrF,jnrG,jnrH;
514 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
515 int jnrlistE,jnrlistF,jnrlistG,jnrlistH;
516 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
517 int j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
518 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
520 real *shiftvec,*fshift,*x,*f;
521 real *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
523 __m256 tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
524 real * vdwioffsetptr0;
525 __m256 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
526 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
527 __m256 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
528 __m256 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
529 __m256 velec,felec,velecsum,facel,crf,krf,krf2;
532 __m256 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
535 __m256 one_sixth = _mm256_set1_ps(1.0/6.0);
536 __m256 one_twelfth = _mm256_set1_ps(1.0/12.0);
538 __m128i ewitab_lo,ewitab_hi;
539 __m256 ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
540 __m256 beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
542 __m256 dummy_mask,cutoff_mask;
543 __m256 signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
544 __m256 one = _mm256_set1_ps(1.0);
545 __m256 two = _mm256_set1_ps(2.0);
551 jindex = nlist->jindex;
553 shiftidx = nlist->shift;
555 shiftvec = fr->shift_vec[0];
556 fshift = fr->fshift[0];
557 facel = _mm256_set1_ps(fr->epsfac);
558 charge = mdatoms->chargeA;
559 nvdwtype = fr->ntype;
561 vdwtype = mdatoms->typeA;
563 sh_ewald = _mm256_set1_ps(fr->ic->sh_ewald);
564 beta = _mm256_set1_ps(fr->ic->ewaldcoeff);
565 beta2 = _mm256_mul_ps(beta,beta);
566 beta3 = _mm256_mul_ps(beta,beta2);
568 ewtab = fr->ic->tabq_coul_F;
569 ewtabscale = _mm256_set1_ps(fr->ic->tabq_scale);
570 ewtabhalfspace = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
572 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
573 rcutoff_scalar = fr->rcoulomb;
574 rcutoff = _mm256_set1_ps(rcutoff_scalar);
575 rcutoff2 = _mm256_mul_ps(rcutoff,rcutoff);
577 sh_vdw_invrcut6 = _mm256_set1_ps(fr->ic->sh_invrc6);
578 rvdw = _mm256_set1_ps(fr->rvdw);
580 /* Avoid stupid compiler warnings */
581 jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
594 for(iidx=0;iidx<4*DIM;iidx++)
599 /* Start outer loop over neighborlists */
600 for(iidx=0; iidx<nri; iidx++)
602 /* Load shift vector for this list */
603 i_shift_offset = DIM*shiftidx[iidx];
605 /* Load limits for loop over neighbors */
606 j_index_start = jindex[iidx];
607 j_index_end = jindex[iidx+1];
609 /* Get outer coordinate index */
611 i_coord_offset = DIM*inr;
613 /* Load i particle coords and add shift vector */
614 gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
616 fix0 = _mm256_setzero_ps();
617 fiy0 = _mm256_setzero_ps();
618 fiz0 = _mm256_setzero_ps();
620 /* Load parameters for i particles */
621 iq0 = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
622 vdwioffsetptr0 = vdwparam+2*nvdwtype*vdwtype[inr+0];
624 /* Start inner kernel loop */
625 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
628 /* Get j neighbor index, and coordinate index */
637 j_coord_offsetA = DIM*jnrA;
638 j_coord_offsetB = DIM*jnrB;
639 j_coord_offsetC = DIM*jnrC;
640 j_coord_offsetD = DIM*jnrD;
641 j_coord_offsetE = DIM*jnrE;
642 j_coord_offsetF = DIM*jnrF;
643 j_coord_offsetG = DIM*jnrG;
644 j_coord_offsetH = DIM*jnrH;
646 /* load j atom coordinates */
647 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
648 x+j_coord_offsetC,x+j_coord_offsetD,
649 x+j_coord_offsetE,x+j_coord_offsetF,
650 x+j_coord_offsetG,x+j_coord_offsetH,
653 /* Calculate displacement vector */
654 dx00 = _mm256_sub_ps(ix0,jx0);
655 dy00 = _mm256_sub_ps(iy0,jy0);
656 dz00 = _mm256_sub_ps(iz0,jz0);
658 /* Calculate squared distance and things based on it */
659 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
661 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
663 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
665 /* Load parameters for j particles */
666 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
667 charge+jnrC+0,charge+jnrD+0,
668 charge+jnrE+0,charge+jnrF+0,
669 charge+jnrG+0,charge+jnrH+0);
670 vdwjidx0A = 2*vdwtype[jnrA+0];
671 vdwjidx0B = 2*vdwtype[jnrB+0];
672 vdwjidx0C = 2*vdwtype[jnrC+0];
673 vdwjidx0D = 2*vdwtype[jnrD+0];
674 vdwjidx0E = 2*vdwtype[jnrE+0];
675 vdwjidx0F = 2*vdwtype[jnrF+0];
676 vdwjidx0G = 2*vdwtype[jnrG+0];
677 vdwjidx0H = 2*vdwtype[jnrH+0];
679 /**************************
680 * CALCULATE INTERACTIONS *
681 **************************/
683 if (gmx_mm256_any_lt(rsq00,rcutoff2))
686 r00 = _mm256_mul_ps(rsq00,rinv00);
688 /* Compute parameters for interactions between i and j atoms */
689 qq00 = _mm256_mul_ps(iq0,jq0);
690 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
691 vdwioffsetptr0+vdwjidx0B,
692 vdwioffsetptr0+vdwjidx0C,
693 vdwioffsetptr0+vdwjidx0D,
694 vdwioffsetptr0+vdwjidx0E,
695 vdwioffsetptr0+vdwjidx0F,
696 vdwioffsetptr0+vdwjidx0G,
697 vdwioffsetptr0+vdwjidx0H,
700 /* EWALD ELECTROSTATICS */
702 /* Analytical PME correction */
703 zeta2 = _mm256_mul_ps(beta2,rsq00);
704 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
705 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
706 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
707 felec = _mm256_mul_ps(qq00,felec);
709 /* LENNARD-JONES DISPERSION/REPULSION */
711 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
712 fvdw = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_00,rinvsix),c6_00),_mm256_mul_ps(rinvsix,rinvsq00));
714 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
716 fscal = _mm256_add_ps(felec,fvdw);
718 fscal = _mm256_and_ps(fscal,cutoff_mask);
720 /* Calculate temporary vectorial force */
721 tx = _mm256_mul_ps(fscal,dx00);
722 ty = _mm256_mul_ps(fscal,dy00);
723 tz = _mm256_mul_ps(fscal,dz00);
725 /* Update vectorial force */
726 fix0 = _mm256_add_ps(fix0,tx);
727 fiy0 = _mm256_add_ps(fiy0,ty);
728 fiz0 = _mm256_add_ps(fiz0,tz);
730 fjptrA = f+j_coord_offsetA;
731 fjptrB = f+j_coord_offsetB;
732 fjptrC = f+j_coord_offsetC;
733 fjptrD = f+j_coord_offsetD;
734 fjptrE = f+j_coord_offsetE;
735 fjptrF = f+j_coord_offsetF;
736 fjptrG = f+j_coord_offsetG;
737 fjptrH = f+j_coord_offsetH;
738 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
742 /* Inner loop uses 66 flops */
748 /* Get j neighbor index, and coordinate index */
749 jnrlistA = jjnr[jidx];
750 jnrlistB = jjnr[jidx+1];
751 jnrlistC = jjnr[jidx+2];
752 jnrlistD = jjnr[jidx+3];
753 jnrlistE = jjnr[jidx+4];
754 jnrlistF = jjnr[jidx+5];
755 jnrlistG = jjnr[jidx+6];
756 jnrlistH = jjnr[jidx+7];
757 /* Sign of each element will be negative for non-real atoms.
758 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
759 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
761 dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
762 gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
764 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
765 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
766 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
767 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
768 jnrE = (jnrlistE>=0) ? jnrlistE : 0;
769 jnrF = (jnrlistF>=0) ? jnrlistF : 0;
770 jnrG = (jnrlistG>=0) ? jnrlistG : 0;
771 jnrH = (jnrlistH>=0) ? jnrlistH : 0;
772 j_coord_offsetA = DIM*jnrA;
773 j_coord_offsetB = DIM*jnrB;
774 j_coord_offsetC = DIM*jnrC;
775 j_coord_offsetD = DIM*jnrD;
776 j_coord_offsetE = DIM*jnrE;
777 j_coord_offsetF = DIM*jnrF;
778 j_coord_offsetG = DIM*jnrG;
779 j_coord_offsetH = DIM*jnrH;
781 /* load j atom coordinates */
782 gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
783 x+j_coord_offsetC,x+j_coord_offsetD,
784 x+j_coord_offsetE,x+j_coord_offsetF,
785 x+j_coord_offsetG,x+j_coord_offsetH,
788 /* Calculate displacement vector */
789 dx00 = _mm256_sub_ps(ix0,jx0);
790 dy00 = _mm256_sub_ps(iy0,jy0);
791 dz00 = _mm256_sub_ps(iz0,jz0);
793 /* Calculate squared distance and things based on it */
794 rsq00 = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
796 rinv00 = gmx_mm256_invsqrt_ps(rsq00);
798 rinvsq00 = _mm256_mul_ps(rinv00,rinv00);
800 /* Load parameters for j particles */
801 jq0 = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
802 charge+jnrC+0,charge+jnrD+0,
803 charge+jnrE+0,charge+jnrF+0,
804 charge+jnrG+0,charge+jnrH+0);
805 vdwjidx0A = 2*vdwtype[jnrA+0];
806 vdwjidx0B = 2*vdwtype[jnrB+0];
807 vdwjidx0C = 2*vdwtype[jnrC+0];
808 vdwjidx0D = 2*vdwtype[jnrD+0];
809 vdwjidx0E = 2*vdwtype[jnrE+0];
810 vdwjidx0F = 2*vdwtype[jnrF+0];
811 vdwjidx0G = 2*vdwtype[jnrG+0];
812 vdwjidx0H = 2*vdwtype[jnrH+0];
814 /**************************
815 * CALCULATE INTERACTIONS *
816 **************************/
818 if (gmx_mm256_any_lt(rsq00,rcutoff2))
821 r00 = _mm256_mul_ps(rsq00,rinv00);
822 r00 = _mm256_andnot_ps(dummy_mask,r00);
824 /* Compute parameters for interactions between i and j atoms */
825 qq00 = _mm256_mul_ps(iq0,jq0);
826 gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
827 vdwioffsetptr0+vdwjidx0B,
828 vdwioffsetptr0+vdwjidx0C,
829 vdwioffsetptr0+vdwjidx0D,
830 vdwioffsetptr0+vdwjidx0E,
831 vdwioffsetptr0+vdwjidx0F,
832 vdwioffsetptr0+vdwjidx0G,
833 vdwioffsetptr0+vdwjidx0H,
836 /* EWALD ELECTROSTATICS */
838 /* Analytical PME correction */
839 zeta2 = _mm256_mul_ps(beta2,rsq00);
840 rinv3 = _mm256_mul_ps(rinvsq00,rinv00);
841 pmecorrF = gmx_mm256_pmecorrF_ps(zeta2);
842 felec = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
843 felec = _mm256_mul_ps(qq00,felec);
845 /* LENNARD-JONES DISPERSION/REPULSION */
847 rinvsix = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
848 fvdw = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_00,rinvsix),c6_00),_mm256_mul_ps(rinvsix,rinvsq00));
850 cutoff_mask = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
852 fscal = _mm256_add_ps(felec,fvdw);
854 fscal = _mm256_and_ps(fscal,cutoff_mask);
856 fscal = _mm256_andnot_ps(dummy_mask,fscal);
858 /* Calculate temporary vectorial force */
859 tx = _mm256_mul_ps(fscal,dx00);
860 ty = _mm256_mul_ps(fscal,dy00);
861 tz = _mm256_mul_ps(fscal,dz00);
863 /* Update vectorial force */
864 fix0 = _mm256_add_ps(fix0,tx);
865 fiy0 = _mm256_add_ps(fiy0,ty);
866 fiz0 = _mm256_add_ps(fiz0,tz);
868 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
869 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
870 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
871 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
872 fjptrE = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
873 fjptrF = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
874 fjptrG = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
875 fjptrH = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
876 gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
880 /* Inner loop uses 67 flops */
883 /* End of innermost loop */
885 gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
886 f+i_coord_offset,fshift+i_shift_offset);
888 /* Increment number of inner iterations */
889 inneriter += j_index_end - j_index_start;
891 /* Outer loop uses 7 flops */
894 /* Increment number of outer iterations */
897 /* Update outer/inner flops */
899 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*67);