2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013, 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_double kernel generator.
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
49 #include "gromacs/simd/math_x86_avx_128_fma_double.h"
50 #include "kernelutil_x86_avx_128_fma_double.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_VF_avx_128_fma_double
54 * Electrostatics interaction: ReactionField
55 * VdW interaction: CubicSplineTable
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_VF_avx_128_fma_double
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 refer to j loop unrolling done with SSE double precision, e.g. for the two 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;
77 int j_coord_offsetA,j_coord_offsetB;
78 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
80 real *shiftvec,*fshift,*x,*f;
81 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
83 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
84 int vdwjidx0A,vdwjidx0B;
85 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
86 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
87 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
90 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
93 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
94 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
96 __m128i ifour = _mm_set1_epi32(4);
97 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
99 __m128d dummy_mask,cutoff_mask;
100 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
101 __m128d one = _mm_set1_pd(1.0);
102 __m128d two = _mm_set1_pd(2.0);
108 jindex = nlist->jindex;
110 shiftidx = nlist->shift;
112 shiftvec = fr->shift_vec[0];
113 fshift = fr->fshift[0];
114 facel = _mm_set1_pd(fr->epsfac);
115 charge = mdatoms->chargeA;
116 krf = _mm_set1_pd(fr->ic->k_rf);
117 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
118 crf = _mm_set1_pd(fr->ic->c_rf);
119 nvdwtype = fr->ntype;
121 vdwtype = mdatoms->typeA;
123 vftab = kernel_data->table_vdw->data;
124 vftabscale = _mm_set1_pd(kernel_data->table_vdw->scale);
126 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
127 rcutoff_scalar = fr->rcoulomb;
128 rcutoff = _mm_set1_pd(rcutoff_scalar);
129 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
131 /* Avoid stupid compiler warnings */
139 /* Start outer loop over neighborlists */
140 for(iidx=0; iidx<nri; iidx++)
142 /* Load shift vector for this list */
143 i_shift_offset = DIM*shiftidx[iidx];
145 /* Load limits for loop over neighbors */
146 j_index_start = jindex[iidx];
147 j_index_end = jindex[iidx+1];
149 /* Get outer coordinate index */
151 i_coord_offset = DIM*inr;
153 /* Load i particle coords and add shift vector */
154 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
156 fix0 = _mm_setzero_pd();
157 fiy0 = _mm_setzero_pd();
158 fiz0 = _mm_setzero_pd();
160 /* Load parameters for i particles */
161 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
162 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
164 /* Reset potential sums */
165 velecsum = _mm_setzero_pd();
166 vvdwsum = _mm_setzero_pd();
168 /* Start inner kernel loop */
169 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
172 /* Get j neighbor index, and coordinate index */
175 j_coord_offsetA = DIM*jnrA;
176 j_coord_offsetB = DIM*jnrB;
178 /* load j atom coordinates */
179 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
182 /* Calculate displacement vector */
183 dx00 = _mm_sub_pd(ix0,jx0);
184 dy00 = _mm_sub_pd(iy0,jy0);
185 dz00 = _mm_sub_pd(iz0,jz0);
187 /* Calculate squared distance and things based on it */
188 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
190 rinv00 = gmx_mm_invsqrt_pd(rsq00);
192 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
194 /* Load parameters for j particles */
195 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
196 vdwjidx0A = 2*vdwtype[jnrA+0];
197 vdwjidx0B = 2*vdwtype[jnrB+0];
199 /**************************
200 * CALCULATE INTERACTIONS *
201 **************************/
203 if (gmx_mm_any_lt(rsq00,rcutoff2))
206 r00 = _mm_mul_pd(rsq00,rinv00);
208 /* Compute parameters for interactions between i and j atoms */
209 qq00 = _mm_mul_pd(iq0,jq0);
210 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
211 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
213 /* Calculate table index by multiplying r with table scale and truncate to integer */
214 rt = _mm_mul_pd(r00,vftabscale);
215 vfitab = _mm_cvttpd_epi32(rt);
217 vfeps = _mm_frcz_pd(rt);
219 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
221 twovfeps = _mm_add_pd(vfeps,vfeps);
222 vfitab = _mm_slli_epi32(vfitab,3);
224 /* REACTION-FIELD ELECTROSTATICS */
225 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_macc_pd(krf,rsq00,rinv00),crf));
226 felec = _mm_mul_pd(qq00,_mm_msub_pd(rinv00,rinvsq00,krf2));
228 /* CUBIC SPLINE TABLE DISPERSION */
229 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
230 F = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
231 GMX_MM_TRANSPOSE2_PD(Y,F);
232 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
233 H = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
234 GMX_MM_TRANSPOSE2_PD(G,H);
235 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
236 VV = _mm_macc_pd(vfeps,Fp,Y);
237 vvdw6 = _mm_mul_pd(c6_00,VV);
238 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
239 fvdw6 = _mm_mul_pd(c6_00,FF);
241 /* CUBIC SPLINE TABLE REPULSION */
242 vfitab = _mm_add_epi32(vfitab,ifour);
243 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
244 F = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
245 GMX_MM_TRANSPOSE2_PD(Y,F);
246 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
247 H = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
248 GMX_MM_TRANSPOSE2_PD(G,H);
249 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
250 VV = _mm_macc_pd(vfeps,Fp,Y);
251 vvdw12 = _mm_mul_pd(c12_00,VV);
252 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
253 fvdw12 = _mm_mul_pd(c12_00,FF);
254 vvdw = _mm_add_pd(vvdw12,vvdw6);
255 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
257 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
259 /* Update potential sum for this i atom from the interaction with this j atom. */
260 velec = _mm_and_pd(velec,cutoff_mask);
261 velecsum = _mm_add_pd(velecsum,velec);
262 vvdw = _mm_and_pd(vvdw,cutoff_mask);
263 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
265 fscal = _mm_add_pd(felec,fvdw);
267 fscal = _mm_and_pd(fscal,cutoff_mask);
269 /* Update vectorial force */
270 fix0 = _mm_macc_pd(dx00,fscal,fix0);
271 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
272 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
274 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
275 _mm_mul_pd(dx00,fscal),
276 _mm_mul_pd(dy00,fscal),
277 _mm_mul_pd(dz00,fscal));
281 /* Inner loop uses 75 flops */
288 j_coord_offsetA = DIM*jnrA;
290 /* load j atom coordinates */
291 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
294 /* Calculate displacement vector */
295 dx00 = _mm_sub_pd(ix0,jx0);
296 dy00 = _mm_sub_pd(iy0,jy0);
297 dz00 = _mm_sub_pd(iz0,jz0);
299 /* Calculate squared distance and things based on it */
300 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
302 rinv00 = gmx_mm_invsqrt_pd(rsq00);
304 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
306 /* Load parameters for j particles */
307 jq0 = _mm_load_sd(charge+jnrA+0);
308 vdwjidx0A = 2*vdwtype[jnrA+0];
310 /**************************
311 * CALCULATE INTERACTIONS *
312 **************************/
314 if (gmx_mm_any_lt(rsq00,rcutoff2))
317 r00 = _mm_mul_pd(rsq00,rinv00);
319 /* Compute parameters for interactions between i and j atoms */
320 qq00 = _mm_mul_pd(iq0,jq0);
321 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
323 /* Calculate table index by multiplying r with table scale and truncate to integer */
324 rt = _mm_mul_pd(r00,vftabscale);
325 vfitab = _mm_cvttpd_epi32(rt);
327 vfeps = _mm_frcz_pd(rt);
329 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
331 twovfeps = _mm_add_pd(vfeps,vfeps);
332 vfitab = _mm_slli_epi32(vfitab,3);
334 /* REACTION-FIELD ELECTROSTATICS */
335 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_macc_pd(krf,rsq00,rinv00),crf));
336 felec = _mm_mul_pd(qq00,_mm_msub_pd(rinv00,rinvsq00,krf2));
338 /* CUBIC SPLINE TABLE DISPERSION */
339 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
340 F = _mm_setzero_pd();
341 GMX_MM_TRANSPOSE2_PD(Y,F);
342 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
343 H = _mm_setzero_pd();
344 GMX_MM_TRANSPOSE2_PD(G,H);
345 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
346 VV = _mm_macc_pd(vfeps,Fp,Y);
347 vvdw6 = _mm_mul_pd(c6_00,VV);
348 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
349 fvdw6 = _mm_mul_pd(c6_00,FF);
351 /* CUBIC SPLINE TABLE REPULSION */
352 vfitab = _mm_add_epi32(vfitab,ifour);
353 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
354 F = _mm_setzero_pd();
355 GMX_MM_TRANSPOSE2_PD(Y,F);
356 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
357 H = _mm_setzero_pd();
358 GMX_MM_TRANSPOSE2_PD(G,H);
359 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
360 VV = _mm_macc_pd(vfeps,Fp,Y);
361 vvdw12 = _mm_mul_pd(c12_00,VV);
362 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
363 fvdw12 = _mm_mul_pd(c12_00,FF);
364 vvdw = _mm_add_pd(vvdw12,vvdw6);
365 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
367 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
369 /* Update potential sum for this i atom from the interaction with this j atom. */
370 velec = _mm_and_pd(velec,cutoff_mask);
371 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
372 velecsum = _mm_add_pd(velecsum,velec);
373 vvdw = _mm_and_pd(vvdw,cutoff_mask);
374 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
375 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
377 fscal = _mm_add_pd(felec,fvdw);
379 fscal = _mm_and_pd(fscal,cutoff_mask);
381 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
383 /* Update vectorial force */
384 fix0 = _mm_macc_pd(dx00,fscal,fix0);
385 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
386 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
388 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
389 _mm_mul_pd(dx00,fscal),
390 _mm_mul_pd(dy00,fscal),
391 _mm_mul_pd(dz00,fscal));
395 /* Inner loop uses 75 flops */
398 /* End of innermost loop */
400 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
401 f+i_coord_offset,fshift+i_shift_offset);
404 /* Update potential energies */
405 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
406 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
408 /* Increment number of inner iterations */
409 inneriter += j_index_end - j_index_start;
411 /* Outer loop uses 9 flops */
414 /* Increment number of outer iterations */
417 /* Update outer/inner flops */
419 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*75);
422 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_F_avx_128_fma_double
423 * Electrostatics interaction: ReactionField
424 * VdW interaction: CubicSplineTable
425 * Geometry: Particle-Particle
426 * Calculate force/pot: Force
429 nb_kernel_ElecRFCut_VdwCSTab_GeomP1P1_F_avx_128_fma_double
430 (t_nblist * gmx_restrict nlist,
431 rvec * gmx_restrict xx,
432 rvec * gmx_restrict ff,
433 t_forcerec * gmx_restrict fr,
434 t_mdatoms * gmx_restrict mdatoms,
435 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
436 t_nrnb * gmx_restrict nrnb)
438 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
439 * just 0 for non-waters.
440 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
441 * jnr indices corresponding to data put in the four positions in the SIMD register.
443 int i_shift_offset,i_coord_offset,outeriter,inneriter;
444 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
446 int j_coord_offsetA,j_coord_offsetB;
447 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
449 real *shiftvec,*fshift,*x,*f;
450 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
452 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
453 int vdwjidx0A,vdwjidx0B;
454 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
455 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
456 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
459 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
462 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
463 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
465 __m128i ifour = _mm_set1_epi32(4);
466 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
468 __m128d dummy_mask,cutoff_mask;
469 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
470 __m128d one = _mm_set1_pd(1.0);
471 __m128d two = _mm_set1_pd(2.0);
477 jindex = nlist->jindex;
479 shiftidx = nlist->shift;
481 shiftvec = fr->shift_vec[0];
482 fshift = fr->fshift[0];
483 facel = _mm_set1_pd(fr->epsfac);
484 charge = mdatoms->chargeA;
485 krf = _mm_set1_pd(fr->ic->k_rf);
486 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
487 crf = _mm_set1_pd(fr->ic->c_rf);
488 nvdwtype = fr->ntype;
490 vdwtype = mdatoms->typeA;
492 vftab = kernel_data->table_vdw->data;
493 vftabscale = _mm_set1_pd(kernel_data->table_vdw->scale);
495 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
496 rcutoff_scalar = fr->rcoulomb;
497 rcutoff = _mm_set1_pd(rcutoff_scalar);
498 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
500 /* Avoid stupid compiler warnings */
508 /* Start outer loop over neighborlists */
509 for(iidx=0; iidx<nri; iidx++)
511 /* Load shift vector for this list */
512 i_shift_offset = DIM*shiftidx[iidx];
514 /* Load limits for loop over neighbors */
515 j_index_start = jindex[iidx];
516 j_index_end = jindex[iidx+1];
518 /* Get outer coordinate index */
520 i_coord_offset = DIM*inr;
522 /* Load i particle coords and add shift vector */
523 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
525 fix0 = _mm_setzero_pd();
526 fiy0 = _mm_setzero_pd();
527 fiz0 = _mm_setzero_pd();
529 /* Load parameters for i particles */
530 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
531 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
533 /* Start inner kernel loop */
534 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
537 /* Get j neighbor index, and coordinate index */
540 j_coord_offsetA = DIM*jnrA;
541 j_coord_offsetB = DIM*jnrB;
543 /* load j atom coordinates */
544 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
547 /* Calculate displacement vector */
548 dx00 = _mm_sub_pd(ix0,jx0);
549 dy00 = _mm_sub_pd(iy0,jy0);
550 dz00 = _mm_sub_pd(iz0,jz0);
552 /* Calculate squared distance and things based on it */
553 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
555 rinv00 = gmx_mm_invsqrt_pd(rsq00);
557 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
559 /* Load parameters for j particles */
560 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
561 vdwjidx0A = 2*vdwtype[jnrA+0];
562 vdwjidx0B = 2*vdwtype[jnrB+0];
564 /**************************
565 * CALCULATE INTERACTIONS *
566 **************************/
568 if (gmx_mm_any_lt(rsq00,rcutoff2))
571 r00 = _mm_mul_pd(rsq00,rinv00);
573 /* Compute parameters for interactions between i and j atoms */
574 qq00 = _mm_mul_pd(iq0,jq0);
575 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
576 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
578 /* Calculate table index by multiplying r with table scale and truncate to integer */
579 rt = _mm_mul_pd(r00,vftabscale);
580 vfitab = _mm_cvttpd_epi32(rt);
582 vfeps = _mm_frcz_pd(rt);
584 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
586 twovfeps = _mm_add_pd(vfeps,vfeps);
587 vfitab = _mm_slli_epi32(vfitab,3);
589 /* REACTION-FIELD ELECTROSTATICS */
590 felec = _mm_mul_pd(qq00,_mm_msub_pd(rinv00,rinvsq00,krf2));
592 /* CUBIC SPLINE TABLE DISPERSION */
593 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
594 F = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
595 GMX_MM_TRANSPOSE2_PD(Y,F);
596 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
597 H = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
598 GMX_MM_TRANSPOSE2_PD(G,H);
599 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
600 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
601 fvdw6 = _mm_mul_pd(c6_00,FF);
603 /* CUBIC SPLINE TABLE REPULSION */
604 vfitab = _mm_add_epi32(vfitab,ifour);
605 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
606 F = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
607 GMX_MM_TRANSPOSE2_PD(Y,F);
608 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
609 H = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
610 GMX_MM_TRANSPOSE2_PD(G,H);
611 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
612 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
613 fvdw12 = _mm_mul_pd(c12_00,FF);
614 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
616 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
618 fscal = _mm_add_pd(felec,fvdw);
620 fscal = _mm_and_pd(fscal,cutoff_mask);
622 /* Update vectorial force */
623 fix0 = _mm_macc_pd(dx00,fscal,fix0);
624 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
625 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
627 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
628 _mm_mul_pd(dx00,fscal),
629 _mm_mul_pd(dy00,fscal),
630 _mm_mul_pd(dz00,fscal));
634 /* Inner loop uses 60 flops */
641 j_coord_offsetA = DIM*jnrA;
643 /* load j atom coordinates */
644 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
647 /* Calculate displacement vector */
648 dx00 = _mm_sub_pd(ix0,jx0);
649 dy00 = _mm_sub_pd(iy0,jy0);
650 dz00 = _mm_sub_pd(iz0,jz0);
652 /* Calculate squared distance and things based on it */
653 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
655 rinv00 = gmx_mm_invsqrt_pd(rsq00);
657 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
659 /* Load parameters for j particles */
660 jq0 = _mm_load_sd(charge+jnrA+0);
661 vdwjidx0A = 2*vdwtype[jnrA+0];
663 /**************************
664 * CALCULATE INTERACTIONS *
665 **************************/
667 if (gmx_mm_any_lt(rsq00,rcutoff2))
670 r00 = _mm_mul_pd(rsq00,rinv00);
672 /* Compute parameters for interactions between i and j atoms */
673 qq00 = _mm_mul_pd(iq0,jq0);
674 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
676 /* Calculate table index by multiplying r with table scale and truncate to integer */
677 rt = _mm_mul_pd(r00,vftabscale);
678 vfitab = _mm_cvttpd_epi32(rt);
680 vfeps = _mm_frcz_pd(rt);
682 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
684 twovfeps = _mm_add_pd(vfeps,vfeps);
685 vfitab = _mm_slli_epi32(vfitab,3);
687 /* REACTION-FIELD ELECTROSTATICS */
688 felec = _mm_mul_pd(qq00,_mm_msub_pd(rinv00,rinvsq00,krf2));
690 /* CUBIC SPLINE TABLE DISPERSION */
691 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
692 F = _mm_setzero_pd();
693 GMX_MM_TRANSPOSE2_PD(Y,F);
694 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
695 H = _mm_setzero_pd();
696 GMX_MM_TRANSPOSE2_PD(G,H);
697 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
698 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
699 fvdw6 = _mm_mul_pd(c6_00,FF);
701 /* CUBIC SPLINE TABLE REPULSION */
702 vfitab = _mm_add_epi32(vfitab,ifour);
703 Y = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
704 F = _mm_setzero_pd();
705 GMX_MM_TRANSPOSE2_PD(Y,F);
706 G = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
707 H = _mm_setzero_pd();
708 GMX_MM_TRANSPOSE2_PD(G,H);
709 Fp = _mm_macc_pd(vfeps,_mm_macc_pd(H,vfeps,G),F);
710 FF = _mm_macc_pd(vfeps,_mm_macc_pd(twovfeps,H,G),Fp);
711 fvdw12 = _mm_mul_pd(c12_00,FF);
712 fvdw = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
714 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
716 fscal = _mm_add_pd(felec,fvdw);
718 fscal = _mm_and_pd(fscal,cutoff_mask);
720 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
722 /* Update vectorial force */
723 fix0 = _mm_macc_pd(dx00,fscal,fix0);
724 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
725 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
727 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
728 _mm_mul_pd(dx00,fscal),
729 _mm_mul_pd(dy00,fscal),
730 _mm_mul_pd(dz00,fscal));
734 /* Inner loop uses 60 flops */
737 /* End of innermost loop */
739 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
740 f+i_coord_offset,fshift+i_shift_offset);
742 /* Increment number of inner iterations */
743 inneriter += j_index_end - j_index_start;
745 /* Outer loop uses 7 flops */
748 /* Increment number of outer iterations */
751 /* Update outer/inner flops */
753 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*60);