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_ElecRFCut_VdwLJSw_GeomP1P1_VF_avx_128_fma_double
38 * Electrostatics interaction: ReactionField
39 * VdW interaction: LennardJones
40 * Geometry: Particle-Particle
41 * Calculate force/pot: PotentialAndForce
44 nb_kernel_ElecRFCut_VdwLJSw_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);
79 __m128d rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
80 real rswitch_scalar,d_scalar;
81 __m128d dummy_mask,cutoff_mask;
82 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
83 __m128d one = _mm_set1_pd(1.0);
84 __m128d two = _mm_set1_pd(2.0);
90 jindex = nlist->jindex;
92 shiftidx = nlist->shift;
94 shiftvec = fr->shift_vec[0];
95 fshift = fr->fshift[0];
96 facel = _mm_set1_pd(fr->epsfac);
97 charge = mdatoms->chargeA;
98 krf = _mm_set1_pd(fr->ic->k_rf);
99 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
100 crf = _mm_set1_pd(fr->ic->c_rf);
101 nvdwtype = fr->ntype;
103 vdwtype = mdatoms->typeA;
105 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
106 rcutoff_scalar = fr->rcoulomb;
107 rcutoff = _mm_set1_pd(rcutoff_scalar);
108 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
110 rswitch_scalar = fr->rvdw_switch;
111 rswitch = _mm_set1_pd(rswitch_scalar);
112 /* Setup switch parameters */
113 d_scalar = rcutoff_scalar-rswitch_scalar;
114 d = _mm_set1_pd(d_scalar);
115 swV3 = _mm_set1_pd(-10.0/(d_scalar*d_scalar*d_scalar));
116 swV4 = _mm_set1_pd( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
117 swV5 = _mm_set1_pd( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
118 swF2 = _mm_set1_pd(-30.0/(d_scalar*d_scalar*d_scalar));
119 swF3 = _mm_set1_pd( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
120 swF4 = _mm_set1_pd(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
122 /* Avoid stupid compiler warnings */
130 /* Start outer loop over neighborlists */
131 for(iidx=0; iidx<nri; iidx++)
133 /* Load shift vector for this list */
134 i_shift_offset = DIM*shiftidx[iidx];
136 /* Load limits for loop over neighbors */
137 j_index_start = jindex[iidx];
138 j_index_end = jindex[iidx+1];
140 /* Get outer coordinate index */
142 i_coord_offset = DIM*inr;
144 /* Load i particle coords and add shift vector */
145 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
147 fix0 = _mm_setzero_pd();
148 fiy0 = _mm_setzero_pd();
149 fiz0 = _mm_setzero_pd();
151 /* Load parameters for i particles */
152 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
153 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
155 /* Reset potential sums */
156 velecsum = _mm_setzero_pd();
157 vvdwsum = _mm_setzero_pd();
159 /* Start inner kernel loop */
160 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
163 /* Get j neighbor index, and coordinate index */
166 j_coord_offsetA = DIM*jnrA;
167 j_coord_offsetB = DIM*jnrB;
169 /* load j atom coordinates */
170 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
173 /* Calculate displacement vector */
174 dx00 = _mm_sub_pd(ix0,jx0);
175 dy00 = _mm_sub_pd(iy0,jy0);
176 dz00 = _mm_sub_pd(iz0,jz0);
178 /* Calculate squared distance and things based on it */
179 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
181 rinv00 = gmx_mm_invsqrt_pd(rsq00);
183 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
185 /* Load parameters for j particles */
186 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
187 vdwjidx0A = 2*vdwtype[jnrA+0];
188 vdwjidx0B = 2*vdwtype[jnrB+0];
190 /**************************
191 * CALCULATE INTERACTIONS *
192 **************************/
194 if (gmx_mm_any_lt(rsq00,rcutoff2))
197 r00 = _mm_mul_pd(rsq00,rinv00);
199 /* Compute parameters for interactions between i and j atoms */
200 qq00 = _mm_mul_pd(iq0,jq0);
201 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
202 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
204 /* REACTION-FIELD ELECTROSTATICS */
205 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_macc_pd(krf,rsq00,rinv00),crf));
206 felec = _mm_mul_pd(qq00,_mm_msub_pd(rinv00,rinvsq00,krf2));
208 /* LENNARD-JONES DISPERSION/REPULSION */
210 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
211 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
212 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
213 vvdw = _mm_msub_pd( vvdw12,one_twelfth, _mm_mul_pd(vvdw6,one_sixth) );
214 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
216 d = _mm_sub_pd(r00,rswitch);
217 d = _mm_max_pd(d,_mm_setzero_pd());
218 d2 = _mm_mul_pd(d,d);
219 sw = _mm_add_pd(one,_mm_mul_pd(d2,_mm_mul_pd(d,_mm_macc_pd(d,_mm_macc_pd(d,swV5,swV4),swV3))));
221 dsw = _mm_mul_pd(d2,_mm_macc_pd(d,_mm_macc_pd(d,swF4,swF3),swF2));
223 /* Evaluate switch function */
224 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
225 fvdw = _mm_msub_pd( fvdw,sw , _mm_mul_pd(rinv00,_mm_mul_pd(vvdw,dsw)) );
226 vvdw = _mm_mul_pd(vvdw,sw);
227 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
229 /* Update potential sum for this i atom from the interaction with this j atom. */
230 velec = _mm_and_pd(velec,cutoff_mask);
231 velecsum = _mm_add_pd(velecsum,velec);
232 vvdw = _mm_and_pd(vvdw,cutoff_mask);
233 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
235 fscal = _mm_add_pd(felec,fvdw);
237 fscal = _mm_and_pd(fscal,cutoff_mask);
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));
251 /* Inner loop uses 73 flops */
258 j_coord_offsetA = DIM*jnrA;
260 /* load j atom coordinates */
261 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
264 /* Calculate displacement vector */
265 dx00 = _mm_sub_pd(ix0,jx0);
266 dy00 = _mm_sub_pd(iy0,jy0);
267 dz00 = _mm_sub_pd(iz0,jz0);
269 /* Calculate squared distance and things based on it */
270 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
272 rinv00 = gmx_mm_invsqrt_pd(rsq00);
274 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
276 /* Load parameters for j particles */
277 jq0 = _mm_load_sd(charge+jnrA+0);
278 vdwjidx0A = 2*vdwtype[jnrA+0];
280 /**************************
281 * CALCULATE INTERACTIONS *
282 **************************/
284 if (gmx_mm_any_lt(rsq00,rcutoff2))
287 r00 = _mm_mul_pd(rsq00,rinv00);
289 /* Compute parameters for interactions between i and j atoms */
290 qq00 = _mm_mul_pd(iq0,jq0);
291 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
293 /* REACTION-FIELD ELECTROSTATICS */
294 velec = _mm_mul_pd(qq00,_mm_sub_pd(_mm_macc_pd(krf,rsq00,rinv00),crf));
295 felec = _mm_mul_pd(qq00,_mm_msub_pd(rinv00,rinvsq00,krf2));
297 /* LENNARD-JONES DISPERSION/REPULSION */
299 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
300 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
301 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
302 vvdw = _mm_msub_pd( vvdw12,one_twelfth, _mm_mul_pd(vvdw6,one_sixth) );
303 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
305 d = _mm_sub_pd(r00,rswitch);
306 d = _mm_max_pd(d,_mm_setzero_pd());
307 d2 = _mm_mul_pd(d,d);
308 sw = _mm_add_pd(one,_mm_mul_pd(d2,_mm_mul_pd(d,_mm_macc_pd(d,_mm_macc_pd(d,swV5,swV4),swV3))));
310 dsw = _mm_mul_pd(d2,_mm_macc_pd(d,_mm_macc_pd(d,swF4,swF3),swF2));
312 /* Evaluate switch function */
313 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
314 fvdw = _mm_msub_pd( fvdw,sw , _mm_mul_pd(rinv00,_mm_mul_pd(vvdw,dsw)) );
315 vvdw = _mm_mul_pd(vvdw,sw);
316 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
318 /* Update potential sum for this i atom from the interaction with this j atom. */
319 velec = _mm_and_pd(velec,cutoff_mask);
320 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
321 velecsum = _mm_add_pd(velecsum,velec);
322 vvdw = _mm_and_pd(vvdw,cutoff_mask);
323 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
324 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
326 fscal = _mm_add_pd(felec,fvdw);
328 fscal = _mm_and_pd(fscal,cutoff_mask);
330 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
332 /* Update vectorial force */
333 fix0 = _mm_macc_pd(dx00,fscal,fix0);
334 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
335 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
337 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
338 _mm_mul_pd(dx00,fscal),
339 _mm_mul_pd(dy00,fscal),
340 _mm_mul_pd(dz00,fscal));
344 /* Inner loop uses 73 flops */
347 /* End of innermost loop */
349 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
350 f+i_coord_offset,fshift+i_shift_offset);
353 /* Update potential energies */
354 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
355 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
357 /* Increment number of inner iterations */
358 inneriter += j_index_end - j_index_start;
360 /* Outer loop uses 9 flops */
363 /* Increment number of outer iterations */
366 /* Update outer/inner flops */
368 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*73);
371 * Gromacs nonbonded kernel: nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_F_avx_128_fma_double
372 * Electrostatics interaction: ReactionField
373 * VdW interaction: LennardJones
374 * Geometry: Particle-Particle
375 * Calculate force/pot: Force
378 nb_kernel_ElecRFCut_VdwLJSw_GeomP1P1_F_avx_128_fma_double
379 (t_nblist * gmx_restrict nlist,
380 rvec * gmx_restrict xx,
381 rvec * gmx_restrict ff,
382 t_forcerec * gmx_restrict fr,
383 t_mdatoms * gmx_restrict mdatoms,
384 nb_kernel_data_t * gmx_restrict kernel_data,
385 t_nrnb * gmx_restrict nrnb)
387 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
388 * just 0 for non-waters.
389 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
390 * jnr indices corresponding to data put in the four positions in the SIMD register.
392 int i_shift_offset,i_coord_offset,outeriter,inneriter;
393 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
395 int j_coord_offsetA,j_coord_offsetB;
396 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
398 real *shiftvec,*fshift,*x,*f;
399 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
401 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
402 int vdwjidx0A,vdwjidx0B;
403 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
404 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
405 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
408 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
411 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
412 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
413 __m128d rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
414 real rswitch_scalar,d_scalar;
415 __m128d dummy_mask,cutoff_mask;
416 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
417 __m128d one = _mm_set1_pd(1.0);
418 __m128d two = _mm_set1_pd(2.0);
424 jindex = nlist->jindex;
426 shiftidx = nlist->shift;
428 shiftvec = fr->shift_vec[0];
429 fshift = fr->fshift[0];
430 facel = _mm_set1_pd(fr->epsfac);
431 charge = mdatoms->chargeA;
432 krf = _mm_set1_pd(fr->ic->k_rf);
433 krf2 = _mm_set1_pd(fr->ic->k_rf*2.0);
434 crf = _mm_set1_pd(fr->ic->c_rf);
435 nvdwtype = fr->ntype;
437 vdwtype = mdatoms->typeA;
439 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
440 rcutoff_scalar = fr->rcoulomb;
441 rcutoff = _mm_set1_pd(rcutoff_scalar);
442 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
444 rswitch_scalar = fr->rvdw_switch;
445 rswitch = _mm_set1_pd(rswitch_scalar);
446 /* Setup switch parameters */
447 d_scalar = rcutoff_scalar-rswitch_scalar;
448 d = _mm_set1_pd(d_scalar);
449 swV3 = _mm_set1_pd(-10.0/(d_scalar*d_scalar*d_scalar));
450 swV4 = _mm_set1_pd( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
451 swV5 = _mm_set1_pd( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
452 swF2 = _mm_set1_pd(-30.0/(d_scalar*d_scalar*d_scalar));
453 swF3 = _mm_set1_pd( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
454 swF4 = _mm_set1_pd(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
456 /* Avoid stupid compiler warnings */
464 /* Start outer loop over neighborlists */
465 for(iidx=0; iidx<nri; iidx++)
467 /* Load shift vector for this list */
468 i_shift_offset = DIM*shiftidx[iidx];
470 /* Load limits for loop over neighbors */
471 j_index_start = jindex[iidx];
472 j_index_end = jindex[iidx+1];
474 /* Get outer coordinate index */
476 i_coord_offset = DIM*inr;
478 /* Load i particle coords and add shift vector */
479 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
481 fix0 = _mm_setzero_pd();
482 fiy0 = _mm_setzero_pd();
483 fiz0 = _mm_setzero_pd();
485 /* Load parameters for i particles */
486 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
487 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
489 /* Start inner kernel loop */
490 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
493 /* Get j neighbor index, and coordinate index */
496 j_coord_offsetA = DIM*jnrA;
497 j_coord_offsetB = DIM*jnrB;
499 /* load j atom coordinates */
500 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
503 /* Calculate displacement vector */
504 dx00 = _mm_sub_pd(ix0,jx0);
505 dy00 = _mm_sub_pd(iy0,jy0);
506 dz00 = _mm_sub_pd(iz0,jz0);
508 /* Calculate squared distance and things based on it */
509 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
511 rinv00 = gmx_mm_invsqrt_pd(rsq00);
513 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
515 /* Load parameters for j particles */
516 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
517 vdwjidx0A = 2*vdwtype[jnrA+0];
518 vdwjidx0B = 2*vdwtype[jnrB+0];
520 /**************************
521 * CALCULATE INTERACTIONS *
522 **************************/
524 if (gmx_mm_any_lt(rsq00,rcutoff2))
527 r00 = _mm_mul_pd(rsq00,rinv00);
529 /* Compute parameters for interactions between i and j atoms */
530 qq00 = _mm_mul_pd(iq0,jq0);
531 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
532 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
534 /* REACTION-FIELD ELECTROSTATICS */
535 felec = _mm_mul_pd(qq00,_mm_msub_pd(rinv00,rinvsq00,krf2));
537 /* LENNARD-JONES DISPERSION/REPULSION */
539 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
540 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
541 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
542 vvdw = _mm_msub_pd( vvdw12,one_twelfth, _mm_mul_pd(vvdw6,one_sixth) );
543 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
545 d = _mm_sub_pd(r00,rswitch);
546 d = _mm_max_pd(d,_mm_setzero_pd());
547 d2 = _mm_mul_pd(d,d);
548 sw = _mm_add_pd(one,_mm_mul_pd(d2,_mm_mul_pd(d,_mm_macc_pd(d,_mm_macc_pd(d,swV5,swV4),swV3))));
550 dsw = _mm_mul_pd(d2,_mm_macc_pd(d,_mm_macc_pd(d,swF4,swF3),swF2));
552 /* Evaluate switch function */
553 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
554 fvdw = _mm_msub_pd( fvdw,sw , _mm_mul_pd(rinv00,_mm_mul_pd(vvdw,dsw)) );
555 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
557 fscal = _mm_add_pd(felec,fvdw);
559 fscal = _mm_and_pd(fscal,cutoff_mask);
561 /* Update vectorial force */
562 fix0 = _mm_macc_pd(dx00,fscal,fix0);
563 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
564 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
566 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
567 _mm_mul_pd(dx00,fscal),
568 _mm_mul_pd(dy00,fscal),
569 _mm_mul_pd(dz00,fscal));
573 /* Inner loop uses 64 flops */
580 j_coord_offsetA = DIM*jnrA;
582 /* load j atom coordinates */
583 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
586 /* Calculate displacement vector */
587 dx00 = _mm_sub_pd(ix0,jx0);
588 dy00 = _mm_sub_pd(iy0,jy0);
589 dz00 = _mm_sub_pd(iz0,jz0);
591 /* Calculate squared distance and things based on it */
592 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
594 rinv00 = gmx_mm_invsqrt_pd(rsq00);
596 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
598 /* Load parameters for j particles */
599 jq0 = _mm_load_sd(charge+jnrA+0);
600 vdwjidx0A = 2*vdwtype[jnrA+0];
602 /**************************
603 * CALCULATE INTERACTIONS *
604 **************************/
606 if (gmx_mm_any_lt(rsq00,rcutoff2))
609 r00 = _mm_mul_pd(rsq00,rinv00);
611 /* Compute parameters for interactions between i and j atoms */
612 qq00 = _mm_mul_pd(iq0,jq0);
613 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
615 /* REACTION-FIELD ELECTROSTATICS */
616 felec = _mm_mul_pd(qq00,_mm_msub_pd(rinv00,rinvsq00,krf2));
618 /* LENNARD-JONES DISPERSION/REPULSION */
620 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
621 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
622 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
623 vvdw = _mm_msub_pd( vvdw12,one_twelfth, _mm_mul_pd(vvdw6,one_sixth) );
624 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
626 d = _mm_sub_pd(r00,rswitch);
627 d = _mm_max_pd(d,_mm_setzero_pd());
628 d2 = _mm_mul_pd(d,d);
629 sw = _mm_add_pd(one,_mm_mul_pd(d2,_mm_mul_pd(d,_mm_macc_pd(d,_mm_macc_pd(d,swV5,swV4),swV3))));
631 dsw = _mm_mul_pd(d2,_mm_macc_pd(d,_mm_macc_pd(d,swF4,swF3),swF2));
633 /* Evaluate switch function */
634 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
635 fvdw = _mm_msub_pd( fvdw,sw , _mm_mul_pd(rinv00,_mm_mul_pd(vvdw,dsw)) );
636 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
638 fscal = _mm_add_pd(felec,fvdw);
640 fscal = _mm_and_pd(fscal,cutoff_mask);
642 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
644 /* Update vectorial force */
645 fix0 = _mm_macc_pd(dx00,fscal,fix0);
646 fiy0 = _mm_macc_pd(dy00,fscal,fiy0);
647 fiz0 = _mm_macc_pd(dz00,fscal,fiz0);
649 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
650 _mm_mul_pd(dx00,fscal),
651 _mm_mul_pd(dy00,fscal),
652 _mm_mul_pd(dz00,fscal));
656 /* Inner loop uses 64 flops */
659 /* End of innermost loop */
661 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
662 f+i_coord_offset,fshift+i_shift_offset);
664 /* Increment number of inner iterations */
665 inneriter += j_index_end - j_index_start;
667 /* Outer loop uses 7 flops */
670 /* Increment number of outer iterations */
673 /* Update outer/inner flops */
675 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*64);