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_ElecRF_VdwCSTab_GeomW4W4_VF_avx_128_fma_single
38 * Electrostatics interaction: ReactionField
39 * VdW interaction: CubicSplineTable
40 * Geometry: Water4-Water4
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecRF_VdwCSTab_GeomW4W4_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;
72 __m128 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
74 __m128 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
76 __m128 ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
77 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
78 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
79 int vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D;
80 __m128 jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
81 int vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D;
82 __m128 jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
83 int vdwjidx3A,vdwjidx3B,vdwjidx3C,vdwjidx3D;
84 __m128 jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
85 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
86 __m128 dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
87 __m128 dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
88 __m128 dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
89 __m128 dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
90 __m128 dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
91 __m128 dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
92 __m128 dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
93 __m128 dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
94 __m128 dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
95 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
98 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
101 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
102 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
104 __m128i ifour = _mm_set1_epi32(4);
105 __m128 rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
107 __m128 dummy_mask,cutoff_mask;
108 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
109 __m128 one = _mm_set1_ps(1.0);
110 __m128 two = _mm_set1_ps(2.0);
116 jindex = nlist->jindex;
118 shiftidx = nlist->shift;
120 shiftvec = fr->shift_vec[0];
121 fshift = fr->fshift[0];
122 facel = _mm_set1_ps(fr->epsfac);
123 charge = mdatoms->chargeA;
124 krf = _mm_set1_ps(fr->ic->k_rf);
125 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
126 crf = _mm_set1_ps(fr->ic->c_rf);
127 nvdwtype = fr->ntype;
129 vdwtype = mdatoms->typeA;
131 vftab = kernel_data->table_vdw->data;
132 vftabscale = _mm_set1_ps(kernel_data->table_vdw->scale);
134 /* Setup water-specific parameters */
135 inr = nlist->iinr[0];
136 iq1 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
137 iq2 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
138 iq3 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
139 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
141 jq1 = _mm_set1_ps(charge[inr+1]);
142 jq2 = _mm_set1_ps(charge[inr+2]);
143 jq3 = _mm_set1_ps(charge[inr+3]);
144 vdwjidx0A = 2*vdwtype[inr+0];
145 c6_00 = _mm_set1_ps(vdwparam[vdwioffset0+vdwjidx0A]);
146 c12_00 = _mm_set1_ps(vdwparam[vdwioffset0+vdwjidx0A+1]);
147 qq11 = _mm_mul_ps(iq1,jq1);
148 qq12 = _mm_mul_ps(iq1,jq2);
149 qq13 = _mm_mul_ps(iq1,jq3);
150 qq21 = _mm_mul_ps(iq2,jq1);
151 qq22 = _mm_mul_ps(iq2,jq2);
152 qq23 = _mm_mul_ps(iq2,jq3);
153 qq31 = _mm_mul_ps(iq3,jq1);
154 qq32 = _mm_mul_ps(iq3,jq2);
155 qq33 = _mm_mul_ps(iq3,jq3);
157 /* Avoid stupid compiler warnings */
158 jnrA = jnrB = jnrC = jnrD = 0;
167 for(iidx=0;iidx<4*DIM;iidx++)
172 /* Start outer loop over neighborlists */
173 for(iidx=0; iidx<nri; iidx++)
175 /* Load shift vector for this list */
176 i_shift_offset = DIM*shiftidx[iidx];
178 /* Load limits for loop over neighbors */
179 j_index_start = jindex[iidx];
180 j_index_end = jindex[iidx+1];
182 /* Get outer coordinate index */
184 i_coord_offset = DIM*inr;
186 /* Load i particle coords and add shift vector */
187 gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
188 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
190 fix0 = _mm_setzero_ps();
191 fiy0 = _mm_setzero_ps();
192 fiz0 = _mm_setzero_ps();
193 fix1 = _mm_setzero_ps();
194 fiy1 = _mm_setzero_ps();
195 fiz1 = _mm_setzero_ps();
196 fix2 = _mm_setzero_ps();
197 fiy2 = _mm_setzero_ps();
198 fiz2 = _mm_setzero_ps();
199 fix3 = _mm_setzero_ps();
200 fiy3 = _mm_setzero_ps();
201 fiz3 = _mm_setzero_ps();
203 /* Reset potential sums */
204 velecsum = _mm_setzero_ps();
205 vvdwsum = _mm_setzero_ps();
207 /* Start inner kernel loop */
208 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
211 /* Get j neighbor index, and coordinate index */
216 j_coord_offsetA = DIM*jnrA;
217 j_coord_offsetB = DIM*jnrB;
218 j_coord_offsetC = DIM*jnrC;
219 j_coord_offsetD = DIM*jnrD;
221 /* load j atom coordinates */
222 gmx_mm_load_4rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
223 x+j_coord_offsetC,x+j_coord_offsetD,
224 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
225 &jy2,&jz2,&jx3,&jy3,&jz3);
227 /* Calculate displacement vector */
228 dx00 = _mm_sub_ps(ix0,jx0);
229 dy00 = _mm_sub_ps(iy0,jy0);
230 dz00 = _mm_sub_ps(iz0,jz0);
231 dx11 = _mm_sub_ps(ix1,jx1);
232 dy11 = _mm_sub_ps(iy1,jy1);
233 dz11 = _mm_sub_ps(iz1,jz1);
234 dx12 = _mm_sub_ps(ix1,jx2);
235 dy12 = _mm_sub_ps(iy1,jy2);
236 dz12 = _mm_sub_ps(iz1,jz2);
237 dx13 = _mm_sub_ps(ix1,jx3);
238 dy13 = _mm_sub_ps(iy1,jy3);
239 dz13 = _mm_sub_ps(iz1,jz3);
240 dx21 = _mm_sub_ps(ix2,jx1);
241 dy21 = _mm_sub_ps(iy2,jy1);
242 dz21 = _mm_sub_ps(iz2,jz1);
243 dx22 = _mm_sub_ps(ix2,jx2);
244 dy22 = _mm_sub_ps(iy2,jy2);
245 dz22 = _mm_sub_ps(iz2,jz2);
246 dx23 = _mm_sub_ps(ix2,jx3);
247 dy23 = _mm_sub_ps(iy2,jy3);
248 dz23 = _mm_sub_ps(iz2,jz3);
249 dx31 = _mm_sub_ps(ix3,jx1);
250 dy31 = _mm_sub_ps(iy3,jy1);
251 dz31 = _mm_sub_ps(iz3,jz1);
252 dx32 = _mm_sub_ps(ix3,jx2);
253 dy32 = _mm_sub_ps(iy3,jy2);
254 dz32 = _mm_sub_ps(iz3,jz2);
255 dx33 = _mm_sub_ps(ix3,jx3);
256 dy33 = _mm_sub_ps(iy3,jy3);
257 dz33 = _mm_sub_ps(iz3,jz3);
259 /* Calculate squared distance and things based on it */
260 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
261 rsq11 = gmx_mm_calc_rsq_ps(dx11,dy11,dz11);
262 rsq12 = gmx_mm_calc_rsq_ps(dx12,dy12,dz12);
263 rsq13 = gmx_mm_calc_rsq_ps(dx13,dy13,dz13);
264 rsq21 = gmx_mm_calc_rsq_ps(dx21,dy21,dz21);
265 rsq22 = gmx_mm_calc_rsq_ps(dx22,dy22,dz22);
266 rsq23 = gmx_mm_calc_rsq_ps(dx23,dy23,dz23);
267 rsq31 = gmx_mm_calc_rsq_ps(dx31,dy31,dz31);
268 rsq32 = gmx_mm_calc_rsq_ps(dx32,dy32,dz32);
269 rsq33 = gmx_mm_calc_rsq_ps(dx33,dy33,dz33);
271 rinv00 = gmx_mm_invsqrt_ps(rsq00);
272 rinv11 = gmx_mm_invsqrt_ps(rsq11);
273 rinv12 = gmx_mm_invsqrt_ps(rsq12);
274 rinv13 = gmx_mm_invsqrt_ps(rsq13);
275 rinv21 = gmx_mm_invsqrt_ps(rsq21);
276 rinv22 = gmx_mm_invsqrt_ps(rsq22);
277 rinv23 = gmx_mm_invsqrt_ps(rsq23);
278 rinv31 = gmx_mm_invsqrt_ps(rsq31);
279 rinv32 = gmx_mm_invsqrt_ps(rsq32);
280 rinv33 = gmx_mm_invsqrt_ps(rsq33);
282 rinvsq11 = _mm_mul_ps(rinv11,rinv11);
283 rinvsq12 = _mm_mul_ps(rinv12,rinv12);
284 rinvsq13 = _mm_mul_ps(rinv13,rinv13);
285 rinvsq21 = _mm_mul_ps(rinv21,rinv21);
286 rinvsq22 = _mm_mul_ps(rinv22,rinv22);
287 rinvsq23 = _mm_mul_ps(rinv23,rinv23);
288 rinvsq31 = _mm_mul_ps(rinv31,rinv31);
289 rinvsq32 = _mm_mul_ps(rinv32,rinv32);
290 rinvsq33 = _mm_mul_ps(rinv33,rinv33);
292 fjx0 = _mm_setzero_ps();
293 fjy0 = _mm_setzero_ps();
294 fjz0 = _mm_setzero_ps();
295 fjx1 = _mm_setzero_ps();
296 fjy1 = _mm_setzero_ps();
297 fjz1 = _mm_setzero_ps();
298 fjx2 = _mm_setzero_ps();
299 fjy2 = _mm_setzero_ps();
300 fjz2 = _mm_setzero_ps();
301 fjx3 = _mm_setzero_ps();
302 fjy3 = _mm_setzero_ps();
303 fjz3 = _mm_setzero_ps();
305 /**************************
306 * CALCULATE INTERACTIONS *
307 **************************/
309 r00 = _mm_mul_ps(rsq00,rinv00);
311 /* Calculate table index by multiplying r with table scale and truncate to integer */
312 rt = _mm_mul_ps(r00,vftabscale);
313 vfitab = _mm_cvttps_epi32(rt);
315 vfeps = _mm_frcz_ps(rt);
317 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
319 twovfeps = _mm_add_ps(vfeps,vfeps);
320 vfitab = _mm_slli_epi32(vfitab,3);
322 /* CUBIC SPLINE TABLE DISPERSION */
323 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
324 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
325 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
326 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
327 _MM_TRANSPOSE4_PS(Y,F,G,H);
328 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
329 VV = _mm_macc_ps(vfeps,Fp,Y);
330 vvdw6 = _mm_mul_ps(c6_00,VV);
331 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
332 fvdw6 = _mm_mul_ps(c6_00,FF);
334 /* CUBIC SPLINE TABLE REPULSION */
335 vfitab = _mm_add_epi32(vfitab,ifour);
336 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
337 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
338 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
339 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
340 _MM_TRANSPOSE4_PS(Y,F,G,H);
341 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
342 VV = _mm_macc_ps(vfeps,Fp,Y);
343 vvdw12 = _mm_mul_ps(c12_00,VV);
344 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
345 fvdw12 = _mm_mul_ps(c12_00,FF);
346 vvdw = _mm_add_ps(vvdw12,vvdw6);
347 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
349 /* Update potential sum for this i atom from the interaction with this j atom. */
350 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
354 /* Update vectorial force */
355 fix0 = _mm_macc_ps(dx00,fscal,fix0);
356 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
357 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
359 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
360 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
361 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
363 /**************************
364 * CALCULATE INTERACTIONS *
365 **************************/
367 /* REACTION-FIELD ELECTROSTATICS */
368 velec = _mm_mul_ps(qq11,_mm_sub_ps(_mm_macc_ps(krf,rsq11,rinv11),crf));
369 felec = _mm_mul_ps(qq11,_mm_msub_ps(rinv11,rinvsq11,krf2));
371 /* Update potential sum for this i atom from the interaction with this j atom. */
372 velecsum = _mm_add_ps(velecsum,velec);
376 /* Update vectorial force */
377 fix1 = _mm_macc_ps(dx11,fscal,fix1);
378 fiy1 = _mm_macc_ps(dy11,fscal,fiy1);
379 fiz1 = _mm_macc_ps(dz11,fscal,fiz1);
381 fjx1 = _mm_macc_ps(dx11,fscal,fjx1);
382 fjy1 = _mm_macc_ps(dy11,fscal,fjy1);
383 fjz1 = _mm_macc_ps(dz11,fscal,fjz1);
385 /**************************
386 * CALCULATE INTERACTIONS *
387 **************************/
389 /* REACTION-FIELD ELECTROSTATICS */
390 velec = _mm_mul_ps(qq12,_mm_sub_ps(_mm_macc_ps(krf,rsq12,rinv12),crf));
391 felec = _mm_mul_ps(qq12,_mm_msub_ps(rinv12,rinvsq12,krf2));
393 /* Update potential sum for this i atom from the interaction with this j atom. */
394 velecsum = _mm_add_ps(velecsum,velec);
398 /* Update vectorial force */
399 fix1 = _mm_macc_ps(dx12,fscal,fix1);
400 fiy1 = _mm_macc_ps(dy12,fscal,fiy1);
401 fiz1 = _mm_macc_ps(dz12,fscal,fiz1);
403 fjx2 = _mm_macc_ps(dx12,fscal,fjx2);
404 fjy2 = _mm_macc_ps(dy12,fscal,fjy2);
405 fjz2 = _mm_macc_ps(dz12,fscal,fjz2);
407 /**************************
408 * CALCULATE INTERACTIONS *
409 **************************/
411 /* REACTION-FIELD ELECTROSTATICS */
412 velec = _mm_mul_ps(qq13,_mm_sub_ps(_mm_macc_ps(krf,rsq13,rinv13),crf));
413 felec = _mm_mul_ps(qq13,_mm_msub_ps(rinv13,rinvsq13,krf2));
415 /* Update potential sum for this i atom from the interaction with this j atom. */
416 velecsum = _mm_add_ps(velecsum,velec);
420 /* Update vectorial force */
421 fix1 = _mm_macc_ps(dx13,fscal,fix1);
422 fiy1 = _mm_macc_ps(dy13,fscal,fiy1);
423 fiz1 = _mm_macc_ps(dz13,fscal,fiz1);
425 fjx3 = _mm_macc_ps(dx13,fscal,fjx3);
426 fjy3 = _mm_macc_ps(dy13,fscal,fjy3);
427 fjz3 = _mm_macc_ps(dz13,fscal,fjz3);
429 /**************************
430 * CALCULATE INTERACTIONS *
431 **************************/
433 /* REACTION-FIELD ELECTROSTATICS */
434 velec = _mm_mul_ps(qq21,_mm_sub_ps(_mm_macc_ps(krf,rsq21,rinv21),crf));
435 felec = _mm_mul_ps(qq21,_mm_msub_ps(rinv21,rinvsq21,krf2));
437 /* Update potential sum for this i atom from the interaction with this j atom. */
438 velecsum = _mm_add_ps(velecsum,velec);
442 /* Update vectorial force */
443 fix2 = _mm_macc_ps(dx21,fscal,fix2);
444 fiy2 = _mm_macc_ps(dy21,fscal,fiy2);
445 fiz2 = _mm_macc_ps(dz21,fscal,fiz2);
447 fjx1 = _mm_macc_ps(dx21,fscal,fjx1);
448 fjy1 = _mm_macc_ps(dy21,fscal,fjy1);
449 fjz1 = _mm_macc_ps(dz21,fscal,fjz1);
451 /**************************
452 * CALCULATE INTERACTIONS *
453 **************************/
455 /* REACTION-FIELD ELECTROSTATICS */
456 velec = _mm_mul_ps(qq22,_mm_sub_ps(_mm_macc_ps(krf,rsq22,rinv22),crf));
457 felec = _mm_mul_ps(qq22,_mm_msub_ps(rinv22,rinvsq22,krf2));
459 /* Update potential sum for this i atom from the interaction with this j atom. */
460 velecsum = _mm_add_ps(velecsum,velec);
464 /* Update vectorial force */
465 fix2 = _mm_macc_ps(dx22,fscal,fix2);
466 fiy2 = _mm_macc_ps(dy22,fscal,fiy2);
467 fiz2 = _mm_macc_ps(dz22,fscal,fiz2);
469 fjx2 = _mm_macc_ps(dx22,fscal,fjx2);
470 fjy2 = _mm_macc_ps(dy22,fscal,fjy2);
471 fjz2 = _mm_macc_ps(dz22,fscal,fjz2);
473 /**************************
474 * CALCULATE INTERACTIONS *
475 **************************/
477 /* REACTION-FIELD ELECTROSTATICS */
478 velec = _mm_mul_ps(qq23,_mm_sub_ps(_mm_macc_ps(krf,rsq23,rinv23),crf));
479 felec = _mm_mul_ps(qq23,_mm_msub_ps(rinv23,rinvsq23,krf2));
481 /* Update potential sum for this i atom from the interaction with this j atom. */
482 velecsum = _mm_add_ps(velecsum,velec);
486 /* Update vectorial force */
487 fix2 = _mm_macc_ps(dx23,fscal,fix2);
488 fiy2 = _mm_macc_ps(dy23,fscal,fiy2);
489 fiz2 = _mm_macc_ps(dz23,fscal,fiz2);
491 fjx3 = _mm_macc_ps(dx23,fscal,fjx3);
492 fjy3 = _mm_macc_ps(dy23,fscal,fjy3);
493 fjz3 = _mm_macc_ps(dz23,fscal,fjz3);
495 /**************************
496 * CALCULATE INTERACTIONS *
497 **************************/
499 /* REACTION-FIELD ELECTROSTATICS */
500 velec = _mm_mul_ps(qq31,_mm_sub_ps(_mm_macc_ps(krf,rsq31,rinv31),crf));
501 felec = _mm_mul_ps(qq31,_mm_msub_ps(rinv31,rinvsq31,krf2));
503 /* Update potential sum for this i atom from the interaction with this j atom. */
504 velecsum = _mm_add_ps(velecsum,velec);
508 /* Update vectorial force */
509 fix3 = _mm_macc_ps(dx31,fscal,fix3);
510 fiy3 = _mm_macc_ps(dy31,fscal,fiy3);
511 fiz3 = _mm_macc_ps(dz31,fscal,fiz3);
513 fjx1 = _mm_macc_ps(dx31,fscal,fjx1);
514 fjy1 = _mm_macc_ps(dy31,fscal,fjy1);
515 fjz1 = _mm_macc_ps(dz31,fscal,fjz1);
517 /**************************
518 * CALCULATE INTERACTIONS *
519 **************************/
521 /* REACTION-FIELD ELECTROSTATICS */
522 velec = _mm_mul_ps(qq32,_mm_sub_ps(_mm_macc_ps(krf,rsq32,rinv32),crf));
523 felec = _mm_mul_ps(qq32,_mm_msub_ps(rinv32,rinvsq32,krf2));
525 /* Update potential sum for this i atom from the interaction with this j atom. */
526 velecsum = _mm_add_ps(velecsum,velec);
530 /* Update vectorial force */
531 fix3 = _mm_macc_ps(dx32,fscal,fix3);
532 fiy3 = _mm_macc_ps(dy32,fscal,fiy3);
533 fiz3 = _mm_macc_ps(dz32,fscal,fiz3);
535 fjx2 = _mm_macc_ps(dx32,fscal,fjx2);
536 fjy2 = _mm_macc_ps(dy32,fscal,fjy2);
537 fjz2 = _mm_macc_ps(dz32,fscal,fjz2);
539 /**************************
540 * CALCULATE INTERACTIONS *
541 **************************/
543 /* REACTION-FIELD ELECTROSTATICS */
544 velec = _mm_mul_ps(qq33,_mm_sub_ps(_mm_macc_ps(krf,rsq33,rinv33),crf));
545 felec = _mm_mul_ps(qq33,_mm_msub_ps(rinv33,rinvsq33,krf2));
547 /* Update potential sum for this i atom from the interaction with this j atom. */
548 velecsum = _mm_add_ps(velecsum,velec);
552 /* Update vectorial force */
553 fix3 = _mm_macc_ps(dx33,fscal,fix3);
554 fiy3 = _mm_macc_ps(dy33,fscal,fiy3);
555 fiz3 = _mm_macc_ps(dz33,fscal,fiz3);
557 fjx3 = _mm_macc_ps(dx33,fscal,fjx3);
558 fjy3 = _mm_macc_ps(dy33,fscal,fjy3);
559 fjz3 = _mm_macc_ps(dz33,fscal,fjz3);
561 fjptrA = f+j_coord_offsetA;
562 fjptrB = f+j_coord_offsetB;
563 fjptrC = f+j_coord_offsetC;
564 fjptrD = f+j_coord_offsetD;
566 gmx_mm_decrement_4rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
567 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
568 fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
570 /* Inner loop uses 377 flops */
576 /* Get j neighbor index, and coordinate index */
577 jnrlistA = jjnr[jidx];
578 jnrlistB = jjnr[jidx+1];
579 jnrlistC = jjnr[jidx+2];
580 jnrlistD = jjnr[jidx+3];
581 /* Sign of each element will be negative for non-real atoms.
582 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
583 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
585 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
586 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
587 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
588 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
589 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
590 j_coord_offsetA = DIM*jnrA;
591 j_coord_offsetB = DIM*jnrB;
592 j_coord_offsetC = DIM*jnrC;
593 j_coord_offsetD = DIM*jnrD;
595 /* load j atom coordinates */
596 gmx_mm_load_4rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
597 x+j_coord_offsetC,x+j_coord_offsetD,
598 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
599 &jy2,&jz2,&jx3,&jy3,&jz3);
601 /* Calculate displacement vector */
602 dx00 = _mm_sub_ps(ix0,jx0);
603 dy00 = _mm_sub_ps(iy0,jy0);
604 dz00 = _mm_sub_ps(iz0,jz0);
605 dx11 = _mm_sub_ps(ix1,jx1);
606 dy11 = _mm_sub_ps(iy1,jy1);
607 dz11 = _mm_sub_ps(iz1,jz1);
608 dx12 = _mm_sub_ps(ix1,jx2);
609 dy12 = _mm_sub_ps(iy1,jy2);
610 dz12 = _mm_sub_ps(iz1,jz2);
611 dx13 = _mm_sub_ps(ix1,jx3);
612 dy13 = _mm_sub_ps(iy1,jy3);
613 dz13 = _mm_sub_ps(iz1,jz3);
614 dx21 = _mm_sub_ps(ix2,jx1);
615 dy21 = _mm_sub_ps(iy2,jy1);
616 dz21 = _mm_sub_ps(iz2,jz1);
617 dx22 = _mm_sub_ps(ix2,jx2);
618 dy22 = _mm_sub_ps(iy2,jy2);
619 dz22 = _mm_sub_ps(iz2,jz2);
620 dx23 = _mm_sub_ps(ix2,jx3);
621 dy23 = _mm_sub_ps(iy2,jy3);
622 dz23 = _mm_sub_ps(iz2,jz3);
623 dx31 = _mm_sub_ps(ix3,jx1);
624 dy31 = _mm_sub_ps(iy3,jy1);
625 dz31 = _mm_sub_ps(iz3,jz1);
626 dx32 = _mm_sub_ps(ix3,jx2);
627 dy32 = _mm_sub_ps(iy3,jy2);
628 dz32 = _mm_sub_ps(iz3,jz2);
629 dx33 = _mm_sub_ps(ix3,jx3);
630 dy33 = _mm_sub_ps(iy3,jy3);
631 dz33 = _mm_sub_ps(iz3,jz3);
633 /* Calculate squared distance and things based on it */
634 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
635 rsq11 = gmx_mm_calc_rsq_ps(dx11,dy11,dz11);
636 rsq12 = gmx_mm_calc_rsq_ps(dx12,dy12,dz12);
637 rsq13 = gmx_mm_calc_rsq_ps(dx13,dy13,dz13);
638 rsq21 = gmx_mm_calc_rsq_ps(dx21,dy21,dz21);
639 rsq22 = gmx_mm_calc_rsq_ps(dx22,dy22,dz22);
640 rsq23 = gmx_mm_calc_rsq_ps(dx23,dy23,dz23);
641 rsq31 = gmx_mm_calc_rsq_ps(dx31,dy31,dz31);
642 rsq32 = gmx_mm_calc_rsq_ps(dx32,dy32,dz32);
643 rsq33 = gmx_mm_calc_rsq_ps(dx33,dy33,dz33);
645 rinv00 = gmx_mm_invsqrt_ps(rsq00);
646 rinv11 = gmx_mm_invsqrt_ps(rsq11);
647 rinv12 = gmx_mm_invsqrt_ps(rsq12);
648 rinv13 = gmx_mm_invsqrt_ps(rsq13);
649 rinv21 = gmx_mm_invsqrt_ps(rsq21);
650 rinv22 = gmx_mm_invsqrt_ps(rsq22);
651 rinv23 = gmx_mm_invsqrt_ps(rsq23);
652 rinv31 = gmx_mm_invsqrt_ps(rsq31);
653 rinv32 = gmx_mm_invsqrt_ps(rsq32);
654 rinv33 = gmx_mm_invsqrt_ps(rsq33);
656 rinvsq11 = _mm_mul_ps(rinv11,rinv11);
657 rinvsq12 = _mm_mul_ps(rinv12,rinv12);
658 rinvsq13 = _mm_mul_ps(rinv13,rinv13);
659 rinvsq21 = _mm_mul_ps(rinv21,rinv21);
660 rinvsq22 = _mm_mul_ps(rinv22,rinv22);
661 rinvsq23 = _mm_mul_ps(rinv23,rinv23);
662 rinvsq31 = _mm_mul_ps(rinv31,rinv31);
663 rinvsq32 = _mm_mul_ps(rinv32,rinv32);
664 rinvsq33 = _mm_mul_ps(rinv33,rinv33);
666 fjx0 = _mm_setzero_ps();
667 fjy0 = _mm_setzero_ps();
668 fjz0 = _mm_setzero_ps();
669 fjx1 = _mm_setzero_ps();
670 fjy1 = _mm_setzero_ps();
671 fjz1 = _mm_setzero_ps();
672 fjx2 = _mm_setzero_ps();
673 fjy2 = _mm_setzero_ps();
674 fjz2 = _mm_setzero_ps();
675 fjx3 = _mm_setzero_ps();
676 fjy3 = _mm_setzero_ps();
677 fjz3 = _mm_setzero_ps();
679 /**************************
680 * CALCULATE INTERACTIONS *
681 **************************/
683 r00 = _mm_mul_ps(rsq00,rinv00);
684 r00 = _mm_andnot_ps(dummy_mask,r00);
686 /* Calculate table index by multiplying r with table scale and truncate to integer */
687 rt = _mm_mul_ps(r00,vftabscale);
688 vfitab = _mm_cvttps_epi32(rt);
690 vfeps = _mm_frcz_ps(rt);
692 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
694 twovfeps = _mm_add_ps(vfeps,vfeps);
695 vfitab = _mm_slli_epi32(vfitab,3);
697 /* CUBIC SPLINE TABLE DISPERSION */
698 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
699 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
700 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
701 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
702 _MM_TRANSPOSE4_PS(Y,F,G,H);
703 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
704 VV = _mm_macc_ps(vfeps,Fp,Y);
705 vvdw6 = _mm_mul_ps(c6_00,VV);
706 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
707 fvdw6 = _mm_mul_ps(c6_00,FF);
709 /* CUBIC SPLINE TABLE REPULSION */
710 vfitab = _mm_add_epi32(vfitab,ifour);
711 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
712 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
713 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
714 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
715 _MM_TRANSPOSE4_PS(Y,F,G,H);
716 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
717 VV = _mm_macc_ps(vfeps,Fp,Y);
718 vvdw12 = _mm_mul_ps(c12_00,VV);
719 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
720 fvdw12 = _mm_mul_ps(c12_00,FF);
721 vvdw = _mm_add_ps(vvdw12,vvdw6);
722 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
724 /* Update potential sum for this i atom from the interaction with this j atom. */
725 vvdw = _mm_andnot_ps(dummy_mask,vvdw);
726 vvdwsum = _mm_add_ps(vvdwsum,vvdw);
730 fscal = _mm_andnot_ps(dummy_mask,fscal);
732 /* Update vectorial force */
733 fix0 = _mm_macc_ps(dx00,fscal,fix0);
734 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
735 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
737 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
738 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
739 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
741 /**************************
742 * CALCULATE INTERACTIONS *
743 **************************/
745 /* REACTION-FIELD ELECTROSTATICS */
746 velec = _mm_mul_ps(qq11,_mm_sub_ps(_mm_macc_ps(krf,rsq11,rinv11),crf));
747 felec = _mm_mul_ps(qq11,_mm_msub_ps(rinv11,rinvsq11,krf2));
749 /* Update potential sum for this i atom from the interaction with this j atom. */
750 velec = _mm_andnot_ps(dummy_mask,velec);
751 velecsum = _mm_add_ps(velecsum,velec);
755 fscal = _mm_andnot_ps(dummy_mask,fscal);
757 /* Update vectorial force */
758 fix1 = _mm_macc_ps(dx11,fscal,fix1);
759 fiy1 = _mm_macc_ps(dy11,fscal,fiy1);
760 fiz1 = _mm_macc_ps(dz11,fscal,fiz1);
762 fjx1 = _mm_macc_ps(dx11,fscal,fjx1);
763 fjy1 = _mm_macc_ps(dy11,fscal,fjy1);
764 fjz1 = _mm_macc_ps(dz11,fscal,fjz1);
766 /**************************
767 * CALCULATE INTERACTIONS *
768 **************************/
770 /* REACTION-FIELD ELECTROSTATICS */
771 velec = _mm_mul_ps(qq12,_mm_sub_ps(_mm_macc_ps(krf,rsq12,rinv12),crf));
772 felec = _mm_mul_ps(qq12,_mm_msub_ps(rinv12,rinvsq12,krf2));
774 /* Update potential sum for this i atom from the interaction with this j atom. */
775 velec = _mm_andnot_ps(dummy_mask,velec);
776 velecsum = _mm_add_ps(velecsum,velec);
780 fscal = _mm_andnot_ps(dummy_mask,fscal);
782 /* Update vectorial force */
783 fix1 = _mm_macc_ps(dx12,fscal,fix1);
784 fiy1 = _mm_macc_ps(dy12,fscal,fiy1);
785 fiz1 = _mm_macc_ps(dz12,fscal,fiz1);
787 fjx2 = _mm_macc_ps(dx12,fscal,fjx2);
788 fjy2 = _mm_macc_ps(dy12,fscal,fjy2);
789 fjz2 = _mm_macc_ps(dz12,fscal,fjz2);
791 /**************************
792 * CALCULATE INTERACTIONS *
793 **************************/
795 /* REACTION-FIELD ELECTROSTATICS */
796 velec = _mm_mul_ps(qq13,_mm_sub_ps(_mm_macc_ps(krf,rsq13,rinv13),crf));
797 felec = _mm_mul_ps(qq13,_mm_msub_ps(rinv13,rinvsq13,krf2));
799 /* Update potential sum for this i atom from the interaction with this j atom. */
800 velec = _mm_andnot_ps(dummy_mask,velec);
801 velecsum = _mm_add_ps(velecsum,velec);
805 fscal = _mm_andnot_ps(dummy_mask,fscal);
807 /* Update vectorial force */
808 fix1 = _mm_macc_ps(dx13,fscal,fix1);
809 fiy1 = _mm_macc_ps(dy13,fscal,fiy1);
810 fiz1 = _mm_macc_ps(dz13,fscal,fiz1);
812 fjx3 = _mm_macc_ps(dx13,fscal,fjx3);
813 fjy3 = _mm_macc_ps(dy13,fscal,fjy3);
814 fjz3 = _mm_macc_ps(dz13,fscal,fjz3);
816 /**************************
817 * CALCULATE INTERACTIONS *
818 **************************/
820 /* REACTION-FIELD ELECTROSTATICS */
821 velec = _mm_mul_ps(qq21,_mm_sub_ps(_mm_macc_ps(krf,rsq21,rinv21),crf));
822 felec = _mm_mul_ps(qq21,_mm_msub_ps(rinv21,rinvsq21,krf2));
824 /* Update potential sum for this i atom from the interaction with this j atom. */
825 velec = _mm_andnot_ps(dummy_mask,velec);
826 velecsum = _mm_add_ps(velecsum,velec);
830 fscal = _mm_andnot_ps(dummy_mask,fscal);
832 /* Update vectorial force */
833 fix2 = _mm_macc_ps(dx21,fscal,fix2);
834 fiy2 = _mm_macc_ps(dy21,fscal,fiy2);
835 fiz2 = _mm_macc_ps(dz21,fscal,fiz2);
837 fjx1 = _mm_macc_ps(dx21,fscal,fjx1);
838 fjy1 = _mm_macc_ps(dy21,fscal,fjy1);
839 fjz1 = _mm_macc_ps(dz21,fscal,fjz1);
841 /**************************
842 * CALCULATE INTERACTIONS *
843 **************************/
845 /* REACTION-FIELD ELECTROSTATICS */
846 velec = _mm_mul_ps(qq22,_mm_sub_ps(_mm_macc_ps(krf,rsq22,rinv22),crf));
847 felec = _mm_mul_ps(qq22,_mm_msub_ps(rinv22,rinvsq22,krf2));
849 /* Update potential sum for this i atom from the interaction with this j atom. */
850 velec = _mm_andnot_ps(dummy_mask,velec);
851 velecsum = _mm_add_ps(velecsum,velec);
855 fscal = _mm_andnot_ps(dummy_mask,fscal);
857 /* Update vectorial force */
858 fix2 = _mm_macc_ps(dx22,fscal,fix2);
859 fiy2 = _mm_macc_ps(dy22,fscal,fiy2);
860 fiz2 = _mm_macc_ps(dz22,fscal,fiz2);
862 fjx2 = _mm_macc_ps(dx22,fscal,fjx2);
863 fjy2 = _mm_macc_ps(dy22,fscal,fjy2);
864 fjz2 = _mm_macc_ps(dz22,fscal,fjz2);
866 /**************************
867 * CALCULATE INTERACTIONS *
868 **************************/
870 /* REACTION-FIELD ELECTROSTATICS */
871 velec = _mm_mul_ps(qq23,_mm_sub_ps(_mm_macc_ps(krf,rsq23,rinv23),crf));
872 felec = _mm_mul_ps(qq23,_mm_msub_ps(rinv23,rinvsq23,krf2));
874 /* Update potential sum for this i atom from the interaction with this j atom. */
875 velec = _mm_andnot_ps(dummy_mask,velec);
876 velecsum = _mm_add_ps(velecsum,velec);
880 fscal = _mm_andnot_ps(dummy_mask,fscal);
882 /* Update vectorial force */
883 fix2 = _mm_macc_ps(dx23,fscal,fix2);
884 fiy2 = _mm_macc_ps(dy23,fscal,fiy2);
885 fiz2 = _mm_macc_ps(dz23,fscal,fiz2);
887 fjx3 = _mm_macc_ps(dx23,fscal,fjx3);
888 fjy3 = _mm_macc_ps(dy23,fscal,fjy3);
889 fjz3 = _mm_macc_ps(dz23,fscal,fjz3);
891 /**************************
892 * CALCULATE INTERACTIONS *
893 **************************/
895 /* REACTION-FIELD ELECTROSTATICS */
896 velec = _mm_mul_ps(qq31,_mm_sub_ps(_mm_macc_ps(krf,rsq31,rinv31),crf));
897 felec = _mm_mul_ps(qq31,_mm_msub_ps(rinv31,rinvsq31,krf2));
899 /* Update potential sum for this i atom from the interaction with this j atom. */
900 velec = _mm_andnot_ps(dummy_mask,velec);
901 velecsum = _mm_add_ps(velecsum,velec);
905 fscal = _mm_andnot_ps(dummy_mask,fscal);
907 /* Update vectorial force */
908 fix3 = _mm_macc_ps(dx31,fscal,fix3);
909 fiy3 = _mm_macc_ps(dy31,fscal,fiy3);
910 fiz3 = _mm_macc_ps(dz31,fscal,fiz3);
912 fjx1 = _mm_macc_ps(dx31,fscal,fjx1);
913 fjy1 = _mm_macc_ps(dy31,fscal,fjy1);
914 fjz1 = _mm_macc_ps(dz31,fscal,fjz1);
916 /**************************
917 * CALCULATE INTERACTIONS *
918 **************************/
920 /* REACTION-FIELD ELECTROSTATICS */
921 velec = _mm_mul_ps(qq32,_mm_sub_ps(_mm_macc_ps(krf,rsq32,rinv32),crf));
922 felec = _mm_mul_ps(qq32,_mm_msub_ps(rinv32,rinvsq32,krf2));
924 /* Update potential sum for this i atom from the interaction with this j atom. */
925 velec = _mm_andnot_ps(dummy_mask,velec);
926 velecsum = _mm_add_ps(velecsum,velec);
930 fscal = _mm_andnot_ps(dummy_mask,fscal);
932 /* Update vectorial force */
933 fix3 = _mm_macc_ps(dx32,fscal,fix3);
934 fiy3 = _mm_macc_ps(dy32,fscal,fiy3);
935 fiz3 = _mm_macc_ps(dz32,fscal,fiz3);
937 fjx2 = _mm_macc_ps(dx32,fscal,fjx2);
938 fjy2 = _mm_macc_ps(dy32,fscal,fjy2);
939 fjz2 = _mm_macc_ps(dz32,fscal,fjz2);
941 /**************************
942 * CALCULATE INTERACTIONS *
943 **************************/
945 /* REACTION-FIELD ELECTROSTATICS */
946 velec = _mm_mul_ps(qq33,_mm_sub_ps(_mm_macc_ps(krf,rsq33,rinv33),crf));
947 felec = _mm_mul_ps(qq33,_mm_msub_ps(rinv33,rinvsq33,krf2));
949 /* Update potential sum for this i atom from the interaction with this j atom. */
950 velec = _mm_andnot_ps(dummy_mask,velec);
951 velecsum = _mm_add_ps(velecsum,velec);
955 fscal = _mm_andnot_ps(dummy_mask,fscal);
957 /* Update vectorial force */
958 fix3 = _mm_macc_ps(dx33,fscal,fix3);
959 fiy3 = _mm_macc_ps(dy33,fscal,fiy3);
960 fiz3 = _mm_macc_ps(dz33,fscal,fiz3);
962 fjx3 = _mm_macc_ps(dx33,fscal,fjx3);
963 fjy3 = _mm_macc_ps(dy33,fscal,fjy3);
964 fjz3 = _mm_macc_ps(dz33,fscal,fjz3);
966 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
967 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
968 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
969 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
971 gmx_mm_decrement_4rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
972 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
973 fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
975 /* Inner loop uses 378 flops */
978 /* End of innermost loop */
980 gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
981 f+i_coord_offset,fshift+i_shift_offset);
984 /* Update potential energies */
985 gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
986 gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
988 /* Increment number of inner iterations */
989 inneriter += j_index_end - j_index_start;
991 /* Outer loop uses 26 flops */
994 /* Increment number of outer iterations */
997 /* Update outer/inner flops */
999 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4W4_VF,outeriter*26 + inneriter*378);
1002 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwCSTab_GeomW4W4_F_avx_128_fma_single
1003 * Electrostatics interaction: ReactionField
1004 * VdW interaction: CubicSplineTable
1005 * Geometry: Water4-Water4
1006 * Calculate force/pot: Force
1009 nb_kernel_ElecRF_VdwCSTab_GeomW4W4_F_avx_128_fma_single
1010 (t_nblist * gmx_restrict nlist,
1011 rvec * gmx_restrict xx,
1012 rvec * gmx_restrict ff,
1013 t_forcerec * gmx_restrict fr,
1014 t_mdatoms * gmx_restrict mdatoms,
1015 nb_kernel_data_t * gmx_restrict kernel_data,
1016 t_nrnb * gmx_restrict nrnb)
1018 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
1019 * just 0 for non-waters.
1020 * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
1021 * jnr indices corresponding to data put in the four positions in the SIMD register.
1023 int i_shift_offset,i_coord_offset,outeriter,inneriter;
1024 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
1025 int jnrA,jnrB,jnrC,jnrD;
1026 int jnrlistA,jnrlistB,jnrlistC,jnrlistD;
1027 int j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
1028 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
1029 real rcutoff_scalar;
1030 real *shiftvec,*fshift,*x,*f;
1031 real *fjptrA,*fjptrB,*fjptrC,*fjptrD;
1032 real scratch[4*DIM];
1033 __m128 fscal,rcutoff,rcutoff2,jidxall;
1035 __m128 ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
1037 __m128 ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
1039 __m128 ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
1041 __m128 ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
1042 int vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
1043 __m128 jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
1044 int vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D;
1045 __m128 jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
1046 int vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D;
1047 __m128 jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
1048 int vdwjidx3A,vdwjidx3B,vdwjidx3C,vdwjidx3D;
1049 __m128 jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
1050 __m128 dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
1051 __m128 dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
1052 __m128 dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
1053 __m128 dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
1054 __m128 dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
1055 __m128 dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
1056 __m128 dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
1057 __m128 dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
1058 __m128 dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
1059 __m128 dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
1060 __m128 velec,felec,velecsum,facel,crf,krf,krf2;
1063 __m128 rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
1066 __m128 one_sixth = _mm_set1_ps(1.0/6.0);
1067 __m128 one_twelfth = _mm_set1_ps(1.0/12.0);
1069 __m128i ifour = _mm_set1_epi32(4);
1070 __m128 rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
1072 __m128 dummy_mask,cutoff_mask;
1073 __m128 signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
1074 __m128 one = _mm_set1_ps(1.0);
1075 __m128 two = _mm_set1_ps(2.0);
1081 jindex = nlist->jindex;
1083 shiftidx = nlist->shift;
1085 shiftvec = fr->shift_vec[0];
1086 fshift = fr->fshift[0];
1087 facel = _mm_set1_ps(fr->epsfac);
1088 charge = mdatoms->chargeA;
1089 krf = _mm_set1_ps(fr->ic->k_rf);
1090 krf2 = _mm_set1_ps(fr->ic->k_rf*2.0);
1091 crf = _mm_set1_ps(fr->ic->c_rf);
1092 nvdwtype = fr->ntype;
1093 vdwparam = fr->nbfp;
1094 vdwtype = mdatoms->typeA;
1096 vftab = kernel_data->table_vdw->data;
1097 vftabscale = _mm_set1_ps(kernel_data->table_vdw->scale);
1099 /* Setup water-specific parameters */
1100 inr = nlist->iinr[0];
1101 iq1 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
1102 iq2 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
1103 iq3 = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
1104 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
1106 jq1 = _mm_set1_ps(charge[inr+1]);
1107 jq2 = _mm_set1_ps(charge[inr+2]);
1108 jq3 = _mm_set1_ps(charge[inr+3]);
1109 vdwjidx0A = 2*vdwtype[inr+0];
1110 c6_00 = _mm_set1_ps(vdwparam[vdwioffset0+vdwjidx0A]);
1111 c12_00 = _mm_set1_ps(vdwparam[vdwioffset0+vdwjidx0A+1]);
1112 qq11 = _mm_mul_ps(iq1,jq1);
1113 qq12 = _mm_mul_ps(iq1,jq2);
1114 qq13 = _mm_mul_ps(iq1,jq3);
1115 qq21 = _mm_mul_ps(iq2,jq1);
1116 qq22 = _mm_mul_ps(iq2,jq2);
1117 qq23 = _mm_mul_ps(iq2,jq3);
1118 qq31 = _mm_mul_ps(iq3,jq1);
1119 qq32 = _mm_mul_ps(iq3,jq2);
1120 qq33 = _mm_mul_ps(iq3,jq3);
1122 /* Avoid stupid compiler warnings */
1123 jnrA = jnrB = jnrC = jnrD = 0;
1124 j_coord_offsetA = 0;
1125 j_coord_offsetB = 0;
1126 j_coord_offsetC = 0;
1127 j_coord_offsetD = 0;
1132 for(iidx=0;iidx<4*DIM;iidx++)
1134 scratch[iidx] = 0.0;
1137 /* Start outer loop over neighborlists */
1138 for(iidx=0; iidx<nri; iidx++)
1140 /* Load shift vector for this list */
1141 i_shift_offset = DIM*shiftidx[iidx];
1143 /* Load limits for loop over neighbors */
1144 j_index_start = jindex[iidx];
1145 j_index_end = jindex[iidx+1];
1147 /* Get outer coordinate index */
1149 i_coord_offset = DIM*inr;
1151 /* Load i particle coords and add shift vector */
1152 gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
1153 &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
1155 fix0 = _mm_setzero_ps();
1156 fiy0 = _mm_setzero_ps();
1157 fiz0 = _mm_setzero_ps();
1158 fix1 = _mm_setzero_ps();
1159 fiy1 = _mm_setzero_ps();
1160 fiz1 = _mm_setzero_ps();
1161 fix2 = _mm_setzero_ps();
1162 fiy2 = _mm_setzero_ps();
1163 fiz2 = _mm_setzero_ps();
1164 fix3 = _mm_setzero_ps();
1165 fiy3 = _mm_setzero_ps();
1166 fiz3 = _mm_setzero_ps();
1168 /* Start inner kernel loop */
1169 for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
1172 /* Get j neighbor index, and coordinate index */
1174 jnrB = jjnr[jidx+1];
1175 jnrC = jjnr[jidx+2];
1176 jnrD = jjnr[jidx+3];
1177 j_coord_offsetA = DIM*jnrA;
1178 j_coord_offsetB = DIM*jnrB;
1179 j_coord_offsetC = DIM*jnrC;
1180 j_coord_offsetD = DIM*jnrD;
1182 /* load j atom coordinates */
1183 gmx_mm_load_4rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1184 x+j_coord_offsetC,x+j_coord_offsetD,
1185 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
1186 &jy2,&jz2,&jx3,&jy3,&jz3);
1188 /* Calculate displacement vector */
1189 dx00 = _mm_sub_ps(ix0,jx0);
1190 dy00 = _mm_sub_ps(iy0,jy0);
1191 dz00 = _mm_sub_ps(iz0,jz0);
1192 dx11 = _mm_sub_ps(ix1,jx1);
1193 dy11 = _mm_sub_ps(iy1,jy1);
1194 dz11 = _mm_sub_ps(iz1,jz1);
1195 dx12 = _mm_sub_ps(ix1,jx2);
1196 dy12 = _mm_sub_ps(iy1,jy2);
1197 dz12 = _mm_sub_ps(iz1,jz2);
1198 dx13 = _mm_sub_ps(ix1,jx3);
1199 dy13 = _mm_sub_ps(iy1,jy3);
1200 dz13 = _mm_sub_ps(iz1,jz3);
1201 dx21 = _mm_sub_ps(ix2,jx1);
1202 dy21 = _mm_sub_ps(iy2,jy1);
1203 dz21 = _mm_sub_ps(iz2,jz1);
1204 dx22 = _mm_sub_ps(ix2,jx2);
1205 dy22 = _mm_sub_ps(iy2,jy2);
1206 dz22 = _mm_sub_ps(iz2,jz2);
1207 dx23 = _mm_sub_ps(ix2,jx3);
1208 dy23 = _mm_sub_ps(iy2,jy3);
1209 dz23 = _mm_sub_ps(iz2,jz3);
1210 dx31 = _mm_sub_ps(ix3,jx1);
1211 dy31 = _mm_sub_ps(iy3,jy1);
1212 dz31 = _mm_sub_ps(iz3,jz1);
1213 dx32 = _mm_sub_ps(ix3,jx2);
1214 dy32 = _mm_sub_ps(iy3,jy2);
1215 dz32 = _mm_sub_ps(iz3,jz2);
1216 dx33 = _mm_sub_ps(ix3,jx3);
1217 dy33 = _mm_sub_ps(iy3,jy3);
1218 dz33 = _mm_sub_ps(iz3,jz3);
1220 /* Calculate squared distance and things based on it */
1221 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
1222 rsq11 = gmx_mm_calc_rsq_ps(dx11,dy11,dz11);
1223 rsq12 = gmx_mm_calc_rsq_ps(dx12,dy12,dz12);
1224 rsq13 = gmx_mm_calc_rsq_ps(dx13,dy13,dz13);
1225 rsq21 = gmx_mm_calc_rsq_ps(dx21,dy21,dz21);
1226 rsq22 = gmx_mm_calc_rsq_ps(dx22,dy22,dz22);
1227 rsq23 = gmx_mm_calc_rsq_ps(dx23,dy23,dz23);
1228 rsq31 = gmx_mm_calc_rsq_ps(dx31,dy31,dz31);
1229 rsq32 = gmx_mm_calc_rsq_ps(dx32,dy32,dz32);
1230 rsq33 = gmx_mm_calc_rsq_ps(dx33,dy33,dz33);
1232 rinv00 = gmx_mm_invsqrt_ps(rsq00);
1233 rinv11 = gmx_mm_invsqrt_ps(rsq11);
1234 rinv12 = gmx_mm_invsqrt_ps(rsq12);
1235 rinv13 = gmx_mm_invsqrt_ps(rsq13);
1236 rinv21 = gmx_mm_invsqrt_ps(rsq21);
1237 rinv22 = gmx_mm_invsqrt_ps(rsq22);
1238 rinv23 = gmx_mm_invsqrt_ps(rsq23);
1239 rinv31 = gmx_mm_invsqrt_ps(rsq31);
1240 rinv32 = gmx_mm_invsqrt_ps(rsq32);
1241 rinv33 = gmx_mm_invsqrt_ps(rsq33);
1243 rinvsq11 = _mm_mul_ps(rinv11,rinv11);
1244 rinvsq12 = _mm_mul_ps(rinv12,rinv12);
1245 rinvsq13 = _mm_mul_ps(rinv13,rinv13);
1246 rinvsq21 = _mm_mul_ps(rinv21,rinv21);
1247 rinvsq22 = _mm_mul_ps(rinv22,rinv22);
1248 rinvsq23 = _mm_mul_ps(rinv23,rinv23);
1249 rinvsq31 = _mm_mul_ps(rinv31,rinv31);
1250 rinvsq32 = _mm_mul_ps(rinv32,rinv32);
1251 rinvsq33 = _mm_mul_ps(rinv33,rinv33);
1253 fjx0 = _mm_setzero_ps();
1254 fjy0 = _mm_setzero_ps();
1255 fjz0 = _mm_setzero_ps();
1256 fjx1 = _mm_setzero_ps();
1257 fjy1 = _mm_setzero_ps();
1258 fjz1 = _mm_setzero_ps();
1259 fjx2 = _mm_setzero_ps();
1260 fjy2 = _mm_setzero_ps();
1261 fjz2 = _mm_setzero_ps();
1262 fjx3 = _mm_setzero_ps();
1263 fjy3 = _mm_setzero_ps();
1264 fjz3 = _mm_setzero_ps();
1266 /**************************
1267 * CALCULATE INTERACTIONS *
1268 **************************/
1270 r00 = _mm_mul_ps(rsq00,rinv00);
1272 /* Calculate table index by multiplying r with table scale and truncate to integer */
1273 rt = _mm_mul_ps(r00,vftabscale);
1274 vfitab = _mm_cvttps_epi32(rt);
1276 vfeps = _mm_frcz_ps(rt);
1278 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
1280 twovfeps = _mm_add_ps(vfeps,vfeps);
1281 vfitab = _mm_slli_epi32(vfitab,3);
1283 /* CUBIC SPLINE TABLE DISPERSION */
1284 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
1285 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
1286 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
1287 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
1288 _MM_TRANSPOSE4_PS(Y,F,G,H);
1289 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
1290 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
1291 fvdw6 = _mm_mul_ps(c6_00,FF);
1293 /* CUBIC SPLINE TABLE REPULSION */
1294 vfitab = _mm_add_epi32(vfitab,ifour);
1295 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
1296 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
1297 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
1298 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
1299 _MM_TRANSPOSE4_PS(Y,F,G,H);
1300 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
1301 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
1302 fvdw12 = _mm_mul_ps(c12_00,FF);
1303 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
1307 /* Update vectorial force */
1308 fix0 = _mm_macc_ps(dx00,fscal,fix0);
1309 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
1310 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
1312 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
1313 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
1314 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
1316 /**************************
1317 * CALCULATE INTERACTIONS *
1318 **************************/
1320 /* REACTION-FIELD ELECTROSTATICS */
1321 felec = _mm_mul_ps(qq11,_mm_msub_ps(rinv11,rinvsq11,krf2));
1325 /* Update vectorial force */
1326 fix1 = _mm_macc_ps(dx11,fscal,fix1);
1327 fiy1 = _mm_macc_ps(dy11,fscal,fiy1);
1328 fiz1 = _mm_macc_ps(dz11,fscal,fiz1);
1330 fjx1 = _mm_macc_ps(dx11,fscal,fjx1);
1331 fjy1 = _mm_macc_ps(dy11,fscal,fjy1);
1332 fjz1 = _mm_macc_ps(dz11,fscal,fjz1);
1334 /**************************
1335 * CALCULATE INTERACTIONS *
1336 **************************/
1338 /* REACTION-FIELD ELECTROSTATICS */
1339 felec = _mm_mul_ps(qq12,_mm_msub_ps(rinv12,rinvsq12,krf2));
1343 /* Update vectorial force */
1344 fix1 = _mm_macc_ps(dx12,fscal,fix1);
1345 fiy1 = _mm_macc_ps(dy12,fscal,fiy1);
1346 fiz1 = _mm_macc_ps(dz12,fscal,fiz1);
1348 fjx2 = _mm_macc_ps(dx12,fscal,fjx2);
1349 fjy2 = _mm_macc_ps(dy12,fscal,fjy2);
1350 fjz2 = _mm_macc_ps(dz12,fscal,fjz2);
1352 /**************************
1353 * CALCULATE INTERACTIONS *
1354 **************************/
1356 /* REACTION-FIELD ELECTROSTATICS */
1357 felec = _mm_mul_ps(qq13,_mm_msub_ps(rinv13,rinvsq13,krf2));
1361 /* Update vectorial force */
1362 fix1 = _mm_macc_ps(dx13,fscal,fix1);
1363 fiy1 = _mm_macc_ps(dy13,fscal,fiy1);
1364 fiz1 = _mm_macc_ps(dz13,fscal,fiz1);
1366 fjx3 = _mm_macc_ps(dx13,fscal,fjx3);
1367 fjy3 = _mm_macc_ps(dy13,fscal,fjy3);
1368 fjz3 = _mm_macc_ps(dz13,fscal,fjz3);
1370 /**************************
1371 * CALCULATE INTERACTIONS *
1372 **************************/
1374 /* REACTION-FIELD ELECTROSTATICS */
1375 felec = _mm_mul_ps(qq21,_mm_msub_ps(rinv21,rinvsq21,krf2));
1379 /* Update vectorial force */
1380 fix2 = _mm_macc_ps(dx21,fscal,fix2);
1381 fiy2 = _mm_macc_ps(dy21,fscal,fiy2);
1382 fiz2 = _mm_macc_ps(dz21,fscal,fiz2);
1384 fjx1 = _mm_macc_ps(dx21,fscal,fjx1);
1385 fjy1 = _mm_macc_ps(dy21,fscal,fjy1);
1386 fjz1 = _mm_macc_ps(dz21,fscal,fjz1);
1388 /**************************
1389 * CALCULATE INTERACTIONS *
1390 **************************/
1392 /* REACTION-FIELD ELECTROSTATICS */
1393 felec = _mm_mul_ps(qq22,_mm_msub_ps(rinv22,rinvsq22,krf2));
1397 /* Update vectorial force */
1398 fix2 = _mm_macc_ps(dx22,fscal,fix2);
1399 fiy2 = _mm_macc_ps(dy22,fscal,fiy2);
1400 fiz2 = _mm_macc_ps(dz22,fscal,fiz2);
1402 fjx2 = _mm_macc_ps(dx22,fscal,fjx2);
1403 fjy2 = _mm_macc_ps(dy22,fscal,fjy2);
1404 fjz2 = _mm_macc_ps(dz22,fscal,fjz2);
1406 /**************************
1407 * CALCULATE INTERACTIONS *
1408 **************************/
1410 /* REACTION-FIELD ELECTROSTATICS */
1411 felec = _mm_mul_ps(qq23,_mm_msub_ps(rinv23,rinvsq23,krf2));
1415 /* Update vectorial force */
1416 fix2 = _mm_macc_ps(dx23,fscal,fix2);
1417 fiy2 = _mm_macc_ps(dy23,fscal,fiy2);
1418 fiz2 = _mm_macc_ps(dz23,fscal,fiz2);
1420 fjx3 = _mm_macc_ps(dx23,fscal,fjx3);
1421 fjy3 = _mm_macc_ps(dy23,fscal,fjy3);
1422 fjz3 = _mm_macc_ps(dz23,fscal,fjz3);
1424 /**************************
1425 * CALCULATE INTERACTIONS *
1426 **************************/
1428 /* REACTION-FIELD ELECTROSTATICS */
1429 felec = _mm_mul_ps(qq31,_mm_msub_ps(rinv31,rinvsq31,krf2));
1433 /* Update vectorial force */
1434 fix3 = _mm_macc_ps(dx31,fscal,fix3);
1435 fiy3 = _mm_macc_ps(dy31,fscal,fiy3);
1436 fiz3 = _mm_macc_ps(dz31,fscal,fiz3);
1438 fjx1 = _mm_macc_ps(dx31,fscal,fjx1);
1439 fjy1 = _mm_macc_ps(dy31,fscal,fjy1);
1440 fjz1 = _mm_macc_ps(dz31,fscal,fjz1);
1442 /**************************
1443 * CALCULATE INTERACTIONS *
1444 **************************/
1446 /* REACTION-FIELD ELECTROSTATICS */
1447 felec = _mm_mul_ps(qq32,_mm_msub_ps(rinv32,rinvsq32,krf2));
1451 /* Update vectorial force */
1452 fix3 = _mm_macc_ps(dx32,fscal,fix3);
1453 fiy3 = _mm_macc_ps(dy32,fscal,fiy3);
1454 fiz3 = _mm_macc_ps(dz32,fscal,fiz3);
1456 fjx2 = _mm_macc_ps(dx32,fscal,fjx2);
1457 fjy2 = _mm_macc_ps(dy32,fscal,fjy2);
1458 fjz2 = _mm_macc_ps(dz32,fscal,fjz2);
1460 /**************************
1461 * CALCULATE INTERACTIONS *
1462 **************************/
1464 /* REACTION-FIELD ELECTROSTATICS */
1465 felec = _mm_mul_ps(qq33,_mm_msub_ps(rinv33,rinvsq33,krf2));
1469 /* Update vectorial force */
1470 fix3 = _mm_macc_ps(dx33,fscal,fix3);
1471 fiy3 = _mm_macc_ps(dy33,fscal,fiy3);
1472 fiz3 = _mm_macc_ps(dz33,fscal,fiz3);
1474 fjx3 = _mm_macc_ps(dx33,fscal,fjx3);
1475 fjy3 = _mm_macc_ps(dy33,fscal,fjy3);
1476 fjz3 = _mm_macc_ps(dz33,fscal,fjz3);
1478 fjptrA = f+j_coord_offsetA;
1479 fjptrB = f+j_coord_offsetB;
1480 fjptrC = f+j_coord_offsetC;
1481 fjptrD = f+j_coord_offsetD;
1483 gmx_mm_decrement_4rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
1484 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
1485 fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1487 /* Inner loop uses 324 flops */
1490 if(jidx<j_index_end)
1493 /* Get j neighbor index, and coordinate index */
1494 jnrlistA = jjnr[jidx];
1495 jnrlistB = jjnr[jidx+1];
1496 jnrlistC = jjnr[jidx+2];
1497 jnrlistD = jjnr[jidx+3];
1498 /* Sign of each element will be negative for non-real atoms.
1499 * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1500 * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1502 dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
1503 jnrA = (jnrlistA>=0) ? jnrlistA : 0;
1504 jnrB = (jnrlistB>=0) ? jnrlistB : 0;
1505 jnrC = (jnrlistC>=0) ? jnrlistC : 0;
1506 jnrD = (jnrlistD>=0) ? jnrlistD : 0;
1507 j_coord_offsetA = DIM*jnrA;
1508 j_coord_offsetB = DIM*jnrB;
1509 j_coord_offsetC = DIM*jnrC;
1510 j_coord_offsetD = DIM*jnrD;
1512 /* load j atom coordinates */
1513 gmx_mm_load_4rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1514 x+j_coord_offsetC,x+j_coord_offsetD,
1515 &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
1516 &jy2,&jz2,&jx3,&jy3,&jz3);
1518 /* Calculate displacement vector */
1519 dx00 = _mm_sub_ps(ix0,jx0);
1520 dy00 = _mm_sub_ps(iy0,jy0);
1521 dz00 = _mm_sub_ps(iz0,jz0);
1522 dx11 = _mm_sub_ps(ix1,jx1);
1523 dy11 = _mm_sub_ps(iy1,jy1);
1524 dz11 = _mm_sub_ps(iz1,jz1);
1525 dx12 = _mm_sub_ps(ix1,jx2);
1526 dy12 = _mm_sub_ps(iy1,jy2);
1527 dz12 = _mm_sub_ps(iz1,jz2);
1528 dx13 = _mm_sub_ps(ix1,jx3);
1529 dy13 = _mm_sub_ps(iy1,jy3);
1530 dz13 = _mm_sub_ps(iz1,jz3);
1531 dx21 = _mm_sub_ps(ix2,jx1);
1532 dy21 = _mm_sub_ps(iy2,jy1);
1533 dz21 = _mm_sub_ps(iz2,jz1);
1534 dx22 = _mm_sub_ps(ix2,jx2);
1535 dy22 = _mm_sub_ps(iy2,jy2);
1536 dz22 = _mm_sub_ps(iz2,jz2);
1537 dx23 = _mm_sub_ps(ix2,jx3);
1538 dy23 = _mm_sub_ps(iy2,jy3);
1539 dz23 = _mm_sub_ps(iz2,jz3);
1540 dx31 = _mm_sub_ps(ix3,jx1);
1541 dy31 = _mm_sub_ps(iy3,jy1);
1542 dz31 = _mm_sub_ps(iz3,jz1);
1543 dx32 = _mm_sub_ps(ix3,jx2);
1544 dy32 = _mm_sub_ps(iy3,jy2);
1545 dz32 = _mm_sub_ps(iz3,jz2);
1546 dx33 = _mm_sub_ps(ix3,jx3);
1547 dy33 = _mm_sub_ps(iy3,jy3);
1548 dz33 = _mm_sub_ps(iz3,jz3);
1550 /* Calculate squared distance and things based on it */
1551 rsq00 = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
1552 rsq11 = gmx_mm_calc_rsq_ps(dx11,dy11,dz11);
1553 rsq12 = gmx_mm_calc_rsq_ps(dx12,dy12,dz12);
1554 rsq13 = gmx_mm_calc_rsq_ps(dx13,dy13,dz13);
1555 rsq21 = gmx_mm_calc_rsq_ps(dx21,dy21,dz21);
1556 rsq22 = gmx_mm_calc_rsq_ps(dx22,dy22,dz22);
1557 rsq23 = gmx_mm_calc_rsq_ps(dx23,dy23,dz23);
1558 rsq31 = gmx_mm_calc_rsq_ps(dx31,dy31,dz31);
1559 rsq32 = gmx_mm_calc_rsq_ps(dx32,dy32,dz32);
1560 rsq33 = gmx_mm_calc_rsq_ps(dx33,dy33,dz33);
1562 rinv00 = gmx_mm_invsqrt_ps(rsq00);
1563 rinv11 = gmx_mm_invsqrt_ps(rsq11);
1564 rinv12 = gmx_mm_invsqrt_ps(rsq12);
1565 rinv13 = gmx_mm_invsqrt_ps(rsq13);
1566 rinv21 = gmx_mm_invsqrt_ps(rsq21);
1567 rinv22 = gmx_mm_invsqrt_ps(rsq22);
1568 rinv23 = gmx_mm_invsqrt_ps(rsq23);
1569 rinv31 = gmx_mm_invsqrt_ps(rsq31);
1570 rinv32 = gmx_mm_invsqrt_ps(rsq32);
1571 rinv33 = gmx_mm_invsqrt_ps(rsq33);
1573 rinvsq11 = _mm_mul_ps(rinv11,rinv11);
1574 rinvsq12 = _mm_mul_ps(rinv12,rinv12);
1575 rinvsq13 = _mm_mul_ps(rinv13,rinv13);
1576 rinvsq21 = _mm_mul_ps(rinv21,rinv21);
1577 rinvsq22 = _mm_mul_ps(rinv22,rinv22);
1578 rinvsq23 = _mm_mul_ps(rinv23,rinv23);
1579 rinvsq31 = _mm_mul_ps(rinv31,rinv31);
1580 rinvsq32 = _mm_mul_ps(rinv32,rinv32);
1581 rinvsq33 = _mm_mul_ps(rinv33,rinv33);
1583 fjx0 = _mm_setzero_ps();
1584 fjy0 = _mm_setzero_ps();
1585 fjz0 = _mm_setzero_ps();
1586 fjx1 = _mm_setzero_ps();
1587 fjy1 = _mm_setzero_ps();
1588 fjz1 = _mm_setzero_ps();
1589 fjx2 = _mm_setzero_ps();
1590 fjy2 = _mm_setzero_ps();
1591 fjz2 = _mm_setzero_ps();
1592 fjx3 = _mm_setzero_ps();
1593 fjy3 = _mm_setzero_ps();
1594 fjz3 = _mm_setzero_ps();
1596 /**************************
1597 * CALCULATE INTERACTIONS *
1598 **************************/
1600 r00 = _mm_mul_ps(rsq00,rinv00);
1601 r00 = _mm_andnot_ps(dummy_mask,r00);
1603 /* Calculate table index by multiplying r with table scale and truncate to integer */
1604 rt = _mm_mul_ps(r00,vftabscale);
1605 vfitab = _mm_cvttps_epi32(rt);
1607 vfeps = _mm_frcz_ps(rt);
1609 vfeps = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
1611 twovfeps = _mm_add_ps(vfeps,vfeps);
1612 vfitab = _mm_slli_epi32(vfitab,3);
1614 /* CUBIC SPLINE TABLE DISPERSION */
1615 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
1616 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
1617 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
1618 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
1619 _MM_TRANSPOSE4_PS(Y,F,G,H);
1620 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
1621 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
1622 fvdw6 = _mm_mul_ps(c6_00,FF);
1624 /* CUBIC SPLINE TABLE REPULSION */
1625 vfitab = _mm_add_epi32(vfitab,ifour);
1626 Y = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
1627 F = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
1628 G = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
1629 H = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
1630 _MM_TRANSPOSE4_PS(Y,F,G,H);
1631 Fp = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
1632 FF = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
1633 fvdw12 = _mm_mul_ps(c12_00,FF);
1634 fvdw = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
1638 fscal = _mm_andnot_ps(dummy_mask,fscal);
1640 /* Update vectorial force */
1641 fix0 = _mm_macc_ps(dx00,fscal,fix0);
1642 fiy0 = _mm_macc_ps(dy00,fscal,fiy0);
1643 fiz0 = _mm_macc_ps(dz00,fscal,fiz0);
1645 fjx0 = _mm_macc_ps(dx00,fscal,fjx0);
1646 fjy0 = _mm_macc_ps(dy00,fscal,fjy0);
1647 fjz0 = _mm_macc_ps(dz00,fscal,fjz0);
1649 /**************************
1650 * CALCULATE INTERACTIONS *
1651 **************************/
1653 /* REACTION-FIELD ELECTROSTATICS */
1654 felec = _mm_mul_ps(qq11,_mm_msub_ps(rinv11,rinvsq11,krf2));
1658 fscal = _mm_andnot_ps(dummy_mask,fscal);
1660 /* Update vectorial force */
1661 fix1 = _mm_macc_ps(dx11,fscal,fix1);
1662 fiy1 = _mm_macc_ps(dy11,fscal,fiy1);
1663 fiz1 = _mm_macc_ps(dz11,fscal,fiz1);
1665 fjx1 = _mm_macc_ps(dx11,fscal,fjx1);
1666 fjy1 = _mm_macc_ps(dy11,fscal,fjy1);
1667 fjz1 = _mm_macc_ps(dz11,fscal,fjz1);
1669 /**************************
1670 * CALCULATE INTERACTIONS *
1671 **************************/
1673 /* REACTION-FIELD ELECTROSTATICS */
1674 felec = _mm_mul_ps(qq12,_mm_msub_ps(rinv12,rinvsq12,krf2));
1678 fscal = _mm_andnot_ps(dummy_mask,fscal);
1680 /* Update vectorial force */
1681 fix1 = _mm_macc_ps(dx12,fscal,fix1);
1682 fiy1 = _mm_macc_ps(dy12,fscal,fiy1);
1683 fiz1 = _mm_macc_ps(dz12,fscal,fiz1);
1685 fjx2 = _mm_macc_ps(dx12,fscal,fjx2);
1686 fjy2 = _mm_macc_ps(dy12,fscal,fjy2);
1687 fjz2 = _mm_macc_ps(dz12,fscal,fjz2);
1689 /**************************
1690 * CALCULATE INTERACTIONS *
1691 **************************/
1693 /* REACTION-FIELD ELECTROSTATICS */
1694 felec = _mm_mul_ps(qq13,_mm_msub_ps(rinv13,rinvsq13,krf2));
1698 fscal = _mm_andnot_ps(dummy_mask,fscal);
1700 /* Update vectorial force */
1701 fix1 = _mm_macc_ps(dx13,fscal,fix1);
1702 fiy1 = _mm_macc_ps(dy13,fscal,fiy1);
1703 fiz1 = _mm_macc_ps(dz13,fscal,fiz1);
1705 fjx3 = _mm_macc_ps(dx13,fscal,fjx3);
1706 fjy3 = _mm_macc_ps(dy13,fscal,fjy3);
1707 fjz3 = _mm_macc_ps(dz13,fscal,fjz3);
1709 /**************************
1710 * CALCULATE INTERACTIONS *
1711 **************************/
1713 /* REACTION-FIELD ELECTROSTATICS */
1714 felec = _mm_mul_ps(qq21,_mm_msub_ps(rinv21,rinvsq21,krf2));
1718 fscal = _mm_andnot_ps(dummy_mask,fscal);
1720 /* Update vectorial force */
1721 fix2 = _mm_macc_ps(dx21,fscal,fix2);
1722 fiy2 = _mm_macc_ps(dy21,fscal,fiy2);
1723 fiz2 = _mm_macc_ps(dz21,fscal,fiz2);
1725 fjx1 = _mm_macc_ps(dx21,fscal,fjx1);
1726 fjy1 = _mm_macc_ps(dy21,fscal,fjy1);
1727 fjz1 = _mm_macc_ps(dz21,fscal,fjz1);
1729 /**************************
1730 * CALCULATE INTERACTIONS *
1731 **************************/
1733 /* REACTION-FIELD ELECTROSTATICS */
1734 felec = _mm_mul_ps(qq22,_mm_msub_ps(rinv22,rinvsq22,krf2));
1738 fscal = _mm_andnot_ps(dummy_mask,fscal);
1740 /* Update vectorial force */
1741 fix2 = _mm_macc_ps(dx22,fscal,fix2);
1742 fiy2 = _mm_macc_ps(dy22,fscal,fiy2);
1743 fiz2 = _mm_macc_ps(dz22,fscal,fiz2);
1745 fjx2 = _mm_macc_ps(dx22,fscal,fjx2);
1746 fjy2 = _mm_macc_ps(dy22,fscal,fjy2);
1747 fjz2 = _mm_macc_ps(dz22,fscal,fjz2);
1749 /**************************
1750 * CALCULATE INTERACTIONS *
1751 **************************/
1753 /* REACTION-FIELD ELECTROSTATICS */
1754 felec = _mm_mul_ps(qq23,_mm_msub_ps(rinv23,rinvsq23,krf2));
1758 fscal = _mm_andnot_ps(dummy_mask,fscal);
1760 /* Update vectorial force */
1761 fix2 = _mm_macc_ps(dx23,fscal,fix2);
1762 fiy2 = _mm_macc_ps(dy23,fscal,fiy2);
1763 fiz2 = _mm_macc_ps(dz23,fscal,fiz2);
1765 fjx3 = _mm_macc_ps(dx23,fscal,fjx3);
1766 fjy3 = _mm_macc_ps(dy23,fscal,fjy3);
1767 fjz3 = _mm_macc_ps(dz23,fscal,fjz3);
1769 /**************************
1770 * CALCULATE INTERACTIONS *
1771 **************************/
1773 /* REACTION-FIELD ELECTROSTATICS */
1774 felec = _mm_mul_ps(qq31,_mm_msub_ps(rinv31,rinvsq31,krf2));
1778 fscal = _mm_andnot_ps(dummy_mask,fscal);
1780 /* Update vectorial force */
1781 fix3 = _mm_macc_ps(dx31,fscal,fix3);
1782 fiy3 = _mm_macc_ps(dy31,fscal,fiy3);
1783 fiz3 = _mm_macc_ps(dz31,fscal,fiz3);
1785 fjx1 = _mm_macc_ps(dx31,fscal,fjx1);
1786 fjy1 = _mm_macc_ps(dy31,fscal,fjy1);
1787 fjz1 = _mm_macc_ps(dz31,fscal,fjz1);
1789 /**************************
1790 * CALCULATE INTERACTIONS *
1791 **************************/
1793 /* REACTION-FIELD ELECTROSTATICS */
1794 felec = _mm_mul_ps(qq32,_mm_msub_ps(rinv32,rinvsq32,krf2));
1798 fscal = _mm_andnot_ps(dummy_mask,fscal);
1800 /* Update vectorial force */
1801 fix3 = _mm_macc_ps(dx32,fscal,fix3);
1802 fiy3 = _mm_macc_ps(dy32,fscal,fiy3);
1803 fiz3 = _mm_macc_ps(dz32,fscal,fiz3);
1805 fjx2 = _mm_macc_ps(dx32,fscal,fjx2);
1806 fjy2 = _mm_macc_ps(dy32,fscal,fjy2);
1807 fjz2 = _mm_macc_ps(dz32,fscal,fjz2);
1809 /**************************
1810 * CALCULATE INTERACTIONS *
1811 **************************/
1813 /* REACTION-FIELD ELECTROSTATICS */
1814 felec = _mm_mul_ps(qq33,_mm_msub_ps(rinv33,rinvsq33,krf2));
1818 fscal = _mm_andnot_ps(dummy_mask,fscal);
1820 /* Update vectorial force */
1821 fix3 = _mm_macc_ps(dx33,fscal,fix3);
1822 fiy3 = _mm_macc_ps(dy33,fscal,fiy3);
1823 fiz3 = _mm_macc_ps(dz33,fscal,fiz3);
1825 fjx3 = _mm_macc_ps(dx33,fscal,fjx3);
1826 fjy3 = _mm_macc_ps(dy33,fscal,fjy3);
1827 fjz3 = _mm_macc_ps(dz33,fscal,fjz3);
1829 fjptrA = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1830 fjptrB = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1831 fjptrC = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1832 fjptrD = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1834 gmx_mm_decrement_4rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
1835 fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
1836 fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1838 /* Inner loop uses 325 flops */
1841 /* End of innermost loop */
1843 gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1844 f+i_coord_offset,fshift+i_shift_offset);
1846 /* Increment number of inner iterations */
1847 inneriter += j_index_end - j_index_start;
1849 /* Outer loop uses 24 flops */
1852 /* Increment number of outer iterations */
1855 /* Update outer/inner flops */
1857 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4W4_F,outeriter*24 + inneriter*325);