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