Compile nonbonded kernels as C++
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_single / nb_kernel_ElecEw_VdwCSTab_GeomW4P1_avx_256_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_256_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_256_single.h"
48
49 /*
50  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwCSTab_GeomW4P1_VF_avx_256_single
51  * Electrostatics interaction: Ewald
52  * VdW interaction:            CubicSplineTable
53  * Geometry:                   Water4-Particle
54  * Calculate force/pot:        PotentialAndForce
55  */
56 void
57 nb_kernel_ElecEw_VdwCSTab_GeomW4P1_VF_avx_256_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,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight 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              jnrE,jnrF,jnrG,jnrH;
75     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
76     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
77     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
78     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
79     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
80     real             rcutoff_scalar;
81     real             *shiftvec,*fshift,*x,*f;
82     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
83     real             scratch[4*DIM];
84     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
85     real *           vdwioffsetptr0;
86     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
87     real *           vdwioffsetptr1;
88     __m256           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
89     real *           vdwioffsetptr2;
90     __m256           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
91     real *           vdwioffsetptr3;
92     __m256           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
93     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
94     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
95     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
96     __m256           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
97     __m256           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
98     __m256           dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
99     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
100     real             *charge;
101     int              nvdwtype;
102     __m256           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
103     int              *vdwtype;
104     real             *vdwparam;
105     __m256           one_sixth   = _mm256_set1_ps(1.0/6.0);
106     __m256           one_twelfth = _mm256_set1_ps(1.0/12.0);
107     __m256i          vfitab;
108     __m128i          vfitab_lo,vfitab_hi;
109     __m128i          ifour       = _mm_set1_epi32(4);
110     __m256           rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
111     real             *vftab;
112     __m256i          ewitab;
113     __m128i          ewitab_lo,ewitab_hi;
114     __m256           ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
115     __m256           beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
116     real             *ewtab;
117     __m256           dummy_mask,cutoff_mask;
118     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
119     __m256           one     = _mm256_set1_ps(1.0);
120     __m256           two     = _mm256_set1_ps(2.0);
121     x                = xx[0];
122     f                = ff[0];
123
124     nri              = nlist->nri;
125     iinr             = nlist->iinr;
126     jindex           = nlist->jindex;
127     jjnr             = nlist->jjnr;
128     shiftidx         = nlist->shift;
129     gid              = nlist->gid;
130     shiftvec         = fr->shift_vec[0];
131     fshift           = fr->fshift[0];
132     facel            = _mm256_set1_ps(fr->ic->epsfac);
133     charge           = mdatoms->chargeA;
134     nvdwtype         = fr->ntype;
135     vdwparam         = fr->nbfp;
136     vdwtype          = mdatoms->typeA;
137
138     vftab            = kernel_data->table_vdw->data;
139     vftabscale       = _mm256_set1_ps(kernel_data->table_vdw->scale);
140
141     sh_ewald         = _mm256_set1_ps(fr->ic->sh_ewald);
142     beta             = _mm256_set1_ps(fr->ic->ewaldcoeff_q);
143     beta2            = _mm256_mul_ps(beta,beta);
144     beta3            = _mm256_mul_ps(beta,beta2);
145
146     ewtab            = fr->ic->tabq_coul_FDV0;
147     ewtabscale       = _mm256_set1_ps(fr->ic->tabq_scale);
148     ewtabhalfspace   = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
149
150     /* Setup water-specific parameters */
151     inr              = nlist->iinr[0];
152     iq1              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
153     iq2              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
154     iq3              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+3]));
155     vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
156
157     /* Avoid stupid compiler warnings */
158     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
159     j_coord_offsetA = 0;
160     j_coord_offsetB = 0;
161     j_coord_offsetC = 0;
162     j_coord_offsetD = 0;
163     j_coord_offsetE = 0;
164     j_coord_offsetF = 0;
165     j_coord_offsetG = 0;
166     j_coord_offsetH = 0;
167
168     outeriter        = 0;
169     inneriter        = 0;
170
171     for(iidx=0;iidx<4*DIM;iidx++)
172     {
173         scratch[iidx] = 0.0;
174     }
175
176     /* Start outer loop over neighborlists */
177     for(iidx=0; iidx<nri; iidx++)
178     {
179         /* Load shift vector for this list */
180         i_shift_offset   = DIM*shiftidx[iidx];
181
182         /* Load limits for loop over neighbors */
183         j_index_start    = jindex[iidx];
184         j_index_end      = jindex[iidx+1];
185
186         /* Get outer coordinate index */
187         inr              = iinr[iidx];
188         i_coord_offset   = DIM*inr;
189
190         /* Load i particle coords and add shift vector */
191         gmx_mm256_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
192                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
193
194         fix0             = _mm256_setzero_ps();
195         fiy0             = _mm256_setzero_ps();
196         fiz0             = _mm256_setzero_ps();
197         fix1             = _mm256_setzero_ps();
198         fiy1             = _mm256_setzero_ps();
199         fiz1             = _mm256_setzero_ps();
200         fix2             = _mm256_setzero_ps();
201         fiy2             = _mm256_setzero_ps();
202         fiz2             = _mm256_setzero_ps();
203         fix3             = _mm256_setzero_ps();
204         fiy3             = _mm256_setzero_ps();
205         fiz3             = _mm256_setzero_ps();
206
207         /* Reset potential sums */
208         velecsum         = _mm256_setzero_ps();
209         vvdwsum          = _mm256_setzero_ps();
210
211         /* Start inner kernel loop */
212         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
213         {
214
215             /* Get j neighbor index, and coordinate index */
216             jnrA             = jjnr[jidx];
217             jnrB             = jjnr[jidx+1];
218             jnrC             = jjnr[jidx+2];
219             jnrD             = jjnr[jidx+3];
220             jnrE             = jjnr[jidx+4];
221             jnrF             = jjnr[jidx+5];
222             jnrG             = jjnr[jidx+6];
223             jnrH             = jjnr[jidx+7];
224             j_coord_offsetA  = DIM*jnrA;
225             j_coord_offsetB  = DIM*jnrB;
226             j_coord_offsetC  = DIM*jnrC;
227             j_coord_offsetD  = DIM*jnrD;
228             j_coord_offsetE  = DIM*jnrE;
229             j_coord_offsetF  = DIM*jnrF;
230             j_coord_offsetG  = DIM*jnrG;
231             j_coord_offsetH  = DIM*jnrH;
232
233             /* load j atom coordinates */
234             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
235                                                  x+j_coord_offsetC,x+j_coord_offsetD,
236                                                  x+j_coord_offsetE,x+j_coord_offsetF,
237                                                  x+j_coord_offsetG,x+j_coord_offsetH,
238                                                  &jx0,&jy0,&jz0);
239
240             /* Calculate displacement vector */
241             dx00             = _mm256_sub_ps(ix0,jx0);
242             dy00             = _mm256_sub_ps(iy0,jy0);
243             dz00             = _mm256_sub_ps(iz0,jz0);
244             dx10             = _mm256_sub_ps(ix1,jx0);
245             dy10             = _mm256_sub_ps(iy1,jy0);
246             dz10             = _mm256_sub_ps(iz1,jz0);
247             dx20             = _mm256_sub_ps(ix2,jx0);
248             dy20             = _mm256_sub_ps(iy2,jy0);
249             dz20             = _mm256_sub_ps(iz2,jz0);
250             dx30             = _mm256_sub_ps(ix3,jx0);
251             dy30             = _mm256_sub_ps(iy3,jy0);
252             dz30             = _mm256_sub_ps(iz3,jz0);
253
254             /* Calculate squared distance and things based on it */
255             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
256             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
257             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
258             rsq30            = gmx_mm256_calc_rsq_ps(dx30,dy30,dz30);
259
260             rinv00           = avx256_invsqrt_f(rsq00);
261             rinv10           = avx256_invsqrt_f(rsq10);
262             rinv20           = avx256_invsqrt_f(rsq20);
263             rinv30           = avx256_invsqrt_f(rsq30);
264
265             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
266             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
267             rinvsq30         = _mm256_mul_ps(rinv30,rinv30);
268
269             /* Load parameters for j particles */
270             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
271                                                                  charge+jnrC+0,charge+jnrD+0,
272                                                                  charge+jnrE+0,charge+jnrF+0,
273                                                                  charge+jnrG+0,charge+jnrH+0);
274             vdwjidx0A        = 2*vdwtype[jnrA+0];
275             vdwjidx0B        = 2*vdwtype[jnrB+0];
276             vdwjidx0C        = 2*vdwtype[jnrC+0];
277             vdwjidx0D        = 2*vdwtype[jnrD+0];
278             vdwjidx0E        = 2*vdwtype[jnrE+0];
279             vdwjidx0F        = 2*vdwtype[jnrF+0];
280             vdwjidx0G        = 2*vdwtype[jnrG+0];
281             vdwjidx0H        = 2*vdwtype[jnrH+0];
282
283             fjx0             = _mm256_setzero_ps();
284             fjy0             = _mm256_setzero_ps();
285             fjz0             = _mm256_setzero_ps();
286
287             /**************************
288              * CALCULATE INTERACTIONS *
289              **************************/
290
291             r00              = _mm256_mul_ps(rsq00,rinv00);
292
293             /* Compute parameters for interactions between i and j atoms */
294             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
295                                             vdwioffsetptr0+vdwjidx0B,
296                                             vdwioffsetptr0+vdwjidx0C,
297                                             vdwioffsetptr0+vdwjidx0D,
298                                             vdwioffsetptr0+vdwjidx0E,
299                                             vdwioffsetptr0+vdwjidx0F,
300                                             vdwioffsetptr0+vdwjidx0G,
301                                             vdwioffsetptr0+vdwjidx0H,
302                                             &c6_00,&c12_00);
303
304             /* Calculate table index by multiplying r with table scale and truncate to integer */
305             rt               = _mm256_mul_ps(r00,vftabscale);
306             vfitab           = _mm256_cvttps_epi32(rt);
307             vfeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
308             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
309             vfitab_lo        = _mm256_extractf128_si256(vfitab,0x0);
310             vfitab_hi        = _mm256_extractf128_si256(vfitab,0x1);
311             vfitab_lo        = _mm_slli_epi32(vfitab_lo,3);
312             vfitab_hi        = _mm_slli_epi32(vfitab_hi,3);
313
314             /* CUBIC SPLINE TABLE DISPERSION */
315             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
316                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
317             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
318                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
319             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
320                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
321             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
322                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
323             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
324             Heps             = _mm256_mul_ps(vfeps,H);
325             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
326             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
327             vvdw6            = _mm256_mul_ps(c6_00,VV);
328             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
329             fvdw6            = _mm256_mul_ps(c6_00,FF);
330
331             /* CUBIC SPLINE TABLE REPULSION */
332             vfitab_lo        = _mm_add_epi32(vfitab_lo,ifour);
333             vfitab_hi        = _mm_add_epi32(vfitab_hi,ifour);
334             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
335                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
336             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
337                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
338             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
339                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
340             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
341                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
342             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
343             Heps             = _mm256_mul_ps(vfeps,H);
344             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
345             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
346             vvdw12           = _mm256_mul_ps(c12_00,VV);
347             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
348             fvdw12           = _mm256_mul_ps(c12_00,FF);
349             vvdw             = _mm256_add_ps(vvdw12,vvdw6);
350             fvdw             = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv00)));
351
352             /* Update potential sum for this i atom from the interaction with this j atom. */
353             vvdwsum          = _mm256_add_ps(vvdwsum,vvdw);
354
355             fscal            = fvdw;
356
357             /* Calculate temporary vectorial force */
358             tx               = _mm256_mul_ps(fscal,dx00);
359             ty               = _mm256_mul_ps(fscal,dy00);
360             tz               = _mm256_mul_ps(fscal,dz00);
361
362             /* Update vectorial force */
363             fix0             = _mm256_add_ps(fix0,tx);
364             fiy0             = _mm256_add_ps(fiy0,ty);
365             fiz0             = _mm256_add_ps(fiz0,tz);
366
367             fjx0             = _mm256_add_ps(fjx0,tx);
368             fjy0             = _mm256_add_ps(fjy0,ty);
369             fjz0             = _mm256_add_ps(fjz0,tz);
370
371             /**************************
372              * CALCULATE INTERACTIONS *
373              **************************/
374
375             r10              = _mm256_mul_ps(rsq10,rinv10);
376
377             /* Compute parameters for interactions between i and j atoms */
378             qq10             = _mm256_mul_ps(iq1,jq0);
379
380             /* EWALD ELECTROSTATICS */
381             
382             /* Analytical PME correction */
383             zeta2            = _mm256_mul_ps(beta2,rsq10);
384             rinv3            = _mm256_mul_ps(rinvsq10,rinv10);
385             pmecorrF         = avx256_pmecorrF_f(zeta2);
386             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
387             felec            = _mm256_mul_ps(qq10,felec);
388             pmecorrV         = avx256_pmecorrV_f(zeta2);
389             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
390             velec            = _mm256_sub_ps(rinv10,pmecorrV);
391             velec            = _mm256_mul_ps(qq10,velec);
392             
393             /* Update potential sum for this i atom from the interaction with this j atom. */
394             velecsum         = _mm256_add_ps(velecsum,velec);
395
396             fscal            = felec;
397
398             /* Calculate temporary vectorial force */
399             tx               = _mm256_mul_ps(fscal,dx10);
400             ty               = _mm256_mul_ps(fscal,dy10);
401             tz               = _mm256_mul_ps(fscal,dz10);
402
403             /* Update vectorial force */
404             fix1             = _mm256_add_ps(fix1,tx);
405             fiy1             = _mm256_add_ps(fiy1,ty);
406             fiz1             = _mm256_add_ps(fiz1,tz);
407
408             fjx0             = _mm256_add_ps(fjx0,tx);
409             fjy0             = _mm256_add_ps(fjy0,ty);
410             fjz0             = _mm256_add_ps(fjz0,tz);
411
412             /**************************
413              * CALCULATE INTERACTIONS *
414              **************************/
415
416             r20              = _mm256_mul_ps(rsq20,rinv20);
417
418             /* Compute parameters for interactions between i and j atoms */
419             qq20             = _mm256_mul_ps(iq2,jq0);
420
421             /* EWALD ELECTROSTATICS */
422             
423             /* Analytical PME correction */
424             zeta2            = _mm256_mul_ps(beta2,rsq20);
425             rinv3            = _mm256_mul_ps(rinvsq20,rinv20);
426             pmecorrF         = avx256_pmecorrF_f(zeta2);
427             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
428             felec            = _mm256_mul_ps(qq20,felec);
429             pmecorrV         = avx256_pmecorrV_f(zeta2);
430             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
431             velec            = _mm256_sub_ps(rinv20,pmecorrV);
432             velec            = _mm256_mul_ps(qq20,velec);
433             
434             /* Update potential sum for this i atom from the interaction with this j atom. */
435             velecsum         = _mm256_add_ps(velecsum,velec);
436
437             fscal            = felec;
438
439             /* Calculate temporary vectorial force */
440             tx               = _mm256_mul_ps(fscal,dx20);
441             ty               = _mm256_mul_ps(fscal,dy20);
442             tz               = _mm256_mul_ps(fscal,dz20);
443
444             /* Update vectorial force */
445             fix2             = _mm256_add_ps(fix2,tx);
446             fiy2             = _mm256_add_ps(fiy2,ty);
447             fiz2             = _mm256_add_ps(fiz2,tz);
448
449             fjx0             = _mm256_add_ps(fjx0,tx);
450             fjy0             = _mm256_add_ps(fjy0,ty);
451             fjz0             = _mm256_add_ps(fjz0,tz);
452
453             /**************************
454              * CALCULATE INTERACTIONS *
455              **************************/
456
457             r30              = _mm256_mul_ps(rsq30,rinv30);
458
459             /* Compute parameters for interactions between i and j atoms */
460             qq30             = _mm256_mul_ps(iq3,jq0);
461
462             /* EWALD ELECTROSTATICS */
463             
464             /* Analytical PME correction */
465             zeta2            = _mm256_mul_ps(beta2,rsq30);
466             rinv3            = _mm256_mul_ps(rinvsq30,rinv30);
467             pmecorrF         = avx256_pmecorrF_f(zeta2);
468             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
469             felec            = _mm256_mul_ps(qq30,felec);
470             pmecorrV         = avx256_pmecorrV_f(zeta2);
471             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
472             velec            = _mm256_sub_ps(rinv30,pmecorrV);
473             velec            = _mm256_mul_ps(qq30,velec);
474             
475             /* Update potential sum for this i atom from the interaction with this j atom. */
476             velecsum         = _mm256_add_ps(velecsum,velec);
477
478             fscal            = felec;
479
480             /* Calculate temporary vectorial force */
481             tx               = _mm256_mul_ps(fscal,dx30);
482             ty               = _mm256_mul_ps(fscal,dy30);
483             tz               = _mm256_mul_ps(fscal,dz30);
484
485             /* Update vectorial force */
486             fix3             = _mm256_add_ps(fix3,tx);
487             fiy3             = _mm256_add_ps(fiy3,ty);
488             fiz3             = _mm256_add_ps(fiz3,tz);
489
490             fjx0             = _mm256_add_ps(fjx0,tx);
491             fjy0             = _mm256_add_ps(fjy0,ty);
492             fjz0             = _mm256_add_ps(fjz0,tz);
493
494             fjptrA             = f+j_coord_offsetA;
495             fjptrB             = f+j_coord_offsetB;
496             fjptrC             = f+j_coord_offsetC;
497             fjptrD             = f+j_coord_offsetD;
498             fjptrE             = f+j_coord_offsetE;
499             fjptrF             = f+j_coord_offsetF;
500             fjptrG             = f+j_coord_offsetG;
501             fjptrH             = f+j_coord_offsetH;
502
503             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
504
505             /* Inner loop uses 311 flops */
506         }
507
508         if(jidx<j_index_end)
509         {
510
511             /* Get j neighbor index, and coordinate index */
512             jnrlistA         = jjnr[jidx];
513             jnrlistB         = jjnr[jidx+1];
514             jnrlistC         = jjnr[jidx+2];
515             jnrlistD         = jjnr[jidx+3];
516             jnrlistE         = jjnr[jidx+4];
517             jnrlistF         = jjnr[jidx+5];
518             jnrlistG         = jjnr[jidx+6];
519             jnrlistH         = jjnr[jidx+7];
520             /* Sign of each element will be negative for non-real atoms.
521              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
522              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
523              */
524             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
525                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
526                                             
527             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
528             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
529             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
530             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
531             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
532             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
533             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
534             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
535             j_coord_offsetA  = DIM*jnrA;
536             j_coord_offsetB  = DIM*jnrB;
537             j_coord_offsetC  = DIM*jnrC;
538             j_coord_offsetD  = DIM*jnrD;
539             j_coord_offsetE  = DIM*jnrE;
540             j_coord_offsetF  = DIM*jnrF;
541             j_coord_offsetG  = DIM*jnrG;
542             j_coord_offsetH  = DIM*jnrH;
543
544             /* load j atom coordinates */
545             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
546                                                  x+j_coord_offsetC,x+j_coord_offsetD,
547                                                  x+j_coord_offsetE,x+j_coord_offsetF,
548                                                  x+j_coord_offsetG,x+j_coord_offsetH,
549                                                  &jx0,&jy0,&jz0);
550
551             /* Calculate displacement vector */
552             dx00             = _mm256_sub_ps(ix0,jx0);
553             dy00             = _mm256_sub_ps(iy0,jy0);
554             dz00             = _mm256_sub_ps(iz0,jz0);
555             dx10             = _mm256_sub_ps(ix1,jx0);
556             dy10             = _mm256_sub_ps(iy1,jy0);
557             dz10             = _mm256_sub_ps(iz1,jz0);
558             dx20             = _mm256_sub_ps(ix2,jx0);
559             dy20             = _mm256_sub_ps(iy2,jy0);
560             dz20             = _mm256_sub_ps(iz2,jz0);
561             dx30             = _mm256_sub_ps(ix3,jx0);
562             dy30             = _mm256_sub_ps(iy3,jy0);
563             dz30             = _mm256_sub_ps(iz3,jz0);
564
565             /* Calculate squared distance and things based on it */
566             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
567             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
568             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
569             rsq30            = gmx_mm256_calc_rsq_ps(dx30,dy30,dz30);
570
571             rinv00           = avx256_invsqrt_f(rsq00);
572             rinv10           = avx256_invsqrt_f(rsq10);
573             rinv20           = avx256_invsqrt_f(rsq20);
574             rinv30           = avx256_invsqrt_f(rsq30);
575
576             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
577             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
578             rinvsq30         = _mm256_mul_ps(rinv30,rinv30);
579
580             /* Load parameters for j particles */
581             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
582                                                                  charge+jnrC+0,charge+jnrD+0,
583                                                                  charge+jnrE+0,charge+jnrF+0,
584                                                                  charge+jnrG+0,charge+jnrH+0);
585             vdwjidx0A        = 2*vdwtype[jnrA+0];
586             vdwjidx0B        = 2*vdwtype[jnrB+0];
587             vdwjidx0C        = 2*vdwtype[jnrC+0];
588             vdwjidx0D        = 2*vdwtype[jnrD+0];
589             vdwjidx0E        = 2*vdwtype[jnrE+0];
590             vdwjidx0F        = 2*vdwtype[jnrF+0];
591             vdwjidx0G        = 2*vdwtype[jnrG+0];
592             vdwjidx0H        = 2*vdwtype[jnrH+0];
593
594             fjx0             = _mm256_setzero_ps();
595             fjy0             = _mm256_setzero_ps();
596             fjz0             = _mm256_setzero_ps();
597
598             /**************************
599              * CALCULATE INTERACTIONS *
600              **************************/
601
602             r00              = _mm256_mul_ps(rsq00,rinv00);
603             r00              = _mm256_andnot_ps(dummy_mask,r00);
604
605             /* Compute parameters for interactions between i and j atoms */
606             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
607                                             vdwioffsetptr0+vdwjidx0B,
608                                             vdwioffsetptr0+vdwjidx0C,
609                                             vdwioffsetptr0+vdwjidx0D,
610                                             vdwioffsetptr0+vdwjidx0E,
611                                             vdwioffsetptr0+vdwjidx0F,
612                                             vdwioffsetptr0+vdwjidx0G,
613                                             vdwioffsetptr0+vdwjidx0H,
614                                             &c6_00,&c12_00);
615
616             /* Calculate table index by multiplying r with table scale and truncate to integer */
617             rt               = _mm256_mul_ps(r00,vftabscale);
618             vfitab           = _mm256_cvttps_epi32(rt);
619             vfeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
620             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
621             vfitab_lo        = _mm256_extractf128_si256(vfitab,0x0);
622             vfitab_hi        = _mm256_extractf128_si256(vfitab,0x1);
623             vfitab_lo        = _mm_slli_epi32(vfitab_lo,3);
624             vfitab_hi        = _mm_slli_epi32(vfitab_hi,3);
625
626             /* CUBIC SPLINE TABLE DISPERSION */
627             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
628                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
629             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
630                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
631             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
632                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
633             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
634                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
635             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
636             Heps             = _mm256_mul_ps(vfeps,H);
637             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
638             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
639             vvdw6            = _mm256_mul_ps(c6_00,VV);
640             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
641             fvdw6            = _mm256_mul_ps(c6_00,FF);
642
643             /* CUBIC SPLINE TABLE REPULSION */
644             vfitab_lo        = _mm_add_epi32(vfitab_lo,ifour);
645             vfitab_hi        = _mm_add_epi32(vfitab_hi,ifour);
646             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
647                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
648             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
649                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
650             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
651                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
652             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
653                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
654             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
655             Heps             = _mm256_mul_ps(vfeps,H);
656             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
657             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
658             vvdw12           = _mm256_mul_ps(c12_00,VV);
659             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
660             fvdw12           = _mm256_mul_ps(c12_00,FF);
661             vvdw             = _mm256_add_ps(vvdw12,vvdw6);
662             fvdw             = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv00)));
663
664             /* Update potential sum for this i atom from the interaction with this j atom. */
665             vvdw             = _mm256_andnot_ps(dummy_mask,vvdw);
666             vvdwsum          = _mm256_add_ps(vvdwsum,vvdw);
667
668             fscal            = fvdw;
669
670             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
671
672             /* Calculate temporary vectorial force */
673             tx               = _mm256_mul_ps(fscal,dx00);
674             ty               = _mm256_mul_ps(fscal,dy00);
675             tz               = _mm256_mul_ps(fscal,dz00);
676
677             /* Update vectorial force */
678             fix0             = _mm256_add_ps(fix0,tx);
679             fiy0             = _mm256_add_ps(fiy0,ty);
680             fiz0             = _mm256_add_ps(fiz0,tz);
681
682             fjx0             = _mm256_add_ps(fjx0,tx);
683             fjy0             = _mm256_add_ps(fjy0,ty);
684             fjz0             = _mm256_add_ps(fjz0,tz);
685
686             /**************************
687              * CALCULATE INTERACTIONS *
688              **************************/
689
690             r10              = _mm256_mul_ps(rsq10,rinv10);
691             r10              = _mm256_andnot_ps(dummy_mask,r10);
692
693             /* Compute parameters for interactions between i and j atoms */
694             qq10             = _mm256_mul_ps(iq1,jq0);
695
696             /* EWALD ELECTROSTATICS */
697             
698             /* Analytical PME correction */
699             zeta2            = _mm256_mul_ps(beta2,rsq10);
700             rinv3            = _mm256_mul_ps(rinvsq10,rinv10);
701             pmecorrF         = avx256_pmecorrF_f(zeta2);
702             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
703             felec            = _mm256_mul_ps(qq10,felec);
704             pmecorrV         = avx256_pmecorrV_f(zeta2);
705             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
706             velec            = _mm256_sub_ps(rinv10,pmecorrV);
707             velec            = _mm256_mul_ps(qq10,velec);
708             
709             /* Update potential sum for this i atom from the interaction with this j atom. */
710             velec            = _mm256_andnot_ps(dummy_mask,velec);
711             velecsum         = _mm256_add_ps(velecsum,velec);
712
713             fscal            = felec;
714
715             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
716
717             /* Calculate temporary vectorial force */
718             tx               = _mm256_mul_ps(fscal,dx10);
719             ty               = _mm256_mul_ps(fscal,dy10);
720             tz               = _mm256_mul_ps(fscal,dz10);
721
722             /* Update vectorial force */
723             fix1             = _mm256_add_ps(fix1,tx);
724             fiy1             = _mm256_add_ps(fiy1,ty);
725             fiz1             = _mm256_add_ps(fiz1,tz);
726
727             fjx0             = _mm256_add_ps(fjx0,tx);
728             fjy0             = _mm256_add_ps(fjy0,ty);
729             fjz0             = _mm256_add_ps(fjz0,tz);
730
731             /**************************
732              * CALCULATE INTERACTIONS *
733              **************************/
734
735             r20              = _mm256_mul_ps(rsq20,rinv20);
736             r20              = _mm256_andnot_ps(dummy_mask,r20);
737
738             /* Compute parameters for interactions between i and j atoms */
739             qq20             = _mm256_mul_ps(iq2,jq0);
740
741             /* EWALD ELECTROSTATICS */
742             
743             /* Analytical PME correction */
744             zeta2            = _mm256_mul_ps(beta2,rsq20);
745             rinv3            = _mm256_mul_ps(rinvsq20,rinv20);
746             pmecorrF         = avx256_pmecorrF_f(zeta2);
747             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
748             felec            = _mm256_mul_ps(qq20,felec);
749             pmecorrV         = avx256_pmecorrV_f(zeta2);
750             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
751             velec            = _mm256_sub_ps(rinv20,pmecorrV);
752             velec            = _mm256_mul_ps(qq20,velec);
753             
754             /* Update potential sum for this i atom from the interaction with this j atom. */
755             velec            = _mm256_andnot_ps(dummy_mask,velec);
756             velecsum         = _mm256_add_ps(velecsum,velec);
757
758             fscal            = felec;
759
760             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
761
762             /* Calculate temporary vectorial force */
763             tx               = _mm256_mul_ps(fscal,dx20);
764             ty               = _mm256_mul_ps(fscal,dy20);
765             tz               = _mm256_mul_ps(fscal,dz20);
766
767             /* Update vectorial force */
768             fix2             = _mm256_add_ps(fix2,tx);
769             fiy2             = _mm256_add_ps(fiy2,ty);
770             fiz2             = _mm256_add_ps(fiz2,tz);
771
772             fjx0             = _mm256_add_ps(fjx0,tx);
773             fjy0             = _mm256_add_ps(fjy0,ty);
774             fjz0             = _mm256_add_ps(fjz0,tz);
775
776             /**************************
777              * CALCULATE INTERACTIONS *
778              **************************/
779
780             r30              = _mm256_mul_ps(rsq30,rinv30);
781             r30              = _mm256_andnot_ps(dummy_mask,r30);
782
783             /* Compute parameters for interactions between i and j atoms */
784             qq30             = _mm256_mul_ps(iq3,jq0);
785
786             /* EWALD ELECTROSTATICS */
787             
788             /* Analytical PME correction */
789             zeta2            = _mm256_mul_ps(beta2,rsq30);
790             rinv3            = _mm256_mul_ps(rinvsq30,rinv30);
791             pmecorrF         = avx256_pmecorrF_f(zeta2);
792             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
793             felec            = _mm256_mul_ps(qq30,felec);
794             pmecorrV         = avx256_pmecorrV_f(zeta2);
795             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
796             velec            = _mm256_sub_ps(rinv30,pmecorrV);
797             velec            = _mm256_mul_ps(qq30,velec);
798             
799             /* Update potential sum for this i atom from the interaction with this j atom. */
800             velec            = _mm256_andnot_ps(dummy_mask,velec);
801             velecsum         = _mm256_add_ps(velecsum,velec);
802
803             fscal            = felec;
804
805             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
806
807             /* Calculate temporary vectorial force */
808             tx               = _mm256_mul_ps(fscal,dx30);
809             ty               = _mm256_mul_ps(fscal,dy30);
810             tz               = _mm256_mul_ps(fscal,dz30);
811
812             /* Update vectorial force */
813             fix3             = _mm256_add_ps(fix3,tx);
814             fiy3             = _mm256_add_ps(fiy3,ty);
815             fiz3             = _mm256_add_ps(fiz3,tz);
816
817             fjx0             = _mm256_add_ps(fjx0,tx);
818             fjy0             = _mm256_add_ps(fjy0,ty);
819             fjz0             = _mm256_add_ps(fjz0,tz);
820
821             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
822             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
823             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
824             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
825             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
826             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
827             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
828             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
829
830             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
831
832             /* Inner loop uses 315 flops */
833         }
834
835         /* End of innermost loop */
836
837         gmx_mm256_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
838                                                  f+i_coord_offset,fshift+i_shift_offset);
839
840         ggid                        = gid[iidx];
841         /* Update potential energies */
842         gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
843         gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
844
845         /* Increment number of inner iterations */
846         inneriter                  += j_index_end - j_index_start;
847
848         /* Outer loop uses 26 flops */
849     }
850
851     /* Increment number of outer iterations */
852     outeriter        += nri;
853
854     /* Update outer/inner flops */
855
856     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*315);
857 }
858 /*
859  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwCSTab_GeomW4P1_F_avx_256_single
860  * Electrostatics interaction: Ewald
861  * VdW interaction:            CubicSplineTable
862  * Geometry:                   Water4-Particle
863  * Calculate force/pot:        Force
864  */
865 void
866 nb_kernel_ElecEw_VdwCSTab_GeomW4P1_F_avx_256_single
867                     (t_nblist                    * gmx_restrict       nlist,
868                      rvec                        * gmx_restrict          xx,
869                      rvec                        * gmx_restrict          ff,
870                      struct t_forcerec           * gmx_restrict          fr,
871                      t_mdatoms                   * gmx_restrict     mdatoms,
872                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
873                      t_nrnb                      * gmx_restrict        nrnb)
874 {
875     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
876      * just 0 for non-waters.
877      * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
878      * jnr indices corresponding to data put in the four positions in the SIMD register.
879      */
880     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
881     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
882     int              jnrA,jnrB,jnrC,jnrD;
883     int              jnrE,jnrF,jnrG,jnrH;
884     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
885     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
886     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
887     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
888     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
889     real             rcutoff_scalar;
890     real             *shiftvec,*fshift,*x,*f;
891     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
892     real             scratch[4*DIM];
893     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
894     real *           vdwioffsetptr0;
895     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
896     real *           vdwioffsetptr1;
897     __m256           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
898     real *           vdwioffsetptr2;
899     __m256           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
900     real *           vdwioffsetptr3;
901     __m256           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
902     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
903     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
904     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
905     __m256           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
906     __m256           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
907     __m256           dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
908     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
909     real             *charge;
910     int              nvdwtype;
911     __m256           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
912     int              *vdwtype;
913     real             *vdwparam;
914     __m256           one_sixth   = _mm256_set1_ps(1.0/6.0);
915     __m256           one_twelfth = _mm256_set1_ps(1.0/12.0);
916     __m256i          vfitab;
917     __m128i          vfitab_lo,vfitab_hi;
918     __m128i          ifour       = _mm_set1_epi32(4);
919     __m256           rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
920     real             *vftab;
921     __m256i          ewitab;
922     __m128i          ewitab_lo,ewitab_hi;
923     __m256           ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
924     __m256           beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
925     real             *ewtab;
926     __m256           dummy_mask,cutoff_mask;
927     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
928     __m256           one     = _mm256_set1_ps(1.0);
929     __m256           two     = _mm256_set1_ps(2.0);
930     x                = xx[0];
931     f                = ff[0];
932
933     nri              = nlist->nri;
934     iinr             = nlist->iinr;
935     jindex           = nlist->jindex;
936     jjnr             = nlist->jjnr;
937     shiftidx         = nlist->shift;
938     gid              = nlist->gid;
939     shiftvec         = fr->shift_vec[0];
940     fshift           = fr->fshift[0];
941     facel            = _mm256_set1_ps(fr->ic->epsfac);
942     charge           = mdatoms->chargeA;
943     nvdwtype         = fr->ntype;
944     vdwparam         = fr->nbfp;
945     vdwtype          = mdatoms->typeA;
946
947     vftab            = kernel_data->table_vdw->data;
948     vftabscale       = _mm256_set1_ps(kernel_data->table_vdw->scale);
949
950     sh_ewald         = _mm256_set1_ps(fr->ic->sh_ewald);
951     beta             = _mm256_set1_ps(fr->ic->ewaldcoeff_q);
952     beta2            = _mm256_mul_ps(beta,beta);
953     beta3            = _mm256_mul_ps(beta,beta2);
954
955     ewtab            = fr->ic->tabq_coul_F;
956     ewtabscale       = _mm256_set1_ps(fr->ic->tabq_scale);
957     ewtabhalfspace   = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
958
959     /* Setup water-specific parameters */
960     inr              = nlist->iinr[0];
961     iq1              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
962     iq2              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
963     iq3              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+3]));
964     vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
965
966     /* Avoid stupid compiler warnings */
967     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
968     j_coord_offsetA = 0;
969     j_coord_offsetB = 0;
970     j_coord_offsetC = 0;
971     j_coord_offsetD = 0;
972     j_coord_offsetE = 0;
973     j_coord_offsetF = 0;
974     j_coord_offsetG = 0;
975     j_coord_offsetH = 0;
976
977     outeriter        = 0;
978     inneriter        = 0;
979
980     for(iidx=0;iidx<4*DIM;iidx++)
981     {
982         scratch[iidx] = 0.0;
983     }
984
985     /* Start outer loop over neighborlists */
986     for(iidx=0; iidx<nri; iidx++)
987     {
988         /* Load shift vector for this list */
989         i_shift_offset   = DIM*shiftidx[iidx];
990
991         /* Load limits for loop over neighbors */
992         j_index_start    = jindex[iidx];
993         j_index_end      = jindex[iidx+1];
994
995         /* Get outer coordinate index */
996         inr              = iinr[iidx];
997         i_coord_offset   = DIM*inr;
998
999         /* Load i particle coords and add shift vector */
1000         gmx_mm256_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
1001                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
1002
1003         fix0             = _mm256_setzero_ps();
1004         fiy0             = _mm256_setzero_ps();
1005         fiz0             = _mm256_setzero_ps();
1006         fix1             = _mm256_setzero_ps();
1007         fiy1             = _mm256_setzero_ps();
1008         fiz1             = _mm256_setzero_ps();
1009         fix2             = _mm256_setzero_ps();
1010         fiy2             = _mm256_setzero_ps();
1011         fiz2             = _mm256_setzero_ps();
1012         fix3             = _mm256_setzero_ps();
1013         fiy3             = _mm256_setzero_ps();
1014         fiz3             = _mm256_setzero_ps();
1015
1016         /* Start inner kernel loop */
1017         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
1018         {
1019
1020             /* Get j neighbor index, and coordinate index */
1021             jnrA             = jjnr[jidx];
1022             jnrB             = jjnr[jidx+1];
1023             jnrC             = jjnr[jidx+2];
1024             jnrD             = jjnr[jidx+3];
1025             jnrE             = jjnr[jidx+4];
1026             jnrF             = jjnr[jidx+5];
1027             jnrG             = jjnr[jidx+6];
1028             jnrH             = jjnr[jidx+7];
1029             j_coord_offsetA  = DIM*jnrA;
1030             j_coord_offsetB  = DIM*jnrB;
1031             j_coord_offsetC  = DIM*jnrC;
1032             j_coord_offsetD  = DIM*jnrD;
1033             j_coord_offsetE  = DIM*jnrE;
1034             j_coord_offsetF  = DIM*jnrF;
1035             j_coord_offsetG  = DIM*jnrG;
1036             j_coord_offsetH  = DIM*jnrH;
1037
1038             /* load j atom coordinates */
1039             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1040                                                  x+j_coord_offsetC,x+j_coord_offsetD,
1041                                                  x+j_coord_offsetE,x+j_coord_offsetF,
1042                                                  x+j_coord_offsetG,x+j_coord_offsetH,
1043                                                  &jx0,&jy0,&jz0);
1044
1045             /* Calculate displacement vector */
1046             dx00             = _mm256_sub_ps(ix0,jx0);
1047             dy00             = _mm256_sub_ps(iy0,jy0);
1048             dz00             = _mm256_sub_ps(iz0,jz0);
1049             dx10             = _mm256_sub_ps(ix1,jx0);
1050             dy10             = _mm256_sub_ps(iy1,jy0);
1051             dz10             = _mm256_sub_ps(iz1,jz0);
1052             dx20             = _mm256_sub_ps(ix2,jx0);
1053             dy20             = _mm256_sub_ps(iy2,jy0);
1054             dz20             = _mm256_sub_ps(iz2,jz0);
1055             dx30             = _mm256_sub_ps(ix3,jx0);
1056             dy30             = _mm256_sub_ps(iy3,jy0);
1057             dz30             = _mm256_sub_ps(iz3,jz0);
1058
1059             /* Calculate squared distance and things based on it */
1060             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
1061             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
1062             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
1063             rsq30            = gmx_mm256_calc_rsq_ps(dx30,dy30,dz30);
1064
1065             rinv00           = avx256_invsqrt_f(rsq00);
1066             rinv10           = avx256_invsqrt_f(rsq10);
1067             rinv20           = avx256_invsqrt_f(rsq20);
1068             rinv30           = avx256_invsqrt_f(rsq30);
1069
1070             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
1071             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
1072             rinvsq30         = _mm256_mul_ps(rinv30,rinv30);
1073
1074             /* Load parameters for j particles */
1075             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1076                                                                  charge+jnrC+0,charge+jnrD+0,
1077                                                                  charge+jnrE+0,charge+jnrF+0,
1078                                                                  charge+jnrG+0,charge+jnrH+0);
1079             vdwjidx0A        = 2*vdwtype[jnrA+0];
1080             vdwjidx0B        = 2*vdwtype[jnrB+0];
1081             vdwjidx0C        = 2*vdwtype[jnrC+0];
1082             vdwjidx0D        = 2*vdwtype[jnrD+0];
1083             vdwjidx0E        = 2*vdwtype[jnrE+0];
1084             vdwjidx0F        = 2*vdwtype[jnrF+0];
1085             vdwjidx0G        = 2*vdwtype[jnrG+0];
1086             vdwjidx0H        = 2*vdwtype[jnrH+0];
1087
1088             fjx0             = _mm256_setzero_ps();
1089             fjy0             = _mm256_setzero_ps();
1090             fjz0             = _mm256_setzero_ps();
1091
1092             /**************************
1093              * CALCULATE INTERACTIONS *
1094              **************************/
1095
1096             r00              = _mm256_mul_ps(rsq00,rinv00);
1097
1098             /* Compute parameters for interactions between i and j atoms */
1099             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
1100                                             vdwioffsetptr0+vdwjidx0B,
1101                                             vdwioffsetptr0+vdwjidx0C,
1102                                             vdwioffsetptr0+vdwjidx0D,
1103                                             vdwioffsetptr0+vdwjidx0E,
1104                                             vdwioffsetptr0+vdwjidx0F,
1105                                             vdwioffsetptr0+vdwjidx0G,
1106                                             vdwioffsetptr0+vdwjidx0H,
1107                                             &c6_00,&c12_00);
1108
1109             /* Calculate table index by multiplying r with table scale and truncate to integer */
1110             rt               = _mm256_mul_ps(r00,vftabscale);
1111             vfitab           = _mm256_cvttps_epi32(rt);
1112             vfeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
1113             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
1114             vfitab_lo        = _mm256_extractf128_si256(vfitab,0x0);
1115             vfitab_hi        = _mm256_extractf128_si256(vfitab,0x1);
1116             vfitab_lo        = _mm_slli_epi32(vfitab_lo,3);
1117             vfitab_hi        = _mm_slli_epi32(vfitab_hi,3);
1118
1119             /* CUBIC SPLINE TABLE DISPERSION */
1120             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
1121                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
1122             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
1123                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
1124             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
1125                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
1126             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
1127                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
1128             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
1129             Heps             = _mm256_mul_ps(vfeps,H);
1130             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
1131             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
1132             fvdw6            = _mm256_mul_ps(c6_00,FF);
1133
1134             /* CUBIC SPLINE TABLE REPULSION */
1135             vfitab_lo        = _mm_add_epi32(vfitab_lo,ifour);
1136             vfitab_hi        = _mm_add_epi32(vfitab_hi,ifour);
1137             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
1138                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
1139             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
1140                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
1141             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
1142                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
1143             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
1144                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
1145             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
1146             Heps             = _mm256_mul_ps(vfeps,H);
1147             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
1148             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
1149             fvdw12           = _mm256_mul_ps(c12_00,FF);
1150             fvdw             = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv00)));
1151
1152             fscal            = fvdw;
1153
1154             /* Calculate temporary vectorial force */
1155             tx               = _mm256_mul_ps(fscal,dx00);
1156             ty               = _mm256_mul_ps(fscal,dy00);
1157             tz               = _mm256_mul_ps(fscal,dz00);
1158
1159             /* Update vectorial force */
1160             fix0             = _mm256_add_ps(fix0,tx);
1161             fiy0             = _mm256_add_ps(fiy0,ty);
1162             fiz0             = _mm256_add_ps(fiz0,tz);
1163
1164             fjx0             = _mm256_add_ps(fjx0,tx);
1165             fjy0             = _mm256_add_ps(fjy0,ty);
1166             fjz0             = _mm256_add_ps(fjz0,tz);
1167
1168             /**************************
1169              * CALCULATE INTERACTIONS *
1170              **************************/
1171
1172             r10              = _mm256_mul_ps(rsq10,rinv10);
1173
1174             /* Compute parameters for interactions between i and j atoms */
1175             qq10             = _mm256_mul_ps(iq1,jq0);
1176
1177             /* EWALD ELECTROSTATICS */
1178             
1179             /* Analytical PME correction */
1180             zeta2            = _mm256_mul_ps(beta2,rsq10);
1181             rinv3            = _mm256_mul_ps(rinvsq10,rinv10);
1182             pmecorrF         = avx256_pmecorrF_f(zeta2);
1183             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1184             felec            = _mm256_mul_ps(qq10,felec);
1185             
1186             fscal            = felec;
1187
1188             /* Calculate temporary vectorial force */
1189             tx               = _mm256_mul_ps(fscal,dx10);
1190             ty               = _mm256_mul_ps(fscal,dy10);
1191             tz               = _mm256_mul_ps(fscal,dz10);
1192
1193             /* Update vectorial force */
1194             fix1             = _mm256_add_ps(fix1,tx);
1195             fiy1             = _mm256_add_ps(fiy1,ty);
1196             fiz1             = _mm256_add_ps(fiz1,tz);
1197
1198             fjx0             = _mm256_add_ps(fjx0,tx);
1199             fjy0             = _mm256_add_ps(fjy0,ty);
1200             fjz0             = _mm256_add_ps(fjz0,tz);
1201
1202             /**************************
1203              * CALCULATE INTERACTIONS *
1204              **************************/
1205
1206             r20              = _mm256_mul_ps(rsq20,rinv20);
1207
1208             /* Compute parameters for interactions between i and j atoms */
1209             qq20             = _mm256_mul_ps(iq2,jq0);
1210
1211             /* EWALD ELECTROSTATICS */
1212             
1213             /* Analytical PME correction */
1214             zeta2            = _mm256_mul_ps(beta2,rsq20);
1215             rinv3            = _mm256_mul_ps(rinvsq20,rinv20);
1216             pmecorrF         = avx256_pmecorrF_f(zeta2);
1217             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1218             felec            = _mm256_mul_ps(qq20,felec);
1219             
1220             fscal            = felec;
1221
1222             /* Calculate temporary vectorial force */
1223             tx               = _mm256_mul_ps(fscal,dx20);
1224             ty               = _mm256_mul_ps(fscal,dy20);
1225             tz               = _mm256_mul_ps(fscal,dz20);
1226
1227             /* Update vectorial force */
1228             fix2             = _mm256_add_ps(fix2,tx);
1229             fiy2             = _mm256_add_ps(fiy2,ty);
1230             fiz2             = _mm256_add_ps(fiz2,tz);
1231
1232             fjx0             = _mm256_add_ps(fjx0,tx);
1233             fjy0             = _mm256_add_ps(fjy0,ty);
1234             fjz0             = _mm256_add_ps(fjz0,tz);
1235
1236             /**************************
1237              * CALCULATE INTERACTIONS *
1238              **************************/
1239
1240             r30              = _mm256_mul_ps(rsq30,rinv30);
1241
1242             /* Compute parameters for interactions between i and j atoms */
1243             qq30             = _mm256_mul_ps(iq3,jq0);
1244
1245             /* EWALD ELECTROSTATICS */
1246             
1247             /* Analytical PME correction */
1248             zeta2            = _mm256_mul_ps(beta2,rsq30);
1249             rinv3            = _mm256_mul_ps(rinvsq30,rinv30);
1250             pmecorrF         = avx256_pmecorrF_f(zeta2);
1251             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1252             felec            = _mm256_mul_ps(qq30,felec);
1253             
1254             fscal            = felec;
1255
1256             /* Calculate temporary vectorial force */
1257             tx               = _mm256_mul_ps(fscal,dx30);
1258             ty               = _mm256_mul_ps(fscal,dy30);
1259             tz               = _mm256_mul_ps(fscal,dz30);
1260
1261             /* Update vectorial force */
1262             fix3             = _mm256_add_ps(fix3,tx);
1263             fiy3             = _mm256_add_ps(fiy3,ty);
1264             fiz3             = _mm256_add_ps(fiz3,tz);
1265
1266             fjx0             = _mm256_add_ps(fjx0,tx);
1267             fjy0             = _mm256_add_ps(fjy0,ty);
1268             fjz0             = _mm256_add_ps(fjz0,tz);
1269
1270             fjptrA             = f+j_coord_offsetA;
1271             fjptrB             = f+j_coord_offsetB;
1272             fjptrC             = f+j_coord_offsetC;
1273             fjptrD             = f+j_coord_offsetD;
1274             fjptrE             = f+j_coord_offsetE;
1275             fjptrF             = f+j_coord_offsetF;
1276             fjptrG             = f+j_coord_offsetG;
1277             fjptrH             = f+j_coord_offsetH;
1278
1279             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
1280
1281             /* Inner loop uses 219 flops */
1282         }
1283
1284         if(jidx<j_index_end)
1285         {
1286
1287             /* Get j neighbor index, and coordinate index */
1288             jnrlistA         = jjnr[jidx];
1289             jnrlistB         = jjnr[jidx+1];
1290             jnrlistC         = jjnr[jidx+2];
1291             jnrlistD         = jjnr[jidx+3];
1292             jnrlistE         = jjnr[jidx+4];
1293             jnrlistF         = jjnr[jidx+5];
1294             jnrlistG         = jjnr[jidx+6];
1295             jnrlistH         = jjnr[jidx+7];
1296             /* Sign of each element will be negative for non-real atoms.
1297              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1298              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1299              */
1300             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
1301                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
1302                                             
1303             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
1304             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
1305             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
1306             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
1307             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
1308             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
1309             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
1310             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
1311             j_coord_offsetA  = DIM*jnrA;
1312             j_coord_offsetB  = DIM*jnrB;
1313             j_coord_offsetC  = DIM*jnrC;
1314             j_coord_offsetD  = DIM*jnrD;
1315             j_coord_offsetE  = DIM*jnrE;
1316             j_coord_offsetF  = DIM*jnrF;
1317             j_coord_offsetG  = DIM*jnrG;
1318             j_coord_offsetH  = DIM*jnrH;
1319
1320             /* load j atom coordinates */
1321             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1322                                                  x+j_coord_offsetC,x+j_coord_offsetD,
1323                                                  x+j_coord_offsetE,x+j_coord_offsetF,
1324                                                  x+j_coord_offsetG,x+j_coord_offsetH,
1325                                                  &jx0,&jy0,&jz0);
1326
1327             /* Calculate displacement vector */
1328             dx00             = _mm256_sub_ps(ix0,jx0);
1329             dy00             = _mm256_sub_ps(iy0,jy0);
1330             dz00             = _mm256_sub_ps(iz0,jz0);
1331             dx10             = _mm256_sub_ps(ix1,jx0);
1332             dy10             = _mm256_sub_ps(iy1,jy0);
1333             dz10             = _mm256_sub_ps(iz1,jz0);
1334             dx20             = _mm256_sub_ps(ix2,jx0);
1335             dy20             = _mm256_sub_ps(iy2,jy0);
1336             dz20             = _mm256_sub_ps(iz2,jz0);
1337             dx30             = _mm256_sub_ps(ix3,jx0);
1338             dy30             = _mm256_sub_ps(iy3,jy0);
1339             dz30             = _mm256_sub_ps(iz3,jz0);
1340
1341             /* Calculate squared distance and things based on it */
1342             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
1343             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
1344             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
1345             rsq30            = gmx_mm256_calc_rsq_ps(dx30,dy30,dz30);
1346
1347             rinv00           = avx256_invsqrt_f(rsq00);
1348             rinv10           = avx256_invsqrt_f(rsq10);
1349             rinv20           = avx256_invsqrt_f(rsq20);
1350             rinv30           = avx256_invsqrt_f(rsq30);
1351
1352             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
1353             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
1354             rinvsq30         = _mm256_mul_ps(rinv30,rinv30);
1355
1356             /* Load parameters for j particles */
1357             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1358                                                                  charge+jnrC+0,charge+jnrD+0,
1359                                                                  charge+jnrE+0,charge+jnrF+0,
1360                                                                  charge+jnrG+0,charge+jnrH+0);
1361             vdwjidx0A        = 2*vdwtype[jnrA+0];
1362             vdwjidx0B        = 2*vdwtype[jnrB+0];
1363             vdwjidx0C        = 2*vdwtype[jnrC+0];
1364             vdwjidx0D        = 2*vdwtype[jnrD+0];
1365             vdwjidx0E        = 2*vdwtype[jnrE+0];
1366             vdwjidx0F        = 2*vdwtype[jnrF+0];
1367             vdwjidx0G        = 2*vdwtype[jnrG+0];
1368             vdwjidx0H        = 2*vdwtype[jnrH+0];
1369
1370             fjx0             = _mm256_setzero_ps();
1371             fjy0             = _mm256_setzero_ps();
1372             fjz0             = _mm256_setzero_ps();
1373
1374             /**************************
1375              * CALCULATE INTERACTIONS *
1376              **************************/
1377
1378             r00              = _mm256_mul_ps(rsq00,rinv00);
1379             r00              = _mm256_andnot_ps(dummy_mask,r00);
1380
1381             /* Compute parameters for interactions between i and j atoms */
1382             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
1383                                             vdwioffsetptr0+vdwjidx0B,
1384                                             vdwioffsetptr0+vdwjidx0C,
1385                                             vdwioffsetptr0+vdwjidx0D,
1386                                             vdwioffsetptr0+vdwjidx0E,
1387                                             vdwioffsetptr0+vdwjidx0F,
1388                                             vdwioffsetptr0+vdwjidx0G,
1389                                             vdwioffsetptr0+vdwjidx0H,
1390                                             &c6_00,&c12_00);
1391
1392             /* Calculate table index by multiplying r with table scale and truncate to integer */
1393             rt               = _mm256_mul_ps(r00,vftabscale);
1394             vfitab           = _mm256_cvttps_epi32(rt);
1395             vfeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
1396             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
1397             vfitab_lo        = _mm256_extractf128_si256(vfitab,0x0);
1398             vfitab_hi        = _mm256_extractf128_si256(vfitab,0x1);
1399             vfitab_lo        = _mm_slli_epi32(vfitab_lo,3);
1400             vfitab_hi        = _mm_slli_epi32(vfitab_hi,3);
1401
1402             /* CUBIC SPLINE TABLE DISPERSION */
1403             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
1404                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
1405             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
1406                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
1407             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
1408                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
1409             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
1410                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
1411             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
1412             Heps             = _mm256_mul_ps(vfeps,H);
1413             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
1414             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
1415             fvdw6            = _mm256_mul_ps(c6_00,FF);
1416
1417             /* CUBIC SPLINE TABLE REPULSION */
1418             vfitab_lo        = _mm_add_epi32(vfitab_lo,ifour);
1419             vfitab_hi        = _mm_add_epi32(vfitab_hi,ifour);
1420             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
1421                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
1422             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
1423                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
1424             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
1425                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
1426             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
1427                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
1428             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
1429             Heps             = _mm256_mul_ps(vfeps,H);
1430             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
1431             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
1432             fvdw12           = _mm256_mul_ps(c12_00,FF);
1433             fvdw             = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv00)));
1434
1435             fscal            = fvdw;
1436
1437             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1438
1439             /* Calculate temporary vectorial force */
1440             tx               = _mm256_mul_ps(fscal,dx00);
1441             ty               = _mm256_mul_ps(fscal,dy00);
1442             tz               = _mm256_mul_ps(fscal,dz00);
1443
1444             /* Update vectorial force */
1445             fix0             = _mm256_add_ps(fix0,tx);
1446             fiy0             = _mm256_add_ps(fiy0,ty);
1447             fiz0             = _mm256_add_ps(fiz0,tz);
1448
1449             fjx0             = _mm256_add_ps(fjx0,tx);
1450             fjy0             = _mm256_add_ps(fjy0,ty);
1451             fjz0             = _mm256_add_ps(fjz0,tz);
1452
1453             /**************************
1454              * CALCULATE INTERACTIONS *
1455              **************************/
1456
1457             r10              = _mm256_mul_ps(rsq10,rinv10);
1458             r10              = _mm256_andnot_ps(dummy_mask,r10);
1459
1460             /* Compute parameters for interactions between i and j atoms */
1461             qq10             = _mm256_mul_ps(iq1,jq0);
1462
1463             /* EWALD ELECTROSTATICS */
1464             
1465             /* Analytical PME correction */
1466             zeta2            = _mm256_mul_ps(beta2,rsq10);
1467             rinv3            = _mm256_mul_ps(rinvsq10,rinv10);
1468             pmecorrF         = avx256_pmecorrF_f(zeta2);
1469             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1470             felec            = _mm256_mul_ps(qq10,felec);
1471             
1472             fscal            = felec;
1473
1474             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1475
1476             /* Calculate temporary vectorial force */
1477             tx               = _mm256_mul_ps(fscal,dx10);
1478             ty               = _mm256_mul_ps(fscal,dy10);
1479             tz               = _mm256_mul_ps(fscal,dz10);
1480
1481             /* Update vectorial force */
1482             fix1             = _mm256_add_ps(fix1,tx);
1483             fiy1             = _mm256_add_ps(fiy1,ty);
1484             fiz1             = _mm256_add_ps(fiz1,tz);
1485
1486             fjx0             = _mm256_add_ps(fjx0,tx);
1487             fjy0             = _mm256_add_ps(fjy0,ty);
1488             fjz0             = _mm256_add_ps(fjz0,tz);
1489
1490             /**************************
1491              * CALCULATE INTERACTIONS *
1492              **************************/
1493
1494             r20              = _mm256_mul_ps(rsq20,rinv20);
1495             r20              = _mm256_andnot_ps(dummy_mask,r20);
1496
1497             /* Compute parameters for interactions between i and j atoms */
1498             qq20             = _mm256_mul_ps(iq2,jq0);
1499
1500             /* EWALD ELECTROSTATICS */
1501             
1502             /* Analytical PME correction */
1503             zeta2            = _mm256_mul_ps(beta2,rsq20);
1504             rinv3            = _mm256_mul_ps(rinvsq20,rinv20);
1505             pmecorrF         = avx256_pmecorrF_f(zeta2);
1506             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1507             felec            = _mm256_mul_ps(qq20,felec);
1508             
1509             fscal            = felec;
1510
1511             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1512
1513             /* Calculate temporary vectorial force */
1514             tx               = _mm256_mul_ps(fscal,dx20);
1515             ty               = _mm256_mul_ps(fscal,dy20);
1516             tz               = _mm256_mul_ps(fscal,dz20);
1517
1518             /* Update vectorial force */
1519             fix2             = _mm256_add_ps(fix2,tx);
1520             fiy2             = _mm256_add_ps(fiy2,ty);
1521             fiz2             = _mm256_add_ps(fiz2,tz);
1522
1523             fjx0             = _mm256_add_ps(fjx0,tx);
1524             fjy0             = _mm256_add_ps(fjy0,ty);
1525             fjz0             = _mm256_add_ps(fjz0,tz);
1526
1527             /**************************
1528              * CALCULATE INTERACTIONS *
1529              **************************/
1530
1531             r30              = _mm256_mul_ps(rsq30,rinv30);
1532             r30              = _mm256_andnot_ps(dummy_mask,r30);
1533
1534             /* Compute parameters for interactions between i and j atoms */
1535             qq30             = _mm256_mul_ps(iq3,jq0);
1536
1537             /* EWALD ELECTROSTATICS */
1538             
1539             /* Analytical PME correction */
1540             zeta2            = _mm256_mul_ps(beta2,rsq30);
1541             rinv3            = _mm256_mul_ps(rinvsq30,rinv30);
1542             pmecorrF         = avx256_pmecorrF_f(zeta2);
1543             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1544             felec            = _mm256_mul_ps(qq30,felec);
1545             
1546             fscal            = felec;
1547
1548             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1549
1550             /* Calculate temporary vectorial force */
1551             tx               = _mm256_mul_ps(fscal,dx30);
1552             ty               = _mm256_mul_ps(fscal,dy30);
1553             tz               = _mm256_mul_ps(fscal,dz30);
1554
1555             /* Update vectorial force */
1556             fix3             = _mm256_add_ps(fix3,tx);
1557             fiy3             = _mm256_add_ps(fiy3,ty);
1558             fiz3             = _mm256_add_ps(fiz3,tz);
1559
1560             fjx0             = _mm256_add_ps(fjx0,tx);
1561             fjy0             = _mm256_add_ps(fjy0,ty);
1562             fjz0             = _mm256_add_ps(fjz0,tz);
1563
1564             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1565             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1566             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1567             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1568             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1569             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1570             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1571             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1572
1573             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
1574
1575             /* Inner loop uses 223 flops */
1576         }
1577
1578         /* End of innermost loop */
1579
1580         gmx_mm256_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1581                                                  f+i_coord_offset,fshift+i_shift_offset);
1582
1583         /* Increment number of inner iterations */
1584         inneriter                  += j_index_end - j_index_start;
1585
1586         /* Outer loop uses 24 flops */
1587     }
1588
1589     /* Increment number of outer iterations */
1590     outeriter        += nri;
1591
1592     /* Update outer/inner flops */
1593
1594     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*223);
1595 }