Remove all unnecessary HAVE_CONFIG_H
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_double / nb_kernel_ElecCSTab_VdwLJ_GeomW4P1_avx_128_fma_double.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 /*
36  * Note: this file was generated by the GROMACS avx_128_fma_double kernel generator.
37  */
38 #include "config.h"
39
40 #include <math.h>
41
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "nrnb.h"
46
47 #include "gromacs/simd/math_x86_avx_128_fma_double.h"
48 #include "kernelutil_x86_avx_128_fma_double.h"
49
50 /*
51  * Gromacs nonbonded kernel:   nb_kernel_ElecCSTab_VdwLJ_GeomW4P1_VF_avx_128_fma_double
52  * Electrostatics interaction: CubicSplineTable
53  * VdW interaction:            LennardJones
54  * Geometry:                   Water4-Particle
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecCSTab_VdwLJ_GeomW4P1_VF_avx_128_fma_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)
66 {
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.
71      */
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              vdwioffset1;
83     __m128d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
84     int              vdwioffset2;
85     __m128d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
86     int              vdwioffset3;
87     __m128d          ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
88     int              vdwjidx0A,vdwjidx0B;
89     __m128d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
90     __m128d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
91     __m128d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
92     __m128d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
93     __m128d          dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
94     __m128d          velec,felec,velecsum,facel,crf,krf,krf2;
95     real             *charge;
96     int              nvdwtype;
97     __m128d          rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
98     int              *vdwtype;
99     real             *vdwparam;
100     __m128d          one_sixth   = _mm_set1_pd(1.0/6.0);
101     __m128d          one_twelfth = _mm_set1_pd(1.0/12.0);
102     __m128i          vfitab;
103     __m128i          ifour       = _mm_set1_epi32(4);
104     __m128d          rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
105     real             *vftab;
106     __m128d          dummy_mask,cutoff_mask;
107     __m128d          signbit   = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
108     __m128d          one     = _mm_set1_pd(1.0);
109     __m128d          two     = _mm_set1_pd(2.0);
110     x                = xx[0];
111     f                = ff[0];
112
113     nri              = nlist->nri;
114     iinr             = nlist->iinr;
115     jindex           = nlist->jindex;
116     jjnr             = nlist->jjnr;
117     shiftidx         = nlist->shift;
118     gid              = nlist->gid;
119     shiftvec         = fr->shift_vec[0];
120     fshift           = fr->fshift[0];
121     facel            = _mm_set1_pd(fr->epsfac);
122     charge           = mdatoms->chargeA;
123     nvdwtype         = fr->ntype;
124     vdwparam         = fr->nbfp;
125     vdwtype          = mdatoms->typeA;
126
127     vftab            = kernel_data->table_elec->data;
128     vftabscale       = _mm_set1_pd(kernel_data->table_elec->scale);
129
130     /* Setup water-specific parameters */
131     inr              = nlist->iinr[0];
132     iq1              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
133     iq2              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
134     iq3              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
135     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
136
137     /* Avoid stupid compiler warnings */
138     jnrA = jnrB = 0;
139     j_coord_offsetA = 0;
140     j_coord_offsetB = 0;
141
142     outeriter        = 0;
143     inneriter        = 0;
144
145     /* Start outer loop over neighborlists */
146     for(iidx=0; iidx<nri; iidx++)
147     {
148         /* Load shift vector for this list */
149         i_shift_offset   = DIM*shiftidx[iidx];
150
151         /* Load limits for loop over neighbors */
152         j_index_start    = jindex[iidx];
153         j_index_end      = jindex[iidx+1];
154
155         /* Get outer coordinate index */
156         inr              = iinr[iidx];
157         i_coord_offset   = DIM*inr;
158
159         /* Load i particle coords and add shift vector */
160         gmx_mm_load_shift_and_4rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
161                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
162
163         fix0             = _mm_setzero_pd();
164         fiy0             = _mm_setzero_pd();
165         fiz0             = _mm_setzero_pd();
166         fix1             = _mm_setzero_pd();
167         fiy1             = _mm_setzero_pd();
168         fiz1             = _mm_setzero_pd();
169         fix2             = _mm_setzero_pd();
170         fiy2             = _mm_setzero_pd();
171         fiz2             = _mm_setzero_pd();
172         fix3             = _mm_setzero_pd();
173         fiy3             = _mm_setzero_pd();
174         fiz3             = _mm_setzero_pd();
175
176         /* Reset potential sums */
177         velecsum         = _mm_setzero_pd();
178         vvdwsum          = _mm_setzero_pd();
179
180         /* Start inner kernel loop */
181         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
182         {
183
184             /* Get j neighbor index, and coordinate index */
185             jnrA             = jjnr[jidx];
186             jnrB             = jjnr[jidx+1];
187             j_coord_offsetA  = DIM*jnrA;
188             j_coord_offsetB  = DIM*jnrB;
189
190             /* load j atom coordinates */
191             gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
192                                               &jx0,&jy0,&jz0);
193
194             /* Calculate displacement vector */
195             dx00             = _mm_sub_pd(ix0,jx0);
196             dy00             = _mm_sub_pd(iy0,jy0);
197             dz00             = _mm_sub_pd(iz0,jz0);
198             dx10             = _mm_sub_pd(ix1,jx0);
199             dy10             = _mm_sub_pd(iy1,jy0);
200             dz10             = _mm_sub_pd(iz1,jz0);
201             dx20             = _mm_sub_pd(ix2,jx0);
202             dy20             = _mm_sub_pd(iy2,jy0);
203             dz20             = _mm_sub_pd(iz2,jz0);
204             dx30             = _mm_sub_pd(ix3,jx0);
205             dy30             = _mm_sub_pd(iy3,jy0);
206             dz30             = _mm_sub_pd(iz3,jz0);
207
208             /* Calculate squared distance and things based on it */
209             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
210             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
211             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
212             rsq30            = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
213
214             rinv10           = gmx_mm_invsqrt_pd(rsq10);
215             rinv20           = gmx_mm_invsqrt_pd(rsq20);
216             rinv30           = gmx_mm_invsqrt_pd(rsq30);
217
218             rinvsq00         = gmx_mm_inv_pd(rsq00);
219
220             /* Load parameters for j particles */
221             jq0              = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
222             vdwjidx0A        = 2*vdwtype[jnrA+0];
223             vdwjidx0B        = 2*vdwtype[jnrB+0];
224
225             fjx0             = _mm_setzero_pd();
226             fjy0             = _mm_setzero_pd();
227             fjz0             = _mm_setzero_pd();
228
229             /**************************
230              * CALCULATE INTERACTIONS *
231              **************************/
232
233             /* Compute parameters for interactions between i and j atoms */
234             gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
235                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
236
237             /* LENNARD-JONES DISPERSION/REPULSION */
238
239             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
240             vvdw6            = _mm_mul_pd(c6_00,rinvsix);
241             vvdw12           = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
242             vvdw             = _mm_msub_pd( vvdw12,one_twelfth, _mm_mul_pd(vvdw6,one_sixth) );
243             fvdw             = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
244
245             /* Update potential sum for this i atom from the interaction with this j atom. */
246             vvdwsum          = _mm_add_pd(vvdwsum,vvdw);
247
248             fscal            = fvdw;
249
250             /* Update vectorial force */
251             fix0             = _mm_macc_pd(dx00,fscal,fix0);
252             fiy0             = _mm_macc_pd(dy00,fscal,fiy0);
253             fiz0             = _mm_macc_pd(dz00,fscal,fiz0);
254             
255             fjx0             = _mm_macc_pd(dx00,fscal,fjx0);
256             fjy0             = _mm_macc_pd(dy00,fscal,fjy0);
257             fjz0             = _mm_macc_pd(dz00,fscal,fjz0);
258
259             /**************************
260              * CALCULATE INTERACTIONS *
261              **************************/
262
263             r10              = _mm_mul_pd(rsq10,rinv10);
264
265             /* Compute parameters for interactions between i and j atoms */
266             qq10             = _mm_mul_pd(iq1,jq0);
267
268             /* Calculate table index by multiplying r with table scale and truncate to integer */
269             rt               = _mm_mul_pd(r10,vftabscale);
270             vfitab           = _mm_cvttpd_epi32(rt);
271 #ifdef __XOP__
272             vfeps            = _mm_frcz_pd(rt);
273 #else
274             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
275 #endif
276             twovfeps         = _mm_add_pd(vfeps,vfeps);
277             vfitab           = _mm_slli_epi32(vfitab,2);
278
279             /* CUBIC SPLINE TABLE ELECTROSTATICS */
280             Y                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
281             F                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
282             GMX_MM_TRANSPOSE2_PD(Y,F);
283             G                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
284             H                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
285             GMX_MM_TRANSPOSE2_PD(G,H);
286             Fp               = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
287             VV               = _mm_macc_pd(vfeps,Fp,Y);
288             velec            = _mm_mul_pd(qq10,VV);
289             FF               = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
290             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq10,FF),_mm_mul_pd(vftabscale,rinv10)));
291
292             /* Update potential sum for this i atom from the interaction with this j atom. */
293             velecsum         = _mm_add_pd(velecsum,velec);
294
295             fscal            = felec;
296
297             /* Update vectorial force */
298             fix1             = _mm_macc_pd(dx10,fscal,fix1);
299             fiy1             = _mm_macc_pd(dy10,fscal,fiy1);
300             fiz1             = _mm_macc_pd(dz10,fscal,fiz1);
301             
302             fjx0             = _mm_macc_pd(dx10,fscal,fjx0);
303             fjy0             = _mm_macc_pd(dy10,fscal,fjy0);
304             fjz0             = _mm_macc_pd(dz10,fscal,fjz0);
305
306             /**************************
307              * CALCULATE INTERACTIONS *
308              **************************/
309
310             r20              = _mm_mul_pd(rsq20,rinv20);
311
312             /* Compute parameters for interactions between i and j atoms */
313             qq20             = _mm_mul_pd(iq2,jq0);
314
315             /* Calculate table index by multiplying r with table scale and truncate to integer */
316             rt               = _mm_mul_pd(r20,vftabscale);
317             vfitab           = _mm_cvttpd_epi32(rt);
318 #ifdef __XOP__
319             vfeps            = _mm_frcz_pd(rt);
320 #else
321             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
322 #endif
323             twovfeps         = _mm_add_pd(vfeps,vfeps);
324             vfitab           = _mm_slli_epi32(vfitab,2);
325
326             /* CUBIC SPLINE TABLE ELECTROSTATICS */
327             Y                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
328             F                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
329             GMX_MM_TRANSPOSE2_PD(Y,F);
330             G                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
331             H                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
332             GMX_MM_TRANSPOSE2_PD(G,H);
333             Fp               = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
334             VV               = _mm_macc_pd(vfeps,Fp,Y);
335             velec            = _mm_mul_pd(qq20,VV);
336             FF               = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
337             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq20,FF),_mm_mul_pd(vftabscale,rinv20)));
338
339             /* Update potential sum for this i atom from the interaction with this j atom. */
340             velecsum         = _mm_add_pd(velecsum,velec);
341
342             fscal            = felec;
343
344             /* Update vectorial force */
345             fix2             = _mm_macc_pd(dx20,fscal,fix2);
346             fiy2             = _mm_macc_pd(dy20,fscal,fiy2);
347             fiz2             = _mm_macc_pd(dz20,fscal,fiz2);
348             
349             fjx0             = _mm_macc_pd(dx20,fscal,fjx0);
350             fjy0             = _mm_macc_pd(dy20,fscal,fjy0);
351             fjz0             = _mm_macc_pd(dz20,fscal,fjz0);
352
353             /**************************
354              * CALCULATE INTERACTIONS *
355              **************************/
356
357             r30              = _mm_mul_pd(rsq30,rinv30);
358
359             /* Compute parameters for interactions between i and j atoms */
360             qq30             = _mm_mul_pd(iq3,jq0);
361
362             /* Calculate table index by multiplying r with table scale and truncate to integer */
363             rt               = _mm_mul_pd(r30,vftabscale);
364             vfitab           = _mm_cvttpd_epi32(rt);
365 #ifdef __XOP__
366             vfeps            = _mm_frcz_pd(rt);
367 #else
368             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
369 #endif
370             twovfeps         = _mm_add_pd(vfeps,vfeps);
371             vfitab           = _mm_slli_epi32(vfitab,2);
372
373             /* CUBIC SPLINE TABLE ELECTROSTATICS */
374             Y                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
375             F                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
376             GMX_MM_TRANSPOSE2_PD(Y,F);
377             G                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
378             H                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
379             GMX_MM_TRANSPOSE2_PD(G,H);
380             Fp               = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
381             VV               = _mm_macc_pd(vfeps,Fp,Y);
382             velec            = _mm_mul_pd(qq30,VV);
383             FF               = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
384             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq30,FF),_mm_mul_pd(vftabscale,rinv30)));
385
386             /* Update potential sum for this i atom from the interaction with this j atom. */
387             velecsum         = _mm_add_pd(velecsum,velec);
388
389             fscal            = felec;
390
391             /* Update vectorial force */
392             fix3             = _mm_macc_pd(dx30,fscal,fix3);
393             fiy3             = _mm_macc_pd(dy30,fscal,fiy3);
394             fiz3             = _mm_macc_pd(dz30,fscal,fiz3);
395             
396             fjx0             = _mm_macc_pd(dx30,fscal,fjx0);
397             fjy0             = _mm_macc_pd(dy30,fscal,fjy0);
398             fjz0             = _mm_macc_pd(dz30,fscal,fjz0);
399
400             gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
401
402             /* Inner loop uses 176 flops */
403         }
404
405         if(jidx<j_index_end)
406         {
407
408             jnrA             = jjnr[jidx];
409             j_coord_offsetA  = DIM*jnrA;
410
411             /* load j atom coordinates */
412             gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
413                                               &jx0,&jy0,&jz0);
414
415             /* Calculate displacement vector */
416             dx00             = _mm_sub_pd(ix0,jx0);
417             dy00             = _mm_sub_pd(iy0,jy0);
418             dz00             = _mm_sub_pd(iz0,jz0);
419             dx10             = _mm_sub_pd(ix1,jx0);
420             dy10             = _mm_sub_pd(iy1,jy0);
421             dz10             = _mm_sub_pd(iz1,jz0);
422             dx20             = _mm_sub_pd(ix2,jx0);
423             dy20             = _mm_sub_pd(iy2,jy0);
424             dz20             = _mm_sub_pd(iz2,jz0);
425             dx30             = _mm_sub_pd(ix3,jx0);
426             dy30             = _mm_sub_pd(iy3,jy0);
427             dz30             = _mm_sub_pd(iz3,jz0);
428
429             /* Calculate squared distance and things based on it */
430             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
431             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
432             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
433             rsq30            = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
434
435             rinv10           = gmx_mm_invsqrt_pd(rsq10);
436             rinv20           = gmx_mm_invsqrt_pd(rsq20);
437             rinv30           = gmx_mm_invsqrt_pd(rsq30);
438
439             rinvsq00         = gmx_mm_inv_pd(rsq00);
440
441             /* Load parameters for j particles */
442             jq0              = _mm_load_sd(charge+jnrA+0);
443             vdwjidx0A        = 2*vdwtype[jnrA+0];
444
445             fjx0             = _mm_setzero_pd();
446             fjy0             = _mm_setzero_pd();
447             fjz0             = _mm_setzero_pd();
448
449             /**************************
450              * CALCULATE INTERACTIONS *
451              **************************/
452
453             /* Compute parameters for interactions between i and j atoms */
454             gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
455
456             /* LENNARD-JONES DISPERSION/REPULSION */
457
458             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
459             vvdw6            = _mm_mul_pd(c6_00,rinvsix);
460             vvdw12           = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
461             vvdw             = _mm_msub_pd( vvdw12,one_twelfth, _mm_mul_pd(vvdw6,one_sixth) );
462             fvdw             = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
463
464             /* Update potential sum for this i atom from the interaction with this j atom. */
465             vvdw             = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
466             vvdwsum          = _mm_add_pd(vvdwsum,vvdw);
467
468             fscal            = fvdw;
469
470             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
471
472             /* Update vectorial force */
473             fix0             = _mm_macc_pd(dx00,fscal,fix0);
474             fiy0             = _mm_macc_pd(dy00,fscal,fiy0);
475             fiz0             = _mm_macc_pd(dz00,fscal,fiz0);
476             
477             fjx0             = _mm_macc_pd(dx00,fscal,fjx0);
478             fjy0             = _mm_macc_pd(dy00,fscal,fjy0);
479             fjz0             = _mm_macc_pd(dz00,fscal,fjz0);
480
481             /**************************
482              * CALCULATE INTERACTIONS *
483              **************************/
484
485             r10              = _mm_mul_pd(rsq10,rinv10);
486
487             /* Compute parameters for interactions between i and j atoms */
488             qq10             = _mm_mul_pd(iq1,jq0);
489
490             /* Calculate table index by multiplying r with table scale and truncate to integer */
491             rt               = _mm_mul_pd(r10,vftabscale);
492             vfitab           = _mm_cvttpd_epi32(rt);
493 #ifdef __XOP__
494             vfeps            = _mm_frcz_pd(rt);
495 #else
496             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
497 #endif
498             twovfeps         = _mm_add_pd(vfeps,vfeps);
499             vfitab           = _mm_slli_epi32(vfitab,2);
500
501             /* CUBIC SPLINE TABLE ELECTROSTATICS */
502             Y                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
503             F                = _mm_setzero_pd();
504             GMX_MM_TRANSPOSE2_PD(Y,F);
505             G                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
506             H                = _mm_setzero_pd();
507             GMX_MM_TRANSPOSE2_PD(G,H);
508             Fp               = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
509             VV               = _mm_macc_pd(vfeps,Fp,Y);
510             velec            = _mm_mul_pd(qq10,VV);
511             FF               = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
512             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq10,FF),_mm_mul_pd(vftabscale,rinv10)));
513
514             /* Update potential sum for this i atom from the interaction with this j atom. */
515             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
516             velecsum         = _mm_add_pd(velecsum,velec);
517
518             fscal            = felec;
519
520             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
521
522             /* Update vectorial force */
523             fix1             = _mm_macc_pd(dx10,fscal,fix1);
524             fiy1             = _mm_macc_pd(dy10,fscal,fiy1);
525             fiz1             = _mm_macc_pd(dz10,fscal,fiz1);
526             
527             fjx0             = _mm_macc_pd(dx10,fscal,fjx0);
528             fjy0             = _mm_macc_pd(dy10,fscal,fjy0);
529             fjz0             = _mm_macc_pd(dz10,fscal,fjz0);
530
531             /**************************
532              * CALCULATE INTERACTIONS *
533              **************************/
534
535             r20              = _mm_mul_pd(rsq20,rinv20);
536
537             /* Compute parameters for interactions between i and j atoms */
538             qq20             = _mm_mul_pd(iq2,jq0);
539
540             /* Calculate table index by multiplying r with table scale and truncate to integer */
541             rt               = _mm_mul_pd(r20,vftabscale);
542             vfitab           = _mm_cvttpd_epi32(rt);
543 #ifdef __XOP__
544             vfeps            = _mm_frcz_pd(rt);
545 #else
546             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
547 #endif
548             twovfeps         = _mm_add_pd(vfeps,vfeps);
549             vfitab           = _mm_slli_epi32(vfitab,2);
550
551             /* CUBIC SPLINE TABLE ELECTROSTATICS */
552             Y                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
553             F                = _mm_setzero_pd();
554             GMX_MM_TRANSPOSE2_PD(Y,F);
555             G                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
556             H                = _mm_setzero_pd();
557             GMX_MM_TRANSPOSE2_PD(G,H);
558             Fp               = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
559             VV               = _mm_macc_pd(vfeps,Fp,Y);
560             velec            = _mm_mul_pd(qq20,VV);
561             FF               = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
562             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq20,FF),_mm_mul_pd(vftabscale,rinv20)));
563
564             /* Update potential sum for this i atom from the interaction with this j atom. */
565             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
566             velecsum         = _mm_add_pd(velecsum,velec);
567
568             fscal            = felec;
569
570             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
571
572             /* Update vectorial force */
573             fix2             = _mm_macc_pd(dx20,fscal,fix2);
574             fiy2             = _mm_macc_pd(dy20,fscal,fiy2);
575             fiz2             = _mm_macc_pd(dz20,fscal,fiz2);
576             
577             fjx0             = _mm_macc_pd(dx20,fscal,fjx0);
578             fjy0             = _mm_macc_pd(dy20,fscal,fjy0);
579             fjz0             = _mm_macc_pd(dz20,fscal,fjz0);
580
581             /**************************
582              * CALCULATE INTERACTIONS *
583              **************************/
584
585             r30              = _mm_mul_pd(rsq30,rinv30);
586
587             /* Compute parameters for interactions between i and j atoms */
588             qq30             = _mm_mul_pd(iq3,jq0);
589
590             /* Calculate table index by multiplying r with table scale and truncate to integer */
591             rt               = _mm_mul_pd(r30,vftabscale);
592             vfitab           = _mm_cvttpd_epi32(rt);
593 #ifdef __XOP__
594             vfeps            = _mm_frcz_pd(rt);
595 #else
596             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
597 #endif
598             twovfeps         = _mm_add_pd(vfeps,vfeps);
599             vfitab           = _mm_slli_epi32(vfitab,2);
600
601             /* CUBIC SPLINE TABLE ELECTROSTATICS */
602             Y                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
603             F                = _mm_setzero_pd();
604             GMX_MM_TRANSPOSE2_PD(Y,F);
605             G                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
606             H                = _mm_setzero_pd();
607             GMX_MM_TRANSPOSE2_PD(G,H);
608             Fp               = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
609             VV               = _mm_macc_pd(vfeps,Fp,Y);
610             velec            = _mm_mul_pd(qq30,VV);
611             FF               = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
612             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq30,FF),_mm_mul_pd(vftabscale,rinv30)));
613
614             /* Update potential sum for this i atom from the interaction with this j atom. */
615             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
616             velecsum         = _mm_add_pd(velecsum,velec);
617
618             fscal            = felec;
619
620             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
621
622             /* Update vectorial force */
623             fix3             = _mm_macc_pd(dx30,fscal,fix3);
624             fiy3             = _mm_macc_pd(dy30,fscal,fiy3);
625             fiz3             = _mm_macc_pd(dz30,fscal,fiz3);
626             
627             fjx0             = _mm_macc_pd(dx30,fscal,fjx0);
628             fjy0             = _mm_macc_pd(dy30,fscal,fjy0);
629             fjz0             = _mm_macc_pd(dz30,fscal,fjz0);
630
631             gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
632
633             /* Inner loop uses 176 flops */
634         }
635
636         /* End of innermost loop */
637
638         gmx_mm_update_iforce_4atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
639                                               f+i_coord_offset,fshift+i_shift_offset);
640
641         ggid                        = gid[iidx];
642         /* Update potential energies */
643         gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
644         gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
645
646         /* Increment number of inner iterations */
647         inneriter                  += j_index_end - j_index_start;
648
649         /* Outer loop uses 26 flops */
650     }
651
652     /* Increment number of outer iterations */
653     outeriter        += nri;
654
655     /* Update outer/inner flops */
656
657     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*176);
658 }
659 /*
660  * Gromacs nonbonded kernel:   nb_kernel_ElecCSTab_VdwLJ_GeomW4P1_F_avx_128_fma_double
661  * Electrostatics interaction: CubicSplineTable
662  * VdW interaction:            LennardJones
663  * Geometry:                   Water4-Particle
664  * Calculate force/pot:        Force
665  */
666 void
667 nb_kernel_ElecCSTab_VdwLJ_GeomW4P1_F_avx_128_fma_double
668                     (t_nblist                    * gmx_restrict       nlist,
669                      rvec                        * gmx_restrict          xx,
670                      rvec                        * gmx_restrict          ff,
671                      t_forcerec                  * gmx_restrict          fr,
672                      t_mdatoms                   * gmx_restrict     mdatoms,
673                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
674                      t_nrnb                      * gmx_restrict        nrnb)
675 {
676     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
677      * just 0 for non-waters.
678      * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
679      * jnr indices corresponding to data put in the four positions in the SIMD register.
680      */
681     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
682     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
683     int              jnrA,jnrB;
684     int              j_coord_offsetA,j_coord_offsetB;
685     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
686     real             rcutoff_scalar;
687     real             *shiftvec,*fshift,*x,*f;
688     __m128d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
689     int              vdwioffset0;
690     __m128d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
691     int              vdwioffset1;
692     __m128d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
693     int              vdwioffset2;
694     __m128d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
695     int              vdwioffset3;
696     __m128d          ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
697     int              vdwjidx0A,vdwjidx0B;
698     __m128d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
699     __m128d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
700     __m128d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
701     __m128d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
702     __m128d          dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
703     __m128d          velec,felec,velecsum,facel,crf,krf,krf2;
704     real             *charge;
705     int              nvdwtype;
706     __m128d          rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
707     int              *vdwtype;
708     real             *vdwparam;
709     __m128d          one_sixth   = _mm_set1_pd(1.0/6.0);
710     __m128d          one_twelfth = _mm_set1_pd(1.0/12.0);
711     __m128i          vfitab;
712     __m128i          ifour       = _mm_set1_epi32(4);
713     __m128d          rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF,twovfeps;
714     real             *vftab;
715     __m128d          dummy_mask,cutoff_mask;
716     __m128d          signbit   = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
717     __m128d          one     = _mm_set1_pd(1.0);
718     __m128d          two     = _mm_set1_pd(2.0);
719     x                = xx[0];
720     f                = ff[0];
721
722     nri              = nlist->nri;
723     iinr             = nlist->iinr;
724     jindex           = nlist->jindex;
725     jjnr             = nlist->jjnr;
726     shiftidx         = nlist->shift;
727     gid              = nlist->gid;
728     shiftvec         = fr->shift_vec[0];
729     fshift           = fr->fshift[0];
730     facel            = _mm_set1_pd(fr->epsfac);
731     charge           = mdatoms->chargeA;
732     nvdwtype         = fr->ntype;
733     vdwparam         = fr->nbfp;
734     vdwtype          = mdatoms->typeA;
735
736     vftab            = kernel_data->table_elec->data;
737     vftabscale       = _mm_set1_pd(kernel_data->table_elec->scale);
738
739     /* Setup water-specific parameters */
740     inr              = nlist->iinr[0];
741     iq1              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
742     iq2              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
743     iq3              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
744     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
745
746     /* Avoid stupid compiler warnings */
747     jnrA = jnrB = 0;
748     j_coord_offsetA = 0;
749     j_coord_offsetB = 0;
750
751     outeriter        = 0;
752     inneriter        = 0;
753
754     /* Start outer loop over neighborlists */
755     for(iidx=0; iidx<nri; iidx++)
756     {
757         /* Load shift vector for this list */
758         i_shift_offset   = DIM*shiftidx[iidx];
759
760         /* Load limits for loop over neighbors */
761         j_index_start    = jindex[iidx];
762         j_index_end      = jindex[iidx+1];
763
764         /* Get outer coordinate index */
765         inr              = iinr[iidx];
766         i_coord_offset   = DIM*inr;
767
768         /* Load i particle coords and add shift vector */
769         gmx_mm_load_shift_and_4rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
770                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
771
772         fix0             = _mm_setzero_pd();
773         fiy0             = _mm_setzero_pd();
774         fiz0             = _mm_setzero_pd();
775         fix1             = _mm_setzero_pd();
776         fiy1             = _mm_setzero_pd();
777         fiz1             = _mm_setzero_pd();
778         fix2             = _mm_setzero_pd();
779         fiy2             = _mm_setzero_pd();
780         fiz2             = _mm_setzero_pd();
781         fix3             = _mm_setzero_pd();
782         fiy3             = _mm_setzero_pd();
783         fiz3             = _mm_setzero_pd();
784
785         /* Start inner kernel loop */
786         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
787         {
788
789             /* Get j neighbor index, and coordinate index */
790             jnrA             = jjnr[jidx];
791             jnrB             = jjnr[jidx+1];
792             j_coord_offsetA  = DIM*jnrA;
793             j_coord_offsetB  = DIM*jnrB;
794
795             /* load j atom coordinates */
796             gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
797                                               &jx0,&jy0,&jz0);
798
799             /* Calculate displacement vector */
800             dx00             = _mm_sub_pd(ix0,jx0);
801             dy00             = _mm_sub_pd(iy0,jy0);
802             dz00             = _mm_sub_pd(iz0,jz0);
803             dx10             = _mm_sub_pd(ix1,jx0);
804             dy10             = _mm_sub_pd(iy1,jy0);
805             dz10             = _mm_sub_pd(iz1,jz0);
806             dx20             = _mm_sub_pd(ix2,jx0);
807             dy20             = _mm_sub_pd(iy2,jy0);
808             dz20             = _mm_sub_pd(iz2,jz0);
809             dx30             = _mm_sub_pd(ix3,jx0);
810             dy30             = _mm_sub_pd(iy3,jy0);
811             dz30             = _mm_sub_pd(iz3,jz0);
812
813             /* Calculate squared distance and things based on it */
814             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
815             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
816             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
817             rsq30            = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
818
819             rinv10           = gmx_mm_invsqrt_pd(rsq10);
820             rinv20           = gmx_mm_invsqrt_pd(rsq20);
821             rinv30           = gmx_mm_invsqrt_pd(rsq30);
822
823             rinvsq00         = gmx_mm_inv_pd(rsq00);
824
825             /* Load parameters for j particles */
826             jq0              = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
827             vdwjidx0A        = 2*vdwtype[jnrA+0];
828             vdwjidx0B        = 2*vdwtype[jnrB+0];
829
830             fjx0             = _mm_setzero_pd();
831             fjy0             = _mm_setzero_pd();
832             fjz0             = _mm_setzero_pd();
833
834             /**************************
835              * CALCULATE INTERACTIONS *
836              **************************/
837
838             /* Compute parameters for interactions between i and j atoms */
839             gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
840                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
841
842             /* LENNARD-JONES DISPERSION/REPULSION */
843
844             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
845             fvdw             = _mm_mul_pd(_mm_msub_pd(c12_00,rinvsix,c6_00),_mm_mul_pd(rinvsix,rinvsq00));
846
847             fscal            = fvdw;
848
849             /* Update vectorial force */
850             fix0             = _mm_macc_pd(dx00,fscal,fix0);
851             fiy0             = _mm_macc_pd(dy00,fscal,fiy0);
852             fiz0             = _mm_macc_pd(dz00,fscal,fiz0);
853             
854             fjx0             = _mm_macc_pd(dx00,fscal,fjx0);
855             fjy0             = _mm_macc_pd(dy00,fscal,fjy0);
856             fjz0             = _mm_macc_pd(dz00,fscal,fjz0);
857
858             /**************************
859              * CALCULATE INTERACTIONS *
860              **************************/
861
862             r10              = _mm_mul_pd(rsq10,rinv10);
863
864             /* Compute parameters for interactions between i and j atoms */
865             qq10             = _mm_mul_pd(iq1,jq0);
866
867             /* Calculate table index by multiplying r with table scale and truncate to integer */
868             rt               = _mm_mul_pd(r10,vftabscale);
869             vfitab           = _mm_cvttpd_epi32(rt);
870 #ifdef __XOP__
871             vfeps            = _mm_frcz_pd(rt);
872 #else
873             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
874 #endif
875             twovfeps         = _mm_add_pd(vfeps,vfeps);
876             vfitab           = _mm_slli_epi32(vfitab,2);
877
878             /* CUBIC SPLINE TABLE ELECTROSTATICS */
879             Y                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
880             F                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
881             GMX_MM_TRANSPOSE2_PD(Y,F);
882             G                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
883             H                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
884             GMX_MM_TRANSPOSE2_PD(G,H);
885             Fp               = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
886             FF               = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
887             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq10,FF),_mm_mul_pd(vftabscale,rinv10)));
888
889             fscal            = felec;
890
891             /* Update vectorial force */
892             fix1             = _mm_macc_pd(dx10,fscal,fix1);
893             fiy1             = _mm_macc_pd(dy10,fscal,fiy1);
894             fiz1             = _mm_macc_pd(dz10,fscal,fiz1);
895             
896             fjx0             = _mm_macc_pd(dx10,fscal,fjx0);
897             fjy0             = _mm_macc_pd(dy10,fscal,fjy0);
898             fjz0             = _mm_macc_pd(dz10,fscal,fjz0);
899
900             /**************************
901              * CALCULATE INTERACTIONS *
902              **************************/
903
904             r20              = _mm_mul_pd(rsq20,rinv20);
905
906             /* Compute parameters for interactions between i and j atoms */
907             qq20             = _mm_mul_pd(iq2,jq0);
908
909             /* Calculate table index by multiplying r with table scale and truncate to integer */
910             rt               = _mm_mul_pd(r20,vftabscale);
911             vfitab           = _mm_cvttpd_epi32(rt);
912 #ifdef __XOP__
913             vfeps            = _mm_frcz_pd(rt);
914 #else
915             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
916 #endif
917             twovfeps         = _mm_add_pd(vfeps,vfeps);
918             vfitab           = _mm_slli_epi32(vfitab,2);
919
920             /* CUBIC SPLINE TABLE ELECTROSTATICS */
921             Y                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
922             F                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
923             GMX_MM_TRANSPOSE2_PD(Y,F);
924             G                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
925             H                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
926             GMX_MM_TRANSPOSE2_PD(G,H);
927             Fp               = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
928             FF               = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
929             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq20,FF),_mm_mul_pd(vftabscale,rinv20)));
930
931             fscal            = felec;
932
933             /* Update vectorial force */
934             fix2             = _mm_macc_pd(dx20,fscal,fix2);
935             fiy2             = _mm_macc_pd(dy20,fscal,fiy2);
936             fiz2             = _mm_macc_pd(dz20,fscal,fiz2);
937             
938             fjx0             = _mm_macc_pd(dx20,fscal,fjx0);
939             fjy0             = _mm_macc_pd(dy20,fscal,fjy0);
940             fjz0             = _mm_macc_pd(dz20,fscal,fjz0);
941
942             /**************************
943              * CALCULATE INTERACTIONS *
944              **************************/
945
946             r30              = _mm_mul_pd(rsq30,rinv30);
947
948             /* Compute parameters for interactions between i and j atoms */
949             qq30             = _mm_mul_pd(iq3,jq0);
950
951             /* Calculate table index by multiplying r with table scale and truncate to integer */
952             rt               = _mm_mul_pd(r30,vftabscale);
953             vfitab           = _mm_cvttpd_epi32(rt);
954 #ifdef __XOP__
955             vfeps            = _mm_frcz_pd(rt);
956 #else
957             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
958 #endif
959             twovfeps         = _mm_add_pd(vfeps,vfeps);
960             vfitab           = _mm_slli_epi32(vfitab,2);
961
962             /* CUBIC SPLINE TABLE ELECTROSTATICS */
963             Y                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
964             F                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
965             GMX_MM_TRANSPOSE2_PD(Y,F);
966             G                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
967             H                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,1) +2);
968             GMX_MM_TRANSPOSE2_PD(G,H);
969             Fp               = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
970             FF               = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
971             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq30,FF),_mm_mul_pd(vftabscale,rinv30)));
972
973             fscal            = felec;
974
975             /* Update vectorial force */
976             fix3             = _mm_macc_pd(dx30,fscal,fix3);
977             fiy3             = _mm_macc_pd(dy30,fscal,fiy3);
978             fiz3             = _mm_macc_pd(dz30,fscal,fiz3);
979             
980             fjx0             = _mm_macc_pd(dx30,fscal,fjx0);
981             fjy0             = _mm_macc_pd(dy30,fscal,fjy0);
982             fjz0             = _mm_macc_pd(dz30,fscal,fjz0);
983
984             gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
985
986             /* Inner loop uses 159 flops */
987         }
988
989         if(jidx<j_index_end)
990         {
991
992             jnrA             = jjnr[jidx];
993             j_coord_offsetA  = DIM*jnrA;
994
995             /* load j atom coordinates */
996             gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
997                                               &jx0,&jy0,&jz0);
998
999             /* Calculate displacement vector */
1000             dx00             = _mm_sub_pd(ix0,jx0);
1001             dy00             = _mm_sub_pd(iy0,jy0);
1002             dz00             = _mm_sub_pd(iz0,jz0);
1003             dx10             = _mm_sub_pd(ix1,jx0);
1004             dy10             = _mm_sub_pd(iy1,jy0);
1005             dz10             = _mm_sub_pd(iz1,jz0);
1006             dx20             = _mm_sub_pd(ix2,jx0);
1007             dy20             = _mm_sub_pd(iy2,jy0);
1008             dz20             = _mm_sub_pd(iz2,jz0);
1009             dx30             = _mm_sub_pd(ix3,jx0);
1010             dy30             = _mm_sub_pd(iy3,jy0);
1011             dz30             = _mm_sub_pd(iz3,jz0);
1012
1013             /* Calculate squared distance and things based on it */
1014             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
1015             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
1016             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
1017             rsq30            = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
1018
1019             rinv10           = gmx_mm_invsqrt_pd(rsq10);
1020             rinv20           = gmx_mm_invsqrt_pd(rsq20);
1021             rinv30           = gmx_mm_invsqrt_pd(rsq30);
1022
1023             rinvsq00         = gmx_mm_inv_pd(rsq00);
1024
1025             /* Load parameters for j particles */
1026             jq0              = _mm_load_sd(charge+jnrA+0);
1027             vdwjidx0A        = 2*vdwtype[jnrA+0];
1028
1029             fjx0             = _mm_setzero_pd();
1030             fjy0             = _mm_setzero_pd();
1031             fjz0             = _mm_setzero_pd();
1032
1033             /**************************
1034              * CALCULATE INTERACTIONS *
1035              **************************/
1036
1037             /* Compute parameters for interactions between i and j atoms */
1038             gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
1039
1040             /* LENNARD-JONES DISPERSION/REPULSION */
1041
1042             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
1043             fvdw             = _mm_mul_pd(_mm_msub_pd(c12_00,rinvsix,c6_00),_mm_mul_pd(rinvsix,rinvsq00));
1044
1045             fscal            = fvdw;
1046
1047             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1048
1049             /* Update vectorial force */
1050             fix0             = _mm_macc_pd(dx00,fscal,fix0);
1051             fiy0             = _mm_macc_pd(dy00,fscal,fiy0);
1052             fiz0             = _mm_macc_pd(dz00,fscal,fiz0);
1053             
1054             fjx0             = _mm_macc_pd(dx00,fscal,fjx0);
1055             fjy0             = _mm_macc_pd(dy00,fscal,fjy0);
1056             fjz0             = _mm_macc_pd(dz00,fscal,fjz0);
1057
1058             /**************************
1059              * CALCULATE INTERACTIONS *
1060              **************************/
1061
1062             r10              = _mm_mul_pd(rsq10,rinv10);
1063
1064             /* Compute parameters for interactions between i and j atoms */
1065             qq10             = _mm_mul_pd(iq1,jq0);
1066
1067             /* Calculate table index by multiplying r with table scale and truncate to integer */
1068             rt               = _mm_mul_pd(r10,vftabscale);
1069             vfitab           = _mm_cvttpd_epi32(rt);
1070 #ifdef __XOP__
1071             vfeps            = _mm_frcz_pd(rt);
1072 #else
1073             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
1074 #endif
1075             twovfeps         = _mm_add_pd(vfeps,vfeps);
1076             vfitab           = _mm_slli_epi32(vfitab,2);
1077
1078             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1079             Y                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
1080             F                = _mm_setzero_pd();
1081             GMX_MM_TRANSPOSE2_PD(Y,F);
1082             G                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
1083             H                = _mm_setzero_pd();
1084             GMX_MM_TRANSPOSE2_PD(G,H);
1085             Fp               = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
1086             FF               = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
1087             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq10,FF),_mm_mul_pd(vftabscale,rinv10)));
1088
1089             fscal            = felec;
1090
1091             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1092
1093             /* Update vectorial force */
1094             fix1             = _mm_macc_pd(dx10,fscal,fix1);
1095             fiy1             = _mm_macc_pd(dy10,fscal,fiy1);
1096             fiz1             = _mm_macc_pd(dz10,fscal,fiz1);
1097             
1098             fjx0             = _mm_macc_pd(dx10,fscal,fjx0);
1099             fjy0             = _mm_macc_pd(dy10,fscal,fjy0);
1100             fjz0             = _mm_macc_pd(dz10,fscal,fjz0);
1101
1102             /**************************
1103              * CALCULATE INTERACTIONS *
1104              **************************/
1105
1106             r20              = _mm_mul_pd(rsq20,rinv20);
1107
1108             /* Compute parameters for interactions between i and j atoms */
1109             qq20             = _mm_mul_pd(iq2,jq0);
1110
1111             /* Calculate table index by multiplying r with table scale and truncate to integer */
1112             rt               = _mm_mul_pd(r20,vftabscale);
1113             vfitab           = _mm_cvttpd_epi32(rt);
1114 #ifdef __XOP__
1115             vfeps            = _mm_frcz_pd(rt);
1116 #else
1117             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
1118 #endif
1119             twovfeps         = _mm_add_pd(vfeps,vfeps);
1120             vfitab           = _mm_slli_epi32(vfitab,2);
1121
1122             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1123             Y                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
1124             F                = _mm_setzero_pd();
1125             GMX_MM_TRANSPOSE2_PD(Y,F);
1126             G                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
1127             H                = _mm_setzero_pd();
1128             GMX_MM_TRANSPOSE2_PD(G,H);
1129             Fp               = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
1130             FF               = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
1131             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq20,FF),_mm_mul_pd(vftabscale,rinv20)));
1132
1133             fscal            = felec;
1134
1135             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1136
1137             /* Update vectorial force */
1138             fix2             = _mm_macc_pd(dx20,fscal,fix2);
1139             fiy2             = _mm_macc_pd(dy20,fscal,fiy2);
1140             fiz2             = _mm_macc_pd(dz20,fscal,fiz2);
1141             
1142             fjx0             = _mm_macc_pd(dx20,fscal,fjx0);
1143             fjy0             = _mm_macc_pd(dy20,fscal,fjy0);
1144             fjz0             = _mm_macc_pd(dz20,fscal,fjz0);
1145
1146             /**************************
1147              * CALCULATE INTERACTIONS *
1148              **************************/
1149
1150             r30              = _mm_mul_pd(rsq30,rinv30);
1151
1152             /* Compute parameters for interactions between i and j atoms */
1153             qq30             = _mm_mul_pd(iq3,jq0);
1154
1155             /* Calculate table index by multiplying r with table scale and truncate to integer */
1156             rt               = _mm_mul_pd(r30,vftabscale);
1157             vfitab           = _mm_cvttpd_epi32(rt);
1158 #ifdef __XOP__
1159             vfeps            = _mm_frcz_pd(rt);
1160 #else
1161             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
1162 #endif
1163             twovfeps         = _mm_add_pd(vfeps,vfeps);
1164             vfitab           = _mm_slli_epi32(vfitab,2);
1165
1166             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1167             Y                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
1168             F                = _mm_setzero_pd();
1169             GMX_MM_TRANSPOSE2_PD(Y,F);
1170             G                = _mm_load_pd( vftab + _mm_extract_epi32(vfitab,0) +2);
1171             H                = _mm_setzero_pd();
1172             GMX_MM_TRANSPOSE2_PD(G,H);
1173             Fp               = _mm_macc_pd(vfeps,_mm_macc_pd(vfeps,H,G),F);
1174             FF               = _mm_macc_pd(_mm_macc_pd(twovfeps,H,G),vfeps,Fp);
1175             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq30,FF),_mm_mul_pd(vftabscale,rinv30)));
1176
1177             fscal            = felec;
1178
1179             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1180
1181             /* Update vectorial force */
1182             fix3             = _mm_macc_pd(dx30,fscal,fix3);
1183             fiy3             = _mm_macc_pd(dy30,fscal,fiy3);
1184             fiz3             = _mm_macc_pd(dz30,fscal,fiz3);
1185             
1186             fjx0             = _mm_macc_pd(dx30,fscal,fjx0);
1187             fjy0             = _mm_macc_pd(dy30,fscal,fjy0);
1188             fjz0             = _mm_macc_pd(dz30,fscal,fjz0);
1189
1190             gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
1191
1192             /* Inner loop uses 159 flops */
1193         }
1194
1195         /* End of innermost loop */
1196
1197         gmx_mm_update_iforce_4atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1198                                               f+i_coord_offset,fshift+i_shift_offset);
1199
1200         /* Increment number of inner iterations */
1201         inneriter                  += j_index_end - j_index_start;
1202
1203         /* Outer loop uses 24 flops */
1204     }
1205
1206     /* Increment number of outer iterations */
1207     outeriter        += nri;
1208
1209     /* Update outer/inner flops */
1210
1211     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*159);
1212 }