0049d02b6ae86db0378ca104a8a97920873a3408
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_single / nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_avx_256_single.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS avx_256_single kernel generator.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
46 #include "vec.h"
47 #include "nrnb.h"
48
49 #include "gmx_math_x86_avx_256_single.h"
50 #include "kernelutil_x86_avx_256_single.h"
51
52 /*
53  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_VF_avx_256_single
54  * Electrostatics interaction: ReactionField
55  * VdW interaction:            CubicSplineTable
56  * Geometry:                   Water3-Particle
57  * Calculate force/pot:        PotentialAndForce
58  */
59 void
60 nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_VF_avx_256_single
61                     (t_nblist                    * gmx_restrict       nlist,
62                      rvec                        * gmx_restrict          xx,
63                      rvec                        * gmx_restrict          ff,
64                      t_forcerec                  * gmx_restrict          fr,
65                      t_mdatoms                   * gmx_restrict     mdatoms,
66                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67                      t_nrnb                      * gmx_restrict        nrnb)
68 {
69     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
70      * just 0 for non-waters.
71      * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
72      * jnr indices corresponding to data put in the four positions in the SIMD register.
73      */
74     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
75     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
76     int              jnrA,jnrB,jnrC,jnrD;
77     int              jnrE,jnrF,jnrG,jnrH;
78     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
79     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
80     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
81     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
82     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
83     real             rcutoff_scalar;
84     real             *shiftvec,*fshift,*x,*f;
85     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
86     real             scratch[4*DIM];
87     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
88     real *           vdwioffsetptr0;
89     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
90     real *           vdwioffsetptr1;
91     __m256           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
92     real *           vdwioffsetptr2;
93     __m256           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
94     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
95     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
96     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
97     __m256           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
98     __m256           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
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     __m256           dummy_mask,cutoff_mask;
113     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
114     __m256           one     = _mm256_set1_ps(1.0);
115     __m256           two     = _mm256_set1_ps(2.0);
116     x                = xx[0];
117     f                = ff[0];
118
119     nri              = nlist->nri;
120     iinr             = nlist->iinr;
121     jindex           = nlist->jindex;
122     jjnr             = nlist->jjnr;
123     shiftidx         = nlist->shift;
124     gid              = nlist->gid;
125     shiftvec         = fr->shift_vec[0];
126     fshift           = fr->fshift[0];
127     facel            = _mm256_set1_ps(fr->epsfac);
128     charge           = mdatoms->chargeA;
129     krf              = _mm256_set1_ps(fr->ic->k_rf);
130     krf2             = _mm256_set1_ps(fr->ic->k_rf*2.0);
131     crf              = _mm256_set1_ps(fr->ic->c_rf);
132     nvdwtype         = fr->ntype;
133     vdwparam         = fr->nbfp;
134     vdwtype          = mdatoms->typeA;
135
136     vftab            = kernel_data->table_vdw->data;
137     vftabscale       = _mm256_set1_ps(kernel_data->table_vdw->scale);
138
139     /* Setup water-specific parameters */
140     inr              = nlist->iinr[0];
141     iq0              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
142     iq1              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
143     iq2              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
144     vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
145
146     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
147     rcutoff_scalar   = fr->rcoulomb;
148     rcutoff          = _mm256_set1_ps(rcutoff_scalar);
149     rcutoff2         = _mm256_mul_ps(rcutoff,rcutoff);
150
151     /* Avoid stupid compiler warnings */
152     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
153     j_coord_offsetA = 0;
154     j_coord_offsetB = 0;
155     j_coord_offsetC = 0;
156     j_coord_offsetD = 0;
157     j_coord_offsetE = 0;
158     j_coord_offsetF = 0;
159     j_coord_offsetG = 0;
160     j_coord_offsetH = 0;
161
162     outeriter        = 0;
163     inneriter        = 0;
164
165     for(iidx=0;iidx<4*DIM;iidx++)
166     {
167         scratch[iidx] = 0.0;
168     }
169
170     /* Start outer loop over neighborlists */
171     for(iidx=0; iidx<nri; iidx++)
172     {
173         /* Load shift vector for this list */
174         i_shift_offset   = DIM*shiftidx[iidx];
175
176         /* Load limits for loop over neighbors */
177         j_index_start    = jindex[iidx];
178         j_index_end      = jindex[iidx+1];
179
180         /* Get outer coordinate index */
181         inr              = iinr[iidx];
182         i_coord_offset   = DIM*inr;
183
184         /* Load i particle coords and add shift vector */
185         gmx_mm256_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
186                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
187
188         fix0             = _mm256_setzero_ps();
189         fiy0             = _mm256_setzero_ps();
190         fiz0             = _mm256_setzero_ps();
191         fix1             = _mm256_setzero_ps();
192         fiy1             = _mm256_setzero_ps();
193         fiz1             = _mm256_setzero_ps();
194         fix2             = _mm256_setzero_ps();
195         fiy2             = _mm256_setzero_ps();
196         fiz2             = _mm256_setzero_ps();
197
198         /* Reset potential sums */
199         velecsum         = _mm256_setzero_ps();
200         vvdwsum          = _mm256_setzero_ps();
201
202         /* Start inner kernel loop */
203         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
204         {
205
206             /* Get j neighbor index, and coordinate index */
207             jnrA             = jjnr[jidx];
208             jnrB             = jjnr[jidx+1];
209             jnrC             = jjnr[jidx+2];
210             jnrD             = jjnr[jidx+3];
211             jnrE             = jjnr[jidx+4];
212             jnrF             = jjnr[jidx+5];
213             jnrG             = jjnr[jidx+6];
214             jnrH             = jjnr[jidx+7];
215             j_coord_offsetA  = DIM*jnrA;
216             j_coord_offsetB  = DIM*jnrB;
217             j_coord_offsetC  = DIM*jnrC;
218             j_coord_offsetD  = DIM*jnrD;
219             j_coord_offsetE  = DIM*jnrE;
220             j_coord_offsetF  = DIM*jnrF;
221             j_coord_offsetG  = DIM*jnrG;
222             j_coord_offsetH  = DIM*jnrH;
223
224             /* load j atom coordinates */
225             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
226                                                  x+j_coord_offsetC,x+j_coord_offsetD,
227                                                  x+j_coord_offsetE,x+j_coord_offsetF,
228                                                  x+j_coord_offsetG,x+j_coord_offsetH,
229                                                  &jx0,&jy0,&jz0);
230
231             /* Calculate displacement vector */
232             dx00             = _mm256_sub_ps(ix0,jx0);
233             dy00             = _mm256_sub_ps(iy0,jy0);
234             dz00             = _mm256_sub_ps(iz0,jz0);
235             dx10             = _mm256_sub_ps(ix1,jx0);
236             dy10             = _mm256_sub_ps(iy1,jy0);
237             dz10             = _mm256_sub_ps(iz1,jz0);
238             dx20             = _mm256_sub_ps(ix2,jx0);
239             dy20             = _mm256_sub_ps(iy2,jy0);
240             dz20             = _mm256_sub_ps(iz2,jz0);
241
242             /* Calculate squared distance and things based on it */
243             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
244             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
245             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
246
247             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
248             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
249             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
250
251             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
252             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
253             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
254
255             /* Load parameters for j particles */
256             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
257                                                                  charge+jnrC+0,charge+jnrD+0,
258                                                                  charge+jnrE+0,charge+jnrF+0,
259                                                                  charge+jnrG+0,charge+jnrH+0);
260             vdwjidx0A        = 2*vdwtype[jnrA+0];
261             vdwjidx0B        = 2*vdwtype[jnrB+0];
262             vdwjidx0C        = 2*vdwtype[jnrC+0];
263             vdwjidx0D        = 2*vdwtype[jnrD+0];
264             vdwjidx0E        = 2*vdwtype[jnrE+0];
265             vdwjidx0F        = 2*vdwtype[jnrF+0];
266             vdwjidx0G        = 2*vdwtype[jnrG+0];
267             vdwjidx0H        = 2*vdwtype[jnrH+0];
268
269             fjx0             = _mm256_setzero_ps();
270             fjy0             = _mm256_setzero_ps();
271             fjz0             = _mm256_setzero_ps();
272
273             /**************************
274              * CALCULATE INTERACTIONS *
275              **************************/
276
277             if (gmx_mm256_any_lt(rsq00,rcutoff2))
278             {
279
280             r00              = _mm256_mul_ps(rsq00,rinv00);
281
282             /* Compute parameters for interactions between i and j atoms */
283             qq00             = _mm256_mul_ps(iq0,jq0);
284             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
285                                             vdwioffsetptr0+vdwjidx0B,
286                                             vdwioffsetptr0+vdwjidx0C,
287                                             vdwioffsetptr0+vdwjidx0D,
288                                             vdwioffsetptr0+vdwjidx0E,
289                                             vdwioffsetptr0+vdwjidx0F,
290                                             vdwioffsetptr0+vdwjidx0G,
291                                             vdwioffsetptr0+vdwjidx0H,
292                                             &c6_00,&c12_00);
293
294             /* Calculate table index by multiplying r with table scale and truncate to integer */
295             rt               = _mm256_mul_ps(r00,vftabscale);
296             vfitab           = _mm256_cvttps_epi32(rt);
297             vfeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
298             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
299             vfitab_lo        = _mm256_extractf128_si256(vfitab,0x0);
300             vfitab_hi        = _mm256_extractf128_si256(vfitab,0x1);
301             vfitab_lo        = _mm_slli_epi32(vfitab_lo,3);
302             vfitab_hi        = _mm_slli_epi32(vfitab_hi,3);
303
304             /* REACTION-FIELD ELECTROSTATICS */
305             velec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
306             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
307
308             /* CUBIC SPLINE TABLE DISPERSION */
309             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
310                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
311             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
312                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
313             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
314                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
315             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
316                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
317             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
318             Heps             = _mm256_mul_ps(vfeps,H);
319             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
320             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
321             vvdw6            = _mm256_mul_ps(c6_00,VV);
322             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
323             fvdw6            = _mm256_mul_ps(c6_00,FF);
324
325             /* CUBIC SPLINE TABLE REPULSION */
326             vfitab_lo        = _mm_add_epi32(vfitab_lo,ifour);
327             vfitab_hi        = _mm_add_epi32(vfitab_hi,ifour);
328             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
329                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
330             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
331                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
332             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
333                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
334             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
335                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
336             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
337             Heps             = _mm256_mul_ps(vfeps,H);
338             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
339             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
340             vvdw12           = _mm256_mul_ps(c12_00,VV);
341             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
342             fvdw12           = _mm256_mul_ps(c12_00,FF);
343             vvdw             = _mm256_add_ps(vvdw12,vvdw6);
344             fvdw             = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv00)));
345
346             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
347
348             /* Update potential sum for this i atom from the interaction with this j atom. */
349             velec            = _mm256_and_ps(velec,cutoff_mask);
350             velecsum         = _mm256_add_ps(velecsum,velec);
351             vvdw             = _mm256_and_ps(vvdw,cutoff_mask);
352             vvdwsum          = _mm256_add_ps(vvdwsum,vvdw);
353
354             fscal            = _mm256_add_ps(felec,fvdw);
355
356             fscal            = _mm256_and_ps(fscal,cutoff_mask);
357
358             /* Calculate temporary vectorial force */
359             tx               = _mm256_mul_ps(fscal,dx00);
360             ty               = _mm256_mul_ps(fscal,dy00);
361             tz               = _mm256_mul_ps(fscal,dz00);
362
363             /* Update vectorial force */
364             fix0             = _mm256_add_ps(fix0,tx);
365             fiy0             = _mm256_add_ps(fiy0,ty);
366             fiz0             = _mm256_add_ps(fiz0,tz);
367
368             fjx0             = _mm256_add_ps(fjx0,tx);
369             fjy0             = _mm256_add_ps(fjy0,ty);
370             fjz0             = _mm256_add_ps(fjz0,tz);
371
372             }
373
374             /**************************
375              * CALCULATE INTERACTIONS *
376              **************************/
377
378             if (gmx_mm256_any_lt(rsq10,rcutoff2))
379             {
380
381             /* Compute parameters for interactions between i and j atoms */
382             qq10             = _mm256_mul_ps(iq1,jq0);
383
384             /* REACTION-FIELD ELECTROSTATICS */
385             velec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_add_ps(rinv10,_mm256_mul_ps(krf,rsq10)),crf));
386             felec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_mul_ps(rinv10,rinvsq10),krf2));
387
388             cutoff_mask      = _mm256_cmp_ps(rsq10,rcutoff2,_CMP_LT_OQ);
389
390             /* Update potential sum for this i atom from the interaction with this j atom. */
391             velec            = _mm256_and_ps(velec,cutoff_mask);
392             velecsum         = _mm256_add_ps(velecsum,velec);
393
394             fscal            = felec;
395
396             fscal            = _mm256_and_ps(fscal,cutoff_mask);
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
414             /**************************
415              * CALCULATE INTERACTIONS *
416              **************************/
417
418             if (gmx_mm256_any_lt(rsq20,rcutoff2))
419             {
420
421             /* Compute parameters for interactions between i and j atoms */
422             qq20             = _mm256_mul_ps(iq2,jq0);
423
424             /* REACTION-FIELD ELECTROSTATICS */
425             velec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_add_ps(rinv20,_mm256_mul_ps(krf,rsq20)),crf));
426             felec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_mul_ps(rinv20,rinvsq20),krf2));
427
428             cutoff_mask      = _mm256_cmp_ps(rsq20,rcutoff2,_CMP_LT_OQ);
429
430             /* Update potential sum for this i atom from the interaction with this j atom. */
431             velec            = _mm256_and_ps(velec,cutoff_mask);
432             velecsum         = _mm256_add_ps(velecsum,velec);
433
434             fscal            = felec;
435
436             fscal            = _mm256_and_ps(fscal,cutoff_mask);
437
438             /* Calculate temporary vectorial force */
439             tx               = _mm256_mul_ps(fscal,dx20);
440             ty               = _mm256_mul_ps(fscal,dy20);
441             tz               = _mm256_mul_ps(fscal,dz20);
442
443             /* Update vectorial force */
444             fix2             = _mm256_add_ps(fix2,tx);
445             fiy2             = _mm256_add_ps(fiy2,ty);
446             fiz2             = _mm256_add_ps(fiz2,tz);
447
448             fjx0             = _mm256_add_ps(fjx0,tx);
449             fjy0             = _mm256_add_ps(fjy0,ty);
450             fjz0             = _mm256_add_ps(fjz0,tz);
451
452             }
453
454             fjptrA             = f+j_coord_offsetA;
455             fjptrB             = f+j_coord_offsetB;
456             fjptrC             = f+j_coord_offsetC;
457             fjptrD             = f+j_coord_offsetD;
458             fjptrE             = f+j_coord_offsetE;
459             fjptrF             = f+j_coord_offsetF;
460             fjptrG             = f+j_coord_offsetG;
461             fjptrH             = f+j_coord_offsetH;
462
463             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
464
465             /* Inner loop uses 147 flops */
466         }
467
468         if(jidx<j_index_end)
469         {
470
471             /* Get j neighbor index, and coordinate index */
472             jnrlistA         = jjnr[jidx];
473             jnrlistB         = jjnr[jidx+1];
474             jnrlistC         = jjnr[jidx+2];
475             jnrlistD         = jjnr[jidx+3];
476             jnrlistE         = jjnr[jidx+4];
477             jnrlistF         = jjnr[jidx+5];
478             jnrlistG         = jjnr[jidx+6];
479             jnrlistH         = jjnr[jidx+7];
480             /* Sign of each element will be negative for non-real atoms.
481              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
482              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
483              */
484             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
485                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
486                                             
487             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
488             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
489             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
490             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
491             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
492             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
493             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
494             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
495             j_coord_offsetA  = DIM*jnrA;
496             j_coord_offsetB  = DIM*jnrB;
497             j_coord_offsetC  = DIM*jnrC;
498             j_coord_offsetD  = DIM*jnrD;
499             j_coord_offsetE  = DIM*jnrE;
500             j_coord_offsetF  = DIM*jnrF;
501             j_coord_offsetG  = DIM*jnrG;
502             j_coord_offsetH  = DIM*jnrH;
503
504             /* load j atom coordinates */
505             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
506                                                  x+j_coord_offsetC,x+j_coord_offsetD,
507                                                  x+j_coord_offsetE,x+j_coord_offsetF,
508                                                  x+j_coord_offsetG,x+j_coord_offsetH,
509                                                  &jx0,&jy0,&jz0);
510
511             /* Calculate displacement vector */
512             dx00             = _mm256_sub_ps(ix0,jx0);
513             dy00             = _mm256_sub_ps(iy0,jy0);
514             dz00             = _mm256_sub_ps(iz0,jz0);
515             dx10             = _mm256_sub_ps(ix1,jx0);
516             dy10             = _mm256_sub_ps(iy1,jy0);
517             dz10             = _mm256_sub_ps(iz1,jz0);
518             dx20             = _mm256_sub_ps(ix2,jx0);
519             dy20             = _mm256_sub_ps(iy2,jy0);
520             dz20             = _mm256_sub_ps(iz2,jz0);
521
522             /* Calculate squared distance and things based on it */
523             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
524             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
525             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
526
527             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
528             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
529             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
530
531             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
532             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
533             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
534
535             /* Load parameters for j particles */
536             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
537                                                                  charge+jnrC+0,charge+jnrD+0,
538                                                                  charge+jnrE+0,charge+jnrF+0,
539                                                                  charge+jnrG+0,charge+jnrH+0);
540             vdwjidx0A        = 2*vdwtype[jnrA+0];
541             vdwjidx0B        = 2*vdwtype[jnrB+0];
542             vdwjidx0C        = 2*vdwtype[jnrC+0];
543             vdwjidx0D        = 2*vdwtype[jnrD+0];
544             vdwjidx0E        = 2*vdwtype[jnrE+0];
545             vdwjidx0F        = 2*vdwtype[jnrF+0];
546             vdwjidx0G        = 2*vdwtype[jnrG+0];
547             vdwjidx0H        = 2*vdwtype[jnrH+0];
548
549             fjx0             = _mm256_setzero_ps();
550             fjy0             = _mm256_setzero_ps();
551             fjz0             = _mm256_setzero_ps();
552
553             /**************************
554              * CALCULATE INTERACTIONS *
555              **************************/
556
557             if (gmx_mm256_any_lt(rsq00,rcutoff2))
558             {
559
560             r00              = _mm256_mul_ps(rsq00,rinv00);
561             r00              = _mm256_andnot_ps(dummy_mask,r00);
562
563             /* Compute parameters for interactions between i and j atoms */
564             qq00             = _mm256_mul_ps(iq0,jq0);
565             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
566                                             vdwioffsetptr0+vdwjidx0B,
567                                             vdwioffsetptr0+vdwjidx0C,
568                                             vdwioffsetptr0+vdwjidx0D,
569                                             vdwioffsetptr0+vdwjidx0E,
570                                             vdwioffsetptr0+vdwjidx0F,
571                                             vdwioffsetptr0+vdwjidx0G,
572                                             vdwioffsetptr0+vdwjidx0H,
573                                             &c6_00,&c12_00);
574
575             /* Calculate table index by multiplying r with table scale and truncate to integer */
576             rt               = _mm256_mul_ps(r00,vftabscale);
577             vfitab           = _mm256_cvttps_epi32(rt);
578             vfeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
579             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
580             vfitab_lo        = _mm256_extractf128_si256(vfitab,0x0);
581             vfitab_hi        = _mm256_extractf128_si256(vfitab,0x1);
582             vfitab_lo        = _mm_slli_epi32(vfitab_lo,3);
583             vfitab_hi        = _mm_slli_epi32(vfitab_hi,3);
584
585             /* REACTION-FIELD ELECTROSTATICS */
586             velec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
587             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
588
589             /* CUBIC SPLINE TABLE DISPERSION */
590             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
591                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
592             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
593                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
594             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
595                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
596             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
597                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
598             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
599             Heps             = _mm256_mul_ps(vfeps,H);
600             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
601             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
602             vvdw6            = _mm256_mul_ps(c6_00,VV);
603             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
604             fvdw6            = _mm256_mul_ps(c6_00,FF);
605
606             /* CUBIC SPLINE TABLE REPULSION */
607             vfitab_lo        = _mm_add_epi32(vfitab_lo,ifour);
608             vfitab_hi        = _mm_add_epi32(vfitab_hi,ifour);
609             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
610                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
611             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
612                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
613             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
614                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
615             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
616                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
617             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
618             Heps             = _mm256_mul_ps(vfeps,H);
619             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
620             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
621             vvdw12           = _mm256_mul_ps(c12_00,VV);
622             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
623             fvdw12           = _mm256_mul_ps(c12_00,FF);
624             vvdw             = _mm256_add_ps(vvdw12,vvdw6);
625             fvdw             = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv00)));
626
627             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
628
629             /* Update potential sum for this i atom from the interaction with this j atom. */
630             velec            = _mm256_and_ps(velec,cutoff_mask);
631             velec            = _mm256_andnot_ps(dummy_mask,velec);
632             velecsum         = _mm256_add_ps(velecsum,velec);
633             vvdw             = _mm256_and_ps(vvdw,cutoff_mask);
634             vvdw             = _mm256_andnot_ps(dummy_mask,vvdw);
635             vvdwsum          = _mm256_add_ps(vvdwsum,vvdw);
636
637             fscal            = _mm256_add_ps(felec,fvdw);
638
639             fscal            = _mm256_and_ps(fscal,cutoff_mask);
640
641             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
642
643             /* Calculate temporary vectorial force */
644             tx               = _mm256_mul_ps(fscal,dx00);
645             ty               = _mm256_mul_ps(fscal,dy00);
646             tz               = _mm256_mul_ps(fscal,dz00);
647
648             /* Update vectorial force */
649             fix0             = _mm256_add_ps(fix0,tx);
650             fiy0             = _mm256_add_ps(fiy0,ty);
651             fiz0             = _mm256_add_ps(fiz0,tz);
652
653             fjx0             = _mm256_add_ps(fjx0,tx);
654             fjy0             = _mm256_add_ps(fjy0,ty);
655             fjz0             = _mm256_add_ps(fjz0,tz);
656
657             }
658
659             /**************************
660              * CALCULATE INTERACTIONS *
661              **************************/
662
663             if (gmx_mm256_any_lt(rsq10,rcutoff2))
664             {
665
666             /* Compute parameters for interactions between i and j atoms */
667             qq10             = _mm256_mul_ps(iq1,jq0);
668
669             /* REACTION-FIELD ELECTROSTATICS */
670             velec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_add_ps(rinv10,_mm256_mul_ps(krf,rsq10)),crf));
671             felec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_mul_ps(rinv10,rinvsq10),krf2));
672
673             cutoff_mask      = _mm256_cmp_ps(rsq10,rcutoff2,_CMP_LT_OQ);
674
675             /* Update potential sum for this i atom from the interaction with this j atom. */
676             velec            = _mm256_and_ps(velec,cutoff_mask);
677             velec            = _mm256_andnot_ps(dummy_mask,velec);
678             velecsum         = _mm256_add_ps(velecsum,velec);
679
680             fscal            = felec;
681
682             fscal            = _mm256_and_ps(fscal,cutoff_mask);
683
684             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
685
686             /* Calculate temporary vectorial force */
687             tx               = _mm256_mul_ps(fscal,dx10);
688             ty               = _mm256_mul_ps(fscal,dy10);
689             tz               = _mm256_mul_ps(fscal,dz10);
690
691             /* Update vectorial force */
692             fix1             = _mm256_add_ps(fix1,tx);
693             fiy1             = _mm256_add_ps(fiy1,ty);
694             fiz1             = _mm256_add_ps(fiz1,tz);
695
696             fjx0             = _mm256_add_ps(fjx0,tx);
697             fjy0             = _mm256_add_ps(fjy0,ty);
698             fjz0             = _mm256_add_ps(fjz0,tz);
699
700             }
701
702             /**************************
703              * CALCULATE INTERACTIONS *
704              **************************/
705
706             if (gmx_mm256_any_lt(rsq20,rcutoff2))
707             {
708
709             /* Compute parameters for interactions between i and j atoms */
710             qq20             = _mm256_mul_ps(iq2,jq0);
711
712             /* REACTION-FIELD ELECTROSTATICS */
713             velec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_add_ps(rinv20,_mm256_mul_ps(krf,rsq20)),crf));
714             felec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_mul_ps(rinv20,rinvsq20),krf2));
715
716             cutoff_mask      = _mm256_cmp_ps(rsq20,rcutoff2,_CMP_LT_OQ);
717
718             /* Update potential sum for this i atom from the interaction with this j atom. */
719             velec            = _mm256_and_ps(velec,cutoff_mask);
720             velec            = _mm256_andnot_ps(dummy_mask,velec);
721             velecsum         = _mm256_add_ps(velecsum,velec);
722
723             fscal            = felec;
724
725             fscal            = _mm256_and_ps(fscal,cutoff_mask);
726
727             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
728
729             /* Calculate temporary vectorial force */
730             tx               = _mm256_mul_ps(fscal,dx20);
731             ty               = _mm256_mul_ps(fscal,dy20);
732             tz               = _mm256_mul_ps(fscal,dz20);
733
734             /* Update vectorial force */
735             fix2             = _mm256_add_ps(fix2,tx);
736             fiy2             = _mm256_add_ps(fiy2,ty);
737             fiz2             = _mm256_add_ps(fiz2,tz);
738
739             fjx0             = _mm256_add_ps(fjx0,tx);
740             fjy0             = _mm256_add_ps(fjy0,ty);
741             fjz0             = _mm256_add_ps(fjz0,tz);
742
743             }
744
745             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
746             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
747             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
748             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
749             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
750             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
751             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
752             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
753
754             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
755
756             /* Inner loop uses 148 flops */
757         }
758
759         /* End of innermost loop */
760
761         gmx_mm256_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
762                                                  f+i_coord_offset,fshift+i_shift_offset);
763
764         ggid                        = gid[iidx];
765         /* Update potential energies */
766         gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
767         gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
768
769         /* Increment number of inner iterations */
770         inneriter                  += j_index_end - j_index_start;
771
772         /* Outer loop uses 20 flops */
773     }
774
775     /* Increment number of outer iterations */
776     outeriter        += nri;
777
778     /* Update outer/inner flops */
779
780     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_VF,outeriter*20 + inneriter*148);
781 }
782 /*
783  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_F_avx_256_single
784  * Electrostatics interaction: ReactionField
785  * VdW interaction:            CubicSplineTable
786  * Geometry:                   Water3-Particle
787  * Calculate force/pot:        Force
788  */
789 void
790 nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_F_avx_256_single
791                     (t_nblist                    * gmx_restrict       nlist,
792                      rvec                        * gmx_restrict          xx,
793                      rvec                        * gmx_restrict          ff,
794                      t_forcerec                  * gmx_restrict          fr,
795                      t_mdatoms                   * gmx_restrict     mdatoms,
796                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
797                      t_nrnb                      * gmx_restrict        nrnb)
798 {
799     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
800      * just 0 for non-waters.
801      * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
802      * jnr indices corresponding to data put in the four positions in the SIMD register.
803      */
804     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
805     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
806     int              jnrA,jnrB,jnrC,jnrD;
807     int              jnrE,jnrF,jnrG,jnrH;
808     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
809     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
810     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
811     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
812     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
813     real             rcutoff_scalar;
814     real             *shiftvec,*fshift,*x,*f;
815     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
816     real             scratch[4*DIM];
817     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
818     real *           vdwioffsetptr0;
819     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
820     real *           vdwioffsetptr1;
821     __m256           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
822     real *           vdwioffsetptr2;
823     __m256           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
824     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
825     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
826     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
827     __m256           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
828     __m256           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
829     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
830     real             *charge;
831     int              nvdwtype;
832     __m256           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
833     int              *vdwtype;
834     real             *vdwparam;
835     __m256           one_sixth   = _mm256_set1_ps(1.0/6.0);
836     __m256           one_twelfth = _mm256_set1_ps(1.0/12.0);
837     __m256i          vfitab;
838     __m128i          vfitab_lo,vfitab_hi;
839     __m128i          ifour       = _mm_set1_epi32(4);
840     __m256           rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
841     real             *vftab;
842     __m256           dummy_mask,cutoff_mask;
843     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
844     __m256           one     = _mm256_set1_ps(1.0);
845     __m256           two     = _mm256_set1_ps(2.0);
846     x                = xx[0];
847     f                = ff[0];
848
849     nri              = nlist->nri;
850     iinr             = nlist->iinr;
851     jindex           = nlist->jindex;
852     jjnr             = nlist->jjnr;
853     shiftidx         = nlist->shift;
854     gid              = nlist->gid;
855     shiftvec         = fr->shift_vec[0];
856     fshift           = fr->fshift[0];
857     facel            = _mm256_set1_ps(fr->epsfac);
858     charge           = mdatoms->chargeA;
859     krf              = _mm256_set1_ps(fr->ic->k_rf);
860     krf2             = _mm256_set1_ps(fr->ic->k_rf*2.0);
861     crf              = _mm256_set1_ps(fr->ic->c_rf);
862     nvdwtype         = fr->ntype;
863     vdwparam         = fr->nbfp;
864     vdwtype          = mdatoms->typeA;
865
866     vftab            = kernel_data->table_vdw->data;
867     vftabscale       = _mm256_set1_ps(kernel_data->table_vdw->scale);
868
869     /* Setup water-specific parameters */
870     inr              = nlist->iinr[0];
871     iq0              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
872     iq1              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
873     iq2              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
874     vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
875
876     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
877     rcutoff_scalar   = fr->rcoulomb;
878     rcutoff          = _mm256_set1_ps(rcutoff_scalar);
879     rcutoff2         = _mm256_mul_ps(rcutoff,rcutoff);
880
881     /* Avoid stupid compiler warnings */
882     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
883     j_coord_offsetA = 0;
884     j_coord_offsetB = 0;
885     j_coord_offsetC = 0;
886     j_coord_offsetD = 0;
887     j_coord_offsetE = 0;
888     j_coord_offsetF = 0;
889     j_coord_offsetG = 0;
890     j_coord_offsetH = 0;
891
892     outeriter        = 0;
893     inneriter        = 0;
894
895     for(iidx=0;iidx<4*DIM;iidx++)
896     {
897         scratch[iidx] = 0.0;
898     }
899
900     /* Start outer loop over neighborlists */
901     for(iidx=0; iidx<nri; iidx++)
902     {
903         /* Load shift vector for this list */
904         i_shift_offset   = DIM*shiftidx[iidx];
905
906         /* Load limits for loop over neighbors */
907         j_index_start    = jindex[iidx];
908         j_index_end      = jindex[iidx+1];
909
910         /* Get outer coordinate index */
911         inr              = iinr[iidx];
912         i_coord_offset   = DIM*inr;
913
914         /* Load i particle coords and add shift vector */
915         gmx_mm256_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
916                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
917
918         fix0             = _mm256_setzero_ps();
919         fiy0             = _mm256_setzero_ps();
920         fiz0             = _mm256_setzero_ps();
921         fix1             = _mm256_setzero_ps();
922         fiy1             = _mm256_setzero_ps();
923         fiz1             = _mm256_setzero_ps();
924         fix2             = _mm256_setzero_ps();
925         fiy2             = _mm256_setzero_ps();
926         fiz2             = _mm256_setzero_ps();
927
928         /* Start inner kernel loop */
929         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
930         {
931
932             /* Get j neighbor index, and coordinate index */
933             jnrA             = jjnr[jidx];
934             jnrB             = jjnr[jidx+1];
935             jnrC             = jjnr[jidx+2];
936             jnrD             = jjnr[jidx+3];
937             jnrE             = jjnr[jidx+4];
938             jnrF             = jjnr[jidx+5];
939             jnrG             = jjnr[jidx+6];
940             jnrH             = jjnr[jidx+7];
941             j_coord_offsetA  = DIM*jnrA;
942             j_coord_offsetB  = DIM*jnrB;
943             j_coord_offsetC  = DIM*jnrC;
944             j_coord_offsetD  = DIM*jnrD;
945             j_coord_offsetE  = DIM*jnrE;
946             j_coord_offsetF  = DIM*jnrF;
947             j_coord_offsetG  = DIM*jnrG;
948             j_coord_offsetH  = DIM*jnrH;
949
950             /* load j atom coordinates */
951             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
952                                                  x+j_coord_offsetC,x+j_coord_offsetD,
953                                                  x+j_coord_offsetE,x+j_coord_offsetF,
954                                                  x+j_coord_offsetG,x+j_coord_offsetH,
955                                                  &jx0,&jy0,&jz0);
956
957             /* Calculate displacement vector */
958             dx00             = _mm256_sub_ps(ix0,jx0);
959             dy00             = _mm256_sub_ps(iy0,jy0);
960             dz00             = _mm256_sub_ps(iz0,jz0);
961             dx10             = _mm256_sub_ps(ix1,jx0);
962             dy10             = _mm256_sub_ps(iy1,jy0);
963             dz10             = _mm256_sub_ps(iz1,jz0);
964             dx20             = _mm256_sub_ps(ix2,jx0);
965             dy20             = _mm256_sub_ps(iy2,jy0);
966             dz20             = _mm256_sub_ps(iz2,jz0);
967
968             /* Calculate squared distance and things based on it */
969             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
970             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
971             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
972
973             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
974             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
975             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
976
977             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
978             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
979             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
980
981             /* Load parameters for j particles */
982             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
983                                                                  charge+jnrC+0,charge+jnrD+0,
984                                                                  charge+jnrE+0,charge+jnrF+0,
985                                                                  charge+jnrG+0,charge+jnrH+0);
986             vdwjidx0A        = 2*vdwtype[jnrA+0];
987             vdwjidx0B        = 2*vdwtype[jnrB+0];
988             vdwjidx0C        = 2*vdwtype[jnrC+0];
989             vdwjidx0D        = 2*vdwtype[jnrD+0];
990             vdwjidx0E        = 2*vdwtype[jnrE+0];
991             vdwjidx0F        = 2*vdwtype[jnrF+0];
992             vdwjidx0G        = 2*vdwtype[jnrG+0];
993             vdwjidx0H        = 2*vdwtype[jnrH+0];
994
995             fjx0             = _mm256_setzero_ps();
996             fjy0             = _mm256_setzero_ps();
997             fjz0             = _mm256_setzero_ps();
998
999             /**************************
1000              * CALCULATE INTERACTIONS *
1001              **************************/
1002
1003             if (gmx_mm256_any_lt(rsq00,rcutoff2))
1004             {
1005
1006             r00              = _mm256_mul_ps(rsq00,rinv00);
1007
1008             /* Compute parameters for interactions between i and j atoms */
1009             qq00             = _mm256_mul_ps(iq0,jq0);
1010             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
1011                                             vdwioffsetptr0+vdwjidx0B,
1012                                             vdwioffsetptr0+vdwjidx0C,
1013                                             vdwioffsetptr0+vdwjidx0D,
1014                                             vdwioffsetptr0+vdwjidx0E,
1015                                             vdwioffsetptr0+vdwjidx0F,
1016                                             vdwioffsetptr0+vdwjidx0G,
1017                                             vdwioffsetptr0+vdwjidx0H,
1018                                             &c6_00,&c12_00);
1019
1020             /* Calculate table index by multiplying r with table scale and truncate to integer */
1021             rt               = _mm256_mul_ps(r00,vftabscale);
1022             vfitab           = _mm256_cvttps_epi32(rt);
1023             vfeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
1024             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
1025             vfitab_lo        = _mm256_extractf128_si256(vfitab,0x0);
1026             vfitab_hi        = _mm256_extractf128_si256(vfitab,0x1);
1027             vfitab_lo        = _mm_slli_epi32(vfitab_lo,3);
1028             vfitab_hi        = _mm_slli_epi32(vfitab_hi,3);
1029
1030             /* REACTION-FIELD ELECTROSTATICS */
1031             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
1032
1033             /* CUBIC SPLINE TABLE DISPERSION */
1034             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
1035                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
1036             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
1037                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
1038             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
1039                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
1040             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
1041                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
1042             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
1043             Heps             = _mm256_mul_ps(vfeps,H);
1044             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
1045             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
1046             fvdw6            = _mm256_mul_ps(c6_00,FF);
1047
1048             /* CUBIC SPLINE TABLE REPULSION */
1049             vfitab_lo        = _mm_add_epi32(vfitab_lo,ifour);
1050             vfitab_hi        = _mm_add_epi32(vfitab_hi,ifour);
1051             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
1052                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
1053             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
1054                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
1055             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
1056                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
1057             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
1058                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
1059             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
1060             Heps             = _mm256_mul_ps(vfeps,H);
1061             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
1062             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
1063             fvdw12           = _mm256_mul_ps(c12_00,FF);
1064             fvdw             = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv00)));
1065
1066             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
1067
1068             fscal            = _mm256_add_ps(felec,fvdw);
1069
1070             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1071
1072             /* Calculate temporary vectorial force */
1073             tx               = _mm256_mul_ps(fscal,dx00);
1074             ty               = _mm256_mul_ps(fscal,dy00);
1075             tz               = _mm256_mul_ps(fscal,dz00);
1076
1077             /* Update vectorial force */
1078             fix0             = _mm256_add_ps(fix0,tx);
1079             fiy0             = _mm256_add_ps(fiy0,ty);
1080             fiz0             = _mm256_add_ps(fiz0,tz);
1081
1082             fjx0             = _mm256_add_ps(fjx0,tx);
1083             fjy0             = _mm256_add_ps(fjy0,ty);
1084             fjz0             = _mm256_add_ps(fjz0,tz);
1085
1086             }
1087
1088             /**************************
1089              * CALCULATE INTERACTIONS *
1090              **************************/
1091
1092             if (gmx_mm256_any_lt(rsq10,rcutoff2))
1093             {
1094
1095             /* Compute parameters for interactions between i and j atoms */
1096             qq10             = _mm256_mul_ps(iq1,jq0);
1097
1098             /* REACTION-FIELD ELECTROSTATICS */
1099             felec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_mul_ps(rinv10,rinvsq10),krf2));
1100
1101             cutoff_mask      = _mm256_cmp_ps(rsq10,rcutoff2,_CMP_LT_OQ);
1102
1103             fscal            = felec;
1104
1105             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1106
1107             /* Calculate temporary vectorial force */
1108             tx               = _mm256_mul_ps(fscal,dx10);
1109             ty               = _mm256_mul_ps(fscal,dy10);
1110             tz               = _mm256_mul_ps(fscal,dz10);
1111
1112             /* Update vectorial force */
1113             fix1             = _mm256_add_ps(fix1,tx);
1114             fiy1             = _mm256_add_ps(fiy1,ty);
1115             fiz1             = _mm256_add_ps(fiz1,tz);
1116
1117             fjx0             = _mm256_add_ps(fjx0,tx);
1118             fjy0             = _mm256_add_ps(fjy0,ty);
1119             fjz0             = _mm256_add_ps(fjz0,tz);
1120
1121             }
1122
1123             /**************************
1124              * CALCULATE INTERACTIONS *
1125              **************************/
1126
1127             if (gmx_mm256_any_lt(rsq20,rcutoff2))
1128             {
1129
1130             /* Compute parameters for interactions between i and j atoms */
1131             qq20             = _mm256_mul_ps(iq2,jq0);
1132
1133             /* REACTION-FIELD ELECTROSTATICS */
1134             felec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_mul_ps(rinv20,rinvsq20),krf2));
1135
1136             cutoff_mask      = _mm256_cmp_ps(rsq20,rcutoff2,_CMP_LT_OQ);
1137
1138             fscal            = felec;
1139
1140             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1141
1142             /* Calculate temporary vectorial force */
1143             tx               = _mm256_mul_ps(fscal,dx20);
1144             ty               = _mm256_mul_ps(fscal,dy20);
1145             tz               = _mm256_mul_ps(fscal,dz20);
1146
1147             /* Update vectorial force */
1148             fix2             = _mm256_add_ps(fix2,tx);
1149             fiy2             = _mm256_add_ps(fiy2,ty);
1150             fiz2             = _mm256_add_ps(fiz2,tz);
1151
1152             fjx0             = _mm256_add_ps(fjx0,tx);
1153             fjy0             = _mm256_add_ps(fjy0,ty);
1154             fjz0             = _mm256_add_ps(fjz0,tz);
1155
1156             }
1157
1158             fjptrA             = f+j_coord_offsetA;
1159             fjptrB             = f+j_coord_offsetB;
1160             fjptrC             = f+j_coord_offsetC;
1161             fjptrD             = f+j_coord_offsetD;
1162             fjptrE             = f+j_coord_offsetE;
1163             fjptrF             = f+j_coord_offsetF;
1164             fjptrG             = f+j_coord_offsetG;
1165             fjptrH             = f+j_coord_offsetH;
1166
1167             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
1168
1169             /* Inner loop uses 120 flops */
1170         }
1171
1172         if(jidx<j_index_end)
1173         {
1174
1175             /* Get j neighbor index, and coordinate index */
1176             jnrlistA         = jjnr[jidx];
1177             jnrlistB         = jjnr[jidx+1];
1178             jnrlistC         = jjnr[jidx+2];
1179             jnrlistD         = jjnr[jidx+3];
1180             jnrlistE         = jjnr[jidx+4];
1181             jnrlistF         = jjnr[jidx+5];
1182             jnrlistG         = jjnr[jidx+6];
1183             jnrlistH         = jjnr[jidx+7];
1184             /* Sign of each element will be negative for non-real atoms.
1185              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1186              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1187              */
1188             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
1189                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
1190                                             
1191             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
1192             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
1193             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
1194             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
1195             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
1196             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
1197             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
1198             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
1199             j_coord_offsetA  = DIM*jnrA;
1200             j_coord_offsetB  = DIM*jnrB;
1201             j_coord_offsetC  = DIM*jnrC;
1202             j_coord_offsetD  = DIM*jnrD;
1203             j_coord_offsetE  = DIM*jnrE;
1204             j_coord_offsetF  = DIM*jnrF;
1205             j_coord_offsetG  = DIM*jnrG;
1206             j_coord_offsetH  = DIM*jnrH;
1207
1208             /* load j atom coordinates */
1209             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1210                                                  x+j_coord_offsetC,x+j_coord_offsetD,
1211                                                  x+j_coord_offsetE,x+j_coord_offsetF,
1212                                                  x+j_coord_offsetG,x+j_coord_offsetH,
1213                                                  &jx0,&jy0,&jz0);
1214
1215             /* Calculate displacement vector */
1216             dx00             = _mm256_sub_ps(ix0,jx0);
1217             dy00             = _mm256_sub_ps(iy0,jy0);
1218             dz00             = _mm256_sub_ps(iz0,jz0);
1219             dx10             = _mm256_sub_ps(ix1,jx0);
1220             dy10             = _mm256_sub_ps(iy1,jy0);
1221             dz10             = _mm256_sub_ps(iz1,jz0);
1222             dx20             = _mm256_sub_ps(ix2,jx0);
1223             dy20             = _mm256_sub_ps(iy2,jy0);
1224             dz20             = _mm256_sub_ps(iz2,jz0);
1225
1226             /* Calculate squared distance and things based on it */
1227             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
1228             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
1229             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
1230
1231             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
1232             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
1233             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
1234
1235             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
1236             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
1237             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
1238
1239             /* Load parameters for j particles */
1240             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1241                                                                  charge+jnrC+0,charge+jnrD+0,
1242                                                                  charge+jnrE+0,charge+jnrF+0,
1243                                                                  charge+jnrG+0,charge+jnrH+0);
1244             vdwjidx0A        = 2*vdwtype[jnrA+0];
1245             vdwjidx0B        = 2*vdwtype[jnrB+0];
1246             vdwjidx0C        = 2*vdwtype[jnrC+0];
1247             vdwjidx0D        = 2*vdwtype[jnrD+0];
1248             vdwjidx0E        = 2*vdwtype[jnrE+0];
1249             vdwjidx0F        = 2*vdwtype[jnrF+0];
1250             vdwjidx0G        = 2*vdwtype[jnrG+0];
1251             vdwjidx0H        = 2*vdwtype[jnrH+0];
1252
1253             fjx0             = _mm256_setzero_ps();
1254             fjy0             = _mm256_setzero_ps();
1255             fjz0             = _mm256_setzero_ps();
1256
1257             /**************************
1258              * CALCULATE INTERACTIONS *
1259              **************************/
1260
1261             if (gmx_mm256_any_lt(rsq00,rcutoff2))
1262             {
1263
1264             r00              = _mm256_mul_ps(rsq00,rinv00);
1265             r00              = _mm256_andnot_ps(dummy_mask,r00);
1266
1267             /* Compute parameters for interactions between i and j atoms */
1268             qq00             = _mm256_mul_ps(iq0,jq0);
1269             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
1270                                             vdwioffsetptr0+vdwjidx0B,
1271                                             vdwioffsetptr0+vdwjidx0C,
1272                                             vdwioffsetptr0+vdwjidx0D,
1273                                             vdwioffsetptr0+vdwjidx0E,
1274                                             vdwioffsetptr0+vdwjidx0F,
1275                                             vdwioffsetptr0+vdwjidx0G,
1276                                             vdwioffsetptr0+vdwjidx0H,
1277                                             &c6_00,&c12_00);
1278
1279             /* Calculate table index by multiplying r with table scale and truncate to integer */
1280             rt               = _mm256_mul_ps(r00,vftabscale);
1281             vfitab           = _mm256_cvttps_epi32(rt);
1282             vfeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
1283             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
1284             vfitab_lo        = _mm256_extractf128_si256(vfitab,0x0);
1285             vfitab_hi        = _mm256_extractf128_si256(vfitab,0x1);
1286             vfitab_lo        = _mm_slli_epi32(vfitab_lo,3);
1287             vfitab_hi        = _mm_slli_epi32(vfitab_hi,3);
1288
1289             /* REACTION-FIELD ELECTROSTATICS */
1290             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
1291
1292             /* CUBIC SPLINE TABLE DISPERSION */
1293             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
1294                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
1295             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
1296                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
1297             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
1298                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
1299             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
1300                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
1301             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
1302             Heps             = _mm256_mul_ps(vfeps,H);
1303             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
1304             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
1305             fvdw6            = _mm256_mul_ps(c6_00,FF);
1306
1307             /* CUBIC SPLINE TABLE REPULSION */
1308             vfitab_lo        = _mm_add_epi32(vfitab_lo,ifour);
1309             vfitab_hi        = _mm_add_epi32(vfitab_hi,ifour);
1310             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
1311                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
1312             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
1313                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
1314             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
1315                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
1316             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
1317                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
1318             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
1319             Heps             = _mm256_mul_ps(vfeps,H);
1320             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
1321             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
1322             fvdw12           = _mm256_mul_ps(c12_00,FF);
1323             fvdw             = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv00)));
1324
1325             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
1326
1327             fscal            = _mm256_add_ps(felec,fvdw);
1328
1329             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1330
1331             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1332
1333             /* Calculate temporary vectorial force */
1334             tx               = _mm256_mul_ps(fscal,dx00);
1335             ty               = _mm256_mul_ps(fscal,dy00);
1336             tz               = _mm256_mul_ps(fscal,dz00);
1337
1338             /* Update vectorial force */
1339             fix0             = _mm256_add_ps(fix0,tx);
1340             fiy0             = _mm256_add_ps(fiy0,ty);
1341             fiz0             = _mm256_add_ps(fiz0,tz);
1342
1343             fjx0             = _mm256_add_ps(fjx0,tx);
1344             fjy0             = _mm256_add_ps(fjy0,ty);
1345             fjz0             = _mm256_add_ps(fjz0,tz);
1346
1347             }
1348
1349             /**************************
1350              * CALCULATE INTERACTIONS *
1351              **************************/
1352
1353             if (gmx_mm256_any_lt(rsq10,rcutoff2))
1354             {
1355
1356             /* Compute parameters for interactions between i and j atoms */
1357             qq10             = _mm256_mul_ps(iq1,jq0);
1358
1359             /* REACTION-FIELD ELECTROSTATICS */
1360             felec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_mul_ps(rinv10,rinvsq10),krf2));
1361
1362             cutoff_mask      = _mm256_cmp_ps(rsq10,rcutoff2,_CMP_LT_OQ);
1363
1364             fscal            = felec;
1365
1366             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1367
1368             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1369
1370             /* Calculate temporary vectorial force */
1371             tx               = _mm256_mul_ps(fscal,dx10);
1372             ty               = _mm256_mul_ps(fscal,dy10);
1373             tz               = _mm256_mul_ps(fscal,dz10);
1374
1375             /* Update vectorial force */
1376             fix1             = _mm256_add_ps(fix1,tx);
1377             fiy1             = _mm256_add_ps(fiy1,ty);
1378             fiz1             = _mm256_add_ps(fiz1,tz);
1379
1380             fjx0             = _mm256_add_ps(fjx0,tx);
1381             fjy0             = _mm256_add_ps(fjy0,ty);
1382             fjz0             = _mm256_add_ps(fjz0,tz);
1383
1384             }
1385
1386             /**************************
1387              * CALCULATE INTERACTIONS *
1388              **************************/
1389
1390             if (gmx_mm256_any_lt(rsq20,rcutoff2))
1391             {
1392
1393             /* Compute parameters for interactions between i and j atoms */
1394             qq20             = _mm256_mul_ps(iq2,jq0);
1395
1396             /* REACTION-FIELD ELECTROSTATICS */
1397             felec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_mul_ps(rinv20,rinvsq20),krf2));
1398
1399             cutoff_mask      = _mm256_cmp_ps(rsq20,rcutoff2,_CMP_LT_OQ);
1400
1401             fscal            = felec;
1402
1403             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1404
1405             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1406
1407             /* Calculate temporary vectorial force */
1408             tx               = _mm256_mul_ps(fscal,dx20);
1409             ty               = _mm256_mul_ps(fscal,dy20);
1410             tz               = _mm256_mul_ps(fscal,dz20);
1411
1412             /* Update vectorial force */
1413             fix2             = _mm256_add_ps(fix2,tx);
1414             fiy2             = _mm256_add_ps(fiy2,ty);
1415             fiz2             = _mm256_add_ps(fiz2,tz);
1416
1417             fjx0             = _mm256_add_ps(fjx0,tx);
1418             fjy0             = _mm256_add_ps(fjy0,ty);
1419             fjz0             = _mm256_add_ps(fjz0,tz);
1420
1421             }
1422
1423             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1424             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1425             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1426             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1427             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1428             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1429             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1430             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1431
1432             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
1433
1434             /* Inner loop uses 121 flops */
1435         }
1436
1437         /* End of innermost loop */
1438
1439         gmx_mm256_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1440                                                  f+i_coord_offset,fshift+i_shift_offset);
1441
1442         /* Increment number of inner iterations */
1443         inneriter                  += j_index_end - j_index_start;
1444
1445         /* Outer loop uses 18 flops */
1446     }
1447
1448     /* Increment number of outer iterations */
1449     outeriter        += nri;
1450
1451     /* Update outer/inner flops */
1452
1453     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_F,outeriter*18 + inneriter*121);
1454 }