2 * Note: this file was generated by the Gromacs avx_128_fma_double 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_double.h"
34 #include "kernelutil_x86_avx_128_fma_double.h"
37 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwCSTab_GeomP1P1_VF_avx_128_fma_double
38 * Electrostatics interaction: ReactionField
39 * VdW interaction: CubicSplineTable
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecRF_VdwCSTab_GeomP1P1_VF_avx_128_fma_double
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 refer to j loop unrolling done with SSE double precision, e.g. for the two 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;
61 int j_coord_offsetA,j_coord_offsetB;
62 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
64 real *shiftvec,*fshift,*x,*f;
65 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
67 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
68 int vdwjidx0A,vdwjidx0B;
69 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
70 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
71 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
74 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
77 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
78 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
80 __m128i ifour = _mm_set1_epi32(4);
81 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
83 __m128d dummy_mask,cutoff_mask;
84 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
85 __m128d one = _mm_set1_pd(1.0);
86 __m128d two = _mm_set1_pd(2.0);
92 jindex = nlist->jindex;
94 shiftidx = nlist->shift;
96 shiftvec = fr->shift_vec[0];
97 fshift = fr->fshift[0];
98 facel = _mm_set1_pd(fr->epsfac);
99 charge = mdatoms->chargeA;
100 krf = _mm_set1_pd(fr->ic->k_rf);
101 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
102 crf = _mm_set1_pd(fr->ic->c_rf);
103 nvdwtype = fr->ntype;
105 vdwtype = mdatoms->typeA;
107 vftab = kernel_data->table_vdw->data;
108 vftabscale = _mm_set1_pd(kernel_data->table_vdw->scale);
110 /* Avoid stupid compiler warnings */
118 /* Start outer loop over neighborlists */
119 for(iidx=0; iidx<nri; iidx++)
121 /* Load shift vector for this list */
122 i_shift_offset = DIM*shiftidx[iidx];
124 /* Load limits for loop over neighbors */
125 j_index_start = jindex[iidx];
126 j_index_end = jindex[iidx+1];
128 /* Get outer coordinate index */
130 i_coord_offset = DIM*inr;
132 /* Load i particle coords and add shift vector */
133 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
135 fix0 = _mm_setzero_pd();
136 fiy0 = _mm_setzero_pd();
137 fiz0 = _mm_setzero_pd();
139 /* Load parameters for i particles */
140 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
141 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
143 /* Reset potential sums */
144 velecsum = _mm_setzero_pd();
145 vvdwsum = _mm_setzero_pd();
147 /* Start inner kernel loop */
148 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
151 /* Get j neighbor index, and coordinate index */
154 j_coord_offsetA = DIM*jnrA;
155 j_coord_offsetB = DIM*jnrB;
157 /* load j atom coordinates */
158 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
161 /* Calculate displacement vector */
162 dx00 = _mm_sub_pd(ix0,jx0);
163 dy00 = _mm_sub_pd(iy0,jy0);
164 dz00 = _mm_sub_pd(iz0,jz0);
166 /* Calculate squared distance and things based on it */
167 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
169 rinv00 = gmx_mm_invsqrt_pd(rsq00);
171 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
173 /* Load parameters for j particles */
174 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
175 vdwjidx0A = 2*vdwtype[jnrA+0];
176 vdwjidx0B = 2*vdwtype[jnrB+0];
178 /**************************
179 * CALCULATE INTERACTIONS *
180 **************************/
182 r00 = _mm_mul_pd(rsq00,rinv00);
184 /* Compute parameters for interactions between i and j atoms */
185 qq00 = _mm_mul_pd(iq0,jq0);
186 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
187 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
189 /* Calculate table index by multiplying r with table scale and truncate to integer */
190 rt = _mm_mul_pd(r00,vftabscale);
191 vfitab = _mm_cvttpd_epi32(rt);
193 vfeps = _mm_frcz_pd(rt);
195 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
197 twovfeps = _mm_add_pd(vfeps,vfeps);
198 vfitab = _mm_slli_epi32(vfitab,3);
200 /* REACTION-FIELD ELECTROSTATICS */
201 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_macc_pd(krf,rsq00,rinv00),crf));
202 felec = _mm_mul_pd(qq00,_mm_msub_pd(rinv00,rinvsq00,krf2));
204 /* CUBIC SPLINE TABLE DISPERSION */
205 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
206 F = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
207 GMX_MM_TRANSPOSE2_PD(Y,F);
208 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
209 H = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
210 GMX_MM_TRANSPOSE2_PD(G,H);
211 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
212 VV = _mm_macc_pd(vfeps,Fp,Y);
213 vvdw6 = _mm_mul_pd(c6_00,VV);
214 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
215 fvdw6 = _mm_mul_pd(c6_00,FF);
217 /* CUBIC SPLINE TABLE REPULSION */
218 vfitab = _mm_add_epi32(vfitab,ifour);
219 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
220 F = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
221 GMX_MM_TRANSPOSE2_PD(Y,F);
222 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
223 H = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
224 GMX_MM_TRANSPOSE2_PD(G,H);
225 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
226 VV = _mm_macc_pd(vfeps,Fp,Y);
227 vvdw12 = _mm_mul_pd(c12_00,VV);
228 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
229 fvdw12 = _mm_mul_pd(c12_00,FF);
230 vvdw = _mm_add_pd(vvdw12,vvdw6);
231 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
233 /* Update potential sum for this i atom from the interaction with this j atom. */
234 velecsum = _mm_add_pd(velecsum,velec);
235 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
237 fscal = _mm_add_pd(felec,fvdw);
239 /* Update vectorial force */
240 fix0 = _mm_macc_pd(dx00,fscal,fix0);
241 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
242 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
244 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
245 _mm_mul_pd(dx00,fscal),
246 _mm_mul_pd(dy00,fscal),
247 _mm_mul_pd(dz00,fscal));
249 /* Inner loop uses 70 flops */
256 j_coord_offsetA = DIM*jnrA;
258 /* load j atom coordinates */
259 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
262 /* Calculate displacement vector */
263 dx00 = _mm_sub_pd(ix0,jx0);
264 dy00 = _mm_sub_pd(iy0,jy0);
265 dz00 = _mm_sub_pd(iz0,jz0);
267 /* Calculate squared distance and things based on it */
268 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
270 rinv00 = gmx_mm_invsqrt_pd(rsq00);
272 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
274 /* Load parameters for j particles */
275 jq0 = _mm_load_sd(charge+jnrA+0);
276 vdwjidx0A = 2*vdwtype[jnrA+0];
278 /**************************
279 * CALCULATE INTERACTIONS *
280 **************************/
282 r00 = _mm_mul_pd(rsq00,rinv00);
284 /* Compute parameters for interactions between i and j atoms */
285 qq00 = _mm_mul_pd(iq0,jq0);
286 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
288 /* Calculate table index by multiplying r with table scale and truncate to integer */
289 rt = _mm_mul_pd(r00,vftabscale);
290 vfitab = _mm_cvttpd_epi32(rt);
292 vfeps = _mm_frcz_pd(rt);
294 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
296 twovfeps = _mm_add_pd(vfeps,vfeps);
297 vfitab = _mm_slli_epi32(vfitab,3);
299 /* REACTION-FIELD ELECTROSTATICS */
300 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_macc_pd(krf,rsq00,rinv00),crf));
301 felec = _mm_mul_pd(qq00,_mm_msub_pd(rinv00,rinvsq00,krf2));
303 /* CUBIC SPLINE TABLE DISPERSION */
304 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
305 F = _mm_setzero_pd();
306 GMX_MM_TRANSPOSE2_PD(Y,F);
307 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
308 H = _mm_setzero_pd();
309 GMX_MM_TRANSPOSE2_PD(G,H);
310 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
311 VV = _mm_macc_pd(vfeps,Fp,Y);
312 vvdw6 = _mm_mul_pd(c6_00,VV);
313 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
314 fvdw6 = _mm_mul_pd(c6_00,FF);
316 /* CUBIC SPLINE TABLE REPULSION */
317 vfitab = _mm_add_epi32(vfitab,ifour);
318 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
319 F = _mm_setzero_pd();
320 GMX_MM_TRANSPOSE2_PD(Y,F);
321 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
322 H = _mm_setzero_pd();
323 GMX_MM_TRANSPOSE2_PD(G,H);
324 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
325 VV = _mm_macc_pd(vfeps,Fp,Y);
326 vvdw12 = _mm_mul_pd(c12_00,VV);
327 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
328 fvdw12 = _mm_mul_pd(c12_00,FF);
329 vvdw = _mm_add_pd(vvdw12,vvdw6);
330 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
332 /* Update potential sum for this i atom from the interaction with this j atom. */
333 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
334 velecsum = _mm_add_pd(velecsum,velec);
335 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
336 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
338 fscal = _mm_add_pd(felec,fvdw);
340 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
342 /* Update vectorial force */
343 fix0 = _mm_macc_pd(dx00,fscal,fix0);
344 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
345 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
347 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
348 _mm_mul_pd(dx00,fscal),
349 _mm_mul_pd(dy00,fscal),
350 _mm_mul_pd(dz00,fscal));
352 /* Inner loop uses 70 flops */
355 /* End of innermost loop */
357 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
358 f+i_coord_offset,fshift+i_shift_offset);
361 /* Update potential energies */
362 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
363 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
365 /* Increment number of inner iterations */
366 inneriter += j_index_end - j_index_start;
368 /* Outer loop uses 9 flops */
371 /* Increment number of outer iterations */
374 /* Update outer/inner flops */
376 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*70);
379 * Gromacs nonbonded kernel: nb_kernel_ElecRF_VdwCSTab_GeomP1P1_F_avx_128_fma_double
380 * Electrostatics interaction: ReactionField
381 * VdW interaction: CubicSplineTable
382 * Geometry: Particle-Particle
383 * Calculate force/pot: Force
386 nb_kernel_ElecRF_VdwCSTab_GeomP1P1_F_avx_128_fma_double
387 (t_nblist * gmx_restrict nlist,
388 rvec * gmx_restrict xx,
389 rvec * gmx_restrict ff,
390 t_forcerec * gmx_restrict fr,
391 t_mdatoms * gmx_restrict mdatoms,
392 nb_kernel_data_t * gmx_restrict kernel_data,
393 t_nrnb * gmx_restrict nrnb)
395 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
396 * just 0 for non-waters.
397 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
398 * jnr indices corresponding to data put in the four positions in the SIMD register.
400 int i_shift_offset,i_coord_offset,outeriter,inneriter;
401 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
403 int j_coord_offsetA,j_coord_offsetB;
404 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
406 real *shiftvec,*fshift,*x,*f;
407 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
409 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
410 int vdwjidx0A,vdwjidx0B;
411 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
412 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
413 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
416 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
419 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
420 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
422 __m128i ifour = _mm_set1_epi32(4);
423 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
425 __m128d dummy_mask,cutoff_mask;
426 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
427 __m128d one = _mm_set1_pd(1.0);
428 __m128d two = _mm_set1_pd(2.0);
434 jindex = nlist->jindex;
436 shiftidx = nlist->shift;
438 shiftvec = fr->shift_vec[0];
439 fshift = fr->fshift[0];
440 facel = _mm_set1_pd(fr->epsfac);
441 charge = mdatoms->chargeA;
442 krf = _mm_set1_pd(fr->ic->k_rf);
443 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
444 crf = _mm_set1_pd(fr->ic->c_rf);
445 nvdwtype = fr->ntype;
447 vdwtype = mdatoms->typeA;
449 vftab = kernel_data->table_vdw->data;
450 vftabscale = _mm_set1_pd(kernel_data->table_vdw->scale);
452 /* Avoid stupid compiler warnings */
460 /* Start outer loop over neighborlists */
461 for(iidx=0; iidx<nri; iidx++)
463 /* Load shift vector for this list */
464 i_shift_offset = DIM*shiftidx[iidx];
466 /* Load limits for loop over neighbors */
467 j_index_start = jindex[iidx];
468 j_index_end = jindex[iidx+1];
470 /* Get outer coordinate index */
472 i_coord_offset = DIM*inr;
474 /* Load i particle coords and add shift vector */
475 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
477 fix0 = _mm_setzero_pd();
478 fiy0 = _mm_setzero_pd();
479 fiz0 = _mm_setzero_pd();
481 /* Load parameters for i particles */
482 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
483 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
485 /* Start inner kernel loop */
486 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
489 /* Get j neighbor index, and coordinate index */
492 j_coord_offsetA = DIM*jnrA;
493 j_coord_offsetB = DIM*jnrB;
495 /* load j atom coordinates */
496 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
499 /* Calculate displacement vector */
500 dx00 = _mm_sub_pd(ix0,jx0);
501 dy00 = _mm_sub_pd(iy0,jy0);
502 dz00 = _mm_sub_pd(iz0,jz0);
504 /* Calculate squared distance and things based on it */
505 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
507 rinv00 = gmx_mm_invsqrt_pd(rsq00);
509 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
511 /* Load parameters for j particles */
512 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
513 vdwjidx0A = 2*vdwtype[jnrA+0];
514 vdwjidx0B = 2*vdwtype[jnrB+0];
516 /**************************
517 * CALCULATE INTERACTIONS *
518 **************************/
520 r00 = _mm_mul_pd(rsq00,rinv00);
522 /* Compute parameters for interactions between i and j atoms */
523 qq00 = _mm_mul_pd(iq0,jq0);
524 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
525 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
527 /* Calculate table index by multiplying r with table scale and truncate to integer */
528 rt = _mm_mul_pd(r00,vftabscale);
529 vfitab = _mm_cvttpd_epi32(rt);
531 vfeps = _mm_frcz_pd(rt);
533 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
535 twovfeps = _mm_add_pd(vfeps,vfeps);
536 vfitab = _mm_slli_epi32(vfitab,3);
538 /* REACTION-FIELD ELECTROSTATICS */
539 felec = _mm_mul_pd(qq00,_mm_msub_pd(rinv00,rinvsq00,krf2));
541 /* CUBIC SPLINE TABLE DISPERSION */
542 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
543 F = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
544 GMX_MM_TRANSPOSE2_PD(Y,F);
545 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
546 H = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
547 GMX_MM_TRANSPOSE2_PD(G,H);
548 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
549 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
550 fvdw6 = _mm_mul_pd(c6_00,FF);
552 /* CUBIC SPLINE TABLE REPULSION */
553 vfitab = _mm_add_epi32(vfitab,ifour);
554 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
555 F = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
556 GMX_MM_TRANSPOSE2_PD(Y,F);
557 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
558 H = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
559 GMX_MM_TRANSPOSE2_PD(G,H);
560 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
561 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
562 fvdw12 = _mm_mul_pd(c12_00,FF);
563 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
565 fscal = _mm_add_pd(felec,fvdw);
567 /* Update vectorial force */
568 fix0 = _mm_macc_pd(dx00,fscal,fix0);
569 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
570 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
572 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
573 _mm_mul_pd(dx00,fscal),
574 _mm_mul_pd(dy00,fscal),
575 _mm_mul_pd(dz00,fscal));
577 /* Inner loop uses 57 flops */
584 j_coord_offsetA = DIM*jnrA;
586 /* load j atom coordinates */
587 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
590 /* Calculate displacement vector */
591 dx00 = _mm_sub_pd(ix0,jx0);
592 dy00 = _mm_sub_pd(iy0,jy0);
593 dz00 = _mm_sub_pd(iz0,jz0);
595 /* Calculate squared distance and things based on it */
596 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
598 rinv00 = gmx_mm_invsqrt_pd(rsq00);
600 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
602 /* Load parameters for j particles */
603 jq0 = _mm_load_sd(charge+jnrA+0);
604 vdwjidx0A = 2*vdwtype[jnrA+0];
606 /**************************
607 * CALCULATE INTERACTIONS *
608 **************************/
610 r00 = _mm_mul_pd(rsq00,rinv00);
612 /* Compute parameters for interactions between i and j atoms */
613 qq00 = _mm_mul_pd(iq0,jq0);
614 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
616 /* Calculate table index by multiplying r with table scale and truncate to integer */
617 rt = _mm_mul_pd(r00,vftabscale);
618 vfitab = _mm_cvttpd_epi32(rt);
620 vfeps = _mm_frcz_pd(rt);
622 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
624 twovfeps = _mm_add_pd(vfeps,vfeps);
625 vfitab = _mm_slli_epi32(vfitab,3);
627 /* REACTION-FIELD ELECTROSTATICS */
628 felec = _mm_mul_pd(qq00,_mm_msub_pd(rinv00,rinvsq00,krf2));
630 /* CUBIC SPLINE TABLE DISPERSION */
631 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
632 F = _mm_setzero_pd();
633 GMX_MM_TRANSPOSE2_PD(Y,F);
634 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
635 H = _mm_setzero_pd();
636 GMX_MM_TRANSPOSE2_PD(G,H);
637 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
638 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
639 fvdw6 = _mm_mul_pd(c6_00,FF);
641 /* CUBIC SPLINE TABLE REPULSION */
642 vfitab = _mm_add_epi32(vfitab,ifour);
643 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
644 F = _mm_setzero_pd();
645 GMX_MM_TRANSPOSE2_PD(Y,F);
646 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
647 H = _mm_setzero_pd();
648 GMX_MM_TRANSPOSE2_PD(G,H);
649 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
650 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
651 fvdw12 = _mm_mul_pd(c12_00,FF);
652 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
654 fscal = _mm_add_pd(felec,fvdw);
656 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
658 /* Update vectorial force */
659 fix0 = _mm_macc_pd(dx00,fscal,fix0);
660 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
661 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
663 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
664 _mm_mul_pd(dx00,fscal),
665 _mm_mul_pd(dy00,fscal),
666 _mm_mul_pd(dz00,fscal));
668 /* Inner loop uses 57 flops */
671 /* End of innermost loop */
673 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
674 f+i_coord_offset,fshift+i_shift_offset);
676 /* Increment number of inner iterations */
677 inneriter += j_index_end - j_index_start;
679 /* Outer loop uses 7 flops */
682 /* Increment number of outer iterations */
685 /* Update outer/inner flops */
687 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*57);