2 * Note: this file was generated by the Gromacs avx_128_fma_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_128_fma_single.h"
34 #include "kernelutil_x86_avx_128_fma_single.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecCSTab_VdwCSTab_GeomP1P1_VF_avx_128_fma_single
38 * Electrostatics interaction: CubicSplineTable
39 * VdW interaction: CubicSplineTable
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecCSTab_VdwCSTab_GeomP1P1_VF_avx_128_fma_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 refer to j loop unrolling done with AVX_128, e.g. for the four 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 jnrlistA,jnrlistB,jnrlistC,jnrlistD;
62 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
63 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
65 real *shiftvec,*fshift,*x,*f;
66 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
68 __m128 fscal,rcutoff,rcutoff2,jidxall;
70 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
71 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
72 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
73 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
74 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
77 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
80 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
81 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
83 __m128i ifour = _mm_set1_epi32(4);
84 __m128 rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
86 __m128 dummy_mask,cutoff_mask;
87 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
88 __m128 one = _mm_set1_ps(1.0);
89 __m128 two = _mm_set1_ps(2.0);
95 jindex = nlist->jindex;
97 shiftidx = nlist->shift;
99 shiftvec = fr->shift_vec[0];
100 fshift = fr->fshift[0];
101 facel = _mm_set1_ps(fr->epsfac);
102 charge = mdatoms->chargeA;
103 nvdwtype = fr->ntype;
105 vdwtype = mdatoms->typeA;
107 vftab = kernel_data->table_elec_vdw->data;
108 vftabscale = _mm_set1_ps(kernel_data->table_elec_vdw->scale);
110 /* Avoid stupid compiler warnings */
111 jnrA = jnrB = jnrC = jnrD = 0;
120 for(iidx=0;iidx<4*DIM;iidx++)
125 /* Start outer loop over neighborlists */
126 for(iidx=0; iidx<nri; iidx++)
128 /* Load shift vector for this list */
129 i_shift_offset = DIM*shiftidx[iidx];
131 /* Load limits for loop over neighbors */
132 j_index_start = jindex[iidx];
133 j_index_end = jindex[iidx+1];
135 /* Get outer coordinate index */
137 i_coord_offset = DIM*inr;
139 /* Load i particle coords and add shift vector */
140 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
142 fix0 = _mm_setzero_ps();
143 fiy0 = _mm_setzero_ps();
144 fiz0 = _mm_setzero_ps();
146 /* Load parameters for i particles */
147 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
148 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
150 /* Reset potential sums */
151 velecsum = _mm_setzero_ps();
152 vvdwsum = _mm_setzero_ps();
154 /* Start inner kernel loop */
155 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
158 /* Get j neighbor index, and coordinate index */
163 j_coord_offsetA = DIM*jnrA;
164 j_coord_offsetB = DIM*jnrB;
165 j_coord_offsetC = DIM*jnrC;
166 j_coord_offsetD = DIM*jnrD;
168 /* load j atom coordinates */
169 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
170 x+j_coord_offsetC,x+j_coord_offsetD,
173 /* Calculate displacement vector */
174 dx00 = _mm_sub_ps(ix0,jx0);
175 dy00 = _mm_sub_ps(iy0,jy0);
176 dz00 = _mm_sub_ps(iz0,jz0);
178 /* Calculate squared distance and things based on it */
179 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
181 rinv00 = gmx_mm_invsqrt_ps(rsq00);
183 /* Load parameters for j particles */
184 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
185 charge+jnrC+0,charge+jnrD+0);
186 vdwjidx0A = 2*vdwtype[jnrA+0];
187 vdwjidx0B = 2*vdwtype[jnrB+0];
188 vdwjidx0C = 2*vdwtype[jnrC+0];
189 vdwjidx0D = 2*vdwtype[jnrD+0];
191 /**************************
192 * CALCULATE INTERACTIONS *
193 **************************/
195 r00 = _mm_mul_ps(rsq00,rinv00);
197 /* Compute parameters for interactions between i and j atoms */
198 qq00 = _mm_mul_ps(iq0,jq0);
199 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
200 vdwparam+vdwioffset0+vdwjidx0B,
201 vdwparam+vdwioffset0+vdwjidx0C,
202 vdwparam+vdwioffset0+vdwjidx0D,
205 /* Calculate table index by multiplying r with table scale and truncate to integer */
206 rt = _mm_mul_ps(r00,vftabscale);
207 vfitab = _mm_cvttps_epi32(rt);
209 vfeps = _mm_frcz_ps(rt);
211 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
213 twovfeps = _mm_add_ps(vfeps,vfeps);
214 vfitab = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
216 /* CUBIC SPLINE TABLE ELECTROSTATICS */
217 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
218 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
219 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
220 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
221 _MM_TRANSPOSE4_PS(Y,F,G,H);
222 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
223 VV = _mm_macc_ps(vfeps,Fp,Y);
224 velec = _mm_mul_ps(qq00,VV);
225 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
226 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
228 /* CUBIC SPLINE TABLE DISPERSION */
229 vfitab = _mm_add_epi32(vfitab,ifour);
230 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
231 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
232 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
233 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
234 _MM_TRANSPOSE4_PS(Y,F,G,H);
235 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
236 VV = _mm_macc_ps(vfeps,Fp,Y);
237 vvdw6 = _mm_mul_ps(c6_00,VV);
238 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
239 fvdw6 = _mm_mul_ps(c6_00,FF);
241 /* CUBIC SPLINE TABLE REPULSION */
242 vfitab = _mm_add_epi32(vfitab,ifour);
243 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
244 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
245 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
246 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
247 _MM_TRANSPOSE4_PS(Y,F,G,H);
248 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
249 VV = _mm_macc_ps(vfeps,Fp,Y);
250 vvdw12 = _mm_mul_ps(c12_00,VV);
251 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
252 fvdw12 = _mm_mul_ps(c12_00,FF);
253 vvdw = _mm_add_ps(vvdw12,vvdw6);
254 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
256 /* Update potential sum for this i atom from the interaction with this j atom. */
257 velecsum = _mm_add_ps(velecsum,velec);
258 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
260 fscal = _mm_add_ps(felec,fvdw);
262 /* Update vectorial force */
263 fix0 = _mm_macc_ps(dx00,fscal,fix0);
264 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
265 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
267 fjptrA = f+j_coord_offsetA;
268 fjptrB = f+j_coord_offsetB;
269 fjptrC = f+j_coord_offsetC;
270 fjptrD = f+j_coord_offsetD;
271 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
272 _mm_mul_ps(dx00,fscal),
273 _mm_mul_ps(dy00,fscal),
274 _mm_mul_ps(dz00,fscal));
276 /* Inner loop uses 76 flops */
282 /* Get j neighbor index, and coordinate index */
283 jnrlistA = jjnr[jidx];
284 jnrlistB = jjnr[jidx+1];
285 jnrlistC = jjnr[jidx+2];
286 jnrlistD = jjnr[jidx+3];
287 /* Sign of each element will be negative for non-real atoms.
288 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
289 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
291 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
292 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
293 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
294 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
295 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
296 j_coord_offsetA = DIM*jnrA;
297 j_coord_offsetB = DIM*jnrB;
298 j_coord_offsetC = DIM*jnrC;
299 j_coord_offsetD = DIM*jnrD;
301 /* load j atom coordinates */
302 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
303 x+j_coord_offsetC,x+j_coord_offsetD,
306 /* Calculate displacement vector */
307 dx00 = _mm_sub_ps(ix0,jx0);
308 dy00 = _mm_sub_ps(iy0,jy0);
309 dz00 = _mm_sub_ps(iz0,jz0);
311 /* Calculate squared distance and things based on it */
312 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
314 rinv00 = gmx_mm_invsqrt_ps(rsq00);
316 /* Load parameters for j particles */
317 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
318 charge+jnrC+0,charge+jnrD+0);
319 vdwjidx0A = 2*vdwtype[jnrA+0];
320 vdwjidx0B = 2*vdwtype[jnrB+0];
321 vdwjidx0C = 2*vdwtype[jnrC+0];
322 vdwjidx0D = 2*vdwtype[jnrD+0];
324 /**************************
325 * CALCULATE INTERACTIONS *
326 **************************/
328 r00 = _mm_mul_ps(rsq00,rinv00);
329 r00 = _mm_andnot_ps(dummy_mask,r00);
331 /* Compute parameters for interactions between i and j atoms */
332 qq00 = _mm_mul_ps(iq0,jq0);
333 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
334 vdwparam+vdwioffset0+vdwjidx0B,
335 vdwparam+vdwioffset0+vdwjidx0C,
336 vdwparam+vdwioffset0+vdwjidx0D,
339 /* Calculate table index by multiplying r with table scale and truncate to integer */
340 rt = _mm_mul_ps(r00,vftabscale);
341 vfitab = _mm_cvttps_epi32(rt);
343 vfeps = _mm_frcz_ps(rt);
345 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
347 twovfeps = _mm_add_ps(vfeps,vfeps);
348 vfitab = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
350 /* CUBIC SPLINE TABLE ELECTROSTATICS */
351 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
352 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
353 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
354 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
355 _MM_TRANSPOSE4_PS(Y,F,G,H);
356 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
357 VV = _mm_macc_ps(vfeps,Fp,Y);
358 velec = _mm_mul_ps(qq00,VV);
359 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
360 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
362 /* CUBIC SPLINE TABLE DISPERSION */
363 vfitab = _mm_add_epi32(vfitab,ifour);
364 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
365 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
366 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
367 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
368 _MM_TRANSPOSE4_PS(Y,F,G,H);
369 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
370 VV = _mm_macc_ps(vfeps,Fp,Y);
371 vvdw6 = _mm_mul_ps(c6_00,VV);
372 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
373 fvdw6 = _mm_mul_ps(c6_00,FF);
375 /* CUBIC SPLINE TABLE REPULSION */
376 vfitab = _mm_add_epi32(vfitab,ifour);
377 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
378 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
379 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
380 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
381 _MM_TRANSPOSE4_PS(Y,F,G,H);
382 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
383 VV = _mm_macc_ps(vfeps,Fp,Y);
384 vvdw12 = _mm_mul_ps(c12_00,VV);
385 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
386 fvdw12 = _mm_mul_ps(c12_00,FF);
387 vvdw = _mm_add_ps(vvdw12,vvdw6);
388 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
390 /* Update potential sum for this i atom from the interaction with this j atom. */
391 velec = _mm_andnot_ps(dummy_mask,velec);
392 velecsum = _mm_add_ps(velecsum,velec);
393 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
394 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
396 fscal = _mm_add_ps(felec,fvdw);
398 fscal = _mm_andnot_ps(dummy_mask,fscal);
400 /* Update vectorial force */
401 fix0 = _mm_macc_ps(dx00,fscal,fix0);
402 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
403 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
405 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
406 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
407 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
408 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
409 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
410 _mm_mul_ps(dx00,fscal),
411 _mm_mul_ps(dy00,fscal),
412 _mm_mul_ps(dz00,fscal));
414 /* Inner loop uses 77 flops */
417 /* End of innermost loop */
419 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
420 f+i_coord_offset,fshift+i_shift_offset);
423 /* Update potential energies */
424 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
425 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
427 /* Increment number of inner iterations */
428 inneriter += j_index_end - j_index_start;
430 /* Outer loop uses 9 flops */
433 /* Increment number of outer iterations */
436 /* Update outer/inner flops */
438 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*77);
441 * Gromacs nonbonded kernel: nb_kernel_ElecCSTab_VdwCSTab_GeomP1P1_F_avx_128_fma_single
442 * Electrostatics interaction: CubicSplineTable
443 * VdW interaction: CubicSplineTable
444 * Geometry: Particle-Particle
445 * Calculate force/pot: Force
448 nb_kernel_ElecCSTab_VdwCSTab_GeomP1P1_F_avx_128_fma_single
449 (t_nblist * gmx_restrict nlist,
450 rvec * gmx_restrict xx,
451 rvec * gmx_restrict ff,
452 t_forcerec * gmx_restrict fr,
453 t_mdatoms * gmx_restrict mdatoms,
454 nb_kernel_data_t * gmx_restrict kernel_data,
455 t_nrnb * gmx_restrict nrnb)
457 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
458 * just 0 for non-waters.
459 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
460 * jnr indices corresponding to data put in the four positions in the SIMD register.
462 int i_shift_offset,i_coord_offset,outeriter,inneriter;
463 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
464 int jnrA,jnrB,jnrC,jnrD;
465 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
466 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
467 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
469 real *shiftvec,*fshift,*x,*f;
470 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
472 __m128 fscal,rcutoff,rcutoff2,jidxall;
474 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
475 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
476 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
477 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
478 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
481 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
484 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
485 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
487 __m128i ifour = _mm_set1_epi32(4);
488 __m128 rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
490 __m128 dummy_mask,cutoff_mask;
491 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
492 __m128 one = _mm_set1_ps(1.0);
493 __m128 two = _mm_set1_ps(2.0);
499 jindex = nlist->jindex;
501 shiftidx = nlist->shift;
503 shiftvec = fr->shift_vec[0];
504 fshift = fr->fshift[0];
505 facel = _mm_set1_ps(fr->epsfac);
506 charge = mdatoms->chargeA;
507 nvdwtype = fr->ntype;
509 vdwtype = mdatoms->typeA;
511 vftab = kernel_data->table_elec_vdw->data;
512 vftabscale = _mm_set1_ps(kernel_data->table_elec_vdw->scale);
514 /* Avoid stupid compiler warnings */
515 jnrA = jnrB = jnrC = jnrD = 0;
524 for(iidx=0;iidx<4*DIM;iidx++)
529 /* Start outer loop over neighborlists */
530 for(iidx=0; iidx<nri; iidx++)
532 /* Load shift vector for this list */
533 i_shift_offset = DIM*shiftidx[iidx];
535 /* Load limits for loop over neighbors */
536 j_index_start = jindex[iidx];
537 j_index_end = jindex[iidx+1];
539 /* Get outer coordinate index */
541 i_coord_offset = DIM*inr;
543 /* Load i particle coords and add shift vector */
544 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
546 fix0 = _mm_setzero_ps();
547 fiy0 = _mm_setzero_ps();
548 fiz0 = _mm_setzero_ps();
550 /* Load parameters for i particles */
551 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
552 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
554 /* Start inner kernel loop */
555 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
558 /* Get j neighbor index, and coordinate index */
563 j_coord_offsetA = DIM*jnrA;
564 j_coord_offsetB = DIM*jnrB;
565 j_coord_offsetC = DIM*jnrC;
566 j_coord_offsetD = DIM*jnrD;
568 /* load j atom coordinates */
569 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
570 x+j_coord_offsetC,x+j_coord_offsetD,
573 /* Calculate displacement vector */
574 dx00 = _mm_sub_ps(ix0,jx0);
575 dy00 = _mm_sub_ps(iy0,jy0);
576 dz00 = _mm_sub_ps(iz0,jz0);
578 /* Calculate squared distance and things based on it */
579 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
581 rinv00 = gmx_mm_invsqrt_ps(rsq00);
583 /* Load parameters for j particles */
584 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
585 charge+jnrC+0,charge+jnrD+0);
586 vdwjidx0A = 2*vdwtype[jnrA+0];
587 vdwjidx0B = 2*vdwtype[jnrB+0];
588 vdwjidx0C = 2*vdwtype[jnrC+0];
589 vdwjidx0D = 2*vdwtype[jnrD+0];
591 /**************************
592 * CALCULATE INTERACTIONS *
593 **************************/
595 r00 = _mm_mul_ps(rsq00,rinv00);
597 /* Compute parameters for interactions between i and j atoms */
598 qq00 = _mm_mul_ps(iq0,jq0);
599 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
600 vdwparam+vdwioffset0+vdwjidx0B,
601 vdwparam+vdwioffset0+vdwjidx0C,
602 vdwparam+vdwioffset0+vdwjidx0D,
605 /* Calculate table index by multiplying r with table scale and truncate to integer */
606 rt = _mm_mul_ps(r00,vftabscale);
607 vfitab = _mm_cvttps_epi32(rt);
609 vfeps = _mm_frcz_ps(rt);
611 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
613 twovfeps = _mm_add_ps(vfeps,vfeps);
614 vfitab = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
616 /* CUBIC SPLINE TABLE ELECTROSTATICS */
617 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
618 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
619 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
620 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
621 _MM_TRANSPOSE4_PS(Y,F,G,H);
622 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
623 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
624 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
626 /* CUBIC SPLINE TABLE DISPERSION */
627 vfitab = _mm_add_epi32(vfitab,ifour);
628 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
629 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
630 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
631 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
632 _MM_TRANSPOSE4_PS(Y,F,G,H);
633 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
634 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
635 fvdw6 = _mm_mul_ps(c6_00,FF);
637 /* CUBIC SPLINE TABLE REPULSION */
638 vfitab = _mm_add_epi32(vfitab,ifour);
639 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
640 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
641 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
642 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
643 _MM_TRANSPOSE4_PS(Y,F,G,H);
644 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
645 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
646 fvdw12 = _mm_mul_ps(c12_00,FF);
647 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
649 fscal = _mm_add_ps(felec,fvdw);
651 /* Update vectorial force */
652 fix0 = _mm_macc_ps(dx00,fscal,fix0);
653 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
654 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
656 fjptrA = f+j_coord_offsetA;
657 fjptrB = f+j_coord_offsetB;
658 fjptrC = f+j_coord_offsetC;
659 fjptrD = f+j_coord_offsetD;
660 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
661 _mm_mul_ps(dx00,fscal),
662 _mm_mul_ps(dy00,fscal),
663 _mm_mul_ps(dz00,fscal));
665 /* Inner loop uses 64 flops */
671 /* Get j neighbor index, and coordinate index */
672 jnrlistA = jjnr[jidx];
673 jnrlistB = jjnr[jidx+1];
674 jnrlistC = jjnr[jidx+2];
675 jnrlistD = jjnr[jidx+3];
676 /* Sign of each element will be negative for non-real atoms.
677 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
678 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
680 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
681 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
682 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
683 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
684 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
685 j_coord_offsetA = DIM*jnrA;
686 j_coord_offsetB = DIM*jnrB;
687 j_coord_offsetC = DIM*jnrC;
688 j_coord_offsetD = DIM*jnrD;
690 /* load j atom coordinates */
691 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
692 x+j_coord_offsetC,x+j_coord_offsetD,
695 /* Calculate displacement vector */
696 dx00 = _mm_sub_ps(ix0,jx0);
697 dy00 = _mm_sub_ps(iy0,jy0);
698 dz00 = _mm_sub_ps(iz0,jz0);
700 /* Calculate squared distance and things based on it */
701 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
703 rinv00 = gmx_mm_invsqrt_ps(rsq00);
705 /* Load parameters for j particles */
706 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
707 charge+jnrC+0,charge+jnrD+0);
708 vdwjidx0A = 2*vdwtype[jnrA+0];
709 vdwjidx0B = 2*vdwtype[jnrB+0];
710 vdwjidx0C = 2*vdwtype[jnrC+0];
711 vdwjidx0D = 2*vdwtype[jnrD+0];
713 /**************************
714 * CALCULATE INTERACTIONS *
715 **************************/
717 r00 = _mm_mul_ps(rsq00,rinv00);
718 r00 = _mm_andnot_ps(dummy_mask,r00);
720 /* Compute parameters for interactions between i and j atoms */
721 qq00 = _mm_mul_ps(iq0,jq0);
722 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
723 vdwparam+vdwioffset0+vdwjidx0B,
724 vdwparam+vdwioffset0+vdwjidx0C,
725 vdwparam+vdwioffset0+vdwjidx0D,
728 /* Calculate table index by multiplying r with table scale and truncate to integer */
729 rt = _mm_mul_ps(r00,vftabscale);
730 vfitab = _mm_cvttps_epi32(rt);
732 vfeps = _mm_frcz_ps(rt);
734 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
736 twovfeps = _mm_add_ps(vfeps,vfeps);
737 vfitab = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
739 /* CUBIC SPLINE TABLE ELECTROSTATICS */
740 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
741 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
742 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
743 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
744 _MM_TRANSPOSE4_PS(Y,F,G,H);
745 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
746 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
747 felec = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq00,FF),_mm_mul_ps(vftabscale,rinv00)));
749 /* CUBIC SPLINE TABLE DISPERSION */
750 vfitab = _mm_add_epi32(vfitab,ifour);
751 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
752 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
753 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
754 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
755 _MM_TRANSPOSE4_PS(Y,F,G,H);
756 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
757 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
758 fvdw6 = _mm_mul_ps(c6_00,FF);
760 /* CUBIC SPLINE TABLE REPULSION */
761 vfitab = _mm_add_epi32(vfitab,ifour);
762 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
763 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
764 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
765 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
766 _MM_TRANSPOSE4_PS(Y,F,G,H);
767 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
768 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
769 fvdw12 = _mm_mul_ps(c12_00,FF);
770 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
772 fscal = _mm_add_ps(felec,fvdw);
774 fscal = _mm_andnot_ps(dummy_mask,fscal);
776 /* Update vectorial force */
777 fix0 = _mm_macc_ps(dx00,fscal,fix0);
778 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
779 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
781 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
782 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
783 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
784 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
785 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
786 _mm_mul_ps(dx00,fscal),
787 _mm_mul_ps(dy00,fscal),
788 _mm_mul_ps(dz00,fscal));
790 /* Inner loop uses 65 flops */
793 /* End of innermost loop */
795 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
796 f+i_coord_offset,fshift+i_shift_offset);
798 /* Increment number of inner iterations */
799 inneriter += j_index_end - j_index_start;
801 /* Outer loop uses 7 flops */
804 /* Increment number of outer iterations */
807 /* Update outer/inner flops */
809 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*65);