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