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