2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
36 * Note: this file was generated by the GROMACS sse2_double kernel generator.
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
49 #include "gromacs/simd/math_x86_sse2_double.h"
50 #include "kernelutil_x86_sse2_double.h"
53 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_sse2_double
54 * Electrostatics interaction: None
55 * VdW interaction: LJEwald
56 * Geometry: Particle-Particle
57 * Calculate force/pot: PotentialAndForce
60 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_sse2_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;
88 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
91 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
92 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
94 __m128d ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
96 __m128d one_half = _mm_set1_pd(0.5);
97 __m128d minus_one = _mm_set1_pd(-1.0);
98 __m128d dummy_mask,cutoff_mask;
99 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
100 __m128d one = _mm_set1_pd(1.0);
101 __m128d two = _mm_set1_pd(2.0);
107 jindex = nlist->jindex;
109 shiftidx = nlist->shift;
111 shiftvec = fr->shift_vec[0];
112 fshift = fr->fshift[0];
113 nvdwtype = fr->ntype;
115 vdwtype = mdatoms->typeA;
116 vdwgridparam = fr->ljpme_c6grid;
117 sh_lj_ewald = _mm_set1_pd(fr->ic->sh_lj_ewald);
118 ewclj = _mm_set1_pd(fr->ewaldcoeff_lj);
119 ewclj2 = _mm_mul_pd(minus_one,_mm_mul_pd(ewclj,ewclj));
121 rcutoff_scalar = fr->rvdw;
122 rcutoff = _mm_set1_pd(rcutoff_scalar);
123 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
125 sh_vdw_invrcut6 = _mm_set1_pd(fr->ic->sh_invrc6);
126 rvdw = _mm_set1_pd(fr->rvdw);
128 /* Avoid stupid compiler warnings */
136 /* Start outer loop over neighborlists */
137 for(iidx=0; iidx<nri; iidx++)
139 /* Load shift vector for this list */
140 i_shift_offset = DIM*shiftidx[iidx];
142 /* Load limits for loop over neighbors */
143 j_index_start = jindex[iidx];
144 j_index_end = jindex[iidx+1];
146 /* Get outer coordinate index */
148 i_coord_offset = DIM*inr;
150 /* Load i particle coords and add shift vector */
151 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
153 fix0 = _mm_setzero_pd();
154 fiy0 = _mm_setzero_pd();
155 fiz0 = _mm_setzero_pd();
157 /* Load parameters for i particles */
158 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
160 /* Reset potential sums */
161 vvdwsum = _mm_setzero_pd();
163 /* Start inner kernel loop */
164 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
167 /* Get j neighbor index, and coordinate index */
170 j_coord_offsetA = DIM*jnrA;
171 j_coord_offsetB = DIM*jnrB;
173 /* load j atom coordinates */
174 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
177 /* Calculate displacement vector */
178 dx00 = _mm_sub_pd(ix0,jx0);
179 dy00 = _mm_sub_pd(iy0,jy0);
180 dz00 = _mm_sub_pd(iz0,jz0);
182 /* Calculate squared distance and things based on it */
183 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
185 rinv00 = gmx_mm_invsqrt_pd(rsq00);
187 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
189 /* Load parameters for j particles */
190 vdwjidx0A = 2*vdwtype[jnrA+0];
191 vdwjidx0B = 2*vdwtype[jnrB+0];
193 /**************************
194 * CALCULATE INTERACTIONS *
195 **************************/
197 if (gmx_mm_any_lt(rsq00,rcutoff2))
200 r00 = _mm_mul_pd(rsq00,rinv00);
202 /* Compute parameters for interactions between i and j atoms */
203 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
204 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
206 c6grid_00 = gmx_mm_load_2real_swizzle_pd(vdwgridparam+vdwioffset0+vdwjidx0A,
207 vdwgridparam+vdwioffset0+vdwjidx0B);
209 /* Analytical LJ-PME */
210 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
211 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
212 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
213 exponent = gmx_simd_exp_d(ewcljrsq);
214 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
215 poly = _mm_mul_pd(exponent,_mm_add_pd(_mm_sub_pd(one,ewcljrsq),_mm_mul_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half)));
216 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
217 vvdw6 = _mm_mul_pd(_mm_sub_pd(c6_00,_mm_mul_pd(c6grid_00,_mm_sub_pd(one,poly))),rinvsix);
218 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
219 vvdw = _mm_sub_pd(_mm_mul_pd( _mm_sub_pd(vvdw12 , _mm_mul_pd(c12_00,_mm_mul_pd(sh_vdw_invrcut6,sh_vdw_invrcut6))),one_twelfth),
220 _mm_mul_pd( _mm_sub_pd(vvdw6,_mm_add_pd(_mm_mul_pd(c6_00,sh_vdw_invrcut6),_mm_mul_pd(c6grid_00,sh_lj_ewald))),one_sixth));
221 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
222 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,_mm_sub_pd(vvdw6,_mm_mul_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6)))),rinvsq00);
224 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
226 /* Update potential sum for this i atom from the interaction with this j atom. */
227 vvdw = _mm_and_pd(vvdw,cutoff_mask);
228 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
232 fscal = _mm_and_pd(fscal,cutoff_mask);
234 /* Calculate temporary vectorial force */
235 tx = _mm_mul_pd(fscal,dx00);
236 ty = _mm_mul_pd(fscal,dy00);
237 tz = _mm_mul_pd(fscal,dz00);
239 /* Update vectorial force */
240 fix0 = _mm_add_pd(fix0,tx);
241 fiy0 = _mm_add_pd(fiy0,ty);
242 fiz0 = _mm_add_pd(fiz0,tz);
244 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
248 /* Inner loop uses 62 flops */
255 j_coord_offsetA = DIM*jnrA;
257 /* load j atom coordinates */
258 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
261 /* Calculate displacement vector */
262 dx00 = _mm_sub_pd(ix0,jx0);
263 dy00 = _mm_sub_pd(iy0,jy0);
264 dz00 = _mm_sub_pd(iz0,jz0);
266 /* Calculate squared distance and things based on it */
267 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
269 rinv00 = gmx_mm_invsqrt_pd(rsq00);
271 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
273 /* Load parameters for j particles */
274 vdwjidx0A = 2*vdwtype[jnrA+0];
276 /**************************
277 * CALCULATE INTERACTIONS *
278 **************************/
280 if (gmx_mm_any_lt(rsq00,rcutoff2))
283 r00 = _mm_mul_pd(rsq00,rinv00);
285 /* Compute parameters for interactions between i and j atoms */
286 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
288 c6grid_00 = gmx_mm_load_1real_pd(vdwgridparam+vdwioffset0+vdwjidx0A);
290 /* Analytical LJ-PME */
291 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
292 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
293 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
294 exponent = gmx_simd_exp_d(ewcljrsq);
295 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
296 poly = _mm_mul_pd(exponent,_mm_add_pd(_mm_sub_pd(one,ewcljrsq),_mm_mul_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half)));
297 /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
298 vvdw6 = _mm_mul_pd(_mm_sub_pd(c6_00,_mm_mul_pd(c6grid_00,_mm_sub_pd(one,poly))),rinvsix);
299 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
300 vvdw = _mm_sub_pd(_mm_mul_pd( _mm_sub_pd(vvdw12 , _mm_mul_pd(c12_00,_mm_mul_pd(sh_vdw_invrcut6,sh_vdw_invrcut6))),one_twelfth),
301 _mm_mul_pd( _mm_sub_pd(vvdw6,_mm_add_pd(_mm_mul_pd(c6_00,sh_vdw_invrcut6),_mm_mul_pd(c6grid_00,sh_lj_ewald))),one_sixth));
302 /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
303 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,_mm_sub_pd(vvdw6,_mm_mul_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6)))),rinvsq00);
305 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
307 /* Update potential sum for this i atom from the interaction with this j atom. */
308 vvdw = _mm_and_pd(vvdw,cutoff_mask);
309 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
310 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
314 fscal = _mm_and_pd(fscal,cutoff_mask);
316 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
318 /* Calculate temporary vectorial force */
319 tx = _mm_mul_pd(fscal,dx00);
320 ty = _mm_mul_pd(fscal,dy00);
321 tz = _mm_mul_pd(fscal,dz00);
323 /* Update vectorial force */
324 fix0 = _mm_add_pd(fix0,tx);
325 fiy0 = _mm_add_pd(fiy0,ty);
326 fiz0 = _mm_add_pd(fiz0,tz);
328 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
332 /* Inner loop uses 62 flops */
335 /* End of innermost loop */
337 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
338 f+i_coord_offset,fshift+i_shift_offset);
341 /* Update potential energies */
342 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
344 /* Increment number of inner iterations */
345 inneriter += j_index_end - j_index_start;
347 /* Outer loop uses 7 flops */
350 /* Increment number of outer iterations */
353 /* Update outer/inner flops */
355 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_VF,outeriter*7 + inneriter*62);
358 * Gromacs nonbonded kernel: nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_sse2_double
359 * Electrostatics interaction: None
360 * VdW interaction: LJEwald
361 * Geometry: Particle-Particle
362 * Calculate force/pot: Force
365 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_sse2_double
366 (t_nblist * gmx_restrict nlist,
367 rvec * gmx_restrict xx,
368 rvec * gmx_restrict ff,
369 t_forcerec * gmx_restrict fr,
370 t_mdatoms * gmx_restrict mdatoms,
371 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
372 t_nrnb * gmx_restrict nrnb)
374 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
375 * just 0 for non-waters.
376 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
377 * jnr indices corresponding to data put in the four positions in the SIMD register.
379 int i_shift_offset,i_coord_offset,outeriter,inneriter;
380 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
382 int j_coord_offsetA,j_coord_offsetB;
383 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
385 real *shiftvec,*fshift,*x,*f;
386 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
388 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
389 int vdwjidx0A,vdwjidx0B;
390 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
391 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
393 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
396 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
397 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
399 __m128d ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
401 __m128d one_half = _mm_set1_pd(0.5);
402 __m128d minus_one = _mm_set1_pd(-1.0);
403 __m128d dummy_mask,cutoff_mask;
404 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
405 __m128d one = _mm_set1_pd(1.0);
406 __m128d two = _mm_set1_pd(2.0);
412 jindex = nlist->jindex;
414 shiftidx = nlist->shift;
416 shiftvec = fr->shift_vec[0];
417 fshift = fr->fshift[0];
418 nvdwtype = fr->ntype;
420 vdwtype = mdatoms->typeA;
421 vdwgridparam = fr->ljpme_c6grid;
422 sh_lj_ewald = _mm_set1_pd(fr->ic->sh_lj_ewald);
423 ewclj = _mm_set1_pd(fr->ewaldcoeff_lj);
424 ewclj2 = _mm_mul_pd(minus_one,_mm_mul_pd(ewclj,ewclj));
426 rcutoff_scalar = fr->rvdw;
427 rcutoff = _mm_set1_pd(rcutoff_scalar);
428 rcutoff2 = _mm_mul_pd(rcutoff,rcutoff);
430 sh_vdw_invrcut6 = _mm_set1_pd(fr->ic->sh_invrc6);
431 rvdw = _mm_set1_pd(fr->rvdw);
433 /* Avoid stupid compiler warnings */
441 /* Start outer loop over neighborlists */
442 for(iidx=0; iidx<nri; iidx++)
444 /* Load shift vector for this list */
445 i_shift_offset = DIM*shiftidx[iidx];
447 /* Load limits for loop over neighbors */
448 j_index_start = jindex[iidx];
449 j_index_end = jindex[iidx+1];
451 /* Get outer coordinate index */
453 i_coord_offset = DIM*inr;
455 /* Load i particle coords and add shift vector */
456 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
458 fix0 = _mm_setzero_pd();
459 fiy0 = _mm_setzero_pd();
460 fiz0 = _mm_setzero_pd();
462 /* Load parameters for i particles */
463 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
465 /* Start inner kernel loop */
466 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
469 /* Get j neighbor index, and coordinate index */
472 j_coord_offsetA = DIM*jnrA;
473 j_coord_offsetB = DIM*jnrB;
475 /* load j atom coordinates */
476 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
479 /* Calculate displacement vector */
480 dx00 = _mm_sub_pd(ix0,jx0);
481 dy00 = _mm_sub_pd(iy0,jy0);
482 dz00 = _mm_sub_pd(iz0,jz0);
484 /* Calculate squared distance and things based on it */
485 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
487 rinv00 = gmx_mm_invsqrt_pd(rsq00);
489 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
491 /* Load parameters for j particles */
492 vdwjidx0A = 2*vdwtype[jnrA+0];
493 vdwjidx0B = 2*vdwtype[jnrB+0];
495 /**************************
496 * CALCULATE INTERACTIONS *
497 **************************/
499 if (gmx_mm_any_lt(rsq00,rcutoff2))
502 r00 = _mm_mul_pd(rsq00,rinv00);
504 /* Compute parameters for interactions between i and j atoms */
505 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
506 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
508 c6grid_00 = gmx_mm_load_2real_swizzle_pd(vdwgridparam+vdwioffset0+vdwjidx0A,
509 vdwgridparam+vdwioffset0+vdwjidx0B);
511 /* Analytical LJ-PME */
512 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
513 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
514 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
515 exponent = gmx_simd_exp_d(ewcljrsq);
516 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
517 poly = _mm_mul_pd(exponent,_mm_add_pd(_mm_sub_pd(one,ewcljrsq),_mm_mul_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half)));
518 /* f6A = 6 * C6grid * (1 - poly) */
519 f6A = _mm_mul_pd(c6grid_00,_mm_sub_pd(one,poly));
520 /* f6B = C6grid * exponent * beta^6 */
521 f6B = _mm_mul_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6));
522 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
523 fvdw = _mm_mul_pd(_mm_add_pd(_mm_mul_pd(_mm_sub_pd(_mm_mul_pd(c12_00,rinvsix),_mm_sub_pd(c6_00,f6A)),rinvsix),f6B),rinvsq00);
525 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
529 fscal = _mm_and_pd(fscal,cutoff_mask);
531 /* Calculate temporary vectorial force */
532 tx = _mm_mul_pd(fscal,dx00);
533 ty = _mm_mul_pd(fscal,dy00);
534 tz = _mm_mul_pd(fscal,dz00);
536 /* Update vectorial force */
537 fix0 = _mm_add_pd(fix0,tx);
538 fiy0 = _mm_add_pd(fiy0,ty);
539 fiz0 = _mm_add_pd(fiz0,tz);
541 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
545 /* Inner loop uses 49 flops */
552 j_coord_offsetA = DIM*jnrA;
554 /* load j atom coordinates */
555 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
558 /* Calculate displacement vector */
559 dx00 = _mm_sub_pd(ix0,jx0);
560 dy00 = _mm_sub_pd(iy0,jy0);
561 dz00 = _mm_sub_pd(iz0,jz0);
563 /* Calculate squared distance and things based on it */
564 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
566 rinv00 = gmx_mm_invsqrt_pd(rsq00);
568 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
570 /* Load parameters for j particles */
571 vdwjidx0A = 2*vdwtype[jnrA+0];
573 /**************************
574 * CALCULATE INTERACTIONS *
575 **************************/
577 if (gmx_mm_any_lt(rsq00,rcutoff2))
580 r00 = _mm_mul_pd(rsq00,rinv00);
582 /* Compute parameters for interactions between i and j atoms */
583 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
585 c6grid_00 = gmx_mm_load_1real_pd(vdwgridparam+vdwioffset0+vdwjidx0A);
587 /* Analytical LJ-PME */
588 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
589 ewcljrsq = _mm_mul_pd(ewclj2,rsq00);
590 ewclj6 = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
591 exponent = gmx_simd_exp_d(ewcljrsq);
592 /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
593 poly = _mm_mul_pd(exponent,_mm_add_pd(_mm_sub_pd(one,ewcljrsq),_mm_mul_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half)));
594 /* f6A = 6 * C6grid * (1 - poly) */
595 f6A = _mm_mul_pd(c6grid_00,_mm_sub_pd(one,poly));
596 /* f6B = C6grid * exponent * beta^6 */
597 f6B = _mm_mul_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6));
598 /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
599 fvdw = _mm_mul_pd(_mm_add_pd(_mm_mul_pd(_mm_sub_pd(_mm_mul_pd(c12_00,rinvsix),_mm_sub_pd(c6_00,f6A)),rinvsix),f6B),rinvsq00);
601 cutoff_mask = _mm_cmplt_pd(rsq00,rcutoff2);
605 fscal = _mm_and_pd(fscal,cutoff_mask);
607 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
609 /* Calculate temporary vectorial force */
610 tx = _mm_mul_pd(fscal,dx00);
611 ty = _mm_mul_pd(fscal,dy00);
612 tz = _mm_mul_pd(fscal,dz00);
614 /* Update vectorial force */
615 fix0 = _mm_add_pd(fix0,tx);
616 fiy0 = _mm_add_pd(fiy0,ty);
617 fiz0 = _mm_add_pd(fiz0,tz);
619 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
623 /* Inner loop uses 49 flops */
626 /* End of innermost loop */
628 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
629 f+i_coord_offset,fshift+i_shift_offset);
631 /* Increment number of inner iterations */
632 inneriter += j_index_end - j_index_start;
634 /* Outer loop uses 6 flops */
637 /* Increment number of outer iterations */
640 /* Update outer/inner flops */
642 inc_nrnb(nrnb,eNR_NBKERNEL_VDW_F,outeriter*6 + inneriter*49);