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