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