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