Remove all unnecessary HAVE_CONFIG_H
[gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sse4_1_double / nb_kernel_ElecCSTab_VdwLJ_GeomP1P1_sse4_1_double.c
blob8afc9da6d936d487382d95c83b1e1e9072fba624
1 /*
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 sse4_1_double kernel generator.
38 #include "config.h"
40 #include <math.h>
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "nrnb.h"
47 #include "gromacs/simd/math_x86_sse4_1_double.h"
48 #include "kernelutil_x86_sse4_1_double.h"
51 * Gromacs nonbonded kernel: nb_kernel_ElecCSTab_VdwLJ_GeomP1P1_VF_sse4_1_double
52 * Electrostatics interaction: CubicSplineTable
53 * VdW interaction: LennardJones
54 * Geometry: Particle-Particle
55 * Calculate force/pot: PotentialAndForce
57 void
58 nb_kernel_ElecCSTab_VdwLJ_GeomP1P1_VF_sse4_1_double
59 (t_nblist * gmx_restrict nlist,
60 rvec * gmx_restrict xx,
61 rvec * gmx_restrict ff,
62 t_forcerec * gmx_restrict fr,
63 t_mdatoms * gmx_restrict mdatoms,
64 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
65 t_nrnb * gmx_restrict nrnb)
67 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
68 * just 0 for non-waters.
69 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
70 * jnr indices corresponding to data put in the four positions in the SIMD register.
72 int i_shift_offset,i_coord_offset,outeriter,inneriter;
73 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
74 int jnrA,jnrB;
75 int j_coord_offsetA,j_coord_offsetB;
76 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
77 real rcutoff_scalar;
78 real *shiftvec,*fshift,*x,*f;
79 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
80 int vdwioffset0;
81 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
82 int vdwjidx0A,vdwjidx0B;
83 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
84 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
85 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
86 real *charge;
87 int nvdwtype;
88 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
89 int *vdwtype;
90 real *vdwparam;
91 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
92 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
93 __m128i vfitab;
94 __m128i ifour = _mm_set1_epi32(4);
95 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
96 real *vftab;
97 __m128d dummy_mask,cutoff_mask;
98 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
99 __m128d one = _mm_set1_pd(1.0);
100 __m128d two = _mm_set1_pd(2.0);
101 x = xx[0];
102 f = ff[0];
104 nri = nlist->nri;
105 iinr = nlist->iinr;
106 jindex = nlist->jindex;
107 jjnr = nlist->jjnr;
108 shiftidx = nlist->shift;
109 gid = nlist->gid;
110 shiftvec = fr->shift_vec[0];
111 fshift = fr->fshift[0];
112 facel = _mm_set1_pd(fr->epsfac);
113 charge = mdatoms->chargeA;
114 nvdwtype = fr->ntype;
115 vdwparam = fr->nbfp;
116 vdwtype = mdatoms->typeA;
118 vftab = kernel_data->table_elec->data;
119 vftabscale = _mm_set1_pd(kernel_data->table_elec->scale);
121 /* Avoid stupid compiler warnings */
122 jnrA = jnrB = 0;
123 j_coord_offsetA = 0;
124 j_coord_offsetB = 0;
126 outeriter = 0;
127 inneriter = 0;
129 /* Start outer loop over neighborlists */
130 for(iidx=0; iidx<nri; iidx++)
132 /* Load shift vector for this list */
133 i_shift_offset = DIM*shiftidx[iidx];
135 /* Load limits for loop over neighbors */
136 j_index_start = jindex[iidx];
137 j_index_end = jindex[iidx+1];
139 /* Get outer coordinate index */
140 inr = iinr[iidx];
141 i_coord_offset = DIM*inr;
143 /* Load i particle coords and add shift vector */
144 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
146 fix0 = _mm_setzero_pd();
147 fiy0 = _mm_setzero_pd();
148 fiz0 = _mm_setzero_pd();
150 /* Load parameters for i particles */
151 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
152 vdwioffset0 = 2*nvdwtype*vdwtype[inr+0];
154 /* Reset potential sums */
155 velecsum = _mm_setzero_pd();
156 vvdwsum = _mm_setzero_pd();
158 /* Start inner kernel loop */
159 for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
162 /* Get j neighbor index, and coordinate index */
163 jnrA = jjnr[jidx];
164 jnrB = jjnr[jidx+1];
165 j_coord_offsetA = DIM*jnrA;
166 j_coord_offsetB = DIM*jnrB;
168 /* load j atom coordinates */
169 gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
170 &jx0,&jy0,&jz0);
172 /* Calculate displacement vector */
173 dx00 = _mm_sub_pd(ix0,jx0);
174 dy00 = _mm_sub_pd(iy0,jy0);
175 dz00 = _mm_sub_pd(iz0,jz0);
177 /* Calculate squared distance and things based on it */
178 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
180 rinv00 = gmx_mm_invsqrt_pd(rsq00);
182 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
184 /* Load parameters for j particles */
185 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
186 vdwjidx0A = 2*vdwtype[jnrA+0];
187 vdwjidx0B = 2*vdwtype[jnrB+0];
189 /**************************
190 * CALCULATE INTERACTIONS *
191 **************************/
193 r00 = _mm_mul_pd(rsq00,rinv00);
195 /* Compute parameters for interactions between i and j atoms */
196 qq00 = _mm_mul_pd(iq0,jq0);
197 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
198 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
200 /* Calculate table index by multiplying r with table scale and truncate to integer */
201 rt = _mm_mul_pd(r00,vftabscale);
202 vfitab = _mm_cvttpd_epi32(rt);
203 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
204 vfitab = _mm_slli_epi32(vfitab,2);
206 /* CUBIC SPLINE TABLE ELECTROSTATICS */
207 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
208 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
209 GMX_MM_TRANSPOSE2_PD(Y,F);
210 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
211 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
212 GMX_MM_TRANSPOSE2_PD(G,H);
213 Heps = _mm_mul_pd(vfeps,H);
214 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
215 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
216 velec = _mm_mul_pd(qq00,VV);
217 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
218 felec = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
220 /* LENNARD-JONES DISPERSION/REPULSION */
222 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
223 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
224 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
225 vvdw = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
226 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
228 /* Update potential sum for this i atom from the interaction with this j atom. */
229 velecsum = _mm_add_pd(velecsum,velec);
230 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
232 fscal = _mm_add_pd(felec,fvdw);
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);
246 /* Inner loop uses 56 flops */
249 if(jidx<j_index_end)
252 jnrA = jjnr[jidx];
253 j_coord_offsetA = DIM*jnrA;
255 /* load j atom coordinates */
256 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
257 &jx0,&jy0,&jz0);
259 /* Calculate displacement vector */
260 dx00 = _mm_sub_pd(ix0,jx0);
261 dy00 = _mm_sub_pd(iy0,jy0);
262 dz00 = _mm_sub_pd(iz0,jz0);
264 /* Calculate squared distance and things based on it */
265 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
267 rinv00 = gmx_mm_invsqrt_pd(rsq00);
269 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
271 /* Load parameters for j particles */
272 jq0 = _mm_load_sd(charge+jnrA+0);
273 vdwjidx0A = 2*vdwtype[jnrA+0];
275 /**************************
276 * CALCULATE INTERACTIONS *
277 **************************/
279 r00 = _mm_mul_pd(rsq00,rinv00);
281 /* Compute parameters for interactions between i and j atoms */
282 qq00 = _mm_mul_pd(iq0,jq0);
283 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
285 /* Calculate table index by multiplying r with table scale and truncate to integer */
286 rt = _mm_mul_pd(r00,vftabscale);
287 vfitab = _mm_cvttpd_epi32(rt);
288 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
289 vfitab = _mm_slli_epi32(vfitab,2);
291 /* CUBIC SPLINE TABLE ELECTROSTATICS */
292 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
293 F = _mm_setzero_pd();
294 GMX_MM_TRANSPOSE2_PD(Y,F);
295 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
296 H = _mm_setzero_pd();
297 GMX_MM_TRANSPOSE2_PD(G,H);
298 Heps = _mm_mul_pd(vfeps,H);
299 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
300 VV = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
301 velec = _mm_mul_pd(qq00,VV);
302 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
303 felec = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
305 /* LENNARD-JONES DISPERSION/REPULSION */
307 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
308 vvdw6 = _mm_mul_pd(c6_00,rinvsix);
309 vvdw12 = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
310 vvdw = _mm_sub_pd( _mm_mul_pd(vvdw12,one_twelfth) , _mm_mul_pd(vvdw6,one_sixth) );
311 fvdw = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
313 /* Update potential sum for this i atom from the interaction with this j atom. */
314 velec = _mm_unpacklo_pd(velec,_mm_setzero_pd());
315 velecsum = _mm_add_pd(velecsum,velec);
316 vvdw = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
317 vvdwsum = _mm_add_pd(vvdwsum,vvdw);
319 fscal = _mm_add_pd(felec,fvdw);
321 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
323 /* Calculate temporary vectorial force */
324 tx = _mm_mul_pd(fscal,dx00);
325 ty = _mm_mul_pd(fscal,dy00);
326 tz = _mm_mul_pd(fscal,dz00);
328 /* Update vectorial force */
329 fix0 = _mm_add_pd(fix0,tx);
330 fiy0 = _mm_add_pd(fiy0,ty);
331 fiz0 = _mm_add_pd(fiz0,tz);
333 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
335 /* Inner loop uses 56 flops */
338 /* End of innermost loop */
340 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
341 f+i_coord_offset,fshift+i_shift_offset);
343 ggid = gid[iidx];
344 /* Update potential energies */
345 gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
346 gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
348 /* Increment number of inner iterations */
349 inneriter += j_index_end - j_index_start;
351 /* Outer loop uses 9 flops */
354 /* Increment number of outer iterations */
355 outeriter += nri;
357 /* Update outer/inner flops */
359 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*56);
362 * Gromacs nonbonded kernel: nb_kernel_ElecCSTab_VdwLJ_GeomP1P1_F_sse4_1_double
363 * Electrostatics interaction: CubicSplineTable
364 * VdW interaction: LennardJones
365 * Geometry: Particle-Particle
366 * Calculate force/pot: Force
368 void
369 nb_kernel_ElecCSTab_VdwLJ_GeomP1P1_F_sse4_1_double
370 (t_nblist * gmx_restrict nlist,
371 rvec * gmx_restrict xx,
372 rvec * gmx_restrict ff,
373 t_forcerec * gmx_restrict fr,
374 t_mdatoms * gmx_restrict mdatoms,
375 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
376 t_nrnb * gmx_restrict nrnb)
378 /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
379 * just 0 for non-waters.
380 * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
381 * jnr indices corresponding to data put in the four positions in the SIMD register.
383 int i_shift_offset,i_coord_offset,outeriter,inneriter;
384 int j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
385 int jnrA,jnrB;
386 int j_coord_offsetA,j_coord_offsetB;
387 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
388 real rcutoff_scalar;
389 real *shiftvec,*fshift,*x,*f;
390 __m128d tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
391 int vdwioffset0;
392 __m128d ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
393 int vdwjidx0A,vdwjidx0B;
394 __m128d jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
395 __m128d dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
396 __m128d velec,felec,velecsum,facel,crf,krf,krf2;
397 real *charge;
398 int nvdwtype;
399 __m128d rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
400 int *vdwtype;
401 real *vdwparam;
402 __m128d one_sixth = _mm_set1_pd(1.0/6.0);
403 __m128d one_twelfth = _mm_set1_pd(1.0/12.0);
404 __m128i vfitab;
405 __m128i ifour = _mm_set1_epi32(4);
406 __m128d rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
407 real *vftab;
408 __m128d dummy_mask,cutoff_mask;
409 __m128d signbit = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
410 __m128d one = _mm_set1_pd(1.0);
411 __m128d two = _mm_set1_pd(2.0);
412 x = xx[0];
413 f = ff[0];
415 nri = nlist->nri;
416 iinr = nlist->iinr;
417 jindex = nlist->jindex;
418 jjnr = nlist->jjnr;
419 shiftidx = nlist->shift;
420 gid = nlist->gid;
421 shiftvec = fr->shift_vec[0];
422 fshift = fr->fshift[0];
423 facel = _mm_set1_pd(fr->epsfac);
424 charge = mdatoms->chargeA;
425 nvdwtype = fr->ntype;
426 vdwparam = fr->nbfp;
427 vdwtype = mdatoms->typeA;
429 vftab = kernel_data->table_elec->data;
430 vftabscale = _mm_set1_pd(kernel_data->table_elec->scale);
432 /* Avoid stupid compiler warnings */
433 jnrA = jnrB = 0;
434 j_coord_offsetA = 0;
435 j_coord_offsetB = 0;
437 outeriter = 0;
438 inneriter = 0;
440 /* Start outer loop over neighborlists */
441 for(iidx=0; iidx<nri; iidx++)
443 /* Load shift vector for this list */
444 i_shift_offset = DIM*shiftidx[iidx];
446 /* Load limits for loop over neighbors */
447 j_index_start = jindex[iidx];
448 j_index_end = jindex[iidx+1];
450 /* Get outer coordinate index */
451 inr = iinr[iidx];
452 i_coord_offset = DIM*inr;
454 /* Load i particle coords and add shift vector */
455 gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
457 fix0 = _mm_setzero_pd();
458 fiy0 = _mm_setzero_pd();
459 fiz0 = _mm_setzero_pd();
461 /* Load parameters for i particles */
462 iq0 = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
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 */
470 jnrA = jjnr[jidx];
471 jnrB = jjnr[jidx+1];
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,
477 &jx0,&jy0,&jz0);
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 jq0 = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
493 vdwjidx0A = 2*vdwtype[jnrA+0];
494 vdwjidx0B = 2*vdwtype[jnrB+0];
496 /**************************
497 * CALCULATE INTERACTIONS *
498 **************************/
500 r00 = _mm_mul_pd(rsq00,rinv00);
502 /* Compute parameters for interactions between i and j atoms */
503 qq00 = _mm_mul_pd(iq0,jq0);
504 gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
505 vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
507 /* Calculate table index by multiplying r with table scale and truncate to integer */
508 rt = _mm_mul_pd(r00,vftabscale);
509 vfitab = _mm_cvttpd_epi32(rt);
510 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
511 vfitab = _mm_slli_epi32(vfitab,2);
513 /* CUBIC SPLINE TABLE ELECTROSTATICS */
514 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
515 F = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
516 GMX_MM_TRANSPOSE2_PD(Y,F);
517 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
518 H = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
519 GMX_MM_TRANSPOSE2_PD(G,H);
520 Heps = _mm_mul_pd(vfeps,H);
521 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
522 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
523 felec = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
525 /* LENNARD-JONES DISPERSION/REPULSION */
527 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
528 fvdw = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(c12_00,rinvsix),c6_00),_mm_mul_pd(rinvsix,rinvsq00));
530 fscal = _mm_add_pd(felec,fvdw);
532 /* Calculate temporary vectorial force */
533 tx = _mm_mul_pd(fscal,dx00);
534 ty = _mm_mul_pd(fscal,dy00);
535 tz = _mm_mul_pd(fscal,dz00);
537 /* Update vectorial force */
538 fix0 = _mm_add_pd(fix0,tx);
539 fiy0 = _mm_add_pd(fiy0,ty);
540 fiz0 = _mm_add_pd(fiz0,tz);
542 gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
544 /* Inner loop uses 47 flops */
547 if(jidx<j_index_end)
550 jnrA = jjnr[jidx];
551 j_coord_offsetA = DIM*jnrA;
553 /* load j atom coordinates */
554 gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
555 &jx0,&jy0,&jz0);
557 /* Calculate displacement vector */
558 dx00 = _mm_sub_pd(ix0,jx0);
559 dy00 = _mm_sub_pd(iy0,jy0);
560 dz00 = _mm_sub_pd(iz0,jz0);
562 /* Calculate squared distance and things based on it */
563 rsq00 = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
565 rinv00 = gmx_mm_invsqrt_pd(rsq00);
567 rinvsq00 = _mm_mul_pd(rinv00,rinv00);
569 /* Load parameters for j particles */
570 jq0 = _mm_load_sd(charge+jnrA+0);
571 vdwjidx0A = 2*vdwtype[jnrA+0];
573 /**************************
574 * CALCULATE INTERACTIONS *
575 **************************/
577 r00 = _mm_mul_pd(rsq00,rinv00);
579 /* Compute parameters for interactions between i and j atoms */
580 qq00 = _mm_mul_pd(iq0,jq0);
581 gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
583 /* Calculate table index by multiplying r with table scale and truncate to integer */
584 rt = _mm_mul_pd(r00,vftabscale);
585 vfitab = _mm_cvttpd_epi32(rt);
586 vfeps = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
587 vfitab = _mm_slli_epi32(vfitab,2);
589 /* CUBIC SPLINE TABLE ELECTROSTATICS */
590 Y = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
591 F = _mm_setzero_pd();
592 GMX_MM_TRANSPOSE2_PD(Y,F);
593 G = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
594 H = _mm_setzero_pd();
595 GMX_MM_TRANSPOSE2_PD(G,H);
596 Heps = _mm_mul_pd(vfeps,H);
597 Fp = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
598 FF = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
599 felec = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
601 /* LENNARD-JONES DISPERSION/REPULSION */
603 rinvsix = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
604 fvdw = _mm_mul_pd(_mm_sub_pd(_mm_mul_pd(c12_00,rinvsix),c6_00),_mm_mul_pd(rinvsix,rinvsq00));
606 fscal = _mm_add_pd(felec,fvdw);
608 fscal = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
610 /* Calculate temporary vectorial force */
611 tx = _mm_mul_pd(fscal,dx00);
612 ty = _mm_mul_pd(fscal,dy00);
613 tz = _mm_mul_pd(fscal,dz00);
615 /* Update vectorial force */
616 fix0 = _mm_add_pd(fix0,tx);
617 fiy0 = _mm_add_pd(fiy0,ty);
618 fiz0 = _mm_add_pd(fiz0,tz);
620 gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
622 /* Inner loop uses 47 flops */
625 /* End of innermost loop */
627 gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
628 f+i_coord_offset,fshift+i_shift_offset);
630 /* Increment number of inner iterations */
631 inneriter += j_index_end - j_index_start;
633 /* Outer loop uses 7 flops */
636 /* Increment number of outer iterations */
637 outeriter += nri;
639 /* Update outer/inner flops */
641 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*47);