Use full path for legacyheaders
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_double / nb_kernel_ElecEw_VdwLJEw_GeomW4W4_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 "gromacs/legacyheaders/types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/legacyheaders/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_ElecEw_VdwLJEw_GeomW4W4_VF_avx_128_fma_double
52  * Electrostatics interaction: Ewald
53  * VdW interaction:            LJEwald
54  * Geometry:                   Water4-Water4
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecEw_VdwLJEw_GeomW4W4_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     int              vdwjidx1A,vdwjidx1B;
91     __m128d          jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
92     int              vdwjidx2A,vdwjidx2B;
93     __m128d          jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
94     int              vdwjidx3A,vdwjidx3B;
95     __m128d          jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
96     __m128d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
97     __m128d          dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
98     __m128d          dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
99     __m128d          dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
100     __m128d          dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
101     __m128d          dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
102     __m128d          dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
103     __m128d          dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
104     __m128d          dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
105     __m128d          dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
106     __m128d          velec,felec,velecsum,facel,crf,krf,krf2;
107     real             *charge;
108     int              nvdwtype;
109     __m128d          rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
110     int              *vdwtype;
111     real             *vdwparam;
112     __m128d          one_sixth   = _mm_set1_pd(1.0/6.0);
113     __m128d          one_twelfth = _mm_set1_pd(1.0/12.0);
114     __m128d           c6grid_00;
115     __m128d           c6grid_11;
116     __m128d           c6grid_12;
117     __m128d           c6grid_13;
118     __m128d           c6grid_21;
119     __m128d           c6grid_22;
120     __m128d           c6grid_23;
121     __m128d           c6grid_31;
122     __m128d           c6grid_32;
123     __m128d           c6grid_33;
124     real             *vdwgridparam;
125     __m128d           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
126     __m128d           one_half  = _mm_set1_pd(0.5);
127     __m128d           minus_one = _mm_set1_pd(-1.0);
128     __m128i          ewitab;
129     __m128d          ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
130     real             *ewtab;
131     __m128d          dummy_mask,cutoff_mask;
132     __m128d          signbit   = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
133     __m128d          one     = _mm_set1_pd(1.0);
134     __m128d          two     = _mm_set1_pd(2.0);
135     x                = xx[0];
136     f                = ff[0];
137
138     nri              = nlist->nri;
139     iinr             = nlist->iinr;
140     jindex           = nlist->jindex;
141     jjnr             = nlist->jjnr;
142     shiftidx         = nlist->shift;
143     gid              = nlist->gid;
144     shiftvec         = fr->shift_vec[0];
145     fshift           = fr->fshift[0];
146     facel            = _mm_set1_pd(fr->epsfac);
147     charge           = mdatoms->chargeA;
148     nvdwtype         = fr->ntype;
149     vdwparam         = fr->nbfp;
150     vdwtype          = mdatoms->typeA;
151     vdwgridparam     = fr->ljpme_c6grid;
152     sh_lj_ewald      = _mm_set1_pd(fr->ic->sh_lj_ewald);
153     ewclj            = _mm_set1_pd(fr->ewaldcoeff_lj);
154     ewclj2           = _mm_mul_pd(minus_one,_mm_mul_pd(ewclj,ewclj));
155
156     sh_ewald         = _mm_set1_pd(fr->ic->sh_ewald);
157     ewtab            = fr->ic->tabq_coul_FDV0;
158     ewtabscale       = _mm_set1_pd(fr->ic->tabq_scale);
159     ewtabhalfspace   = _mm_set1_pd(0.5/fr->ic->tabq_scale);
160
161     /* Setup water-specific parameters */
162     inr              = nlist->iinr[0];
163     iq1              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
164     iq2              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
165     iq3              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
166     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
167
168     jq1              = _mm_set1_pd(charge[inr+1]);
169     jq2              = _mm_set1_pd(charge[inr+2]);
170     jq3              = _mm_set1_pd(charge[inr+3]);
171     vdwjidx0A        = 2*vdwtype[inr+0];
172     c6_00            = _mm_set1_pd(vdwparam[vdwioffset0+vdwjidx0A]);
173     c12_00           = _mm_set1_pd(vdwparam[vdwioffset0+vdwjidx0A+1]);
174     c6grid_00        = _mm_set1_pd(vdwgridparam[vdwioffset0+vdwjidx0A]);
175     qq11             = _mm_mul_pd(iq1,jq1);
176     qq12             = _mm_mul_pd(iq1,jq2);
177     qq13             = _mm_mul_pd(iq1,jq3);
178     qq21             = _mm_mul_pd(iq2,jq1);
179     qq22             = _mm_mul_pd(iq2,jq2);
180     qq23             = _mm_mul_pd(iq2,jq3);
181     qq31             = _mm_mul_pd(iq3,jq1);
182     qq32             = _mm_mul_pd(iq3,jq2);
183     qq33             = _mm_mul_pd(iq3,jq3);
184
185     /* Avoid stupid compiler warnings */
186     jnrA = jnrB = 0;
187     j_coord_offsetA = 0;
188     j_coord_offsetB = 0;
189
190     outeriter        = 0;
191     inneriter        = 0;
192
193     /* Start outer loop over neighborlists */
194     for(iidx=0; iidx<nri; iidx++)
195     {
196         /* Load shift vector for this list */
197         i_shift_offset   = DIM*shiftidx[iidx];
198
199         /* Load limits for loop over neighbors */
200         j_index_start    = jindex[iidx];
201         j_index_end      = jindex[iidx+1];
202
203         /* Get outer coordinate index */
204         inr              = iinr[iidx];
205         i_coord_offset   = DIM*inr;
206
207         /* Load i particle coords and add shift vector */
208         gmx_mm_load_shift_and_4rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
209                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
210
211         fix0             = _mm_setzero_pd();
212         fiy0             = _mm_setzero_pd();
213         fiz0             = _mm_setzero_pd();
214         fix1             = _mm_setzero_pd();
215         fiy1             = _mm_setzero_pd();
216         fiz1             = _mm_setzero_pd();
217         fix2             = _mm_setzero_pd();
218         fiy2             = _mm_setzero_pd();
219         fiz2             = _mm_setzero_pd();
220         fix3             = _mm_setzero_pd();
221         fiy3             = _mm_setzero_pd();
222         fiz3             = _mm_setzero_pd();
223
224         /* Reset potential sums */
225         velecsum         = _mm_setzero_pd();
226         vvdwsum          = _mm_setzero_pd();
227
228         /* Start inner kernel loop */
229         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
230         {
231
232             /* Get j neighbor index, and coordinate index */
233             jnrA             = jjnr[jidx];
234             jnrB             = jjnr[jidx+1];
235             j_coord_offsetA  = DIM*jnrA;
236             j_coord_offsetB  = DIM*jnrB;
237
238             /* load j atom coordinates */
239             gmx_mm_load_4rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
240                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
241                                               &jy2,&jz2,&jx3,&jy3,&jz3);
242
243             /* Calculate displacement vector */
244             dx00             = _mm_sub_pd(ix0,jx0);
245             dy00             = _mm_sub_pd(iy0,jy0);
246             dz00             = _mm_sub_pd(iz0,jz0);
247             dx11             = _mm_sub_pd(ix1,jx1);
248             dy11             = _mm_sub_pd(iy1,jy1);
249             dz11             = _mm_sub_pd(iz1,jz1);
250             dx12             = _mm_sub_pd(ix1,jx2);
251             dy12             = _mm_sub_pd(iy1,jy2);
252             dz12             = _mm_sub_pd(iz1,jz2);
253             dx13             = _mm_sub_pd(ix1,jx3);
254             dy13             = _mm_sub_pd(iy1,jy3);
255             dz13             = _mm_sub_pd(iz1,jz3);
256             dx21             = _mm_sub_pd(ix2,jx1);
257             dy21             = _mm_sub_pd(iy2,jy1);
258             dz21             = _mm_sub_pd(iz2,jz1);
259             dx22             = _mm_sub_pd(ix2,jx2);
260             dy22             = _mm_sub_pd(iy2,jy2);
261             dz22             = _mm_sub_pd(iz2,jz2);
262             dx23             = _mm_sub_pd(ix2,jx3);
263             dy23             = _mm_sub_pd(iy2,jy3);
264             dz23             = _mm_sub_pd(iz2,jz3);
265             dx31             = _mm_sub_pd(ix3,jx1);
266             dy31             = _mm_sub_pd(iy3,jy1);
267             dz31             = _mm_sub_pd(iz3,jz1);
268             dx32             = _mm_sub_pd(ix3,jx2);
269             dy32             = _mm_sub_pd(iy3,jy2);
270             dz32             = _mm_sub_pd(iz3,jz2);
271             dx33             = _mm_sub_pd(ix3,jx3);
272             dy33             = _mm_sub_pd(iy3,jy3);
273             dz33             = _mm_sub_pd(iz3,jz3);
274
275             /* Calculate squared distance and things based on it */
276             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
277             rsq11            = gmx_mm_calc_rsq_pd(dx11,dy11,dz11);
278             rsq12            = gmx_mm_calc_rsq_pd(dx12,dy12,dz12);
279             rsq13            = gmx_mm_calc_rsq_pd(dx13,dy13,dz13);
280             rsq21            = gmx_mm_calc_rsq_pd(dx21,dy21,dz21);
281             rsq22            = gmx_mm_calc_rsq_pd(dx22,dy22,dz22);
282             rsq23            = gmx_mm_calc_rsq_pd(dx23,dy23,dz23);
283             rsq31            = gmx_mm_calc_rsq_pd(dx31,dy31,dz31);
284             rsq32            = gmx_mm_calc_rsq_pd(dx32,dy32,dz32);
285             rsq33            = gmx_mm_calc_rsq_pd(dx33,dy33,dz33);
286
287             rinv00           = gmx_mm_invsqrt_pd(rsq00);
288             rinv11           = gmx_mm_invsqrt_pd(rsq11);
289             rinv12           = gmx_mm_invsqrt_pd(rsq12);
290             rinv13           = gmx_mm_invsqrt_pd(rsq13);
291             rinv21           = gmx_mm_invsqrt_pd(rsq21);
292             rinv22           = gmx_mm_invsqrt_pd(rsq22);
293             rinv23           = gmx_mm_invsqrt_pd(rsq23);
294             rinv31           = gmx_mm_invsqrt_pd(rsq31);
295             rinv32           = gmx_mm_invsqrt_pd(rsq32);
296             rinv33           = gmx_mm_invsqrt_pd(rsq33);
297
298             rinvsq00         = _mm_mul_pd(rinv00,rinv00);
299             rinvsq11         = _mm_mul_pd(rinv11,rinv11);
300             rinvsq12         = _mm_mul_pd(rinv12,rinv12);
301             rinvsq13         = _mm_mul_pd(rinv13,rinv13);
302             rinvsq21         = _mm_mul_pd(rinv21,rinv21);
303             rinvsq22         = _mm_mul_pd(rinv22,rinv22);
304             rinvsq23         = _mm_mul_pd(rinv23,rinv23);
305             rinvsq31         = _mm_mul_pd(rinv31,rinv31);
306             rinvsq32         = _mm_mul_pd(rinv32,rinv32);
307             rinvsq33         = _mm_mul_pd(rinv33,rinv33);
308
309             fjx0             = _mm_setzero_pd();
310             fjy0             = _mm_setzero_pd();
311             fjz0             = _mm_setzero_pd();
312             fjx1             = _mm_setzero_pd();
313             fjy1             = _mm_setzero_pd();
314             fjz1             = _mm_setzero_pd();
315             fjx2             = _mm_setzero_pd();
316             fjy2             = _mm_setzero_pd();
317             fjz2             = _mm_setzero_pd();
318             fjx3             = _mm_setzero_pd();
319             fjy3             = _mm_setzero_pd();
320             fjz3             = _mm_setzero_pd();
321
322             /**************************
323              * CALCULATE INTERACTIONS *
324              **************************/
325
326             r00              = _mm_mul_pd(rsq00,rinv00);
327
328             /* Analytical LJ-PME */
329             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
330             ewcljrsq         = _mm_mul_pd(ewclj2,rsq00);
331             ewclj6           = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
332             exponent         = gmx_simd_exp_d(ewcljrsq);
333             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
334             poly             = _mm_mul_pd(exponent,_mm_macc_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half,_mm_sub_pd(one,ewcljrsq)));
335             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
336             vvdw6            = _mm_mul_pd(_mm_macc_pd(-c6grid_00,_mm_sub_pd(one,poly),c6_00),rinvsix);
337             vvdw12           = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
338             vvdw             = _mm_msub_pd(vvdw12,one_twelfth,_mm_mul_pd(vvdw6,one_sixth));
339             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
340             fvdw             = _mm_mul_pd(_mm_add_pd(vvdw12,_mm_msub_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6),vvdw6)),rinvsq00);
341
342             /* Update potential sum for this i atom from the interaction with this j atom. */
343             vvdwsum          = _mm_add_pd(vvdwsum,vvdw);
344
345             fscal            = fvdw;
346
347             /* Update vectorial force */
348             fix0             = _mm_macc_pd(dx00,fscal,fix0);
349             fiy0             = _mm_macc_pd(dy00,fscal,fiy0);
350             fiz0             = _mm_macc_pd(dz00,fscal,fiz0);
351             
352             fjx0             = _mm_macc_pd(dx00,fscal,fjx0);
353             fjy0             = _mm_macc_pd(dy00,fscal,fjy0);
354             fjz0             = _mm_macc_pd(dz00,fscal,fjz0);
355
356             /**************************
357              * CALCULATE INTERACTIONS *
358              **************************/
359
360             r11              = _mm_mul_pd(rsq11,rinv11);
361
362             /* EWALD ELECTROSTATICS */
363
364             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
365             ewrt             = _mm_mul_pd(r11,ewtabscale);
366             ewitab           = _mm_cvttpd_epi32(ewrt);
367 #ifdef __XOP__
368             eweps            = _mm_frcz_pd(ewrt);
369 #else
370             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
371 #endif
372             twoeweps         = _mm_add_pd(eweps,eweps);
373             ewitab           = _mm_slli_epi32(ewitab,2);
374             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
375             ewtabD           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
376             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
377             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
378             ewtabFn          = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,1) +2);
379             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
380             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
381             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
382             velec            = _mm_mul_pd(qq11,_mm_sub_pd(rinv11,velec));
383             felec            = _mm_mul_pd(_mm_mul_pd(qq11,rinv11),_mm_sub_pd(rinvsq11,felec));
384
385             /* Update potential sum for this i atom from the interaction with this j atom. */
386             velecsum         = _mm_add_pd(velecsum,velec);
387
388             fscal            = felec;
389
390             /* Update vectorial force */
391             fix1             = _mm_macc_pd(dx11,fscal,fix1);
392             fiy1             = _mm_macc_pd(dy11,fscal,fiy1);
393             fiz1             = _mm_macc_pd(dz11,fscal,fiz1);
394             
395             fjx1             = _mm_macc_pd(dx11,fscal,fjx1);
396             fjy1             = _mm_macc_pd(dy11,fscal,fjy1);
397             fjz1             = _mm_macc_pd(dz11,fscal,fjz1);
398
399             /**************************
400              * CALCULATE INTERACTIONS *
401              **************************/
402
403             r12              = _mm_mul_pd(rsq12,rinv12);
404
405             /* EWALD ELECTROSTATICS */
406
407             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
408             ewrt             = _mm_mul_pd(r12,ewtabscale);
409             ewitab           = _mm_cvttpd_epi32(ewrt);
410 #ifdef __XOP__
411             eweps            = _mm_frcz_pd(ewrt);
412 #else
413             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
414 #endif
415             twoeweps         = _mm_add_pd(eweps,eweps);
416             ewitab           = _mm_slli_epi32(ewitab,2);
417             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
418             ewtabD           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
419             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
420             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
421             ewtabFn          = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,1) +2);
422             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
423             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
424             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
425             velec            = _mm_mul_pd(qq12,_mm_sub_pd(rinv12,velec));
426             felec            = _mm_mul_pd(_mm_mul_pd(qq12,rinv12),_mm_sub_pd(rinvsq12,felec));
427
428             /* Update potential sum for this i atom from the interaction with this j atom. */
429             velecsum         = _mm_add_pd(velecsum,velec);
430
431             fscal            = felec;
432
433             /* Update vectorial force */
434             fix1             = _mm_macc_pd(dx12,fscal,fix1);
435             fiy1             = _mm_macc_pd(dy12,fscal,fiy1);
436             fiz1             = _mm_macc_pd(dz12,fscal,fiz1);
437             
438             fjx2             = _mm_macc_pd(dx12,fscal,fjx2);
439             fjy2             = _mm_macc_pd(dy12,fscal,fjy2);
440             fjz2             = _mm_macc_pd(dz12,fscal,fjz2);
441
442             /**************************
443              * CALCULATE INTERACTIONS *
444              **************************/
445
446             r13              = _mm_mul_pd(rsq13,rinv13);
447
448             /* EWALD ELECTROSTATICS */
449
450             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
451             ewrt             = _mm_mul_pd(r13,ewtabscale);
452             ewitab           = _mm_cvttpd_epi32(ewrt);
453 #ifdef __XOP__
454             eweps            = _mm_frcz_pd(ewrt);
455 #else
456             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
457 #endif
458             twoeweps         = _mm_add_pd(eweps,eweps);
459             ewitab           = _mm_slli_epi32(ewitab,2);
460             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
461             ewtabD           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
462             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
463             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
464             ewtabFn          = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,1) +2);
465             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
466             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
467             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
468             velec            = _mm_mul_pd(qq13,_mm_sub_pd(rinv13,velec));
469             felec            = _mm_mul_pd(_mm_mul_pd(qq13,rinv13),_mm_sub_pd(rinvsq13,felec));
470
471             /* Update potential sum for this i atom from the interaction with this j atom. */
472             velecsum         = _mm_add_pd(velecsum,velec);
473
474             fscal            = felec;
475
476             /* Update vectorial force */
477             fix1             = _mm_macc_pd(dx13,fscal,fix1);
478             fiy1             = _mm_macc_pd(dy13,fscal,fiy1);
479             fiz1             = _mm_macc_pd(dz13,fscal,fiz1);
480             
481             fjx3             = _mm_macc_pd(dx13,fscal,fjx3);
482             fjy3             = _mm_macc_pd(dy13,fscal,fjy3);
483             fjz3             = _mm_macc_pd(dz13,fscal,fjz3);
484
485             /**************************
486              * CALCULATE INTERACTIONS *
487              **************************/
488
489             r21              = _mm_mul_pd(rsq21,rinv21);
490
491             /* EWALD ELECTROSTATICS */
492
493             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
494             ewrt             = _mm_mul_pd(r21,ewtabscale);
495             ewitab           = _mm_cvttpd_epi32(ewrt);
496 #ifdef __XOP__
497             eweps            = _mm_frcz_pd(ewrt);
498 #else
499             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
500 #endif
501             twoeweps         = _mm_add_pd(eweps,eweps);
502             ewitab           = _mm_slli_epi32(ewitab,2);
503             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
504             ewtabD           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
505             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
506             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
507             ewtabFn          = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,1) +2);
508             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
509             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
510             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
511             velec            = _mm_mul_pd(qq21,_mm_sub_pd(rinv21,velec));
512             felec            = _mm_mul_pd(_mm_mul_pd(qq21,rinv21),_mm_sub_pd(rinvsq21,felec));
513
514             /* Update potential sum for this i atom from the interaction with this j atom. */
515             velecsum         = _mm_add_pd(velecsum,velec);
516
517             fscal            = felec;
518
519             /* Update vectorial force */
520             fix2             = _mm_macc_pd(dx21,fscal,fix2);
521             fiy2             = _mm_macc_pd(dy21,fscal,fiy2);
522             fiz2             = _mm_macc_pd(dz21,fscal,fiz2);
523             
524             fjx1             = _mm_macc_pd(dx21,fscal,fjx1);
525             fjy1             = _mm_macc_pd(dy21,fscal,fjy1);
526             fjz1             = _mm_macc_pd(dz21,fscal,fjz1);
527
528             /**************************
529              * CALCULATE INTERACTIONS *
530              **************************/
531
532             r22              = _mm_mul_pd(rsq22,rinv22);
533
534             /* EWALD ELECTROSTATICS */
535
536             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
537             ewrt             = _mm_mul_pd(r22,ewtabscale);
538             ewitab           = _mm_cvttpd_epi32(ewrt);
539 #ifdef __XOP__
540             eweps            = _mm_frcz_pd(ewrt);
541 #else
542             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
543 #endif
544             twoeweps         = _mm_add_pd(eweps,eweps);
545             ewitab           = _mm_slli_epi32(ewitab,2);
546             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
547             ewtabD           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
548             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
549             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
550             ewtabFn          = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,1) +2);
551             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
552             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
553             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
554             velec            = _mm_mul_pd(qq22,_mm_sub_pd(rinv22,velec));
555             felec            = _mm_mul_pd(_mm_mul_pd(qq22,rinv22),_mm_sub_pd(rinvsq22,felec));
556
557             /* Update potential sum for this i atom from the interaction with this j atom. */
558             velecsum         = _mm_add_pd(velecsum,velec);
559
560             fscal            = felec;
561
562             /* Update vectorial force */
563             fix2             = _mm_macc_pd(dx22,fscal,fix2);
564             fiy2             = _mm_macc_pd(dy22,fscal,fiy2);
565             fiz2             = _mm_macc_pd(dz22,fscal,fiz2);
566             
567             fjx2             = _mm_macc_pd(dx22,fscal,fjx2);
568             fjy2             = _mm_macc_pd(dy22,fscal,fjy2);
569             fjz2             = _mm_macc_pd(dz22,fscal,fjz2);
570
571             /**************************
572              * CALCULATE INTERACTIONS *
573              **************************/
574
575             r23              = _mm_mul_pd(rsq23,rinv23);
576
577             /* EWALD ELECTROSTATICS */
578
579             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
580             ewrt             = _mm_mul_pd(r23,ewtabscale);
581             ewitab           = _mm_cvttpd_epi32(ewrt);
582 #ifdef __XOP__
583             eweps            = _mm_frcz_pd(ewrt);
584 #else
585             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
586 #endif
587             twoeweps         = _mm_add_pd(eweps,eweps);
588             ewitab           = _mm_slli_epi32(ewitab,2);
589             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
590             ewtabD           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
591             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
592             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
593             ewtabFn          = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,1) +2);
594             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
595             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
596             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
597             velec            = _mm_mul_pd(qq23,_mm_sub_pd(rinv23,velec));
598             felec            = _mm_mul_pd(_mm_mul_pd(qq23,rinv23),_mm_sub_pd(rinvsq23,felec));
599
600             /* Update potential sum for this i atom from the interaction with this j atom. */
601             velecsum         = _mm_add_pd(velecsum,velec);
602
603             fscal            = felec;
604
605             /* Update vectorial force */
606             fix2             = _mm_macc_pd(dx23,fscal,fix2);
607             fiy2             = _mm_macc_pd(dy23,fscal,fiy2);
608             fiz2             = _mm_macc_pd(dz23,fscal,fiz2);
609             
610             fjx3             = _mm_macc_pd(dx23,fscal,fjx3);
611             fjy3             = _mm_macc_pd(dy23,fscal,fjy3);
612             fjz3             = _mm_macc_pd(dz23,fscal,fjz3);
613
614             /**************************
615              * CALCULATE INTERACTIONS *
616              **************************/
617
618             r31              = _mm_mul_pd(rsq31,rinv31);
619
620             /* EWALD ELECTROSTATICS */
621
622             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
623             ewrt             = _mm_mul_pd(r31,ewtabscale);
624             ewitab           = _mm_cvttpd_epi32(ewrt);
625 #ifdef __XOP__
626             eweps            = _mm_frcz_pd(ewrt);
627 #else
628             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
629 #endif
630             twoeweps         = _mm_add_pd(eweps,eweps);
631             ewitab           = _mm_slli_epi32(ewitab,2);
632             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
633             ewtabD           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
634             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
635             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
636             ewtabFn          = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,1) +2);
637             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
638             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
639             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
640             velec            = _mm_mul_pd(qq31,_mm_sub_pd(rinv31,velec));
641             felec            = _mm_mul_pd(_mm_mul_pd(qq31,rinv31),_mm_sub_pd(rinvsq31,felec));
642
643             /* Update potential sum for this i atom from the interaction with this j atom. */
644             velecsum         = _mm_add_pd(velecsum,velec);
645
646             fscal            = felec;
647
648             /* Update vectorial force */
649             fix3             = _mm_macc_pd(dx31,fscal,fix3);
650             fiy3             = _mm_macc_pd(dy31,fscal,fiy3);
651             fiz3             = _mm_macc_pd(dz31,fscal,fiz3);
652             
653             fjx1             = _mm_macc_pd(dx31,fscal,fjx1);
654             fjy1             = _mm_macc_pd(dy31,fscal,fjy1);
655             fjz1             = _mm_macc_pd(dz31,fscal,fjz1);
656
657             /**************************
658              * CALCULATE INTERACTIONS *
659              **************************/
660
661             r32              = _mm_mul_pd(rsq32,rinv32);
662
663             /* EWALD ELECTROSTATICS */
664
665             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
666             ewrt             = _mm_mul_pd(r32,ewtabscale);
667             ewitab           = _mm_cvttpd_epi32(ewrt);
668 #ifdef __XOP__
669             eweps            = _mm_frcz_pd(ewrt);
670 #else
671             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
672 #endif
673             twoeweps         = _mm_add_pd(eweps,eweps);
674             ewitab           = _mm_slli_epi32(ewitab,2);
675             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
676             ewtabD           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
677             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
678             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
679             ewtabFn          = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,1) +2);
680             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
681             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
682             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
683             velec            = _mm_mul_pd(qq32,_mm_sub_pd(rinv32,velec));
684             felec            = _mm_mul_pd(_mm_mul_pd(qq32,rinv32),_mm_sub_pd(rinvsq32,felec));
685
686             /* Update potential sum for this i atom from the interaction with this j atom. */
687             velecsum         = _mm_add_pd(velecsum,velec);
688
689             fscal            = felec;
690
691             /* Update vectorial force */
692             fix3             = _mm_macc_pd(dx32,fscal,fix3);
693             fiy3             = _mm_macc_pd(dy32,fscal,fiy3);
694             fiz3             = _mm_macc_pd(dz32,fscal,fiz3);
695             
696             fjx2             = _mm_macc_pd(dx32,fscal,fjx2);
697             fjy2             = _mm_macc_pd(dy32,fscal,fjy2);
698             fjz2             = _mm_macc_pd(dz32,fscal,fjz2);
699
700             /**************************
701              * CALCULATE INTERACTIONS *
702              **************************/
703
704             r33              = _mm_mul_pd(rsq33,rinv33);
705
706             /* EWALD ELECTROSTATICS */
707
708             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
709             ewrt             = _mm_mul_pd(r33,ewtabscale);
710             ewitab           = _mm_cvttpd_epi32(ewrt);
711 #ifdef __XOP__
712             eweps            = _mm_frcz_pd(ewrt);
713 #else
714             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
715 #endif
716             twoeweps         = _mm_add_pd(eweps,eweps);
717             ewitab           = _mm_slli_epi32(ewitab,2);
718             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
719             ewtabD           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
720             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
721             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
722             ewtabFn          = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,1) +2);
723             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
724             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
725             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
726             velec            = _mm_mul_pd(qq33,_mm_sub_pd(rinv33,velec));
727             felec            = _mm_mul_pd(_mm_mul_pd(qq33,rinv33),_mm_sub_pd(rinvsq33,felec));
728
729             /* Update potential sum for this i atom from the interaction with this j atom. */
730             velecsum         = _mm_add_pd(velecsum,velec);
731
732             fscal            = felec;
733
734             /* Update vectorial force */
735             fix3             = _mm_macc_pd(dx33,fscal,fix3);
736             fiy3             = _mm_macc_pd(dy33,fscal,fiy3);
737             fiz3             = _mm_macc_pd(dz33,fscal,fiz3);
738             
739             fjx3             = _mm_macc_pd(dx33,fscal,fjx3);
740             fjy3             = _mm_macc_pd(dy33,fscal,fjy3);
741             fjz3             = _mm_macc_pd(dz33,fscal,fjz3);
742
743             gmx_mm_decrement_4rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
744
745             /* Inner loop uses 449 flops */
746         }
747
748         if(jidx<j_index_end)
749         {
750
751             jnrA             = jjnr[jidx];
752             j_coord_offsetA  = DIM*jnrA;
753
754             /* load j atom coordinates */
755             gmx_mm_load_4rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
756                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
757                                               &jy2,&jz2,&jx3,&jy3,&jz3);
758
759             /* Calculate displacement vector */
760             dx00             = _mm_sub_pd(ix0,jx0);
761             dy00             = _mm_sub_pd(iy0,jy0);
762             dz00             = _mm_sub_pd(iz0,jz0);
763             dx11             = _mm_sub_pd(ix1,jx1);
764             dy11             = _mm_sub_pd(iy1,jy1);
765             dz11             = _mm_sub_pd(iz1,jz1);
766             dx12             = _mm_sub_pd(ix1,jx2);
767             dy12             = _mm_sub_pd(iy1,jy2);
768             dz12             = _mm_sub_pd(iz1,jz2);
769             dx13             = _mm_sub_pd(ix1,jx3);
770             dy13             = _mm_sub_pd(iy1,jy3);
771             dz13             = _mm_sub_pd(iz1,jz3);
772             dx21             = _mm_sub_pd(ix2,jx1);
773             dy21             = _mm_sub_pd(iy2,jy1);
774             dz21             = _mm_sub_pd(iz2,jz1);
775             dx22             = _mm_sub_pd(ix2,jx2);
776             dy22             = _mm_sub_pd(iy2,jy2);
777             dz22             = _mm_sub_pd(iz2,jz2);
778             dx23             = _mm_sub_pd(ix2,jx3);
779             dy23             = _mm_sub_pd(iy2,jy3);
780             dz23             = _mm_sub_pd(iz2,jz3);
781             dx31             = _mm_sub_pd(ix3,jx1);
782             dy31             = _mm_sub_pd(iy3,jy1);
783             dz31             = _mm_sub_pd(iz3,jz1);
784             dx32             = _mm_sub_pd(ix3,jx2);
785             dy32             = _mm_sub_pd(iy3,jy2);
786             dz32             = _mm_sub_pd(iz3,jz2);
787             dx33             = _mm_sub_pd(ix3,jx3);
788             dy33             = _mm_sub_pd(iy3,jy3);
789             dz33             = _mm_sub_pd(iz3,jz3);
790
791             /* Calculate squared distance and things based on it */
792             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
793             rsq11            = gmx_mm_calc_rsq_pd(dx11,dy11,dz11);
794             rsq12            = gmx_mm_calc_rsq_pd(dx12,dy12,dz12);
795             rsq13            = gmx_mm_calc_rsq_pd(dx13,dy13,dz13);
796             rsq21            = gmx_mm_calc_rsq_pd(dx21,dy21,dz21);
797             rsq22            = gmx_mm_calc_rsq_pd(dx22,dy22,dz22);
798             rsq23            = gmx_mm_calc_rsq_pd(dx23,dy23,dz23);
799             rsq31            = gmx_mm_calc_rsq_pd(dx31,dy31,dz31);
800             rsq32            = gmx_mm_calc_rsq_pd(dx32,dy32,dz32);
801             rsq33            = gmx_mm_calc_rsq_pd(dx33,dy33,dz33);
802
803             rinv00           = gmx_mm_invsqrt_pd(rsq00);
804             rinv11           = gmx_mm_invsqrt_pd(rsq11);
805             rinv12           = gmx_mm_invsqrt_pd(rsq12);
806             rinv13           = gmx_mm_invsqrt_pd(rsq13);
807             rinv21           = gmx_mm_invsqrt_pd(rsq21);
808             rinv22           = gmx_mm_invsqrt_pd(rsq22);
809             rinv23           = gmx_mm_invsqrt_pd(rsq23);
810             rinv31           = gmx_mm_invsqrt_pd(rsq31);
811             rinv32           = gmx_mm_invsqrt_pd(rsq32);
812             rinv33           = gmx_mm_invsqrt_pd(rsq33);
813
814             rinvsq00         = _mm_mul_pd(rinv00,rinv00);
815             rinvsq11         = _mm_mul_pd(rinv11,rinv11);
816             rinvsq12         = _mm_mul_pd(rinv12,rinv12);
817             rinvsq13         = _mm_mul_pd(rinv13,rinv13);
818             rinvsq21         = _mm_mul_pd(rinv21,rinv21);
819             rinvsq22         = _mm_mul_pd(rinv22,rinv22);
820             rinvsq23         = _mm_mul_pd(rinv23,rinv23);
821             rinvsq31         = _mm_mul_pd(rinv31,rinv31);
822             rinvsq32         = _mm_mul_pd(rinv32,rinv32);
823             rinvsq33         = _mm_mul_pd(rinv33,rinv33);
824
825             fjx0             = _mm_setzero_pd();
826             fjy0             = _mm_setzero_pd();
827             fjz0             = _mm_setzero_pd();
828             fjx1             = _mm_setzero_pd();
829             fjy1             = _mm_setzero_pd();
830             fjz1             = _mm_setzero_pd();
831             fjx2             = _mm_setzero_pd();
832             fjy2             = _mm_setzero_pd();
833             fjz2             = _mm_setzero_pd();
834             fjx3             = _mm_setzero_pd();
835             fjy3             = _mm_setzero_pd();
836             fjz3             = _mm_setzero_pd();
837
838             /**************************
839              * CALCULATE INTERACTIONS *
840              **************************/
841
842             r00              = _mm_mul_pd(rsq00,rinv00);
843
844             /* Analytical LJ-PME */
845             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
846             ewcljrsq         = _mm_mul_pd(ewclj2,rsq00);
847             ewclj6           = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
848             exponent         = gmx_simd_exp_d(ewcljrsq);
849             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
850             poly             = _mm_mul_pd(exponent,_mm_macc_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half,_mm_sub_pd(one,ewcljrsq)));
851             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
852             vvdw6            = _mm_mul_pd(_mm_macc_pd(-c6grid_00,_mm_sub_pd(one,poly),c6_00),rinvsix);
853             vvdw12           = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
854             vvdw             = _mm_msub_pd(vvdw12,one_twelfth,_mm_mul_pd(vvdw6,one_sixth));
855             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
856             fvdw             = _mm_mul_pd(_mm_add_pd(vvdw12,_mm_msub_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6),vvdw6)),rinvsq00);
857
858             /* Update potential sum for this i atom from the interaction with this j atom. */
859             vvdw             = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
860             vvdwsum          = _mm_add_pd(vvdwsum,vvdw);
861
862             fscal            = fvdw;
863
864             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
865
866             /* Update vectorial force */
867             fix0             = _mm_macc_pd(dx00,fscal,fix0);
868             fiy0             = _mm_macc_pd(dy00,fscal,fiy0);
869             fiz0             = _mm_macc_pd(dz00,fscal,fiz0);
870             
871             fjx0             = _mm_macc_pd(dx00,fscal,fjx0);
872             fjy0             = _mm_macc_pd(dy00,fscal,fjy0);
873             fjz0             = _mm_macc_pd(dz00,fscal,fjz0);
874
875             /**************************
876              * CALCULATE INTERACTIONS *
877              **************************/
878
879             r11              = _mm_mul_pd(rsq11,rinv11);
880
881             /* EWALD ELECTROSTATICS */
882
883             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
884             ewrt             = _mm_mul_pd(r11,ewtabscale);
885             ewitab           = _mm_cvttpd_epi32(ewrt);
886 #ifdef __XOP__
887             eweps            = _mm_frcz_pd(ewrt);
888 #else
889             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
890 #endif
891             twoeweps         = _mm_add_pd(eweps,eweps);
892             ewitab           = _mm_slli_epi32(ewitab,2);
893             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
894             ewtabD           = _mm_setzero_pd();
895             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
896             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
897             ewtabFn          = _mm_setzero_pd();
898             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
899             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
900             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
901             velec            = _mm_mul_pd(qq11,_mm_sub_pd(rinv11,velec));
902             felec            = _mm_mul_pd(_mm_mul_pd(qq11,rinv11),_mm_sub_pd(rinvsq11,felec));
903
904             /* Update potential sum for this i atom from the interaction with this j atom. */
905             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
906             velecsum         = _mm_add_pd(velecsum,velec);
907
908             fscal            = felec;
909
910             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
911
912             /* Update vectorial force */
913             fix1             = _mm_macc_pd(dx11,fscal,fix1);
914             fiy1             = _mm_macc_pd(dy11,fscal,fiy1);
915             fiz1             = _mm_macc_pd(dz11,fscal,fiz1);
916             
917             fjx1             = _mm_macc_pd(dx11,fscal,fjx1);
918             fjy1             = _mm_macc_pd(dy11,fscal,fjy1);
919             fjz1             = _mm_macc_pd(dz11,fscal,fjz1);
920
921             /**************************
922              * CALCULATE INTERACTIONS *
923              **************************/
924
925             r12              = _mm_mul_pd(rsq12,rinv12);
926
927             /* EWALD ELECTROSTATICS */
928
929             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
930             ewrt             = _mm_mul_pd(r12,ewtabscale);
931             ewitab           = _mm_cvttpd_epi32(ewrt);
932 #ifdef __XOP__
933             eweps            = _mm_frcz_pd(ewrt);
934 #else
935             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
936 #endif
937             twoeweps         = _mm_add_pd(eweps,eweps);
938             ewitab           = _mm_slli_epi32(ewitab,2);
939             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
940             ewtabD           = _mm_setzero_pd();
941             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
942             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
943             ewtabFn          = _mm_setzero_pd();
944             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
945             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
946             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
947             velec            = _mm_mul_pd(qq12,_mm_sub_pd(rinv12,velec));
948             felec            = _mm_mul_pd(_mm_mul_pd(qq12,rinv12),_mm_sub_pd(rinvsq12,felec));
949
950             /* Update potential sum for this i atom from the interaction with this j atom. */
951             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
952             velecsum         = _mm_add_pd(velecsum,velec);
953
954             fscal            = felec;
955
956             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
957
958             /* Update vectorial force */
959             fix1             = _mm_macc_pd(dx12,fscal,fix1);
960             fiy1             = _mm_macc_pd(dy12,fscal,fiy1);
961             fiz1             = _mm_macc_pd(dz12,fscal,fiz1);
962             
963             fjx2             = _mm_macc_pd(dx12,fscal,fjx2);
964             fjy2             = _mm_macc_pd(dy12,fscal,fjy2);
965             fjz2             = _mm_macc_pd(dz12,fscal,fjz2);
966
967             /**************************
968              * CALCULATE INTERACTIONS *
969              **************************/
970
971             r13              = _mm_mul_pd(rsq13,rinv13);
972
973             /* EWALD ELECTROSTATICS */
974
975             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
976             ewrt             = _mm_mul_pd(r13,ewtabscale);
977             ewitab           = _mm_cvttpd_epi32(ewrt);
978 #ifdef __XOP__
979             eweps            = _mm_frcz_pd(ewrt);
980 #else
981             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
982 #endif
983             twoeweps         = _mm_add_pd(eweps,eweps);
984             ewitab           = _mm_slli_epi32(ewitab,2);
985             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
986             ewtabD           = _mm_setzero_pd();
987             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
988             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
989             ewtabFn          = _mm_setzero_pd();
990             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
991             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
992             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
993             velec            = _mm_mul_pd(qq13,_mm_sub_pd(rinv13,velec));
994             felec            = _mm_mul_pd(_mm_mul_pd(qq13,rinv13),_mm_sub_pd(rinvsq13,felec));
995
996             /* Update potential sum for this i atom from the interaction with this j atom. */
997             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
998             velecsum         = _mm_add_pd(velecsum,velec);
999
1000             fscal            = felec;
1001
1002             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1003
1004             /* Update vectorial force */
1005             fix1             = _mm_macc_pd(dx13,fscal,fix1);
1006             fiy1             = _mm_macc_pd(dy13,fscal,fiy1);
1007             fiz1             = _mm_macc_pd(dz13,fscal,fiz1);
1008             
1009             fjx3             = _mm_macc_pd(dx13,fscal,fjx3);
1010             fjy3             = _mm_macc_pd(dy13,fscal,fjy3);
1011             fjz3             = _mm_macc_pd(dz13,fscal,fjz3);
1012
1013             /**************************
1014              * CALCULATE INTERACTIONS *
1015              **************************/
1016
1017             r21              = _mm_mul_pd(rsq21,rinv21);
1018
1019             /* EWALD ELECTROSTATICS */
1020
1021             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1022             ewrt             = _mm_mul_pd(r21,ewtabscale);
1023             ewitab           = _mm_cvttpd_epi32(ewrt);
1024 #ifdef __XOP__
1025             eweps            = _mm_frcz_pd(ewrt);
1026 #else
1027             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1028 #endif
1029             twoeweps         = _mm_add_pd(eweps,eweps);
1030             ewitab           = _mm_slli_epi32(ewitab,2);
1031             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
1032             ewtabD           = _mm_setzero_pd();
1033             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
1034             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
1035             ewtabFn          = _mm_setzero_pd();
1036             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
1037             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
1038             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
1039             velec            = _mm_mul_pd(qq21,_mm_sub_pd(rinv21,velec));
1040             felec            = _mm_mul_pd(_mm_mul_pd(qq21,rinv21),_mm_sub_pd(rinvsq21,felec));
1041
1042             /* Update potential sum for this i atom from the interaction with this j atom. */
1043             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
1044             velecsum         = _mm_add_pd(velecsum,velec);
1045
1046             fscal            = felec;
1047
1048             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1049
1050             /* Update vectorial force */
1051             fix2             = _mm_macc_pd(dx21,fscal,fix2);
1052             fiy2             = _mm_macc_pd(dy21,fscal,fiy2);
1053             fiz2             = _mm_macc_pd(dz21,fscal,fiz2);
1054             
1055             fjx1             = _mm_macc_pd(dx21,fscal,fjx1);
1056             fjy1             = _mm_macc_pd(dy21,fscal,fjy1);
1057             fjz1             = _mm_macc_pd(dz21,fscal,fjz1);
1058
1059             /**************************
1060              * CALCULATE INTERACTIONS *
1061              **************************/
1062
1063             r22              = _mm_mul_pd(rsq22,rinv22);
1064
1065             /* EWALD ELECTROSTATICS */
1066
1067             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1068             ewrt             = _mm_mul_pd(r22,ewtabscale);
1069             ewitab           = _mm_cvttpd_epi32(ewrt);
1070 #ifdef __XOP__
1071             eweps            = _mm_frcz_pd(ewrt);
1072 #else
1073             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1074 #endif
1075             twoeweps         = _mm_add_pd(eweps,eweps);
1076             ewitab           = _mm_slli_epi32(ewitab,2);
1077             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
1078             ewtabD           = _mm_setzero_pd();
1079             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
1080             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
1081             ewtabFn          = _mm_setzero_pd();
1082             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
1083             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
1084             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
1085             velec            = _mm_mul_pd(qq22,_mm_sub_pd(rinv22,velec));
1086             felec            = _mm_mul_pd(_mm_mul_pd(qq22,rinv22),_mm_sub_pd(rinvsq22,felec));
1087
1088             /* Update potential sum for this i atom from the interaction with this j atom. */
1089             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
1090             velecsum         = _mm_add_pd(velecsum,velec);
1091
1092             fscal            = felec;
1093
1094             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1095
1096             /* Update vectorial force */
1097             fix2             = _mm_macc_pd(dx22,fscal,fix2);
1098             fiy2             = _mm_macc_pd(dy22,fscal,fiy2);
1099             fiz2             = _mm_macc_pd(dz22,fscal,fiz2);
1100             
1101             fjx2             = _mm_macc_pd(dx22,fscal,fjx2);
1102             fjy2             = _mm_macc_pd(dy22,fscal,fjy2);
1103             fjz2             = _mm_macc_pd(dz22,fscal,fjz2);
1104
1105             /**************************
1106              * CALCULATE INTERACTIONS *
1107              **************************/
1108
1109             r23              = _mm_mul_pd(rsq23,rinv23);
1110
1111             /* EWALD ELECTROSTATICS */
1112
1113             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1114             ewrt             = _mm_mul_pd(r23,ewtabscale);
1115             ewitab           = _mm_cvttpd_epi32(ewrt);
1116 #ifdef __XOP__
1117             eweps            = _mm_frcz_pd(ewrt);
1118 #else
1119             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1120 #endif
1121             twoeweps         = _mm_add_pd(eweps,eweps);
1122             ewitab           = _mm_slli_epi32(ewitab,2);
1123             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
1124             ewtabD           = _mm_setzero_pd();
1125             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
1126             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
1127             ewtabFn          = _mm_setzero_pd();
1128             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
1129             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
1130             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
1131             velec            = _mm_mul_pd(qq23,_mm_sub_pd(rinv23,velec));
1132             felec            = _mm_mul_pd(_mm_mul_pd(qq23,rinv23),_mm_sub_pd(rinvsq23,felec));
1133
1134             /* Update potential sum for this i atom from the interaction with this j atom. */
1135             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
1136             velecsum         = _mm_add_pd(velecsum,velec);
1137
1138             fscal            = felec;
1139
1140             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1141
1142             /* Update vectorial force */
1143             fix2             = _mm_macc_pd(dx23,fscal,fix2);
1144             fiy2             = _mm_macc_pd(dy23,fscal,fiy2);
1145             fiz2             = _mm_macc_pd(dz23,fscal,fiz2);
1146             
1147             fjx3             = _mm_macc_pd(dx23,fscal,fjx3);
1148             fjy3             = _mm_macc_pd(dy23,fscal,fjy3);
1149             fjz3             = _mm_macc_pd(dz23,fscal,fjz3);
1150
1151             /**************************
1152              * CALCULATE INTERACTIONS *
1153              **************************/
1154
1155             r31              = _mm_mul_pd(rsq31,rinv31);
1156
1157             /* EWALD ELECTROSTATICS */
1158
1159             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1160             ewrt             = _mm_mul_pd(r31,ewtabscale);
1161             ewitab           = _mm_cvttpd_epi32(ewrt);
1162 #ifdef __XOP__
1163             eweps            = _mm_frcz_pd(ewrt);
1164 #else
1165             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1166 #endif
1167             twoeweps         = _mm_add_pd(eweps,eweps);
1168             ewitab           = _mm_slli_epi32(ewitab,2);
1169             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
1170             ewtabD           = _mm_setzero_pd();
1171             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
1172             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
1173             ewtabFn          = _mm_setzero_pd();
1174             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
1175             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
1176             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
1177             velec            = _mm_mul_pd(qq31,_mm_sub_pd(rinv31,velec));
1178             felec            = _mm_mul_pd(_mm_mul_pd(qq31,rinv31),_mm_sub_pd(rinvsq31,felec));
1179
1180             /* Update potential sum for this i atom from the interaction with this j atom. */
1181             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
1182             velecsum         = _mm_add_pd(velecsum,velec);
1183
1184             fscal            = felec;
1185
1186             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1187
1188             /* Update vectorial force */
1189             fix3             = _mm_macc_pd(dx31,fscal,fix3);
1190             fiy3             = _mm_macc_pd(dy31,fscal,fiy3);
1191             fiz3             = _mm_macc_pd(dz31,fscal,fiz3);
1192             
1193             fjx1             = _mm_macc_pd(dx31,fscal,fjx1);
1194             fjy1             = _mm_macc_pd(dy31,fscal,fjy1);
1195             fjz1             = _mm_macc_pd(dz31,fscal,fjz1);
1196
1197             /**************************
1198              * CALCULATE INTERACTIONS *
1199              **************************/
1200
1201             r32              = _mm_mul_pd(rsq32,rinv32);
1202
1203             /* EWALD ELECTROSTATICS */
1204
1205             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1206             ewrt             = _mm_mul_pd(r32,ewtabscale);
1207             ewitab           = _mm_cvttpd_epi32(ewrt);
1208 #ifdef __XOP__
1209             eweps            = _mm_frcz_pd(ewrt);
1210 #else
1211             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1212 #endif
1213             twoeweps         = _mm_add_pd(eweps,eweps);
1214             ewitab           = _mm_slli_epi32(ewitab,2);
1215             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
1216             ewtabD           = _mm_setzero_pd();
1217             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
1218             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
1219             ewtabFn          = _mm_setzero_pd();
1220             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
1221             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
1222             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
1223             velec            = _mm_mul_pd(qq32,_mm_sub_pd(rinv32,velec));
1224             felec            = _mm_mul_pd(_mm_mul_pd(qq32,rinv32),_mm_sub_pd(rinvsq32,felec));
1225
1226             /* Update potential sum for this i atom from the interaction with this j atom. */
1227             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
1228             velecsum         = _mm_add_pd(velecsum,velec);
1229
1230             fscal            = felec;
1231
1232             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1233
1234             /* Update vectorial force */
1235             fix3             = _mm_macc_pd(dx32,fscal,fix3);
1236             fiy3             = _mm_macc_pd(dy32,fscal,fiy3);
1237             fiz3             = _mm_macc_pd(dz32,fscal,fiz3);
1238             
1239             fjx2             = _mm_macc_pd(dx32,fscal,fjx2);
1240             fjy2             = _mm_macc_pd(dy32,fscal,fjy2);
1241             fjz2             = _mm_macc_pd(dz32,fscal,fjz2);
1242
1243             /**************************
1244              * CALCULATE INTERACTIONS *
1245              **************************/
1246
1247             r33              = _mm_mul_pd(rsq33,rinv33);
1248
1249             /* EWALD ELECTROSTATICS */
1250
1251             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1252             ewrt             = _mm_mul_pd(r33,ewtabscale);
1253             ewitab           = _mm_cvttpd_epi32(ewrt);
1254 #ifdef __XOP__
1255             eweps            = _mm_frcz_pd(ewrt);
1256 #else
1257             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1258 #endif
1259             twoeweps         = _mm_add_pd(eweps,eweps);
1260             ewitab           = _mm_slli_epi32(ewitab,2);
1261             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
1262             ewtabD           = _mm_setzero_pd();
1263             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
1264             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
1265             ewtabFn          = _mm_setzero_pd();
1266             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
1267             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
1268             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
1269             velec            = _mm_mul_pd(qq33,_mm_sub_pd(rinv33,velec));
1270             felec            = _mm_mul_pd(_mm_mul_pd(qq33,rinv33),_mm_sub_pd(rinvsq33,felec));
1271
1272             /* Update potential sum for this i atom from the interaction with this j atom. */
1273             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
1274             velecsum         = _mm_add_pd(velecsum,velec);
1275
1276             fscal            = felec;
1277
1278             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
1279
1280             /* Update vectorial force */
1281             fix3             = _mm_macc_pd(dx33,fscal,fix3);
1282             fiy3             = _mm_macc_pd(dy33,fscal,fiy3);
1283             fiz3             = _mm_macc_pd(dz33,fscal,fiz3);
1284             
1285             fjx3             = _mm_macc_pd(dx33,fscal,fjx3);
1286             fjy3             = _mm_macc_pd(dy33,fscal,fjy3);
1287             fjz3             = _mm_macc_pd(dz33,fscal,fjz3);
1288
1289             gmx_mm_decrement_4rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1290
1291             /* Inner loop uses 449 flops */
1292         }
1293
1294         /* End of innermost loop */
1295
1296         gmx_mm_update_iforce_4atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1297                                               f+i_coord_offset,fshift+i_shift_offset);
1298
1299         ggid                        = gid[iidx];
1300         /* Update potential energies */
1301         gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
1302         gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
1303
1304         /* Increment number of inner iterations */
1305         inneriter                  += j_index_end - j_index_start;
1306
1307         /* Outer loop uses 26 flops */
1308     }
1309
1310     /* Increment number of outer iterations */
1311     outeriter        += nri;
1312
1313     /* Update outer/inner flops */
1314
1315     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4W4_VF,outeriter*26 + inneriter*449);
1316 }
1317 /*
1318  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwLJEw_GeomW4W4_F_avx_128_fma_double
1319  * Electrostatics interaction: Ewald
1320  * VdW interaction:            LJEwald
1321  * Geometry:                   Water4-Water4
1322  * Calculate force/pot:        Force
1323  */
1324 void
1325 nb_kernel_ElecEw_VdwLJEw_GeomW4W4_F_avx_128_fma_double
1326                     (t_nblist                    * gmx_restrict       nlist,
1327                      rvec                        * gmx_restrict          xx,
1328                      rvec                        * gmx_restrict          ff,
1329                      t_forcerec                  * gmx_restrict          fr,
1330                      t_mdatoms                   * gmx_restrict     mdatoms,
1331                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
1332                      t_nrnb                      * gmx_restrict        nrnb)
1333 {
1334     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
1335      * just 0 for non-waters.
1336      * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
1337      * jnr indices corresponding to data put in the four positions in the SIMD register.
1338      */
1339     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
1340     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
1341     int              jnrA,jnrB;
1342     int              j_coord_offsetA,j_coord_offsetB;
1343     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
1344     real             rcutoff_scalar;
1345     real             *shiftvec,*fshift,*x,*f;
1346     __m128d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
1347     int              vdwioffset0;
1348     __m128d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
1349     int              vdwioffset1;
1350     __m128d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
1351     int              vdwioffset2;
1352     __m128d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
1353     int              vdwioffset3;
1354     __m128d          ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
1355     int              vdwjidx0A,vdwjidx0B;
1356     __m128d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
1357     int              vdwjidx1A,vdwjidx1B;
1358     __m128d          jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
1359     int              vdwjidx2A,vdwjidx2B;
1360     __m128d          jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
1361     int              vdwjidx3A,vdwjidx3B;
1362     __m128d          jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
1363     __m128d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
1364     __m128d          dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
1365     __m128d          dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
1366     __m128d          dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
1367     __m128d          dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
1368     __m128d          dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
1369     __m128d          dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
1370     __m128d          dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
1371     __m128d          dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
1372     __m128d          dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
1373     __m128d          velec,felec,velecsum,facel,crf,krf,krf2;
1374     real             *charge;
1375     int              nvdwtype;
1376     __m128d          rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
1377     int              *vdwtype;
1378     real             *vdwparam;
1379     __m128d          one_sixth   = _mm_set1_pd(1.0/6.0);
1380     __m128d          one_twelfth = _mm_set1_pd(1.0/12.0);
1381     __m128d           c6grid_00;
1382     __m128d           c6grid_11;
1383     __m128d           c6grid_12;
1384     __m128d           c6grid_13;
1385     __m128d           c6grid_21;
1386     __m128d           c6grid_22;
1387     __m128d           c6grid_23;
1388     __m128d           c6grid_31;
1389     __m128d           c6grid_32;
1390     __m128d           c6grid_33;
1391     real             *vdwgridparam;
1392     __m128d           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
1393     __m128d           one_half  = _mm_set1_pd(0.5);
1394     __m128d           minus_one = _mm_set1_pd(-1.0);
1395     __m128i          ewitab;
1396     __m128d          ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
1397     real             *ewtab;
1398     __m128d          dummy_mask,cutoff_mask;
1399     __m128d          signbit   = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
1400     __m128d          one     = _mm_set1_pd(1.0);
1401     __m128d          two     = _mm_set1_pd(2.0);
1402     x                = xx[0];
1403     f                = ff[0];
1404
1405     nri              = nlist->nri;
1406     iinr             = nlist->iinr;
1407     jindex           = nlist->jindex;
1408     jjnr             = nlist->jjnr;
1409     shiftidx         = nlist->shift;
1410     gid              = nlist->gid;
1411     shiftvec         = fr->shift_vec[0];
1412     fshift           = fr->fshift[0];
1413     facel            = _mm_set1_pd(fr->epsfac);
1414     charge           = mdatoms->chargeA;
1415     nvdwtype         = fr->ntype;
1416     vdwparam         = fr->nbfp;
1417     vdwtype          = mdatoms->typeA;
1418     vdwgridparam     = fr->ljpme_c6grid;
1419     sh_lj_ewald      = _mm_set1_pd(fr->ic->sh_lj_ewald);
1420     ewclj            = _mm_set1_pd(fr->ewaldcoeff_lj);
1421     ewclj2           = _mm_mul_pd(minus_one,_mm_mul_pd(ewclj,ewclj));
1422
1423     sh_ewald         = _mm_set1_pd(fr->ic->sh_ewald);
1424     ewtab            = fr->ic->tabq_coul_F;
1425     ewtabscale       = _mm_set1_pd(fr->ic->tabq_scale);
1426     ewtabhalfspace   = _mm_set1_pd(0.5/fr->ic->tabq_scale);
1427
1428     /* Setup water-specific parameters */
1429     inr              = nlist->iinr[0];
1430     iq1              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
1431     iq2              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
1432     iq3              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
1433     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
1434
1435     jq1              = _mm_set1_pd(charge[inr+1]);
1436     jq2              = _mm_set1_pd(charge[inr+2]);
1437     jq3              = _mm_set1_pd(charge[inr+3]);
1438     vdwjidx0A        = 2*vdwtype[inr+0];
1439     c6_00            = _mm_set1_pd(vdwparam[vdwioffset0+vdwjidx0A]);
1440     c12_00           = _mm_set1_pd(vdwparam[vdwioffset0+vdwjidx0A+1]);
1441     c6grid_00        = _mm_set1_pd(vdwgridparam[vdwioffset0+vdwjidx0A]);
1442     qq11             = _mm_mul_pd(iq1,jq1);
1443     qq12             = _mm_mul_pd(iq1,jq2);
1444     qq13             = _mm_mul_pd(iq1,jq3);
1445     qq21             = _mm_mul_pd(iq2,jq1);
1446     qq22             = _mm_mul_pd(iq2,jq2);
1447     qq23             = _mm_mul_pd(iq2,jq3);
1448     qq31             = _mm_mul_pd(iq3,jq1);
1449     qq32             = _mm_mul_pd(iq3,jq2);
1450     qq33             = _mm_mul_pd(iq3,jq3);
1451
1452     /* Avoid stupid compiler warnings */
1453     jnrA = jnrB = 0;
1454     j_coord_offsetA = 0;
1455     j_coord_offsetB = 0;
1456
1457     outeriter        = 0;
1458     inneriter        = 0;
1459
1460     /* Start outer loop over neighborlists */
1461     for(iidx=0; iidx<nri; iidx++)
1462     {
1463         /* Load shift vector for this list */
1464         i_shift_offset   = DIM*shiftidx[iidx];
1465
1466         /* Load limits for loop over neighbors */
1467         j_index_start    = jindex[iidx];
1468         j_index_end      = jindex[iidx+1];
1469
1470         /* Get outer coordinate index */
1471         inr              = iinr[iidx];
1472         i_coord_offset   = DIM*inr;
1473
1474         /* Load i particle coords and add shift vector */
1475         gmx_mm_load_shift_and_4rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
1476                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
1477
1478         fix0             = _mm_setzero_pd();
1479         fiy0             = _mm_setzero_pd();
1480         fiz0             = _mm_setzero_pd();
1481         fix1             = _mm_setzero_pd();
1482         fiy1             = _mm_setzero_pd();
1483         fiz1             = _mm_setzero_pd();
1484         fix2             = _mm_setzero_pd();
1485         fiy2             = _mm_setzero_pd();
1486         fiz2             = _mm_setzero_pd();
1487         fix3             = _mm_setzero_pd();
1488         fiy3             = _mm_setzero_pd();
1489         fiz3             = _mm_setzero_pd();
1490
1491         /* Start inner kernel loop */
1492         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
1493         {
1494
1495             /* Get j neighbor index, and coordinate index */
1496             jnrA             = jjnr[jidx];
1497             jnrB             = jjnr[jidx+1];
1498             j_coord_offsetA  = DIM*jnrA;
1499             j_coord_offsetB  = DIM*jnrB;
1500
1501             /* load j atom coordinates */
1502             gmx_mm_load_4rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
1503                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
1504                                               &jy2,&jz2,&jx3,&jy3,&jz3);
1505
1506             /* Calculate displacement vector */
1507             dx00             = _mm_sub_pd(ix0,jx0);
1508             dy00             = _mm_sub_pd(iy0,jy0);
1509             dz00             = _mm_sub_pd(iz0,jz0);
1510             dx11             = _mm_sub_pd(ix1,jx1);
1511             dy11             = _mm_sub_pd(iy1,jy1);
1512             dz11             = _mm_sub_pd(iz1,jz1);
1513             dx12             = _mm_sub_pd(ix1,jx2);
1514             dy12             = _mm_sub_pd(iy1,jy2);
1515             dz12             = _mm_sub_pd(iz1,jz2);
1516             dx13             = _mm_sub_pd(ix1,jx3);
1517             dy13             = _mm_sub_pd(iy1,jy3);
1518             dz13             = _mm_sub_pd(iz1,jz3);
1519             dx21             = _mm_sub_pd(ix2,jx1);
1520             dy21             = _mm_sub_pd(iy2,jy1);
1521             dz21             = _mm_sub_pd(iz2,jz1);
1522             dx22             = _mm_sub_pd(ix2,jx2);
1523             dy22             = _mm_sub_pd(iy2,jy2);
1524             dz22             = _mm_sub_pd(iz2,jz2);
1525             dx23             = _mm_sub_pd(ix2,jx3);
1526             dy23             = _mm_sub_pd(iy2,jy3);
1527             dz23             = _mm_sub_pd(iz2,jz3);
1528             dx31             = _mm_sub_pd(ix3,jx1);
1529             dy31             = _mm_sub_pd(iy3,jy1);
1530             dz31             = _mm_sub_pd(iz3,jz1);
1531             dx32             = _mm_sub_pd(ix3,jx2);
1532             dy32             = _mm_sub_pd(iy3,jy2);
1533             dz32             = _mm_sub_pd(iz3,jz2);
1534             dx33             = _mm_sub_pd(ix3,jx3);
1535             dy33             = _mm_sub_pd(iy3,jy3);
1536             dz33             = _mm_sub_pd(iz3,jz3);
1537
1538             /* Calculate squared distance and things based on it */
1539             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
1540             rsq11            = gmx_mm_calc_rsq_pd(dx11,dy11,dz11);
1541             rsq12            = gmx_mm_calc_rsq_pd(dx12,dy12,dz12);
1542             rsq13            = gmx_mm_calc_rsq_pd(dx13,dy13,dz13);
1543             rsq21            = gmx_mm_calc_rsq_pd(dx21,dy21,dz21);
1544             rsq22            = gmx_mm_calc_rsq_pd(dx22,dy22,dz22);
1545             rsq23            = gmx_mm_calc_rsq_pd(dx23,dy23,dz23);
1546             rsq31            = gmx_mm_calc_rsq_pd(dx31,dy31,dz31);
1547             rsq32            = gmx_mm_calc_rsq_pd(dx32,dy32,dz32);
1548             rsq33            = gmx_mm_calc_rsq_pd(dx33,dy33,dz33);
1549
1550             rinv00           = gmx_mm_invsqrt_pd(rsq00);
1551             rinv11           = gmx_mm_invsqrt_pd(rsq11);
1552             rinv12           = gmx_mm_invsqrt_pd(rsq12);
1553             rinv13           = gmx_mm_invsqrt_pd(rsq13);
1554             rinv21           = gmx_mm_invsqrt_pd(rsq21);
1555             rinv22           = gmx_mm_invsqrt_pd(rsq22);
1556             rinv23           = gmx_mm_invsqrt_pd(rsq23);
1557             rinv31           = gmx_mm_invsqrt_pd(rsq31);
1558             rinv32           = gmx_mm_invsqrt_pd(rsq32);
1559             rinv33           = gmx_mm_invsqrt_pd(rsq33);
1560
1561             rinvsq00         = _mm_mul_pd(rinv00,rinv00);
1562             rinvsq11         = _mm_mul_pd(rinv11,rinv11);
1563             rinvsq12         = _mm_mul_pd(rinv12,rinv12);
1564             rinvsq13         = _mm_mul_pd(rinv13,rinv13);
1565             rinvsq21         = _mm_mul_pd(rinv21,rinv21);
1566             rinvsq22         = _mm_mul_pd(rinv22,rinv22);
1567             rinvsq23         = _mm_mul_pd(rinv23,rinv23);
1568             rinvsq31         = _mm_mul_pd(rinv31,rinv31);
1569             rinvsq32         = _mm_mul_pd(rinv32,rinv32);
1570             rinvsq33         = _mm_mul_pd(rinv33,rinv33);
1571
1572             fjx0             = _mm_setzero_pd();
1573             fjy0             = _mm_setzero_pd();
1574             fjz0             = _mm_setzero_pd();
1575             fjx1             = _mm_setzero_pd();
1576             fjy1             = _mm_setzero_pd();
1577             fjz1             = _mm_setzero_pd();
1578             fjx2             = _mm_setzero_pd();
1579             fjy2             = _mm_setzero_pd();
1580             fjz2             = _mm_setzero_pd();
1581             fjx3             = _mm_setzero_pd();
1582             fjy3             = _mm_setzero_pd();
1583             fjz3             = _mm_setzero_pd();
1584
1585             /**************************
1586              * CALCULATE INTERACTIONS *
1587              **************************/
1588
1589             r00              = _mm_mul_pd(rsq00,rinv00);
1590
1591             /* Analytical LJ-PME */
1592             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
1593             ewcljrsq         = _mm_mul_pd(ewclj2,rsq00);
1594             ewclj6           = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
1595             exponent         = gmx_simd_exp_d(ewcljrsq);
1596             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
1597             poly             = _mm_mul_pd(exponent,_mm_macc_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half,_mm_sub_pd(one,ewcljrsq)));
1598             /* f6A = 6 * C6grid * (1 - poly) */
1599             f6A              = _mm_mul_pd(c6grid_00,_mm_sub_pd(one,poly));
1600             /* f6B = C6grid * exponent * beta^6 */
1601             f6B              = _mm_mul_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6));
1602             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
1603             fvdw              = _mm_mul_pd(_mm_macc_pd(_mm_msub_pd(c12_00,rinvsix,_mm_sub_pd(c6_00,f6A)),rinvsix,f6B),rinvsq00);
1604
1605             fscal            = fvdw;
1606
1607             /* Update vectorial force */
1608             fix0             = _mm_macc_pd(dx00,fscal,fix0);
1609             fiy0             = _mm_macc_pd(dy00,fscal,fiy0);
1610             fiz0             = _mm_macc_pd(dz00,fscal,fiz0);
1611             
1612             fjx0             = _mm_macc_pd(dx00,fscal,fjx0);
1613             fjy0             = _mm_macc_pd(dy00,fscal,fjy0);
1614             fjz0             = _mm_macc_pd(dz00,fscal,fjz0);
1615
1616             /**************************
1617              * CALCULATE INTERACTIONS *
1618              **************************/
1619
1620             r11              = _mm_mul_pd(rsq11,rinv11);
1621
1622             /* EWALD ELECTROSTATICS */
1623
1624             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1625             ewrt             = _mm_mul_pd(r11,ewtabscale);
1626             ewitab           = _mm_cvttpd_epi32(ewrt);
1627 #ifdef __XOP__
1628             eweps            = _mm_frcz_pd(ewrt);
1629 #else
1630             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1631 #endif
1632             twoeweps         = _mm_add_pd(eweps,eweps);
1633             gmx_mm_load_2pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),ewtab+_mm_extract_epi32(ewitab,1),
1634                                          &ewtabF,&ewtabFn);
1635             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
1636             felec            = _mm_mul_pd(_mm_mul_pd(qq11,rinv11),_mm_sub_pd(rinvsq11,felec));
1637
1638             fscal            = felec;
1639
1640             /* Update vectorial force */
1641             fix1             = _mm_macc_pd(dx11,fscal,fix1);
1642             fiy1             = _mm_macc_pd(dy11,fscal,fiy1);
1643             fiz1             = _mm_macc_pd(dz11,fscal,fiz1);
1644             
1645             fjx1             = _mm_macc_pd(dx11,fscal,fjx1);
1646             fjy1             = _mm_macc_pd(dy11,fscal,fjy1);
1647             fjz1             = _mm_macc_pd(dz11,fscal,fjz1);
1648
1649             /**************************
1650              * CALCULATE INTERACTIONS *
1651              **************************/
1652
1653             r12              = _mm_mul_pd(rsq12,rinv12);
1654
1655             /* EWALD ELECTROSTATICS */
1656
1657             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1658             ewrt             = _mm_mul_pd(r12,ewtabscale);
1659             ewitab           = _mm_cvttpd_epi32(ewrt);
1660 #ifdef __XOP__
1661             eweps            = _mm_frcz_pd(ewrt);
1662 #else
1663             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1664 #endif
1665             twoeweps         = _mm_add_pd(eweps,eweps);
1666             gmx_mm_load_2pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),ewtab+_mm_extract_epi32(ewitab,1),
1667                                          &ewtabF,&ewtabFn);
1668             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
1669             felec            = _mm_mul_pd(_mm_mul_pd(qq12,rinv12),_mm_sub_pd(rinvsq12,felec));
1670
1671             fscal            = felec;
1672
1673             /* Update vectorial force */
1674             fix1             = _mm_macc_pd(dx12,fscal,fix1);
1675             fiy1             = _mm_macc_pd(dy12,fscal,fiy1);
1676             fiz1             = _mm_macc_pd(dz12,fscal,fiz1);
1677             
1678             fjx2             = _mm_macc_pd(dx12,fscal,fjx2);
1679             fjy2             = _mm_macc_pd(dy12,fscal,fjy2);
1680             fjz2             = _mm_macc_pd(dz12,fscal,fjz2);
1681
1682             /**************************
1683              * CALCULATE INTERACTIONS *
1684              **************************/
1685
1686             r13              = _mm_mul_pd(rsq13,rinv13);
1687
1688             /* EWALD ELECTROSTATICS */
1689
1690             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1691             ewrt             = _mm_mul_pd(r13,ewtabscale);
1692             ewitab           = _mm_cvttpd_epi32(ewrt);
1693 #ifdef __XOP__
1694             eweps            = _mm_frcz_pd(ewrt);
1695 #else
1696             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1697 #endif
1698             twoeweps         = _mm_add_pd(eweps,eweps);
1699             gmx_mm_load_2pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),ewtab+_mm_extract_epi32(ewitab,1),
1700                                          &ewtabF,&ewtabFn);
1701             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
1702             felec            = _mm_mul_pd(_mm_mul_pd(qq13,rinv13),_mm_sub_pd(rinvsq13,felec));
1703
1704             fscal            = felec;
1705
1706             /* Update vectorial force */
1707             fix1             = _mm_macc_pd(dx13,fscal,fix1);
1708             fiy1             = _mm_macc_pd(dy13,fscal,fiy1);
1709             fiz1             = _mm_macc_pd(dz13,fscal,fiz1);
1710             
1711             fjx3             = _mm_macc_pd(dx13,fscal,fjx3);
1712             fjy3             = _mm_macc_pd(dy13,fscal,fjy3);
1713             fjz3             = _mm_macc_pd(dz13,fscal,fjz3);
1714
1715             /**************************
1716              * CALCULATE INTERACTIONS *
1717              **************************/
1718
1719             r21              = _mm_mul_pd(rsq21,rinv21);
1720
1721             /* EWALD ELECTROSTATICS */
1722
1723             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1724             ewrt             = _mm_mul_pd(r21,ewtabscale);
1725             ewitab           = _mm_cvttpd_epi32(ewrt);
1726 #ifdef __XOP__
1727             eweps            = _mm_frcz_pd(ewrt);
1728 #else
1729             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1730 #endif
1731             twoeweps         = _mm_add_pd(eweps,eweps);
1732             gmx_mm_load_2pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),ewtab+_mm_extract_epi32(ewitab,1),
1733                                          &ewtabF,&ewtabFn);
1734             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
1735             felec            = _mm_mul_pd(_mm_mul_pd(qq21,rinv21),_mm_sub_pd(rinvsq21,felec));
1736
1737             fscal            = felec;
1738
1739             /* Update vectorial force */
1740             fix2             = _mm_macc_pd(dx21,fscal,fix2);
1741             fiy2             = _mm_macc_pd(dy21,fscal,fiy2);
1742             fiz2             = _mm_macc_pd(dz21,fscal,fiz2);
1743             
1744             fjx1             = _mm_macc_pd(dx21,fscal,fjx1);
1745             fjy1             = _mm_macc_pd(dy21,fscal,fjy1);
1746             fjz1             = _mm_macc_pd(dz21,fscal,fjz1);
1747
1748             /**************************
1749              * CALCULATE INTERACTIONS *
1750              **************************/
1751
1752             r22              = _mm_mul_pd(rsq22,rinv22);
1753
1754             /* EWALD ELECTROSTATICS */
1755
1756             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1757             ewrt             = _mm_mul_pd(r22,ewtabscale);
1758             ewitab           = _mm_cvttpd_epi32(ewrt);
1759 #ifdef __XOP__
1760             eweps            = _mm_frcz_pd(ewrt);
1761 #else
1762             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1763 #endif
1764             twoeweps         = _mm_add_pd(eweps,eweps);
1765             gmx_mm_load_2pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),ewtab+_mm_extract_epi32(ewitab,1),
1766                                          &ewtabF,&ewtabFn);
1767             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
1768             felec            = _mm_mul_pd(_mm_mul_pd(qq22,rinv22),_mm_sub_pd(rinvsq22,felec));
1769
1770             fscal            = felec;
1771
1772             /* Update vectorial force */
1773             fix2             = _mm_macc_pd(dx22,fscal,fix2);
1774             fiy2             = _mm_macc_pd(dy22,fscal,fiy2);
1775             fiz2             = _mm_macc_pd(dz22,fscal,fiz2);
1776             
1777             fjx2             = _mm_macc_pd(dx22,fscal,fjx2);
1778             fjy2             = _mm_macc_pd(dy22,fscal,fjy2);
1779             fjz2             = _mm_macc_pd(dz22,fscal,fjz2);
1780
1781             /**************************
1782              * CALCULATE INTERACTIONS *
1783              **************************/
1784
1785             r23              = _mm_mul_pd(rsq23,rinv23);
1786
1787             /* EWALD ELECTROSTATICS */
1788
1789             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1790             ewrt             = _mm_mul_pd(r23,ewtabscale);
1791             ewitab           = _mm_cvttpd_epi32(ewrt);
1792 #ifdef __XOP__
1793             eweps            = _mm_frcz_pd(ewrt);
1794 #else
1795             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1796 #endif
1797             twoeweps         = _mm_add_pd(eweps,eweps);
1798             gmx_mm_load_2pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),ewtab+_mm_extract_epi32(ewitab,1),
1799                                          &ewtabF,&ewtabFn);
1800             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
1801             felec            = _mm_mul_pd(_mm_mul_pd(qq23,rinv23),_mm_sub_pd(rinvsq23,felec));
1802
1803             fscal            = felec;
1804
1805             /* Update vectorial force */
1806             fix2             = _mm_macc_pd(dx23,fscal,fix2);
1807             fiy2             = _mm_macc_pd(dy23,fscal,fiy2);
1808             fiz2             = _mm_macc_pd(dz23,fscal,fiz2);
1809             
1810             fjx3             = _mm_macc_pd(dx23,fscal,fjx3);
1811             fjy3             = _mm_macc_pd(dy23,fscal,fjy3);
1812             fjz3             = _mm_macc_pd(dz23,fscal,fjz3);
1813
1814             /**************************
1815              * CALCULATE INTERACTIONS *
1816              **************************/
1817
1818             r31              = _mm_mul_pd(rsq31,rinv31);
1819
1820             /* EWALD ELECTROSTATICS */
1821
1822             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1823             ewrt             = _mm_mul_pd(r31,ewtabscale);
1824             ewitab           = _mm_cvttpd_epi32(ewrt);
1825 #ifdef __XOP__
1826             eweps            = _mm_frcz_pd(ewrt);
1827 #else
1828             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1829 #endif
1830             twoeweps         = _mm_add_pd(eweps,eweps);
1831             gmx_mm_load_2pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),ewtab+_mm_extract_epi32(ewitab,1),
1832                                          &ewtabF,&ewtabFn);
1833             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
1834             felec            = _mm_mul_pd(_mm_mul_pd(qq31,rinv31),_mm_sub_pd(rinvsq31,felec));
1835
1836             fscal            = felec;
1837
1838             /* Update vectorial force */
1839             fix3             = _mm_macc_pd(dx31,fscal,fix3);
1840             fiy3             = _mm_macc_pd(dy31,fscal,fiy3);
1841             fiz3             = _mm_macc_pd(dz31,fscal,fiz3);
1842             
1843             fjx1             = _mm_macc_pd(dx31,fscal,fjx1);
1844             fjy1             = _mm_macc_pd(dy31,fscal,fjy1);
1845             fjz1             = _mm_macc_pd(dz31,fscal,fjz1);
1846
1847             /**************************
1848              * CALCULATE INTERACTIONS *
1849              **************************/
1850
1851             r32              = _mm_mul_pd(rsq32,rinv32);
1852
1853             /* EWALD ELECTROSTATICS */
1854
1855             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1856             ewrt             = _mm_mul_pd(r32,ewtabscale);
1857             ewitab           = _mm_cvttpd_epi32(ewrt);
1858 #ifdef __XOP__
1859             eweps            = _mm_frcz_pd(ewrt);
1860 #else
1861             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1862 #endif
1863             twoeweps         = _mm_add_pd(eweps,eweps);
1864             gmx_mm_load_2pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),ewtab+_mm_extract_epi32(ewitab,1),
1865                                          &ewtabF,&ewtabFn);
1866             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
1867             felec            = _mm_mul_pd(_mm_mul_pd(qq32,rinv32),_mm_sub_pd(rinvsq32,felec));
1868
1869             fscal            = felec;
1870
1871             /* Update vectorial force */
1872             fix3             = _mm_macc_pd(dx32,fscal,fix3);
1873             fiy3             = _mm_macc_pd(dy32,fscal,fiy3);
1874             fiz3             = _mm_macc_pd(dz32,fscal,fiz3);
1875             
1876             fjx2             = _mm_macc_pd(dx32,fscal,fjx2);
1877             fjy2             = _mm_macc_pd(dy32,fscal,fjy2);
1878             fjz2             = _mm_macc_pd(dz32,fscal,fjz2);
1879
1880             /**************************
1881              * CALCULATE INTERACTIONS *
1882              **************************/
1883
1884             r33              = _mm_mul_pd(rsq33,rinv33);
1885
1886             /* EWALD ELECTROSTATICS */
1887
1888             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1889             ewrt             = _mm_mul_pd(r33,ewtabscale);
1890             ewitab           = _mm_cvttpd_epi32(ewrt);
1891 #ifdef __XOP__
1892             eweps            = _mm_frcz_pd(ewrt);
1893 #else
1894             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
1895 #endif
1896             twoeweps         = _mm_add_pd(eweps,eweps);
1897             gmx_mm_load_2pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),ewtab+_mm_extract_epi32(ewitab,1),
1898                                          &ewtabF,&ewtabFn);
1899             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
1900             felec            = _mm_mul_pd(_mm_mul_pd(qq33,rinv33),_mm_sub_pd(rinvsq33,felec));
1901
1902             fscal            = felec;
1903
1904             /* Update vectorial force */
1905             fix3             = _mm_macc_pd(dx33,fscal,fix3);
1906             fiy3             = _mm_macc_pd(dy33,fscal,fiy3);
1907             fiz3             = _mm_macc_pd(dz33,fscal,fiz3);
1908             
1909             fjx3             = _mm_macc_pd(dx33,fscal,fjx3);
1910             fjy3             = _mm_macc_pd(dy33,fscal,fjy3);
1911             fjz3             = _mm_macc_pd(dz33,fscal,fjz3);
1912
1913             gmx_mm_decrement_4rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1914
1915             /* Inner loop uses 401 flops */
1916         }
1917
1918         if(jidx<j_index_end)
1919         {
1920
1921             jnrA             = jjnr[jidx];
1922             j_coord_offsetA  = DIM*jnrA;
1923
1924             /* load j atom coordinates */
1925             gmx_mm_load_4rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
1926                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
1927                                               &jy2,&jz2,&jx3,&jy3,&jz3);
1928
1929             /* Calculate displacement vector */
1930             dx00             = _mm_sub_pd(ix0,jx0);
1931             dy00             = _mm_sub_pd(iy0,jy0);
1932             dz00             = _mm_sub_pd(iz0,jz0);
1933             dx11             = _mm_sub_pd(ix1,jx1);
1934             dy11             = _mm_sub_pd(iy1,jy1);
1935             dz11             = _mm_sub_pd(iz1,jz1);
1936             dx12             = _mm_sub_pd(ix1,jx2);
1937             dy12             = _mm_sub_pd(iy1,jy2);
1938             dz12             = _mm_sub_pd(iz1,jz2);
1939             dx13             = _mm_sub_pd(ix1,jx3);
1940             dy13             = _mm_sub_pd(iy1,jy3);
1941             dz13             = _mm_sub_pd(iz1,jz3);
1942             dx21             = _mm_sub_pd(ix2,jx1);
1943             dy21             = _mm_sub_pd(iy2,jy1);
1944             dz21             = _mm_sub_pd(iz2,jz1);
1945             dx22             = _mm_sub_pd(ix2,jx2);
1946             dy22             = _mm_sub_pd(iy2,jy2);
1947             dz22             = _mm_sub_pd(iz2,jz2);
1948             dx23             = _mm_sub_pd(ix2,jx3);
1949             dy23             = _mm_sub_pd(iy2,jy3);
1950             dz23             = _mm_sub_pd(iz2,jz3);
1951             dx31             = _mm_sub_pd(ix3,jx1);
1952             dy31             = _mm_sub_pd(iy3,jy1);
1953             dz31             = _mm_sub_pd(iz3,jz1);
1954             dx32             = _mm_sub_pd(ix3,jx2);
1955             dy32             = _mm_sub_pd(iy3,jy2);
1956             dz32             = _mm_sub_pd(iz3,jz2);
1957             dx33             = _mm_sub_pd(ix3,jx3);
1958             dy33             = _mm_sub_pd(iy3,jy3);
1959             dz33             = _mm_sub_pd(iz3,jz3);
1960
1961             /* Calculate squared distance and things based on it */
1962             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
1963             rsq11            = gmx_mm_calc_rsq_pd(dx11,dy11,dz11);
1964             rsq12            = gmx_mm_calc_rsq_pd(dx12,dy12,dz12);
1965             rsq13            = gmx_mm_calc_rsq_pd(dx13,dy13,dz13);
1966             rsq21            = gmx_mm_calc_rsq_pd(dx21,dy21,dz21);
1967             rsq22            = gmx_mm_calc_rsq_pd(dx22,dy22,dz22);
1968             rsq23            = gmx_mm_calc_rsq_pd(dx23,dy23,dz23);
1969             rsq31            = gmx_mm_calc_rsq_pd(dx31,dy31,dz31);
1970             rsq32            = gmx_mm_calc_rsq_pd(dx32,dy32,dz32);
1971             rsq33            = gmx_mm_calc_rsq_pd(dx33,dy33,dz33);
1972
1973             rinv00           = gmx_mm_invsqrt_pd(rsq00);
1974             rinv11           = gmx_mm_invsqrt_pd(rsq11);
1975             rinv12           = gmx_mm_invsqrt_pd(rsq12);
1976             rinv13           = gmx_mm_invsqrt_pd(rsq13);
1977             rinv21           = gmx_mm_invsqrt_pd(rsq21);
1978             rinv22           = gmx_mm_invsqrt_pd(rsq22);
1979             rinv23           = gmx_mm_invsqrt_pd(rsq23);
1980             rinv31           = gmx_mm_invsqrt_pd(rsq31);
1981             rinv32           = gmx_mm_invsqrt_pd(rsq32);
1982             rinv33           = gmx_mm_invsqrt_pd(rsq33);
1983
1984             rinvsq00         = _mm_mul_pd(rinv00,rinv00);
1985             rinvsq11         = _mm_mul_pd(rinv11,rinv11);
1986             rinvsq12         = _mm_mul_pd(rinv12,rinv12);
1987             rinvsq13         = _mm_mul_pd(rinv13,rinv13);
1988             rinvsq21         = _mm_mul_pd(rinv21,rinv21);
1989             rinvsq22         = _mm_mul_pd(rinv22,rinv22);
1990             rinvsq23         = _mm_mul_pd(rinv23,rinv23);
1991             rinvsq31         = _mm_mul_pd(rinv31,rinv31);
1992             rinvsq32         = _mm_mul_pd(rinv32,rinv32);
1993             rinvsq33         = _mm_mul_pd(rinv33,rinv33);
1994
1995             fjx0             = _mm_setzero_pd();
1996             fjy0             = _mm_setzero_pd();
1997             fjz0             = _mm_setzero_pd();
1998             fjx1             = _mm_setzero_pd();
1999             fjy1             = _mm_setzero_pd();
2000             fjz1             = _mm_setzero_pd();
2001             fjx2             = _mm_setzero_pd();
2002             fjy2             = _mm_setzero_pd();
2003             fjz2             = _mm_setzero_pd();
2004             fjx3             = _mm_setzero_pd();
2005             fjy3             = _mm_setzero_pd();
2006             fjz3             = _mm_setzero_pd();
2007
2008             /**************************
2009              * CALCULATE INTERACTIONS *
2010              **************************/
2011
2012             r00              = _mm_mul_pd(rsq00,rinv00);
2013
2014             /* Analytical LJ-PME */
2015             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
2016             ewcljrsq         = _mm_mul_pd(ewclj2,rsq00);
2017             ewclj6           = _mm_mul_pd(ewclj2,_mm_mul_pd(ewclj2,ewclj2));
2018             exponent         = gmx_simd_exp_d(ewcljrsq);
2019             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
2020             poly             = _mm_mul_pd(exponent,_mm_macc_pd(_mm_mul_pd(ewcljrsq,ewcljrsq),one_half,_mm_sub_pd(one,ewcljrsq)));
2021             /* f6A = 6 * C6grid * (1 - poly) */
2022             f6A              = _mm_mul_pd(c6grid_00,_mm_sub_pd(one,poly));
2023             /* f6B = C6grid * exponent * beta^6 */
2024             f6B              = _mm_mul_pd(_mm_mul_pd(c6grid_00,one_sixth),_mm_mul_pd(exponent,ewclj6));
2025             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
2026             fvdw              = _mm_mul_pd(_mm_macc_pd(_mm_msub_pd(c12_00,rinvsix,_mm_sub_pd(c6_00,f6A)),rinvsix,f6B),rinvsq00);
2027
2028             fscal            = fvdw;
2029
2030             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
2031
2032             /* Update vectorial force */
2033             fix0             = _mm_macc_pd(dx00,fscal,fix0);
2034             fiy0             = _mm_macc_pd(dy00,fscal,fiy0);
2035             fiz0             = _mm_macc_pd(dz00,fscal,fiz0);
2036             
2037             fjx0             = _mm_macc_pd(dx00,fscal,fjx0);
2038             fjy0             = _mm_macc_pd(dy00,fscal,fjy0);
2039             fjz0             = _mm_macc_pd(dz00,fscal,fjz0);
2040
2041             /**************************
2042              * CALCULATE INTERACTIONS *
2043              **************************/
2044
2045             r11              = _mm_mul_pd(rsq11,rinv11);
2046
2047             /* EWALD ELECTROSTATICS */
2048
2049             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2050             ewrt             = _mm_mul_pd(r11,ewtabscale);
2051             ewitab           = _mm_cvttpd_epi32(ewrt);
2052 #ifdef __XOP__
2053             eweps            = _mm_frcz_pd(ewrt);
2054 #else
2055             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
2056 #endif
2057             twoeweps         = _mm_add_pd(eweps,eweps);
2058             gmx_mm_load_1pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
2059             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
2060             felec            = _mm_mul_pd(_mm_mul_pd(qq11,rinv11),_mm_sub_pd(rinvsq11,felec));
2061
2062             fscal            = felec;
2063
2064             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
2065
2066             /* Update vectorial force */
2067             fix1             = _mm_macc_pd(dx11,fscal,fix1);
2068             fiy1             = _mm_macc_pd(dy11,fscal,fiy1);
2069             fiz1             = _mm_macc_pd(dz11,fscal,fiz1);
2070             
2071             fjx1             = _mm_macc_pd(dx11,fscal,fjx1);
2072             fjy1             = _mm_macc_pd(dy11,fscal,fjy1);
2073             fjz1             = _mm_macc_pd(dz11,fscal,fjz1);
2074
2075             /**************************
2076              * CALCULATE INTERACTIONS *
2077              **************************/
2078
2079             r12              = _mm_mul_pd(rsq12,rinv12);
2080
2081             /* EWALD ELECTROSTATICS */
2082
2083             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2084             ewrt             = _mm_mul_pd(r12,ewtabscale);
2085             ewitab           = _mm_cvttpd_epi32(ewrt);
2086 #ifdef __XOP__
2087             eweps            = _mm_frcz_pd(ewrt);
2088 #else
2089             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
2090 #endif
2091             twoeweps         = _mm_add_pd(eweps,eweps);
2092             gmx_mm_load_1pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
2093             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
2094             felec            = _mm_mul_pd(_mm_mul_pd(qq12,rinv12),_mm_sub_pd(rinvsq12,felec));
2095
2096             fscal            = felec;
2097
2098             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
2099
2100             /* Update vectorial force */
2101             fix1             = _mm_macc_pd(dx12,fscal,fix1);
2102             fiy1             = _mm_macc_pd(dy12,fscal,fiy1);
2103             fiz1             = _mm_macc_pd(dz12,fscal,fiz1);
2104             
2105             fjx2             = _mm_macc_pd(dx12,fscal,fjx2);
2106             fjy2             = _mm_macc_pd(dy12,fscal,fjy2);
2107             fjz2             = _mm_macc_pd(dz12,fscal,fjz2);
2108
2109             /**************************
2110              * CALCULATE INTERACTIONS *
2111              **************************/
2112
2113             r13              = _mm_mul_pd(rsq13,rinv13);
2114
2115             /* EWALD ELECTROSTATICS */
2116
2117             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2118             ewrt             = _mm_mul_pd(r13,ewtabscale);
2119             ewitab           = _mm_cvttpd_epi32(ewrt);
2120 #ifdef __XOP__
2121             eweps            = _mm_frcz_pd(ewrt);
2122 #else
2123             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
2124 #endif
2125             twoeweps         = _mm_add_pd(eweps,eweps);
2126             gmx_mm_load_1pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
2127             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
2128             felec            = _mm_mul_pd(_mm_mul_pd(qq13,rinv13),_mm_sub_pd(rinvsq13,felec));
2129
2130             fscal            = felec;
2131
2132             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
2133
2134             /* Update vectorial force */
2135             fix1             = _mm_macc_pd(dx13,fscal,fix1);
2136             fiy1             = _mm_macc_pd(dy13,fscal,fiy1);
2137             fiz1             = _mm_macc_pd(dz13,fscal,fiz1);
2138             
2139             fjx3             = _mm_macc_pd(dx13,fscal,fjx3);
2140             fjy3             = _mm_macc_pd(dy13,fscal,fjy3);
2141             fjz3             = _mm_macc_pd(dz13,fscal,fjz3);
2142
2143             /**************************
2144              * CALCULATE INTERACTIONS *
2145              **************************/
2146
2147             r21              = _mm_mul_pd(rsq21,rinv21);
2148
2149             /* EWALD ELECTROSTATICS */
2150
2151             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2152             ewrt             = _mm_mul_pd(r21,ewtabscale);
2153             ewitab           = _mm_cvttpd_epi32(ewrt);
2154 #ifdef __XOP__
2155             eweps            = _mm_frcz_pd(ewrt);
2156 #else
2157             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
2158 #endif
2159             twoeweps         = _mm_add_pd(eweps,eweps);
2160             gmx_mm_load_1pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
2161             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
2162             felec            = _mm_mul_pd(_mm_mul_pd(qq21,rinv21),_mm_sub_pd(rinvsq21,felec));
2163
2164             fscal            = felec;
2165
2166             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
2167
2168             /* Update vectorial force */
2169             fix2             = _mm_macc_pd(dx21,fscal,fix2);
2170             fiy2             = _mm_macc_pd(dy21,fscal,fiy2);
2171             fiz2             = _mm_macc_pd(dz21,fscal,fiz2);
2172             
2173             fjx1             = _mm_macc_pd(dx21,fscal,fjx1);
2174             fjy1             = _mm_macc_pd(dy21,fscal,fjy1);
2175             fjz1             = _mm_macc_pd(dz21,fscal,fjz1);
2176
2177             /**************************
2178              * CALCULATE INTERACTIONS *
2179              **************************/
2180
2181             r22              = _mm_mul_pd(rsq22,rinv22);
2182
2183             /* EWALD ELECTROSTATICS */
2184
2185             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2186             ewrt             = _mm_mul_pd(r22,ewtabscale);
2187             ewitab           = _mm_cvttpd_epi32(ewrt);
2188 #ifdef __XOP__
2189             eweps            = _mm_frcz_pd(ewrt);
2190 #else
2191             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
2192 #endif
2193             twoeweps         = _mm_add_pd(eweps,eweps);
2194             gmx_mm_load_1pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
2195             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
2196             felec            = _mm_mul_pd(_mm_mul_pd(qq22,rinv22),_mm_sub_pd(rinvsq22,felec));
2197
2198             fscal            = felec;
2199
2200             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
2201
2202             /* Update vectorial force */
2203             fix2             = _mm_macc_pd(dx22,fscal,fix2);
2204             fiy2             = _mm_macc_pd(dy22,fscal,fiy2);
2205             fiz2             = _mm_macc_pd(dz22,fscal,fiz2);
2206             
2207             fjx2             = _mm_macc_pd(dx22,fscal,fjx2);
2208             fjy2             = _mm_macc_pd(dy22,fscal,fjy2);
2209             fjz2             = _mm_macc_pd(dz22,fscal,fjz2);
2210
2211             /**************************
2212              * CALCULATE INTERACTIONS *
2213              **************************/
2214
2215             r23              = _mm_mul_pd(rsq23,rinv23);
2216
2217             /* EWALD ELECTROSTATICS */
2218
2219             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2220             ewrt             = _mm_mul_pd(r23,ewtabscale);
2221             ewitab           = _mm_cvttpd_epi32(ewrt);
2222 #ifdef __XOP__
2223             eweps            = _mm_frcz_pd(ewrt);
2224 #else
2225             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
2226 #endif
2227             twoeweps         = _mm_add_pd(eweps,eweps);
2228             gmx_mm_load_1pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
2229             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
2230             felec            = _mm_mul_pd(_mm_mul_pd(qq23,rinv23),_mm_sub_pd(rinvsq23,felec));
2231
2232             fscal            = felec;
2233
2234             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
2235
2236             /* Update vectorial force */
2237             fix2             = _mm_macc_pd(dx23,fscal,fix2);
2238             fiy2             = _mm_macc_pd(dy23,fscal,fiy2);
2239             fiz2             = _mm_macc_pd(dz23,fscal,fiz2);
2240             
2241             fjx3             = _mm_macc_pd(dx23,fscal,fjx3);
2242             fjy3             = _mm_macc_pd(dy23,fscal,fjy3);
2243             fjz3             = _mm_macc_pd(dz23,fscal,fjz3);
2244
2245             /**************************
2246              * CALCULATE INTERACTIONS *
2247              **************************/
2248
2249             r31              = _mm_mul_pd(rsq31,rinv31);
2250
2251             /* EWALD ELECTROSTATICS */
2252
2253             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2254             ewrt             = _mm_mul_pd(r31,ewtabscale);
2255             ewitab           = _mm_cvttpd_epi32(ewrt);
2256 #ifdef __XOP__
2257             eweps            = _mm_frcz_pd(ewrt);
2258 #else
2259             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
2260 #endif
2261             twoeweps         = _mm_add_pd(eweps,eweps);
2262             gmx_mm_load_1pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
2263             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
2264             felec            = _mm_mul_pd(_mm_mul_pd(qq31,rinv31),_mm_sub_pd(rinvsq31,felec));
2265
2266             fscal            = felec;
2267
2268             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
2269
2270             /* Update vectorial force */
2271             fix3             = _mm_macc_pd(dx31,fscal,fix3);
2272             fiy3             = _mm_macc_pd(dy31,fscal,fiy3);
2273             fiz3             = _mm_macc_pd(dz31,fscal,fiz3);
2274             
2275             fjx1             = _mm_macc_pd(dx31,fscal,fjx1);
2276             fjy1             = _mm_macc_pd(dy31,fscal,fjy1);
2277             fjz1             = _mm_macc_pd(dz31,fscal,fjz1);
2278
2279             /**************************
2280              * CALCULATE INTERACTIONS *
2281              **************************/
2282
2283             r32              = _mm_mul_pd(rsq32,rinv32);
2284
2285             /* EWALD ELECTROSTATICS */
2286
2287             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2288             ewrt             = _mm_mul_pd(r32,ewtabscale);
2289             ewitab           = _mm_cvttpd_epi32(ewrt);
2290 #ifdef __XOP__
2291             eweps            = _mm_frcz_pd(ewrt);
2292 #else
2293             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
2294 #endif
2295             twoeweps         = _mm_add_pd(eweps,eweps);
2296             gmx_mm_load_1pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
2297             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
2298             felec            = _mm_mul_pd(_mm_mul_pd(qq32,rinv32),_mm_sub_pd(rinvsq32,felec));
2299
2300             fscal            = felec;
2301
2302             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
2303
2304             /* Update vectorial force */
2305             fix3             = _mm_macc_pd(dx32,fscal,fix3);
2306             fiy3             = _mm_macc_pd(dy32,fscal,fiy3);
2307             fiz3             = _mm_macc_pd(dz32,fscal,fiz3);
2308             
2309             fjx2             = _mm_macc_pd(dx32,fscal,fjx2);
2310             fjy2             = _mm_macc_pd(dy32,fscal,fjy2);
2311             fjz2             = _mm_macc_pd(dz32,fscal,fjz2);
2312
2313             /**************************
2314              * CALCULATE INTERACTIONS *
2315              **************************/
2316
2317             r33              = _mm_mul_pd(rsq33,rinv33);
2318
2319             /* EWALD ELECTROSTATICS */
2320
2321             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
2322             ewrt             = _mm_mul_pd(r33,ewtabscale);
2323             ewitab           = _mm_cvttpd_epi32(ewrt);
2324 #ifdef __XOP__
2325             eweps            = _mm_frcz_pd(ewrt);
2326 #else
2327             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
2328 #endif
2329             twoeweps         = _mm_add_pd(eweps,eweps);
2330             gmx_mm_load_1pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
2331             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
2332             felec            = _mm_mul_pd(_mm_mul_pd(qq33,rinv33),_mm_sub_pd(rinvsq33,felec));
2333
2334             fscal            = felec;
2335
2336             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
2337
2338             /* Update vectorial force */
2339             fix3             = _mm_macc_pd(dx33,fscal,fix3);
2340             fiy3             = _mm_macc_pd(dy33,fscal,fiy3);
2341             fiz3             = _mm_macc_pd(dz33,fscal,fiz3);
2342             
2343             fjx3             = _mm_macc_pd(dx33,fscal,fjx3);
2344             fjy3             = _mm_macc_pd(dy33,fscal,fjy3);
2345             fjz3             = _mm_macc_pd(dz33,fscal,fjz3);
2346
2347             gmx_mm_decrement_4rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
2348
2349             /* Inner loop uses 401 flops */
2350         }
2351
2352         /* End of innermost loop */
2353
2354         gmx_mm_update_iforce_4atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
2355                                               f+i_coord_offset,fshift+i_shift_offset);
2356
2357         /* Increment number of inner iterations */
2358         inneriter                  += j_index_end - j_index_start;
2359
2360         /* Outer loop uses 24 flops */
2361     }
2362
2363     /* Increment number of outer iterations */
2364     outeriter        += nri;
2365
2366     /* Update outer/inner flops */
2367
2368     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4W4_F,outeriter*24 + inneriter*401);
2369 }