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