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