2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 * Note: this file was generated by the GROMACS avx_128_fma_single kernel generator.
44 #include "../nb_kernel.h"
45 #include "gromacs/legacyheaders/types/simple.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/legacyheaders/nrnb.h"
49 #include "gromacs/simd/math_x86_avx_128_fma_single.h"
50 #include "kernelutil_x86_avx_128_fma_single.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecCoul_VdwCSTab_GeomP1P1_VF_avx_128_fma_single
54 * Electrostatics interaction: Coulomb
55 * VdW interaction: CubicSplineTable
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecCoul_VdwCSTab_GeomP1P1_VF_avx_128_fma_single
61 (t_nblist * gmx_restrict nlist,
62 rvec * gmx_restrict xx,
63 rvec * gmx_restrict ff,
64 t_forcerec * gmx_restrict fr,
65 t_mdatoms * gmx_restrict mdatoms,
66 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67 t_nrnb * gmx_restrict nrnb)
69 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
70 * just 0 for non-waters.
71 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
72 * jnr indices corresponding to data put in the four positions in the SIMD register.
74 int i_shift_offset,i_coord_offset,outeriter,inneriter;
75 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
76 int jnrA,jnrB,jnrC,jnrD;
77 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
78 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
79 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
81 real *shiftvec,*fshift,*x,*f;
82 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
84 __m128 fscal,rcutoff,rcutoff2,jidxall;
86 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
87 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
88 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
89 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
90 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
93 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
96 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
97 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
99 __m128i ifour = _mm_set1_epi32(4);
100 __m128 rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
102 __m128 dummy_mask,cutoff_mask;
103 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
104 __m128 one = _mm_set1_ps(1.0);
105 __m128 two = _mm_set1_ps(2.0);
111 jindex = nlist->jindex;
113 shiftidx = nlist->shift;
115 shiftvec = fr->shift_vec[0];
116 fshift = fr->fshift[0];
117 facel = _mm_set1_ps(fr->epsfac);
118 charge = mdatoms->chargeA;
119 nvdwtype = fr->ntype;
121 vdwtype = mdatoms->typeA;
123 vftab = kernel_data->table_vdw->data;
124 vftabscale = _mm_set1_ps(kernel_data->table_vdw->scale);
126 /* Avoid stupid compiler warnings */
127 jnrA = jnrB = jnrC = jnrD = 0;
136 for(iidx=0;iidx<4*DIM;iidx++)
141 /* Start outer loop over neighborlists */
142 for(iidx=0; iidx<nri; iidx++)
144 /* Load shift vector for this list */
145 i_shift_offset = DIM*shiftidx[iidx];
147 /* Load limits for loop over neighbors */
148 j_index_start = jindex[iidx];
149 j_index_end = jindex[iidx+1];
151 /* Get outer coordinate index */
153 i_coord_offset = DIM*inr;
155 /* Load i particle coords and add shift vector */
156 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
158 fix0 = _mm_setzero_ps();
159 fiy0 = _mm_setzero_ps();
160 fiz0 = _mm_setzero_ps();
162 /* Load parameters for i particles */
163 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
164 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
166 /* Reset potential sums */
167 velecsum = _mm_setzero_ps();
168 vvdwsum = _mm_setzero_ps();
170 /* Start inner kernel loop */
171 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
174 /* Get j neighbor index, and coordinate index */
179 j_coord_offsetA = DIM*jnrA;
180 j_coord_offsetB = DIM*jnrB;
181 j_coord_offsetC = DIM*jnrC;
182 j_coord_offsetD = DIM*jnrD;
184 /* load j atom coordinates */
185 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
186 x+j_coord_offsetC,x+j_coord_offsetD,
189 /* Calculate displacement vector */
190 dx00 = _mm_sub_ps(ix0,jx0);
191 dy00 = _mm_sub_ps(iy0,jy0);
192 dz00 = _mm_sub_ps(iz0,jz0);
194 /* Calculate squared distance and things based on it */
195 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
197 rinv00 = gmx_mm_invsqrt_ps(rsq00);
199 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
201 /* Load parameters for j particles */
202 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
203 charge+jnrC+0,charge+jnrD+0);
204 vdwjidx0A = 2*vdwtype[jnrA+0];
205 vdwjidx0B = 2*vdwtype[jnrB+0];
206 vdwjidx0C = 2*vdwtype[jnrC+0];
207 vdwjidx0D = 2*vdwtype[jnrD+0];
209 /**************************
210 * CALCULATE INTERACTIONS *
211 **************************/
213 r00 = _mm_mul_ps(rsq00,rinv00);
215 /* Compute parameters for interactions between i and j atoms */
216 qq00 = _mm_mul_ps(iq0,jq0);
217 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
218 vdwparam+vdwioffset0+vdwjidx0B,
219 vdwparam+vdwioffset0+vdwjidx0C,
220 vdwparam+vdwioffset0+vdwjidx0D,
223 /* Calculate table index by multiplying r with table scale and truncate to integer */
224 rt = _mm_mul_ps(r00,vftabscale);
225 vfitab = _mm_cvttps_epi32(rt);
227 vfeps = _mm_frcz_ps(rt);
229 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
231 twovfeps = _mm_add_ps(vfeps,vfeps);
232 vfitab = _mm_slli_epi32(vfitab,3);
234 /* COULOMB ELECTROSTATICS */
235 velec = _mm_mul_ps(qq00,rinv00);
236 felec = _mm_mul_ps(velec,rinvsq00);
238 /* CUBIC SPLINE TABLE DISPERSION */
239 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
240 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
241 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
242 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
243 _MM_TRANSPOSE4_PS(Y,F,G,H);
244 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
245 VV = _mm_macc_ps(vfeps,Fp,Y);
246 vvdw6 = _mm_mul_ps(c6_00,VV);
247 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
248 fvdw6 = _mm_mul_ps(c6_00,FF);
250 /* CUBIC SPLINE TABLE REPULSION */
251 vfitab = _mm_add_epi32(vfitab,ifour);
252 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
253 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
254 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
255 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
256 _MM_TRANSPOSE4_PS(Y,F,G,H);
257 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
258 VV = _mm_macc_ps(vfeps,Fp,Y);
259 vvdw12 = _mm_mul_ps(c12_00,VV);
260 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
261 fvdw12 = _mm_mul_ps(c12_00,FF);
262 vvdw = _mm_add_ps(vvdw12,vvdw6);
263 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
265 /* Update potential sum for this i atom from the interaction with this j atom. */
266 velecsum = _mm_add_ps(velecsum,velec);
267 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
269 fscal = _mm_add_ps(felec,fvdw);
271 /* Update vectorial force */
272 fix0 = _mm_macc_ps(dx00,fscal,fix0);
273 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
274 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
276 fjptrA = f+j_coord_offsetA;
277 fjptrB = f+j_coord_offsetB;
278 fjptrC = f+j_coord_offsetC;
279 fjptrD = f+j_coord_offsetD;
280 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
281 _mm_mul_ps(dx00,fscal),
282 _mm_mul_ps(dy00,fscal),
283 _mm_mul_ps(dz00,fscal));
285 /* Inner loop uses 66 flops */
291 /* Get j neighbor index, and coordinate index */
292 jnrlistA = jjnr[jidx];
293 jnrlistB = jjnr[jidx+1];
294 jnrlistC = jjnr[jidx+2];
295 jnrlistD = jjnr[jidx+3];
296 /* Sign of each element will be negative for non-real atoms.
297 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
298 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
300 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
301 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
302 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
303 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
304 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
305 j_coord_offsetA = DIM*jnrA;
306 j_coord_offsetB = DIM*jnrB;
307 j_coord_offsetC = DIM*jnrC;
308 j_coord_offsetD = DIM*jnrD;
310 /* load j atom coordinates */
311 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
312 x+j_coord_offsetC,x+j_coord_offsetD,
315 /* Calculate displacement vector */
316 dx00 = _mm_sub_ps(ix0,jx0);
317 dy00 = _mm_sub_ps(iy0,jy0);
318 dz00 = _mm_sub_ps(iz0,jz0);
320 /* Calculate squared distance and things based on it */
321 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
323 rinv00 = gmx_mm_invsqrt_ps(rsq00);
325 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
327 /* Load parameters for j particles */
328 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
329 charge+jnrC+0,charge+jnrD+0);
330 vdwjidx0A = 2*vdwtype[jnrA+0];
331 vdwjidx0B = 2*vdwtype[jnrB+0];
332 vdwjidx0C = 2*vdwtype[jnrC+0];
333 vdwjidx0D = 2*vdwtype[jnrD+0];
335 /**************************
336 * CALCULATE INTERACTIONS *
337 **************************/
339 r00 = _mm_mul_ps(rsq00,rinv00);
340 r00 = _mm_andnot_ps(dummy_mask,r00);
342 /* Compute parameters for interactions between i and j atoms */
343 qq00 = _mm_mul_ps(iq0,jq0);
344 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
345 vdwparam+vdwioffset0+vdwjidx0B,
346 vdwparam+vdwioffset0+vdwjidx0C,
347 vdwparam+vdwioffset0+vdwjidx0D,
350 /* Calculate table index by multiplying r with table scale and truncate to integer */
351 rt = _mm_mul_ps(r00,vftabscale);
352 vfitab = _mm_cvttps_epi32(rt);
354 vfeps = _mm_frcz_ps(rt);
356 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
358 twovfeps = _mm_add_ps(vfeps,vfeps);
359 vfitab = _mm_slli_epi32(vfitab,3);
361 /* COULOMB ELECTROSTATICS */
362 velec = _mm_mul_ps(qq00,rinv00);
363 felec = _mm_mul_ps(velec,rinvsq00);
365 /* CUBIC SPLINE TABLE DISPERSION */
366 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
367 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
368 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
369 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
370 _MM_TRANSPOSE4_PS(Y,F,G,H);
371 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
372 VV = _mm_macc_ps(vfeps,Fp,Y);
373 vvdw6 = _mm_mul_ps(c6_00,VV);
374 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
375 fvdw6 = _mm_mul_ps(c6_00,FF);
377 /* CUBIC SPLINE TABLE REPULSION */
378 vfitab = _mm_add_epi32(vfitab,ifour);
379 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
380 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
381 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
382 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
383 _MM_TRANSPOSE4_PS(Y,F,G,H);
384 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
385 VV = _mm_macc_ps(vfeps,Fp,Y);
386 vvdw12 = _mm_mul_ps(c12_00,VV);
387 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
388 fvdw12 = _mm_mul_ps(c12_00,FF);
389 vvdw = _mm_add_ps(vvdw12,vvdw6);
390 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
392 /* Update potential sum for this i atom from the interaction with this j atom. */
393 velec = _mm_andnot_ps(dummy_mask,velec);
394 velecsum = _mm_add_ps(velecsum,velec);
395 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
396 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
398 fscal = _mm_add_ps(felec,fvdw);
400 fscal = _mm_andnot_ps(dummy_mask,fscal);
402 /* Update vectorial force */
403 fix0 = _mm_macc_ps(dx00,fscal,fix0);
404 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
405 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
407 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
408 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
409 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
410 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
411 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
412 _mm_mul_ps(dx00,fscal),
413 _mm_mul_ps(dy00,fscal),
414 _mm_mul_ps(dz00,fscal));
416 /* Inner loop uses 67 flops */
419 /* End of innermost loop */
421 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
422 f+i_coord_offset,fshift+i_shift_offset);
425 /* Update potential energies */
426 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
427 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
429 /* Increment number of inner iterations */
430 inneriter += j_index_end - j_index_start;
432 /* Outer loop uses 9 flops */
435 /* Increment number of outer iterations */
438 /* Update outer/inner flops */
440 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*67);
443 * Gromacs nonbonded kernel: nb_kernel_ElecCoul_VdwCSTab_GeomP1P1_F_avx_128_fma_single
444 * Electrostatics interaction: Coulomb
445 * VdW interaction: CubicSplineTable
446 * Geometry: Particle-Particle
447 * Calculate force/pot: Force
450 nb_kernel_ElecCoul_VdwCSTab_GeomP1P1_F_avx_128_fma_single
451 (t_nblist * gmx_restrict nlist,
452 rvec * gmx_restrict xx,
453 rvec * gmx_restrict ff,
454 t_forcerec * gmx_restrict fr,
455 t_mdatoms * gmx_restrict mdatoms,
456 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
457 t_nrnb * gmx_restrict nrnb)
459 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
460 * just 0 for non-waters.
461 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
462 * jnr indices corresponding to data put in the four positions in the SIMD register.
464 int i_shift_offset,i_coord_offset,outeriter,inneriter;
465 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
466 int jnrA,jnrB,jnrC,jnrD;
467 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
468 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
469 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
471 real *shiftvec,*fshift,*x,*f;
472 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
474 __m128 fscal,rcutoff,rcutoff2,jidxall;
476 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
477 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
478 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
479 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
480 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
483 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
486 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
487 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
489 __m128i ifour = _mm_set1_epi32(4);
490 __m128 rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
492 __m128 dummy_mask,cutoff_mask;
493 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
494 __m128 one = _mm_set1_ps(1.0);
495 __m128 two = _mm_set1_ps(2.0);
501 jindex = nlist->jindex;
503 shiftidx = nlist->shift;
505 shiftvec = fr->shift_vec[0];
506 fshift = fr->fshift[0];
507 facel = _mm_set1_ps(fr->epsfac);
508 charge = mdatoms->chargeA;
509 nvdwtype = fr->ntype;
511 vdwtype = mdatoms->typeA;
513 vftab = kernel_data->table_vdw->data;
514 vftabscale = _mm_set1_ps(kernel_data->table_vdw->scale);
516 /* Avoid stupid compiler warnings */
517 jnrA = jnrB = jnrC = jnrD = 0;
526 for(iidx=0;iidx<4*DIM;iidx++)
531 /* Start outer loop over neighborlists */
532 for(iidx=0; iidx<nri; iidx++)
534 /* Load shift vector for this list */
535 i_shift_offset = DIM*shiftidx[iidx];
537 /* Load limits for loop over neighbors */
538 j_index_start = jindex[iidx];
539 j_index_end = jindex[iidx+1];
541 /* Get outer coordinate index */
543 i_coord_offset = DIM*inr;
545 /* Load i particle coords and add shift vector */
546 gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
548 fix0 = _mm_setzero_ps();
549 fiy0 = _mm_setzero_ps();
550 fiz0 = _mm_setzero_ps();
552 /* Load parameters for i particles */
553 iq0 = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
554 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
556 /* Start inner kernel loop */
557 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
560 /* Get j neighbor index, and coordinate index */
565 j_coord_offsetA = DIM*jnrA;
566 j_coord_offsetB = DIM*jnrB;
567 j_coord_offsetC = DIM*jnrC;
568 j_coord_offsetD = DIM*jnrD;
570 /* load j atom coordinates */
571 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
572 x+j_coord_offsetC,x+j_coord_offsetD,
575 /* Calculate displacement vector */
576 dx00 = _mm_sub_ps(ix0,jx0);
577 dy00 = _mm_sub_ps(iy0,jy0);
578 dz00 = _mm_sub_ps(iz0,jz0);
580 /* Calculate squared distance and things based on it */
581 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
583 rinv00 = gmx_mm_invsqrt_ps(rsq00);
585 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
587 /* Load parameters for j particles */
588 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
589 charge+jnrC+0,charge+jnrD+0);
590 vdwjidx0A = 2*vdwtype[jnrA+0];
591 vdwjidx0B = 2*vdwtype[jnrB+0];
592 vdwjidx0C = 2*vdwtype[jnrC+0];
593 vdwjidx0D = 2*vdwtype[jnrD+0];
595 /**************************
596 * CALCULATE INTERACTIONS *
597 **************************/
599 r00 = _mm_mul_ps(rsq00,rinv00);
601 /* Compute parameters for interactions between i and j atoms */
602 qq00 = _mm_mul_ps(iq0,jq0);
603 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
604 vdwparam+vdwioffset0+vdwjidx0B,
605 vdwparam+vdwioffset0+vdwjidx0C,
606 vdwparam+vdwioffset0+vdwjidx0D,
609 /* Calculate table index by multiplying r with table scale and truncate to integer */
610 rt = _mm_mul_ps(r00,vftabscale);
611 vfitab = _mm_cvttps_epi32(rt);
613 vfeps = _mm_frcz_ps(rt);
615 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
617 twovfeps = _mm_add_ps(vfeps,vfeps);
618 vfitab = _mm_slli_epi32(vfitab,3);
620 /* COULOMB ELECTROSTATICS */
621 velec = _mm_mul_ps(qq00,rinv00);
622 felec = _mm_mul_ps(velec,rinvsq00);
624 /* CUBIC SPLINE TABLE DISPERSION */
625 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
626 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
627 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
628 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
629 _MM_TRANSPOSE4_PS(Y,F,G,H);
630 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
631 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
632 fvdw6 = _mm_mul_ps(c6_00,FF);
634 /* CUBIC SPLINE TABLE REPULSION */
635 vfitab = _mm_add_epi32(vfitab,ifour);
636 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
637 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
638 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
639 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
640 _MM_TRANSPOSE4_PS(Y,F,G,H);
641 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
642 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
643 fvdw12 = _mm_mul_ps(c12_00,FF);
644 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
646 fscal = _mm_add_ps(felec,fvdw);
648 /* Update vectorial force */
649 fix0 = _mm_macc_ps(dx00,fscal,fix0);
650 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
651 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
653 fjptrA = f+j_coord_offsetA;
654 fjptrB = f+j_coord_offsetB;
655 fjptrC = f+j_coord_offsetC;
656 fjptrD = f+j_coord_offsetD;
657 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
658 _mm_mul_ps(dx00,fscal),
659 _mm_mul_ps(dy00,fscal),
660 _mm_mul_ps(dz00,fscal));
662 /* Inner loop uses 57 flops */
668 /* Get j neighbor index, and coordinate index */
669 jnrlistA = jjnr[jidx];
670 jnrlistB = jjnr[jidx+1];
671 jnrlistC = jjnr[jidx+2];
672 jnrlistD = jjnr[jidx+3];
673 /* Sign of each element will be negative for non-real atoms.
674 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
675 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
677 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
678 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
679 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
680 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
681 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
682 j_coord_offsetA = DIM*jnrA;
683 j_coord_offsetB = DIM*jnrB;
684 j_coord_offsetC = DIM*jnrC;
685 j_coord_offsetD = DIM*jnrD;
687 /* load j atom coordinates */
688 gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
689 x+j_coord_offsetC,x+j_coord_offsetD,
692 /* Calculate displacement vector */
693 dx00 = _mm_sub_ps(ix0,jx0);
694 dy00 = _mm_sub_ps(iy0,jy0);
695 dz00 = _mm_sub_ps(iz0,jz0);
697 /* Calculate squared distance and things based on it */
698 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
700 rinv00 = gmx_mm_invsqrt_ps(rsq00);
702 rinvsq00 = _mm_mul_ps(rinv00,rinv00);
704 /* Load parameters for j particles */
705 jq0 = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
706 charge+jnrC+0,charge+jnrD+0);
707 vdwjidx0A = 2*vdwtype[jnrA+0];
708 vdwjidx0B = 2*vdwtype[jnrB+0];
709 vdwjidx0C = 2*vdwtype[jnrC+0];
710 vdwjidx0D = 2*vdwtype[jnrD+0];
712 /**************************
713 * CALCULATE INTERACTIONS *
714 **************************/
716 r00 = _mm_mul_ps(rsq00,rinv00);
717 r00 = _mm_andnot_ps(dummy_mask,r00);
719 /* Compute parameters for interactions between i and j atoms */
720 qq00 = _mm_mul_ps(iq0,jq0);
721 gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
722 vdwparam+vdwioffset0+vdwjidx0B,
723 vdwparam+vdwioffset0+vdwjidx0C,
724 vdwparam+vdwioffset0+vdwjidx0D,
727 /* Calculate table index by multiplying r with table scale and truncate to integer */
728 rt = _mm_mul_ps(r00,vftabscale);
729 vfitab = _mm_cvttps_epi32(rt);
731 vfeps = _mm_frcz_ps(rt);
733 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
735 twovfeps = _mm_add_ps(vfeps,vfeps);
736 vfitab = _mm_slli_epi32(vfitab,3);
738 /* COULOMB ELECTROSTATICS */
739 velec = _mm_mul_ps(qq00,rinv00);
740 felec = _mm_mul_ps(velec,rinvsq00);
742 /* CUBIC SPLINE TABLE DISPERSION */
743 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
744 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
745 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
746 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
747 _MM_TRANSPOSE4_PS(Y,F,G,H);
748 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
749 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
750 fvdw6 = _mm_mul_ps(c6_00,FF);
752 /* CUBIC SPLINE TABLE REPULSION */
753 vfitab = _mm_add_epi32(vfitab,ifour);
754 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
755 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
756 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
757 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
758 _MM_TRANSPOSE4_PS(Y,F,G,H);
759 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
760 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
761 fvdw12 = _mm_mul_ps(c12_00,FF);
762 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
764 fscal = _mm_add_ps(felec,fvdw);
766 fscal = _mm_andnot_ps(dummy_mask,fscal);
768 /* Update vectorial force */
769 fix0 = _mm_macc_ps(dx00,fscal,fix0);
770 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
771 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
773 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
774 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
775 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
776 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
777 gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
778 _mm_mul_ps(dx00,fscal),
779 _mm_mul_ps(dy00,fscal),
780 _mm_mul_ps(dz00,fscal));
782 /* Inner loop uses 58 flops */
785 /* End of innermost loop */
787 gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
788 f+i_coord_offset,fshift+i_shift_offset);
790 /* Increment number of inner iterations */
791 inneriter += j_index_end - j_index_start;
793 /* Outer loop uses 7 flops */
796 /* Increment number of outer iterations */
799 /* Update outer/inner flops */
801 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*58);