Remove all unnecessary HAVE_CONFIG_H
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_single / nb_kernel_ElecRF_VdwCSTab_GeomP1P1_avx_256_single.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS avx_256_single kernel generator.
37  */
38 #include "config.h"
39
40 #include <math.h>
41
42 #include "../nb_kernel.h"
43 #include "types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "nrnb.h"
46
47 #include "gromacs/simd/math_x86_avx_256_single.h"
48 #include "kernelutil_x86_avx_256_single.h"
49
50 /*
51  * Gromacs nonbonded kernel:   nb_kernel_ElecRF_VdwCSTab_GeomP1P1_VF_avx_256_single
52  * Electrostatics interaction: ReactionField
53  * VdW interaction:            CubicSplineTable
54  * Geometry:                   Particle-Particle
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecRF_VdwCSTab_GeomP1P1_VF_avx_256_single
59                     (t_nblist                    * gmx_restrict       nlist,
60                      rvec                        * gmx_restrict          xx,
61                      rvec                        * gmx_restrict          ff,
62                      t_forcerec                  * gmx_restrict          fr,
63                      t_mdatoms                   * gmx_restrict     mdatoms,
64                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
65                      t_nrnb                      * gmx_restrict        nrnb)
66 {
67     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
68      * just 0 for non-waters.
69      * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
70      * jnr indices corresponding to data put in the four positions in the SIMD register.
71      */
72     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
73     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
74     int              jnrA,jnrB,jnrC,jnrD;
75     int              jnrE,jnrF,jnrG,jnrH;
76     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
77     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
78     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
79     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
80     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
81     real             rcutoff_scalar;
82     real             *shiftvec,*fshift,*x,*f;
83     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
84     real             scratch[4*DIM];
85     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
86     real *           vdwioffsetptr0;
87     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
88     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
89     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
90     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
91     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
92     real             *charge;
93     int              nvdwtype;
94     __m256           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
95     int              *vdwtype;
96     real             *vdwparam;
97     __m256           one_sixth   = _mm256_set1_ps(1.0/6.0);
98     __m256           one_twelfth = _mm256_set1_ps(1.0/12.0);
99     __m256i          vfitab;
100     __m128i          vfitab_lo,vfitab_hi;
101     __m128i          ifour       = _mm_set1_epi32(4);
102     __m256           rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
103     real             *vftab;
104     __m256           dummy_mask,cutoff_mask;
105     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
106     __m256           one     = _mm256_set1_ps(1.0);
107     __m256           two     = _mm256_set1_ps(2.0);
108     x                = xx[0];
109     f                = ff[0];
110
111     nri              = nlist->nri;
112     iinr             = nlist->iinr;
113     jindex           = nlist->jindex;
114     jjnr             = nlist->jjnr;
115     shiftidx         = nlist->shift;
116     gid              = nlist->gid;
117     shiftvec         = fr->shift_vec[0];
118     fshift           = fr->fshift[0];
119     facel            = _mm256_set1_ps(fr->epsfac);
120     charge           = mdatoms->chargeA;
121     krf              = _mm256_set1_ps(fr->ic->k_rf);
122     krf2             = _mm256_set1_ps(fr->ic->k_rf*2.0);
123     crf              = _mm256_set1_ps(fr->ic->c_rf);
124     nvdwtype         = fr->ntype;
125     vdwparam         = fr->nbfp;
126     vdwtype          = mdatoms->typeA;
127
128     vftab            = kernel_data->table_vdw->data;
129     vftabscale       = _mm256_set1_ps(kernel_data->table_vdw->scale);
130
131     /* Avoid stupid compiler warnings */
132     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
133     j_coord_offsetA = 0;
134     j_coord_offsetB = 0;
135     j_coord_offsetC = 0;
136     j_coord_offsetD = 0;
137     j_coord_offsetE = 0;
138     j_coord_offsetF = 0;
139     j_coord_offsetG = 0;
140     j_coord_offsetH = 0;
141
142     outeriter        = 0;
143     inneriter        = 0;
144
145     for(iidx=0;iidx<4*DIM;iidx++)
146     {
147         scratch[iidx] = 0.0;
148     }
149
150     /* Start outer loop over neighborlists */
151     for(iidx=0; iidx<nri; iidx++)
152     {
153         /* Load shift vector for this list */
154         i_shift_offset   = DIM*shiftidx[iidx];
155
156         /* Load limits for loop over neighbors */
157         j_index_start    = jindex[iidx];
158         j_index_end      = jindex[iidx+1];
159
160         /* Get outer coordinate index */
161         inr              = iinr[iidx];
162         i_coord_offset   = DIM*inr;
163
164         /* Load i particle coords and add shift vector */
165         gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
166
167         fix0             = _mm256_setzero_ps();
168         fiy0             = _mm256_setzero_ps();
169         fiz0             = _mm256_setzero_ps();
170
171         /* Load parameters for i particles */
172         iq0              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
173         vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
174
175         /* Reset potential sums */
176         velecsum         = _mm256_setzero_ps();
177         vvdwsum          = _mm256_setzero_ps();
178
179         /* Start inner kernel loop */
180         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
181         {
182
183             /* Get j neighbor index, and coordinate index */
184             jnrA             = jjnr[jidx];
185             jnrB             = jjnr[jidx+1];
186             jnrC             = jjnr[jidx+2];
187             jnrD             = jjnr[jidx+3];
188             jnrE             = jjnr[jidx+4];
189             jnrF             = jjnr[jidx+5];
190             jnrG             = jjnr[jidx+6];
191             jnrH             = jjnr[jidx+7];
192             j_coord_offsetA  = DIM*jnrA;
193             j_coord_offsetB  = DIM*jnrB;
194             j_coord_offsetC  = DIM*jnrC;
195             j_coord_offsetD  = DIM*jnrD;
196             j_coord_offsetE  = DIM*jnrE;
197             j_coord_offsetF  = DIM*jnrF;
198             j_coord_offsetG  = DIM*jnrG;
199             j_coord_offsetH  = DIM*jnrH;
200
201             /* load j atom coordinates */
202             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
203                                                  x+j_coord_offsetC,x+j_coord_offsetD,
204                                                  x+j_coord_offsetE,x+j_coord_offsetF,
205                                                  x+j_coord_offsetG,x+j_coord_offsetH,
206                                                  &jx0,&jy0,&jz0);
207
208             /* Calculate displacement vector */
209             dx00             = _mm256_sub_ps(ix0,jx0);
210             dy00             = _mm256_sub_ps(iy0,jy0);
211             dz00             = _mm256_sub_ps(iz0,jz0);
212
213             /* Calculate squared distance and things based on it */
214             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
215
216             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
217
218             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
219
220             /* Load parameters for j particles */
221             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
222                                                                  charge+jnrC+0,charge+jnrD+0,
223                                                                  charge+jnrE+0,charge+jnrF+0,
224                                                                  charge+jnrG+0,charge+jnrH+0);
225             vdwjidx0A        = 2*vdwtype[jnrA+0];
226             vdwjidx0B        = 2*vdwtype[jnrB+0];
227             vdwjidx0C        = 2*vdwtype[jnrC+0];
228             vdwjidx0D        = 2*vdwtype[jnrD+0];
229             vdwjidx0E        = 2*vdwtype[jnrE+0];
230             vdwjidx0F        = 2*vdwtype[jnrF+0];
231             vdwjidx0G        = 2*vdwtype[jnrG+0];
232             vdwjidx0H        = 2*vdwtype[jnrH+0];
233
234             /**************************
235              * CALCULATE INTERACTIONS *
236              **************************/
237
238             r00              = _mm256_mul_ps(rsq00,rinv00);
239
240             /* Compute parameters for interactions between i and j atoms */
241             qq00             = _mm256_mul_ps(iq0,jq0);
242             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
243                                             vdwioffsetptr0+vdwjidx0B,
244                                             vdwioffsetptr0+vdwjidx0C,
245                                             vdwioffsetptr0+vdwjidx0D,
246                                             vdwioffsetptr0+vdwjidx0E,
247                                             vdwioffsetptr0+vdwjidx0F,
248                                             vdwioffsetptr0+vdwjidx0G,
249                                             vdwioffsetptr0+vdwjidx0H,
250                                             &c6_00,&c12_00);
251
252             /* Calculate table index by multiplying r with table scale and truncate to integer */
253             rt               = _mm256_mul_ps(r00,vftabscale);
254             vfitab           = _mm256_cvttps_epi32(rt);
255             vfeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
256             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
257             vfitab_lo        = _mm256_extractf128_si256(vfitab,0x0);
258             vfitab_hi        = _mm256_extractf128_si256(vfitab,0x1);
259             vfitab_lo        = _mm_slli_epi32(vfitab_lo,3);
260             vfitab_hi        = _mm_slli_epi32(vfitab_hi,3);
261
262             /* REACTION-FIELD ELECTROSTATICS */
263             velec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
264             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
265
266             /* CUBIC SPLINE TABLE DISPERSION */
267             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
268                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
269             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
270                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
271             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
272                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
273             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
274                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
275             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
276             Heps             = _mm256_mul_ps(vfeps,H);
277             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
278             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
279             vvdw6            = _mm256_mul_ps(c6_00,VV);
280             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
281             fvdw6            = _mm256_mul_ps(c6_00,FF);
282
283             /* CUBIC SPLINE TABLE REPULSION */
284             vfitab_lo        = _mm_add_epi32(vfitab_lo,ifour);
285             vfitab_hi        = _mm_add_epi32(vfitab_hi,ifour);
286             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
287                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
288             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
289                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
290             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
291                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
292             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
293                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
294             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
295             Heps             = _mm256_mul_ps(vfeps,H);
296             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
297             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
298             vvdw12           = _mm256_mul_ps(c12_00,VV);
299             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
300             fvdw12           = _mm256_mul_ps(c12_00,FF);
301             vvdw             = _mm256_add_ps(vvdw12,vvdw6);
302             fvdw             = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv00)));
303
304             /* Update potential sum for this i atom from the interaction with this j atom. */
305             velecsum         = _mm256_add_ps(velecsum,velec);
306             vvdwsum          = _mm256_add_ps(vvdwsum,vvdw);
307
308             fscal            = _mm256_add_ps(felec,fvdw);
309
310             /* Calculate temporary vectorial force */
311             tx               = _mm256_mul_ps(fscal,dx00);
312             ty               = _mm256_mul_ps(fscal,dy00);
313             tz               = _mm256_mul_ps(fscal,dz00);
314
315             /* Update vectorial force */
316             fix0             = _mm256_add_ps(fix0,tx);
317             fiy0             = _mm256_add_ps(fiy0,ty);
318             fiz0             = _mm256_add_ps(fiz0,tz);
319
320             fjptrA             = f+j_coord_offsetA;
321             fjptrB             = f+j_coord_offsetB;
322             fjptrC             = f+j_coord_offsetC;
323             fjptrD             = f+j_coord_offsetD;
324             fjptrE             = f+j_coord_offsetE;
325             fjptrF             = f+j_coord_offsetF;
326             fjptrG             = f+j_coord_offsetG;
327             fjptrH             = f+j_coord_offsetH;
328             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
329
330             /* Inner loop uses 67 flops */
331         }
332
333         if(jidx<j_index_end)
334         {
335
336             /* Get j neighbor index, and coordinate index */
337             jnrlistA         = jjnr[jidx];
338             jnrlistB         = jjnr[jidx+1];
339             jnrlistC         = jjnr[jidx+2];
340             jnrlistD         = jjnr[jidx+3];
341             jnrlistE         = jjnr[jidx+4];
342             jnrlistF         = jjnr[jidx+5];
343             jnrlistG         = jjnr[jidx+6];
344             jnrlistH         = jjnr[jidx+7];
345             /* Sign of each element will be negative for non-real atoms.
346              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
347              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
348              */
349             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
350                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
351                                             
352             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
353             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
354             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
355             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
356             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
357             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
358             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
359             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
360             j_coord_offsetA  = DIM*jnrA;
361             j_coord_offsetB  = DIM*jnrB;
362             j_coord_offsetC  = DIM*jnrC;
363             j_coord_offsetD  = DIM*jnrD;
364             j_coord_offsetE  = DIM*jnrE;
365             j_coord_offsetF  = DIM*jnrF;
366             j_coord_offsetG  = DIM*jnrG;
367             j_coord_offsetH  = DIM*jnrH;
368
369             /* load j atom coordinates */
370             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
371                                                  x+j_coord_offsetC,x+j_coord_offsetD,
372                                                  x+j_coord_offsetE,x+j_coord_offsetF,
373                                                  x+j_coord_offsetG,x+j_coord_offsetH,
374                                                  &jx0,&jy0,&jz0);
375
376             /* Calculate displacement vector */
377             dx00             = _mm256_sub_ps(ix0,jx0);
378             dy00             = _mm256_sub_ps(iy0,jy0);
379             dz00             = _mm256_sub_ps(iz0,jz0);
380
381             /* Calculate squared distance and things based on it */
382             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
383
384             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
385
386             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
387
388             /* Load parameters for j particles */
389             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
390                                                                  charge+jnrC+0,charge+jnrD+0,
391                                                                  charge+jnrE+0,charge+jnrF+0,
392                                                                  charge+jnrG+0,charge+jnrH+0);
393             vdwjidx0A        = 2*vdwtype[jnrA+0];
394             vdwjidx0B        = 2*vdwtype[jnrB+0];
395             vdwjidx0C        = 2*vdwtype[jnrC+0];
396             vdwjidx0D        = 2*vdwtype[jnrD+0];
397             vdwjidx0E        = 2*vdwtype[jnrE+0];
398             vdwjidx0F        = 2*vdwtype[jnrF+0];
399             vdwjidx0G        = 2*vdwtype[jnrG+0];
400             vdwjidx0H        = 2*vdwtype[jnrH+0];
401
402             /**************************
403              * CALCULATE INTERACTIONS *
404              **************************/
405
406             r00              = _mm256_mul_ps(rsq00,rinv00);
407             r00              = _mm256_andnot_ps(dummy_mask,r00);
408
409             /* Compute parameters for interactions between i and j atoms */
410             qq00             = _mm256_mul_ps(iq0,jq0);
411             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
412                                             vdwioffsetptr0+vdwjidx0B,
413                                             vdwioffsetptr0+vdwjidx0C,
414                                             vdwioffsetptr0+vdwjidx0D,
415                                             vdwioffsetptr0+vdwjidx0E,
416                                             vdwioffsetptr0+vdwjidx0F,
417                                             vdwioffsetptr0+vdwjidx0G,
418                                             vdwioffsetptr0+vdwjidx0H,
419                                             &c6_00,&c12_00);
420
421             /* Calculate table index by multiplying r with table scale and truncate to integer */
422             rt               = _mm256_mul_ps(r00,vftabscale);
423             vfitab           = _mm256_cvttps_epi32(rt);
424             vfeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
425             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
426             vfitab_lo        = _mm256_extractf128_si256(vfitab,0x0);
427             vfitab_hi        = _mm256_extractf128_si256(vfitab,0x1);
428             vfitab_lo        = _mm_slli_epi32(vfitab_lo,3);
429             vfitab_hi        = _mm_slli_epi32(vfitab_hi,3);
430
431             /* REACTION-FIELD ELECTROSTATICS */
432             velec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
433             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
434
435             /* CUBIC SPLINE TABLE DISPERSION */
436             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
437                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
438             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
439                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
440             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
441                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
442             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
443                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
444             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
445             Heps             = _mm256_mul_ps(vfeps,H);
446             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
447             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
448             vvdw6            = _mm256_mul_ps(c6_00,VV);
449             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
450             fvdw6            = _mm256_mul_ps(c6_00,FF);
451
452             /* CUBIC SPLINE TABLE REPULSION */
453             vfitab_lo        = _mm_add_epi32(vfitab_lo,ifour);
454             vfitab_hi        = _mm_add_epi32(vfitab_hi,ifour);
455             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
456                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
457             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
458                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
459             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
460                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
461             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
462                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
463             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
464             Heps             = _mm256_mul_ps(vfeps,H);
465             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
466             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
467             vvdw12           = _mm256_mul_ps(c12_00,VV);
468             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
469             fvdw12           = _mm256_mul_ps(c12_00,FF);
470             vvdw             = _mm256_add_ps(vvdw12,vvdw6);
471             fvdw             = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv00)));
472
473             /* Update potential sum for this i atom from the interaction with this j atom. */
474             velec            = _mm256_andnot_ps(dummy_mask,velec);
475             velecsum         = _mm256_add_ps(velecsum,velec);
476             vvdw             = _mm256_andnot_ps(dummy_mask,vvdw);
477             vvdwsum          = _mm256_add_ps(vvdwsum,vvdw);
478
479             fscal            = _mm256_add_ps(felec,fvdw);
480
481             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
482
483             /* Calculate temporary vectorial force */
484             tx               = _mm256_mul_ps(fscal,dx00);
485             ty               = _mm256_mul_ps(fscal,dy00);
486             tz               = _mm256_mul_ps(fscal,dz00);
487
488             /* Update vectorial force */
489             fix0             = _mm256_add_ps(fix0,tx);
490             fiy0             = _mm256_add_ps(fiy0,ty);
491             fiz0             = _mm256_add_ps(fiz0,tz);
492
493             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
494             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
495             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
496             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
497             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
498             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
499             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
500             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
501             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
502
503             /* Inner loop uses 68 flops */
504         }
505
506         /* End of innermost loop */
507
508         gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
509                                                  f+i_coord_offset,fshift+i_shift_offset);
510
511         ggid                        = gid[iidx];
512         /* Update potential energies */
513         gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
514         gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
515
516         /* Increment number of inner iterations */
517         inneriter                  += j_index_end - j_index_start;
518
519         /* Outer loop uses 9 flops */
520     }
521
522     /* Increment number of outer iterations */
523     outeriter        += nri;
524
525     /* Update outer/inner flops */
526
527     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*68);
528 }
529 /*
530  * Gromacs nonbonded kernel:   nb_kernel_ElecRF_VdwCSTab_GeomP1P1_F_avx_256_single
531  * Electrostatics interaction: ReactionField
532  * VdW interaction:            CubicSplineTable
533  * Geometry:                   Particle-Particle
534  * Calculate force/pot:        Force
535  */
536 void
537 nb_kernel_ElecRF_VdwCSTab_GeomP1P1_F_avx_256_single
538                     (t_nblist                    * gmx_restrict       nlist,
539                      rvec                        * gmx_restrict          xx,
540                      rvec                        * gmx_restrict          ff,
541                      t_forcerec                  * gmx_restrict          fr,
542                      t_mdatoms                   * gmx_restrict     mdatoms,
543                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
544                      t_nrnb                      * gmx_restrict        nrnb)
545 {
546     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
547      * just 0 for non-waters.
548      * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
549      * jnr indices corresponding to data put in the four positions in the SIMD register.
550      */
551     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
552     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
553     int              jnrA,jnrB,jnrC,jnrD;
554     int              jnrE,jnrF,jnrG,jnrH;
555     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
556     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
557     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
558     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
559     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
560     real             rcutoff_scalar;
561     real             *shiftvec,*fshift,*x,*f;
562     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
563     real             scratch[4*DIM];
564     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
565     real *           vdwioffsetptr0;
566     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
567     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
568     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
569     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
570     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
571     real             *charge;
572     int              nvdwtype;
573     __m256           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
574     int              *vdwtype;
575     real             *vdwparam;
576     __m256           one_sixth   = _mm256_set1_ps(1.0/6.0);
577     __m256           one_twelfth = _mm256_set1_ps(1.0/12.0);
578     __m256i          vfitab;
579     __m128i          vfitab_lo,vfitab_hi;
580     __m128i          ifour       = _mm_set1_epi32(4);
581     __m256           rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
582     real             *vftab;
583     __m256           dummy_mask,cutoff_mask;
584     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
585     __m256           one     = _mm256_set1_ps(1.0);
586     __m256           two     = _mm256_set1_ps(2.0);
587     x                = xx[0];
588     f                = ff[0];
589
590     nri              = nlist->nri;
591     iinr             = nlist->iinr;
592     jindex           = nlist->jindex;
593     jjnr             = nlist->jjnr;
594     shiftidx         = nlist->shift;
595     gid              = nlist->gid;
596     shiftvec         = fr->shift_vec[0];
597     fshift           = fr->fshift[0];
598     facel            = _mm256_set1_ps(fr->epsfac);
599     charge           = mdatoms->chargeA;
600     krf              = _mm256_set1_ps(fr->ic->k_rf);
601     krf2             = _mm256_set1_ps(fr->ic->k_rf*2.0);
602     crf              = _mm256_set1_ps(fr->ic->c_rf);
603     nvdwtype         = fr->ntype;
604     vdwparam         = fr->nbfp;
605     vdwtype          = mdatoms->typeA;
606
607     vftab            = kernel_data->table_vdw->data;
608     vftabscale       = _mm256_set1_ps(kernel_data->table_vdw->scale);
609
610     /* Avoid stupid compiler warnings */
611     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
612     j_coord_offsetA = 0;
613     j_coord_offsetB = 0;
614     j_coord_offsetC = 0;
615     j_coord_offsetD = 0;
616     j_coord_offsetE = 0;
617     j_coord_offsetF = 0;
618     j_coord_offsetG = 0;
619     j_coord_offsetH = 0;
620
621     outeriter        = 0;
622     inneriter        = 0;
623
624     for(iidx=0;iidx<4*DIM;iidx++)
625     {
626         scratch[iidx] = 0.0;
627     }
628
629     /* Start outer loop over neighborlists */
630     for(iidx=0; iidx<nri; iidx++)
631     {
632         /* Load shift vector for this list */
633         i_shift_offset   = DIM*shiftidx[iidx];
634
635         /* Load limits for loop over neighbors */
636         j_index_start    = jindex[iidx];
637         j_index_end      = jindex[iidx+1];
638
639         /* Get outer coordinate index */
640         inr              = iinr[iidx];
641         i_coord_offset   = DIM*inr;
642
643         /* Load i particle coords and add shift vector */
644         gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
645
646         fix0             = _mm256_setzero_ps();
647         fiy0             = _mm256_setzero_ps();
648         fiz0             = _mm256_setzero_ps();
649
650         /* Load parameters for i particles */
651         iq0              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
652         vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
653
654         /* Start inner kernel loop */
655         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
656         {
657
658             /* Get j neighbor index, and coordinate index */
659             jnrA             = jjnr[jidx];
660             jnrB             = jjnr[jidx+1];
661             jnrC             = jjnr[jidx+2];
662             jnrD             = jjnr[jidx+3];
663             jnrE             = jjnr[jidx+4];
664             jnrF             = jjnr[jidx+5];
665             jnrG             = jjnr[jidx+6];
666             jnrH             = jjnr[jidx+7];
667             j_coord_offsetA  = DIM*jnrA;
668             j_coord_offsetB  = DIM*jnrB;
669             j_coord_offsetC  = DIM*jnrC;
670             j_coord_offsetD  = DIM*jnrD;
671             j_coord_offsetE  = DIM*jnrE;
672             j_coord_offsetF  = DIM*jnrF;
673             j_coord_offsetG  = DIM*jnrG;
674             j_coord_offsetH  = DIM*jnrH;
675
676             /* load j atom coordinates */
677             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
678                                                  x+j_coord_offsetC,x+j_coord_offsetD,
679                                                  x+j_coord_offsetE,x+j_coord_offsetF,
680                                                  x+j_coord_offsetG,x+j_coord_offsetH,
681                                                  &jx0,&jy0,&jz0);
682
683             /* Calculate displacement vector */
684             dx00             = _mm256_sub_ps(ix0,jx0);
685             dy00             = _mm256_sub_ps(iy0,jy0);
686             dz00             = _mm256_sub_ps(iz0,jz0);
687
688             /* Calculate squared distance and things based on it */
689             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
690
691             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
692
693             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
694
695             /* Load parameters for j particles */
696             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
697                                                                  charge+jnrC+0,charge+jnrD+0,
698                                                                  charge+jnrE+0,charge+jnrF+0,
699                                                                  charge+jnrG+0,charge+jnrH+0);
700             vdwjidx0A        = 2*vdwtype[jnrA+0];
701             vdwjidx0B        = 2*vdwtype[jnrB+0];
702             vdwjidx0C        = 2*vdwtype[jnrC+0];
703             vdwjidx0D        = 2*vdwtype[jnrD+0];
704             vdwjidx0E        = 2*vdwtype[jnrE+0];
705             vdwjidx0F        = 2*vdwtype[jnrF+0];
706             vdwjidx0G        = 2*vdwtype[jnrG+0];
707             vdwjidx0H        = 2*vdwtype[jnrH+0];
708
709             /**************************
710              * CALCULATE INTERACTIONS *
711              **************************/
712
713             r00              = _mm256_mul_ps(rsq00,rinv00);
714
715             /* Compute parameters for interactions between i and j atoms */
716             qq00             = _mm256_mul_ps(iq0,jq0);
717             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
718                                             vdwioffsetptr0+vdwjidx0B,
719                                             vdwioffsetptr0+vdwjidx0C,
720                                             vdwioffsetptr0+vdwjidx0D,
721                                             vdwioffsetptr0+vdwjidx0E,
722                                             vdwioffsetptr0+vdwjidx0F,
723                                             vdwioffsetptr0+vdwjidx0G,
724                                             vdwioffsetptr0+vdwjidx0H,
725                                             &c6_00,&c12_00);
726
727             /* Calculate table index by multiplying r with table scale and truncate to integer */
728             rt               = _mm256_mul_ps(r00,vftabscale);
729             vfitab           = _mm256_cvttps_epi32(rt);
730             vfeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
731             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
732             vfitab_lo        = _mm256_extractf128_si256(vfitab,0x0);
733             vfitab_hi        = _mm256_extractf128_si256(vfitab,0x1);
734             vfitab_lo        = _mm_slli_epi32(vfitab_lo,3);
735             vfitab_hi        = _mm_slli_epi32(vfitab_hi,3);
736
737             /* REACTION-FIELD ELECTROSTATICS */
738             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
739
740             /* CUBIC SPLINE TABLE DISPERSION */
741             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
742                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
743             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
744                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
745             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
746                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
747             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
748                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
749             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
750             Heps             = _mm256_mul_ps(vfeps,H);
751             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
752             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
753             fvdw6            = _mm256_mul_ps(c6_00,FF);
754
755             /* CUBIC SPLINE TABLE REPULSION */
756             vfitab_lo        = _mm_add_epi32(vfitab_lo,ifour);
757             vfitab_hi        = _mm_add_epi32(vfitab_hi,ifour);
758             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
759                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
760             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
761                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
762             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
763                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
764             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
765                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
766             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
767             Heps             = _mm256_mul_ps(vfeps,H);
768             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
769             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
770             fvdw12           = _mm256_mul_ps(c12_00,FF);
771             fvdw             = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv00)));
772
773             fscal            = _mm256_add_ps(felec,fvdw);
774
775             /* Calculate temporary vectorial force */
776             tx               = _mm256_mul_ps(fscal,dx00);
777             ty               = _mm256_mul_ps(fscal,dy00);
778             tz               = _mm256_mul_ps(fscal,dz00);
779
780             /* Update vectorial force */
781             fix0             = _mm256_add_ps(fix0,tx);
782             fiy0             = _mm256_add_ps(fiy0,ty);
783             fiz0             = _mm256_add_ps(fiz0,tz);
784
785             fjptrA             = f+j_coord_offsetA;
786             fjptrB             = f+j_coord_offsetB;
787             fjptrC             = f+j_coord_offsetC;
788             fjptrD             = f+j_coord_offsetD;
789             fjptrE             = f+j_coord_offsetE;
790             fjptrF             = f+j_coord_offsetF;
791             fjptrG             = f+j_coord_offsetG;
792             fjptrH             = f+j_coord_offsetH;
793             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
794
795             /* Inner loop uses 54 flops */
796         }
797
798         if(jidx<j_index_end)
799         {
800
801             /* Get j neighbor index, and coordinate index */
802             jnrlistA         = jjnr[jidx];
803             jnrlistB         = jjnr[jidx+1];
804             jnrlistC         = jjnr[jidx+2];
805             jnrlistD         = jjnr[jidx+3];
806             jnrlistE         = jjnr[jidx+4];
807             jnrlistF         = jjnr[jidx+5];
808             jnrlistG         = jjnr[jidx+6];
809             jnrlistH         = jjnr[jidx+7];
810             /* Sign of each element will be negative for non-real atoms.
811              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
812              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
813              */
814             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
815                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
816                                             
817             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
818             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
819             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
820             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
821             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
822             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
823             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
824             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
825             j_coord_offsetA  = DIM*jnrA;
826             j_coord_offsetB  = DIM*jnrB;
827             j_coord_offsetC  = DIM*jnrC;
828             j_coord_offsetD  = DIM*jnrD;
829             j_coord_offsetE  = DIM*jnrE;
830             j_coord_offsetF  = DIM*jnrF;
831             j_coord_offsetG  = DIM*jnrG;
832             j_coord_offsetH  = DIM*jnrH;
833
834             /* load j atom coordinates */
835             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
836                                                  x+j_coord_offsetC,x+j_coord_offsetD,
837                                                  x+j_coord_offsetE,x+j_coord_offsetF,
838                                                  x+j_coord_offsetG,x+j_coord_offsetH,
839                                                  &jx0,&jy0,&jz0);
840
841             /* Calculate displacement vector */
842             dx00             = _mm256_sub_ps(ix0,jx0);
843             dy00             = _mm256_sub_ps(iy0,jy0);
844             dz00             = _mm256_sub_ps(iz0,jz0);
845
846             /* Calculate squared distance and things based on it */
847             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
848
849             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
850
851             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
852
853             /* Load parameters for j particles */
854             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
855                                                                  charge+jnrC+0,charge+jnrD+0,
856                                                                  charge+jnrE+0,charge+jnrF+0,
857                                                                  charge+jnrG+0,charge+jnrH+0);
858             vdwjidx0A        = 2*vdwtype[jnrA+0];
859             vdwjidx0B        = 2*vdwtype[jnrB+0];
860             vdwjidx0C        = 2*vdwtype[jnrC+0];
861             vdwjidx0D        = 2*vdwtype[jnrD+0];
862             vdwjidx0E        = 2*vdwtype[jnrE+0];
863             vdwjidx0F        = 2*vdwtype[jnrF+0];
864             vdwjidx0G        = 2*vdwtype[jnrG+0];
865             vdwjidx0H        = 2*vdwtype[jnrH+0];
866
867             /**************************
868              * CALCULATE INTERACTIONS *
869              **************************/
870
871             r00              = _mm256_mul_ps(rsq00,rinv00);
872             r00              = _mm256_andnot_ps(dummy_mask,r00);
873
874             /* Compute parameters for interactions between i and j atoms */
875             qq00             = _mm256_mul_ps(iq0,jq0);
876             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
877                                             vdwioffsetptr0+vdwjidx0B,
878                                             vdwioffsetptr0+vdwjidx0C,
879                                             vdwioffsetptr0+vdwjidx0D,
880                                             vdwioffsetptr0+vdwjidx0E,
881                                             vdwioffsetptr0+vdwjidx0F,
882                                             vdwioffsetptr0+vdwjidx0G,
883                                             vdwioffsetptr0+vdwjidx0H,
884                                             &c6_00,&c12_00);
885
886             /* Calculate table index by multiplying r with table scale and truncate to integer */
887             rt               = _mm256_mul_ps(r00,vftabscale);
888             vfitab           = _mm256_cvttps_epi32(rt);
889             vfeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
890             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
891             vfitab_lo        = _mm256_extractf128_si256(vfitab,0x0);
892             vfitab_hi        = _mm256_extractf128_si256(vfitab,0x1);
893             vfitab_lo        = _mm_slli_epi32(vfitab_lo,3);
894             vfitab_hi        = _mm_slli_epi32(vfitab_hi,3);
895
896             /* REACTION-FIELD ELECTROSTATICS */
897             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
898
899             /* CUBIC SPLINE TABLE DISPERSION */
900             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
901                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
902             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
903                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
904             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
905                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
906             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
907                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
908             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
909             Heps             = _mm256_mul_ps(vfeps,H);
910             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
911             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
912             fvdw6            = _mm256_mul_ps(c6_00,FF);
913
914             /* CUBIC SPLINE TABLE REPULSION */
915             vfitab_lo        = _mm_add_epi32(vfitab_lo,ifour);
916             vfitab_hi        = _mm_add_epi32(vfitab_hi,ifour);
917             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
918                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
919             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
920                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
921             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
922                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
923             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
924                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
925             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
926             Heps             = _mm256_mul_ps(vfeps,H);
927             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
928             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
929             fvdw12           = _mm256_mul_ps(c12_00,FF);
930             fvdw             = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv00)));
931
932             fscal            = _mm256_add_ps(felec,fvdw);
933
934             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
935
936             /* Calculate temporary vectorial force */
937             tx               = _mm256_mul_ps(fscal,dx00);
938             ty               = _mm256_mul_ps(fscal,dy00);
939             tz               = _mm256_mul_ps(fscal,dz00);
940
941             /* Update vectorial force */
942             fix0             = _mm256_add_ps(fix0,tx);
943             fiy0             = _mm256_add_ps(fiy0,ty);
944             fiz0             = _mm256_add_ps(fiz0,tz);
945
946             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
947             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
948             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
949             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
950             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
951             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
952             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
953             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
954             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
955
956             /* Inner loop uses 55 flops */
957         }
958
959         /* End of innermost loop */
960
961         gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
962                                                  f+i_coord_offset,fshift+i_shift_offset);
963
964         /* Increment number of inner iterations */
965         inneriter                  += j_index_end - j_index_start;
966
967         /* Outer loop uses 7 flops */
968     }
969
970     /* Increment number of outer iterations */
971     outeriter        += nri;
972
973     /* Update outer/inner flops */
974
975     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*55);
976 }