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