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