Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_single / nb_kernel_ElecEw_VdwLJEw_GeomW4P1_avx_256_single.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS avx_256_single kernel generator.
37  */
38 #include "gmxpre.h"
39
40 #include "config.h"
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "gromacs/legacyheaders/types/simple.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/legacyheaders/nrnb.h"
48
49 #include "gromacs/simd/math_x86_avx_256_single.h"
50 #include "kernelutil_x86_avx_256_single.h"
51
52 /*
53  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwLJEw_GeomW4P1_VF_avx_256_single
54  * Electrostatics interaction: Ewald
55  * VdW interaction:            LJEwald
56  * Geometry:                   Water4-Particle
57  * Calculate force/pot:        PotentialAndForce
58  */
59 void
60 nb_kernel_ElecEw_VdwLJEw_GeomW4P1_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     real *           vdwioffsetptr1;
92     real *           vdwgridioffsetptr1;
93     __m256           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
94     real *           vdwioffsetptr2;
95     real *           vdwgridioffsetptr2;
96     __m256           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
97     real *           vdwioffsetptr3;
98     real *           vdwgridioffsetptr3;
99     __m256           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
100     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
101     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
102     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
103     __m256           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
104     __m256           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
105     __m256           dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
106     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
107     real             *charge;
108     int              nvdwtype;
109     __m256           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
110     int              *vdwtype;
111     real             *vdwparam;
112     __m256           one_sixth   = _mm256_set1_ps(1.0/6.0);
113     __m256           one_twelfth = _mm256_set1_ps(1.0/12.0);
114     __m256           c6grid_00;
115     __m256           c6grid_10;
116     __m256           c6grid_20;
117     __m256           c6grid_30;
118     real             *vdwgridparam;
119     __m256           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
120     __m256           one_half  = _mm256_set1_ps(0.5);
121     __m256           minus_one = _mm256_set1_ps(-1.0);
122     __m256i          ewitab;
123     __m128i          ewitab_lo,ewitab_hi;
124     __m256           ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
125     __m256           beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
126     real             *ewtab;
127     __m256           dummy_mask,cutoff_mask;
128     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
129     __m256           one     = _mm256_set1_ps(1.0);
130     __m256           two     = _mm256_set1_ps(2.0);
131     x                = xx[0];
132     f                = ff[0];
133
134     nri              = nlist->nri;
135     iinr             = nlist->iinr;
136     jindex           = nlist->jindex;
137     jjnr             = nlist->jjnr;
138     shiftidx         = nlist->shift;
139     gid              = nlist->gid;
140     shiftvec         = fr->shift_vec[0];
141     fshift           = fr->fshift[0];
142     facel            = _mm256_set1_ps(fr->epsfac);
143     charge           = mdatoms->chargeA;
144     nvdwtype         = fr->ntype;
145     vdwparam         = fr->nbfp;
146     vdwtype          = mdatoms->typeA;
147     vdwgridparam     = fr->ljpme_c6grid;
148     sh_lj_ewald      = _mm256_set1_ps(fr->ic->sh_lj_ewald);
149     ewclj            = _mm256_set1_ps(fr->ewaldcoeff_lj);
150     ewclj2           = _mm256_mul_ps(minus_one,_mm256_mul_ps(ewclj,ewclj));
151
152     sh_ewald         = _mm256_set1_ps(fr->ic->sh_ewald);
153     beta             = _mm256_set1_ps(fr->ic->ewaldcoeff_q);
154     beta2            = _mm256_mul_ps(beta,beta);
155     beta3            = _mm256_mul_ps(beta,beta2);
156
157     ewtab            = fr->ic->tabq_coul_FDV0;
158     ewtabscale       = _mm256_set1_ps(fr->ic->tabq_scale);
159     ewtabhalfspace   = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
160
161     /* Setup water-specific parameters */
162     inr              = nlist->iinr[0];
163     iq1              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
164     iq2              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
165     iq3              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+3]));
166     vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
167     vdwgridioffsetptr0 = vdwgridparam+2*nvdwtype*vdwtype[inr+0];
168
169     /* Avoid stupid compiler warnings */
170     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
171     j_coord_offsetA = 0;
172     j_coord_offsetB = 0;
173     j_coord_offsetC = 0;
174     j_coord_offsetD = 0;
175     j_coord_offsetE = 0;
176     j_coord_offsetF = 0;
177     j_coord_offsetG = 0;
178     j_coord_offsetH = 0;
179
180     outeriter        = 0;
181     inneriter        = 0;
182
183     for(iidx=0;iidx<4*DIM;iidx++)
184     {
185         scratch[iidx] = 0.0;
186     }
187
188     /* Start outer loop over neighborlists */
189     for(iidx=0; iidx<nri; iidx++)
190     {
191         /* Load shift vector for this list */
192         i_shift_offset   = DIM*shiftidx[iidx];
193
194         /* Load limits for loop over neighbors */
195         j_index_start    = jindex[iidx];
196         j_index_end      = jindex[iidx+1];
197
198         /* Get outer coordinate index */
199         inr              = iinr[iidx];
200         i_coord_offset   = DIM*inr;
201
202         /* Load i particle coords and add shift vector */
203         gmx_mm256_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
204                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
205
206         fix0             = _mm256_setzero_ps();
207         fiy0             = _mm256_setzero_ps();
208         fiz0             = _mm256_setzero_ps();
209         fix1             = _mm256_setzero_ps();
210         fiy1             = _mm256_setzero_ps();
211         fiz1             = _mm256_setzero_ps();
212         fix2             = _mm256_setzero_ps();
213         fiy2             = _mm256_setzero_ps();
214         fiz2             = _mm256_setzero_ps();
215         fix3             = _mm256_setzero_ps();
216         fiy3             = _mm256_setzero_ps();
217         fiz3             = _mm256_setzero_ps();
218
219         /* Reset potential sums */
220         velecsum         = _mm256_setzero_ps();
221         vvdwsum          = _mm256_setzero_ps();
222
223         /* Start inner kernel loop */
224         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
225         {
226
227             /* Get j neighbor index, and coordinate index */
228             jnrA             = jjnr[jidx];
229             jnrB             = jjnr[jidx+1];
230             jnrC             = jjnr[jidx+2];
231             jnrD             = jjnr[jidx+3];
232             jnrE             = jjnr[jidx+4];
233             jnrF             = jjnr[jidx+5];
234             jnrG             = jjnr[jidx+6];
235             jnrH             = jjnr[jidx+7];
236             j_coord_offsetA  = DIM*jnrA;
237             j_coord_offsetB  = DIM*jnrB;
238             j_coord_offsetC  = DIM*jnrC;
239             j_coord_offsetD  = DIM*jnrD;
240             j_coord_offsetE  = DIM*jnrE;
241             j_coord_offsetF  = DIM*jnrF;
242             j_coord_offsetG  = DIM*jnrG;
243             j_coord_offsetH  = DIM*jnrH;
244
245             /* load j atom coordinates */
246             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
247                                                  x+j_coord_offsetC,x+j_coord_offsetD,
248                                                  x+j_coord_offsetE,x+j_coord_offsetF,
249                                                  x+j_coord_offsetG,x+j_coord_offsetH,
250                                                  &jx0,&jy0,&jz0);
251
252             /* Calculate displacement vector */
253             dx00             = _mm256_sub_ps(ix0,jx0);
254             dy00             = _mm256_sub_ps(iy0,jy0);
255             dz00             = _mm256_sub_ps(iz0,jz0);
256             dx10             = _mm256_sub_ps(ix1,jx0);
257             dy10             = _mm256_sub_ps(iy1,jy0);
258             dz10             = _mm256_sub_ps(iz1,jz0);
259             dx20             = _mm256_sub_ps(ix2,jx0);
260             dy20             = _mm256_sub_ps(iy2,jy0);
261             dz20             = _mm256_sub_ps(iz2,jz0);
262             dx30             = _mm256_sub_ps(ix3,jx0);
263             dy30             = _mm256_sub_ps(iy3,jy0);
264             dz30             = _mm256_sub_ps(iz3,jz0);
265
266             /* Calculate squared distance and things based on it */
267             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
268             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
269             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
270             rsq30            = gmx_mm256_calc_rsq_ps(dx30,dy30,dz30);
271
272             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
273             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
274             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
275             rinv30           = gmx_mm256_invsqrt_ps(rsq30);
276
277             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
278             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
279             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
280             rinvsq30         = _mm256_mul_ps(rinv30,rinv30);
281
282             /* Load parameters for j particles */
283             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
284                                                                  charge+jnrC+0,charge+jnrD+0,
285                                                                  charge+jnrE+0,charge+jnrF+0,
286                                                                  charge+jnrG+0,charge+jnrH+0);
287             vdwjidx0A        = 2*vdwtype[jnrA+0];
288             vdwjidx0B        = 2*vdwtype[jnrB+0];
289             vdwjidx0C        = 2*vdwtype[jnrC+0];
290             vdwjidx0D        = 2*vdwtype[jnrD+0];
291             vdwjidx0E        = 2*vdwtype[jnrE+0];
292             vdwjidx0F        = 2*vdwtype[jnrF+0];
293             vdwjidx0G        = 2*vdwtype[jnrG+0];
294             vdwjidx0H        = 2*vdwtype[jnrH+0];
295
296             fjx0             = _mm256_setzero_ps();
297             fjy0             = _mm256_setzero_ps();
298             fjz0             = _mm256_setzero_ps();
299
300             /**************************
301              * CALCULATE INTERACTIONS *
302              **************************/
303
304             r00              = _mm256_mul_ps(rsq00,rinv00);
305
306             /* Compute parameters for interactions between i and j atoms */
307             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
308                                             vdwioffsetptr0+vdwjidx0B,
309                                             vdwioffsetptr0+vdwjidx0C,
310                                             vdwioffsetptr0+vdwjidx0D,
311                                             vdwioffsetptr0+vdwjidx0E,
312                                             vdwioffsetptr0+vdwjidx0F,
313                                             vdwioffsetptr0+vdwjidx0G,
314                                             vdwioffsetptr0+vdwjidx0H,
315                                             &c6_00,&c12_00);
316
317             c6grid_00       = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
318                                                                   vdwgridioffsetptr0+vdwjidx0B,
319                                                                   vdwgridioffsetptr0+vdwjidx0C,
320                                                                   vdwgridioffsetptr0+vdwjidx0D,
321                                                                   vdwgridioffsetptr0+vdwjidx0E,
322                                                                   vdwgridioffsetptr0+vdwjidx0F,
323                                                                   vdwgridioffsetptr0+vdwjidx0G,
324                                                                   vdwgridioffsetptr0+vdwjidx0H);
325
326             /* Analytical LJ-PME */
327             rinvsix          = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
328             ewcljrsq         = _mm256_mul_ps(ewclj2,rsq00);
329             ewclj6           = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
330             exponent         = gmx_simd_exp_r(ewcljrsq);
331             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
332             poly             = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
333             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
334             vvdw6            = _mm256_mul_ps(_mm256_sub_ps(c6_00,_mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly))),rinvsix);
335             vvdw12           = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
336             vvdw             = _mm256_sub_ps(_mm256_mul_ps(vvdw12,one_twelfth),_mm256_mul_ps(vvdw6,one_sixth));
337             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
338             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);
339
340             /* Update potential sum for this i atom from the interaction with this j atom. */
341             vvdwsum          = _mm256_add_ps(vvdwsum,vvdw);
342
343             fscal            = fvdw;
344
345             /* Calculate temporary vectorial force */
346             tx               = _mm256_mul_ps(fscal,dx00);
347             ty               = _mm256_mul_ps(fscal,dy00);
348             tz               = _mm256_mul_ps(fscal,dz00);
349
350             /* Update vectorial force */
351             fix0             = _mm256_add_ps(fix0,tx);
352             fiy0             = _mm256_add_ps(fiy0,ty);
353             fiz0             = _mm256_add_ps(fiz0,tz);
354
355             fjx0             = _mm256_add_ps(fjx0,tx);
356             fjy0             = _mm256_add_ps(fjy0,ty);
357             fjz0             = _mm256_add_ps(fjz0,tz);
358
359             /**************************
360              * CALCULATE INTERACTIONS *
361              **************************/
362
363             r10              = _mm256_mul_ps(rsq10,rinv10);
364
365             /* Compute parameters for interactions between i and j atoms */
366             qq10             = _mm256_mul_ps(iq1,jq0);
367
368             /* EWALD ELECTROSTATICS */
369             
370             /* Analytical PME correction */
371             zeta2            = _mm256_mul_ps(beta2,rsq10);
372             rinv3            = _mm256_mul_ps(rinvsq10,rinv10);
373             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
374             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
375             felec            = _mm256_mul_ps(qq10,felec);
376             pmecorrV         = gmx_mm256_pmecorrV_ps(zeta2);
377             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
378             velec            = _mm256_sub_ps(rinv10,pmecorrV);
379             velec            = _mm256_mul_ps(qq10,velec);
380             
381             /* Update potential sum for this i atom from the interaction with this j atom. */
382             velecsum         = _mm256_add_ps(velecsum,velec);
383
384             fscal            = felec;
385
386             /* Calculate temporary vectorial force */
387             tx               = _mm256_mul_ps(fscal,dx10);
388             ty               = _mm256_mul_ps(fscal,dy10);
389             tz               = _mm256_mul_ps(fscal,dz10);
390
391             /* Update vectorial force */
392             fix1             = _mm256_add_ps(fix1,tx);
393             fiy1             = _mm256_add_ps(fiy1,ty);
394             fiz1             = _mm256_add_ps(fiz1,tz);
395
396             fjx0             = _mm256_add_ps(fjx0,tx);
397             fjy0             = _mm256_add_ps(fjy0,ty);
398             fjz0             = _mm256_add_ps(fjz0,tz);
399
400             /**************************
401              * CALCULATE INTERACTIONS *
402              **************************/
403
404             r20              = _mm256_mul_ps(rsq20,rinv20);
405
406             /* Compute parameters for interactions between i and j atoms */
407             qq20             = _mm256_mul_ps(iq2,jq0);
408
409             /* EWALD ELECTROSTATICS */
410             
411             /* Analytical PME correction */
412             zeta2            = _mm256_mul_ps(beta2,rsq20);
413             rinv3            = _mm256_mul_ps(rinvsq20,rinv20);
414             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
415             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
416             felec            = _mm256_mul_ps(qq20,felec);
417             pmecorrV         = gmx_mm256_pmecorrV_ps(zeta2);
418             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
419             velec            = _mm256_sub_ps(rinv20,pmecorrV);
420             velec            = _mm256_mul_ps(qq20,velec);
421             
422             /* Update potential sum for this i atom from the interaction with this j atom. */
423             velecsum         = _mm256_add_ps(velecsum,velec);
424
425             fscal            = felec;
426
427             /* Calculate temporary vectorial force */
428             tx               = _mm256_mul_ps(fscal,dx20);
429             ty               = _mm256_mul_ps(fscal,dy20);
430             tz               = _mm256_mul_ps(fscal,dz20);
431
432             /* Update vectorial force */
433             fix2             = _mm256_add_ps(fix2,tx);
434             fiy2             = _mm256_add_ps(fiy2,ty);
435             fiz2             = _mm256_add_ps(fiz2,tz);
436
437             fjx0             = _mm256_add_ps(fjx0,tx);
438             fjy0             = _mm256_add_ps(fjy0,ty);
439             fjz0             = _mm256_add_ps(fjz0,tz);
440
441             /**************************
442              * CALCULATE INTERACTIONS *
443              **************************/
444
445             r30              = _mm256_mul_ps(rsq30,rinv30);
446
447             /* Compute parameters for interactions between i and j atoms */
448             qq30             = _mm256_mul_ps(iq3,jq0);
449
450             /* EWALD ELECTROSTATICS */
451             
452             /* Analytical PME correction */
453             zeta2            = _mm256_mul_ps(beta2,rsq30);
454             rinv3            = _mm256_mul_ps(rinvsq30,rinv30);
455             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
456             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
457             felec            = _mm256_mul_ps(qq30,felec);
458             pmecorrV         = gmx_mm256_pmecorrV_ps(zeta2);
459             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
460             velec            = _mm256_sub_ps(rinv30,pmecorrV);
461             velec            = _mm256_mul_ps(qq30,velec);
462             
463             /* Update potential sum for this i atom from the interaction with this j atom. */
464             velecsum         = _mm256_add_ps(velecsum,velec);
465
466             fscal            = felec;
467
468             /* Calculate temporary vectorial force */
469             tx               = _mm256_mul_ps(fscal,dx30);
470             ty               = _mm256_mul_ps(fscal,dy30);
471             tz               = _mm256_mul_ps(fscal,dz30);
472
473             /* Update vectorial force */
474             fix3             = _mm256_add_ps(fix3,tx);
475             fiy3             = _mm256_add_ps(fiy3,ty);
476             fiz3             = _mm256_add_ps(fiz3,tz);
477
478             fjx0             = _mm256_add_ps(fjx0,tx);
479             fjy0             = _mm256_add_ps(fjy0,ty);
480             fjz0             = _mm256_add_ps(fjz0,tz);
481
482             fjptrA             = f+j_coord_offsetA;
483             fjptrB             = f+j_coord_offsetB;
484             fjptrC             = f+j_coord_offsetC;
485             fjptrD             = f+j_coord_offsetD;
486             fjptrE             = f+j_coord_offsetE;
487             fjptrF             = f+j_coord_offsetF;
488             fjptrG             = f+j_coord_offsetG;
489             fjptrH             = f+j_coord_offsetH;
490
491             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
492
493             /* Inner loop uses 306 flops */
494         }
495
496         if(jidx<j_index_end)
497         {
498
499             /* Get j neighbor index, and coordinate index */
500             jnrlistA         = jjnr[jidx];
501             jnrlistB         = jjnr[jidx+1];
502             jnrlistC         = jjnr[jidx+2];
503             jnrlistD         = jjnr[jidx+3];
504             jnrlistE         = jjnr[jidx+4];
505             jnrlistF         = jjnr[jidx+5];
506             jnrlistG         = jjnr[jidx+6];
507             jnrlistH         = jjnr[jidx+7];
508             /* Sign of each element will be negative for non-real atoms.
509              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
510              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
511              */
512             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
513                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
514                                             
515             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
516             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
517             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
518             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
519             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
520             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
521             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
522             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
523             j_coord_offsetA  = DIM*jnrA;
524             j_coord_offsetB  = DIM*jnrB;
525             j_coord_offsetC  = DIM*jnrC;
526             j_coord_offsetD  = DIM*jnrD;
527             j_coord_offsetE  = DIM*jnrE;
528             j_coord_offsetF  = DIM*jnrF;
529             j_coord_offsetG  = DIM*jnrG;
530             j_coord_offsetH  = DIM*jnrH;
531
532             /* load j atom coordinates */
533             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
534                                                  x+j_coord_offsetC,x+j_coord_offsetD,
535                                                  x+j_coord_offsetE,x+j_coord_offsetF,
536                                                  x+j_coord_offsetG,x+j_coord_offsetH,
537                                                  &jx0,&jy0,&jz0);
538
539             /* Calculate displacement vector */
540             dx00             = _mm256_sub_ps(ix0,jx0);
541             dy00             = _mm256_sub_ps(iy0,jy0);
542             dz00             = _mm256_sub_ps(iz0,jz0);
543             dx10             = _mm256_sub_ps(ix1,jx0);
544             dy10             = _mm256_sub_ps(iy1,jy0);
545             dz10             = _mm256_sub_ps(iz1,jz0);
546             dx20             = _mm256_sub_ps(ix2,jx0);
547             dy20             = _mm256_sub_ps(iy2,jy0);
548             dz20             = _mm256_sub_ps(iz2,jz0);
549             dx30             = _mm256_sub_ps(ix3,jx0);
550             dy30             = _mm256_sub_ps(iy3,jy0);
551             dz30             = _mm256_sub_ps(iz3,jz0);
552
553             /* Calculate squared distance and things based on it */
554             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
555             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
556             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
557             rsq30            = gmx_mm256_calc_rsq_ps(dx30,dy30,dz30);
558
559             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
560             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
561             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
562             rinv30           = gmx_mm256_invsqrt_ps(rsq30);
563
564             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
565             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
566             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
567             rinvsq30         = _mm256_mul_ps(rinv30,rinv30);
568
569             /* Load parameters for j particles */
570             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
571                                                                  charge+jnrC+0,charge+jnrD+0,
572                                                                  charge+jnrE+0,charge+jnrF+0,
573                                                                  charge+jnrG+0,charge+jnrH+0);
574             vdwjidx0A        = 2*vdwtype[jnrA+0];
575             vdwjidx0B        = 2*vdwtype[jnrB+0];
576             vdwjidx0C        = 2*vdwtype[jnrC+0];
577             vdwjidx0D        = 2*vdwtype[jnrD+0];
578             vdwjidx0E        = 2*vdwtype[jnrE+0];
579             vdwjidx0F        = 2*vdwtype[jnrF+0];
580             vdwjidx0G        = 2*vdwtype[jnrG+0];
581             vdwjidx0H        = 2*vdwtype[jnrH+0];
582
583             fjx0             = _mm256_setzero_ps();
584             fjy0             = _mm256_setzero_ps();
585             fjz0             = _mm256_setzero_ps();
586
587             /**************************
588              * CALCULATE INTERACTIONS *
589              **************************/
590
591             r00              = _mm256_mul_ps(rsq00,rinv00);
592             r00              = _mm256_andnot_ps(dummy_mask,r00);
593
594             /* Compute parameters for interactions between i and j atoms */
595             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
596                                             vdwioffsetptr0+vdwjidx0B,
597                                             vdwioffsetptr0+vdwjidx0C,
598                                             vdwioffsetptr0+vdwjidx0D,
599                                             vdwioffsetptr0+vdwjidx0E,
600                                             vdwioffsetptr0+vdwjidx0F,
601                                             vdwioffsetptr0+vdwjidx0G,
602                                             vdwioffsetptr0+vdwjidx0H,
603                                             &c6_00,&c12_00);
604
605             c6grid_00       = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
606                                                                   vdwgridioffsetptr0+vdwjidx0B,
607                                                                   vdwgridioffsetptr0+vdwjidx0C,
608                                                                   vdwgridioffsetptr0+vdwjidx0D,
609                                                                   vdwgridioffsetptr0+vdwjidx0E,
610                                                                   vdwgridioffsetptr0+vdwjidx0F,
611                                                                   vdwgridioffsetptr0+vdwjidx0G,
612                                                                   vdwgridioffsetptr0+vdwjidx0H);
613
614             /* Analytical LJ-PME */
615             rinvsix          = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
616             ewcljrsq         = _mm256_mul_ps(ewclj2,rsq00);
617             ewclj6           = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
618             exponent         = gmx_simd_exp_r(ewcljrsq);
619             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
620             poly             = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
621             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
622             vvdw6            = _mm256_mul_ps(_mm256_sub_ps(c6_00,_mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly))),rinvsix);
623             vvdw12           = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
624             vvdw             = _mm256_sub_ps(_mm256_mul_ps(vvdw12,one_twelfth),_mm256_mul_ps(vvdw6,one_sixth));
625             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
626             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);
627
628             /* Update potential sum for this i atom from the interaction with this j atom. */
629             vvdw             = _mm256_andnot_ps(dummy_mask,vvdw);
630             vvdwsum          = _mm256_add_ps(vvdwsum,vvdw);
631
632             fscal            = fvdw;
633
634             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
635
636             /* Calculate temporary vectorial force */
637             tx               = _mm256_mul_ps(fscal,dx00);
638             ty               = _mm256_mul_ps(fscal,dy00);
639             tz               = _mm256_mul_ps(fscal,dz00);
640
641             /* Update vectorial force */
642             fix0             = _mm256_add_ps(fix0,tx);
643             fiy0             = _mm256_add_ps(fiy0,ty);
644             fiz0             = _mm256_add_ps(fiz0,tz);
645
646             fjx0             = _mm256_add_ps(fjx0,tx);
647             fjy0             = _mm256_add_ps(fjy0,ty);
648             fjz0             = _mm256_add_ps(fjz0,tz);
649
650             /**************************
651              * CALCULATE INTERACTIONS *
652              **************************/
653
654             r10              = _mm256_mul_ps(rsq10,rinv10);
655             r10              = _mm256_andnot_ps(dummy_mask,r10);
656
657             /* Compute parameters for interactions between i and j atoms */
658             qq10             = _mm256_mul_ps(iq1,jq0);
659
660             /* EWALD ELECTROSTATICS */
661             
662             /* Analytical PME correction */
663             zeta2            = _mm256_mul_ps(beta2,rsq10);
664             rinv3            = _mm256_mul_ps(rinvsq10,rinv10);
665             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
666             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
667             felec            = _mm256_mul_ps(qq10,felec);
668             pmecorrV         = gmx_mm256_pmecorrV_ps(zeta2);
669             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
670             velec            = _mm256_sub_ps(rinv10,pmecorrV);
671             velec            = _mm256_mul_ps(qq10,velec);
672             
673             /* Update potential sum for this i atom from the interaction with this j atom. */
674             velec            = _mm256_andnot_ps(dummy_mask,velec);
675             velecsum         = _mm256_add_ps(velecsum,velec);
676
677             fscal            = felec;
678
679             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
680
681             /* Calculate temporary vectorial force */
682             tx               = _mm256_mul_ps(fscal,dx10);
683             ty               = _mm256_mul_ps(fscal,dy10);
684             tz               = _mm256_mul_ps(fscal,dz10);
685
686             /* Update vectorial force */
687             fix1             = _mm256_add_ps(fix1,tx);
688             fiy1             = _mm256_add_ps(fiy1,ty);
689             fiz1             = _mm256_add_ps(fiz1,tz);
690
691             fjx0             = _mm256_add_ps(fjx0,tx);
692             fjy0             = _mm256_add_ps(fjy0,ty);
693             fjz0             = _mm256_add_ps(fjz0,tz);
694
695             /**************************
696              * CALCULATE INTERACTIONS *
697              **************************/
698
699             r20              = _mm256_mul_ps(rsq20,rinv20);
700             r20              = _mm256_andnot_ps(dummy_mask,r20);
701
702             /* Compute parameters for interactions between i and j atoms */
703             qq20             = _mm256_mul_ps(iq2,jq0);
704
705             /* EWALD ELECTROSTATICS */
706             
707             /* Analytical PME correction */
708             zeta2            = _mm256_mul_ps(beta2,rsq20);
709             rinv3            = _mm256_mul_ps(rinvsq20,rinv20);
710             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
711             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
712             felec            = _mm256_mul_ps(qq20,felec);
713             pmecorrV         = gmx_mm256_pmecorrV_ps(zeta2);
714             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
715             velec            = _mm256_sub_ps(rinv20,pmecorrV);
716             velec            = _mm256_mul_ps(qq20,velec);
717             
718             /* Update potential sum for this i atom from the interaction with this j atom. */
719             velec            = _mm256_andnot_ps(dummy_mask,velec);
720             velecsum         = _mm256_add_ps(velecsum,velec);
721
722             fscal            = felec;
723
724             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
725
726             /* Calculate temporary vectorial force */
727             tx               = _mm256_mul_ps(fscal,dx20);
728             ty               = _mm256_mul_ps(fscal,dy20);
729             tz               = _mm256_mul_ps(fscal,dz20);
730
731             /* Update vectorial force */
732             fix2             = _mm256_add_ps(fix2,tx);
733             fiy2             = _mm256_add_ps(fiy2,ty);
734             fiz2             = _mm256_add_ps(fiz2,tz);
735
736             fjx0             = _mm256_add_ps(fjx0,tx);
737             fjy0             = _mm256_add_ps(fjy0,ty);
738             fjz0             = _mm256_add_ps(fjz0,tz);
739
740             /**************************
741              * CALCULATE INTERACTIONS *
742              **************************/
743
744             r30              = _mm256_mul_ps(rsq30,rinv30);
745             r30              = _mm256_andnot_ps(dummy_mask,r30);
746
747             /* Compute parameters for interactions between i and j atoms */
748             qq30             = _mm256_mul_ps(iq3,jq0);
749
750             /* EWALD ELECTROSTATICS */
751             
752             /* Analytical PME correction */
753             zeta2            = _mm256_mul_ps(beta2,rsq30);
754             rinv3            = _mm256_mul_ps(rinvsq30,rinv30);
755             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
756             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
757             felec            = _mm256_mul_ps(qq30,felec);
758             pmecorrV         = gmx_mm256_pmecorrV_ps(zeta2);
759             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
760             velec            = _mm256_sub_ps(rinv30,pmecorrV);
761             velec            = _mm256_mul_ps(qq30,velec);
762             
763             /* Update potential sum for this i atom from the interaction with this j atom. */
764             velec            = _mm256_andnot_ps(dummy_mask,velec);
765             velecsum         = _mm256_add_ps(velecsum,velec);
766
767             fscal            = felec;
768
769             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
770
771             /* Calculate temporary vectorial force */
772             tx               = _mm256_mul_ps(fscal,dx30);
773             ty               = _mm256_mul_ps(fscal,dy30);
774             tz               = _mm256_mul_ps(fscal,dz30);
775
776             /* Update vectorial force */
777             fix3             = _mm256_add_ps(fix3,tx);
778             fiy3             = _mm256_add_ps(fiy3,ty);
779             fiz3             = _mm256_add_ps(fiz3,tz);
780
781             fjx0             = _mm256_add_ps(fjx0,tx);
782             fjy0             = _mm256_add_ps(fjy0,ty);
783             fjz0             = _mm256_add_ps(fjz0,tz);
784
785             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
786             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
787             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
788             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
789             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
790             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
791             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
792             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
793
794             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
795
796             /* Inner loop uses 310 flops */
797         }
798
799         /* End of innermost loop */
800
801         gmx_mm256_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
802                                                  f+i_coord_offset,fshift+i_shift_offset);
803
804         ggid                        = gid[iidx];
805         /* Update potential energies */
806         gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
807         gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
808
809         /* Increment number of inner iterations */
810         inneriter                  += j_index_end - j_index_start;
811
812         /* Outer loop uses 26 flops */
813     }
814
815     /* Increment number of outer iterations */
816     outeriter        += nri;
817
818     /* Update outer/inner flops */
819
820     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*310);
821 }
822 /*
823  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwLJEw_GeomW4P1_F_avx_256_single
824  * Electrostatics interaction: Ewald
825  * VdW interaction:            LJEwald
826  * Geometry:                   Water4-Particle
827  * Calculate force/pot:        Force
828  */
829 void
830 nb_kernel_ElecEw_VdwLJEw_GeomW4P1_F_avx_256_single
831                     (t_nblist                    * gmx_restrict       nlist,
832                      rvec                        * gmx_restrict          xx,
833                      rvec                        * gmx_restrict          ff,
834                      t_forcerec                  * gmx_restrict          fr,
835                      t_mdatoms                   * gmx_restrict     mdatoms,
836                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
837                      t_nrnb                      * gmx_restrict        nrnb)
838 {
839     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
840      * just 0 for non-waters.
841      * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
842      * jnr indices corresponding to data put in the four positions in the SIMD register.
843      */
844     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
845     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
846     int              jnrA,jnrB,jnrC,jnrD;
847     int              jnrE,jnrF,jnrG,jnrH;
848     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
849     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
850     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
851     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
852     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
853     real             rcutoff_scalar;
854     real             *shiftvec,*fshift,*x,*f;
855     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
856     real             scratch[4*DIM];
857     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
858     real *           vdwioffsetptr0;
859     real *           vdwgridioffsetptr0;
860     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
861     real *           vdwioffsetptr1;
862     real *           vdwgridioffsetptr1;
863     __m256           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
864     real *           vdwioffsetptr2;
865     real *           vdwgridioffsetptr2;
866     __m256           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
867     real *           vdwioffsetptr3;
868     real *           vdwgridioffsetptr3;
869     __m256           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
870     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
871     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
872     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
873     __m256           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
874     __m256           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
875     __m256           dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
876     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
877     real             *charge;
878     int              nvdwtype;
879     __m256           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
880     int              *vdwtype;
881     real             *vdwparam;
882     __m256           one_sixth   = _mm256_set1_ps(1.0/6.0);
883     __m256           one_twelfth = _mm256_set1_ps(1.0/12.0);
884     __m256           c6grid_00;
885     __m256           c6grid_10;
886     __m256           c6grid_20;
887     __m256           c6grid_30;
888     real             *vdwgridparam;
889     __m256           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
890     __m256           one_half  = _mm256_set1_ps(0.5);
891     __m256           minus_one = _mm256_set1_ps(-1.0);
892     __m256i          ewitab;
893     __m128i          ewitab_lo,ewitab_hi;
894     __m256           ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
895     __m256           beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
896     real             *ewtab;
897     __m256           dummy_mask,cutoff_mask;
898     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
899     __m256           one     = _mm256_set1_ps(1.0);
900     __m256           two     = _mm256_set1_ps(2.0);
901     x                = xx[0];
902     f                = ff[0];
903
904     nri              = nlist->nri;
905     iinr             = nlist->iinr;
906     jindex           = nlist->jindex;
907     jjnr             = nlist->jjnr;
908     shiftidx         = nlist->shift;
909     gid              = nlist->gid;
910     shiftvec         = fr->shift_vec[0];
911     fshift           = fr->fshift[0];
912     facel            = _mm256_set1_ps(fr->epsfac);
913     charge           = mdatoms->chargeA;
914     nvdwtype         = fr->ntype;
915     vdwparam         = fr->nbfp;
916     vdwtype          = mdatoms->typeA;
917     vdwgridparam     = fr->ljpme_c6grid;
918     sh_lj_ewald      = _mm256_set1_ps(fr->ic->sh_lj_ewald);
919     ewclj            = _mm256_set1_ps(fr->ewaldcoeff_lj);
920     ewclj2           = _mm256_mul_ps(minus_one,_mm256_mul_ps(ewclj,ewclj));
921
922     sh_ewald         = _mm256_set1_ps(fr->ic->sh_ewald);
923     beta             = _mm256_set1_ps(fr->ic->ewaldcoeff_q);
924     beta2            = _mm256_mul_ps(beta,beta);
925     beta3            = _mm256_mul_ps(beta,beta2);
926
927     ewtab            = fr->ic->tabq_coul_F;
928     ewtabscale       = _mm256_set1_ps(fr->ic->tabq_scale);
929     ewtabhalfspace   = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
930
931     /* Setup water-specific parameters */
932     inr              = nlist->iinr[0];
933     iq1              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
934     iq2              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
935     iq3              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+3]));
936     vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
937     vdwgridioffsetptr0 = vdwgridparam+2*nvdwtype*vdwtype[inr+0];
938
939     /* Avoid stupid compiler warnings */
940     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
941     j_coord_offsetA = 0;
942     j_coord_offsetB = 0;
943     j_coord_offsetC = 0;
944     j_coord_offsetD = 0;
945     j_coord_offsetE = 0;
946     j_coord_offsetF = 0;
947     j_coord_offsetG = 0;
948     j_coord_offsetH = 0;
949
950     outeriter        = 0;
951     inneriter        = 0;
952
953     for(iidx=0;iidx<4*DIM;iidx++)
954     {
955         scratch[iidx] = 0.0;
956     }
957
958     /* Start outer loop over neighborlists */
959     for(iidx=0; iidx<nri; iidx++)
960     {
961         /* Load shift vector for this list */
962         i_shift_offset   = DIM*shiftidx[iidx];
963
964         /* Load limits for loop over neighbors */
965         j_index_start    = jindex[iidx];
966         j_index_end      = jindex[iidx+1];
967
968         /* Get outer coordinate index */
969         inr              = iinr[iidx];
970         i_coord_offset   = DIM*inr;
971
972         /* Load i particle coords and add shift vector */
973         gmx_mm256_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
974                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
975
976         fix0             = _mm256_setzero_ps();
977         fiy0             = _mm256_setzero_ps();
978         fiz0             = _mm256_setzero_ps();
979         fix1             = _mm256_setzero_ps();
980         fiy1             = _mm256_setzero_ps();
981         fiz1             = _mm256_setzero_ps();
982         fix2             = _mm256_setzero_ps();
983         fiy2             = _mm256_setzero_ps();
984         fiz2             = _mm256_setzero_ps();
985         fix3             = _mm256_setzero_ps();
986         fiy3             = _mm256_setzero_ps();
987         fiz3             = _mm256_setzero_ps();
988
989         /* Start inner kernel loop */
990         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
991         {
992
993             /* Get j neighbor index, and coordinate index */
994             jnrA             = jjnr[jidx];
995             jnrB             = jjnr[jidx+1];
996             jnrC             = jjnr[jidx+2];
997             jnrD             = jjnr[jidx+3];
998             jnrE             = jjnr[jidx+4];
999             jnrF             = jjnr[jidx+5];
1000             jnrG             = jjnr[jidx+6];
1001             jnrH             = jjnr[jidx+7];
1002             j_coord_offsetA  = DIM*jnrA;
1003             j_coord_offsetB  = DIM*jnrB;
1004             j_coord_offsetC  = DIM*jnrC;
1005             j_coord_offsetD  = DIM*jnrD;
1006             j_coord_offsetE  = DIM*jnrE;
1007             j_coord_offsetF  = DIM*jnrF;
1008             j_coord_offsetG  = DIM*jnrG;
1009             j_coord_offsetH  = DIM*jnrH;
1010
1011             /* load j atom coordinates */
1012             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1013                                                  x+j_coord_offsetC,x+j_coord_offsetD,
1014                                                  x+j_coord_offsetE,x+j_coord_offsetF,
1015                                                  x+j_coord_offsetG,x+j_coord_offsetH,
1016                                                  &jx0,&jy0,&jz0);
1017
1018             /* Calculate displacement vector */
1019             dx00             = _mm256_sub_ps(ix0,jx0);
1020             dy00             = _mm256_sub_ps(iy0,jy0);
1021             dz00             = _mm256_sub_ps(iz0,jz0);
1022             dx10             = _mm256_sub_ps(ix1,jx0);
1023             dy10             = _mm256_sub_ps(iy1,jy0);
1024             dz10             = _mm256_sub_ps(iz1,jz0);
1025             dx20             = _mm256_sub_ps(ix2,jx0);
1026             dy20             = _mm256_sub_ps(iy2,jy0);
1027             dz20             = _mm256_sub_ps(iz2,jz0);
1028             dx30             = _mm256_sub_ps(ix3,jx0);
1029             dy30             = _mm256_sub_ps(iy3,jy0);
1030             dz30             = _mm256_sub_ps(iz3,jz0);
1031
1032             /* Calculate squared distance and things based on it */
1033             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
1034             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
1035             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
1036             rsq30            = gmx_mm256_calc_rsq_ps(dx30,dy30,dz30);
1037
1038             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
1039             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
1040             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
1041             rinv30           = gmx_mm256_invsqrt_ps(rsq30);
1042
1043             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
1044             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
1045             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
1046             rinvsq30         = _mm256_mul_ps(rinv30,rinv30);
1047
1048             /* Load parameters for j particles */
1049             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1050                                                                  charge+jnrC+0,charge+jnrD+0,
1051                                                                  charge+jnrE+0,charge+jnrF+0,
1052                                                                  charge+jnrG+0,charge+jnrH+0);
1053             vdwjidx0A        = 2*vdwtype[jnrA+0];
1054             vdwjidx0B        = 2*vdwtype[jnrB+0];
1055             vdwjidx0C        = 2*vdwtype[jnrC+0];
1056             vdwjidx0D        = 2*vdwtype[jnrD+0];
1057             vdwjidx0E        = 2*vdwtype[jnrE+0];
1058             vdwjidx0F        = 2*vdwtype[jnrF+0];
1059             vdwjidx0G        = 2*vdwtype[jnrG+0];
1060             vdwjidx0H        = 2*vdwtype[jnrH+0];
1061
1062             fjx0             = _mm256_setzero_ps();
1063             fjy0             = _mm256_setzero_ps();
1064             fjz0             = _mm256_setzero_ps();
1065
1066             /**************************
1067              * CALCULATE INTERACTIONS *
1068              **************************/
1069
1070             r00              = _mm256_mul_ps(rsq00,rinv00);
1071
1072             /* Compute parameters for interactions between i and j atoms */
1073             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
1074                                             vdwioffsetptr0+vdwjidx0B,
1075                                             vdwioffsetptr0+vdwjidx0C,
1076                                             vdwioffsetptr0+vdwjidx0D,
1077                                             vdwioffsetptr0+vdwjidx0E,
1078                                             vdwioffsetptr0+vdwjidx0F,
1079                                             vdwioffsetptr0+vdwjidx0G,
1080                                             vdwioffsetptr0+vdwjidx0H,
1081                                             &c6_00,&c12_00);
1082
1083             c6grid_00       = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
1084                                                                   vdwgridioffsetptr0+vdwjidx0B,
1085                                                                   vdwgridioffsetptr0+vdwjidx0C,
1086                                                                   vdwgridioffsetptr0+vdwjidx0D,
1087                                                                   vdwgridioffsetptr0+vdwjidx0E,
1088                                                                   vdwgridioffsetptr0+vdwjidx0F,
1089                                                                   vdwgridioffsetptr0+vdwjidx0G,
1090                                                                   vdwgridioffsetptr0+vdwjidx0H);
1091
1092             /* Analytical LJ-PME */
1093             rinvsix          = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
1094             ewcljrsq         = _mm256_mul_ps(ewclj2,rsq00);
1095             ewclj6           = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
1096             exponent         = gmx_simd_exp_r(ewcljrsq);
1097             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
1098             poly             = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
1099             /* f6A = 6 * C6grid * (1 - poly) */
1100             f6A              = _mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly));
1101             /* f6B = C6grid * exponent * beta^6 */
1102             f6B              = _mm256_mul_ps(_mm256_mul_ps(c6grid_00,one_sixth),_mm256_mul_ps(exponent,ewclj6));
1103             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
1104             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);
1105
1106             fscal            = fvdw;
1107
1108             /* Calculate temporary vectorial force */
1109             tx               = _mm256_mul_ps(fscal,dx00);
1110             ty               = _mm256_mul_ps(fscal,dy00);
1111             tz               = _mm256_mul_ps(fscal,dz00);
1112
1113             /* Update vectorial force */
1114             fix0             = _mm256_add_ps(fix0,tx);
1115             fiy0             = _mm256_add_ps(fiy0,ty);
1116             fiz0             = _mm256_add_ps(fiz0,tz);
1117
1118             fjx0             = _mm256_add_ps(fjx0,tx);
1119             fjy0             = _mm256_add_ps(fjy0,ty);
1120             fjz0             = _mm256_add_ps(fjz0,tz);
1121
1122             /**************************
1123              * CALCULATE INTERACTIONS *
1124              **************************/
1125
1126             r10              = _mm256_mul_ps(rsq10,rinv10);
1127
1128             /* Compute parameters for interactions between i and j atoms */
1129             qq10             = _mm256_mul_ps(iq1,jq0);
1130
1131             /* EWALD ELECTROSTATICS */
1132             
1133             /* Analytical PME correction */
1134             zeta2            = _mm256_mul_ps(beta2,rsq10);
1135             rinv3            = _mm256_mul_ps(rinvsq10,rinv10);
1136             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
1137             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1138             felec            = _mm256_mul_ps(qq10,felec);
1139             
1140             fscal            = felec;
1141
1142             /* Calculate temporary vectorial force */
1143             tx               = _mm256_mul_ps(fscal,dx10);
1144             ty               = _mm256_mul_ps(fscal,dy10);
1145             tz               = _mm256_mul_ps(fscal,dz10);
1146
1147             /* Update vectorial force */
1148             fix1             = _mm256_add_ps(fix1,tx);
1149             fiy1             = _mm256_add_ps(fiy1,ty);
1150             fiz1             = _mm256_add_ps(fiz1,tz);
1151
1152             fjx0             = _mm256_add_ps(fjx0,tx);
1153             fjy0             = _mm256_add_ps(fjy0,ty);
1154             fjz0             = _mm256_add_ps(fjz0,tz);
1155
1156             /**************************
1157              * CALCULATE INTERACTIONS *
1158              **************************/
1159
1160             r20              = _mm256_mul_ps(rsq20,rinv20);
1161
1162             /* Compute parameters for interactions between i and j atoms */
1163             qq20             = _mm256_mul_ps(iq2,jq0);
1164
1165             /* EWALD ELECTROSTATICS */
1166             
1167             /* Analytical PME correction */
1168             zeta2            = _mm256_mul_ps(beta2,rsq20);
1169             rinv3            = _mm256_mul_ps(rinvsq20,rinv20);
1170             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
1171             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1172             felec            = _mm256_mul_ps(qq20,felec);
1173             
1174             fscal            = felec;
1175
1176             /* Calculate temporary vectorial force */
1177             tx               = _mm256_mul_ps(fscal,dx20);
1178             ty               = _mm256_mul_ps(fscal,dy20);
1179             tz               = _mm256_mul_ps(fscal,dz20);
1180
1181             /* Update vectorial force */
1182             fix2             = _mm256_add_ps(fix2,tx);
1183             fiy2             = _mm256_add_ps(fiy2,ty);
1184             fiz2             = _mm256_add_ps(fiz2,tz);
1185
1186             fjx0             = _mm256_add_ps(fjx0,tx);
1187             fjy0             = _mm256_add_ps(fjy0,ty);
1188             fjz0             = _mm256_add_ps(fjz0,tz);
1189
1190             /**************************
1191              * CALCULATE INTERACTIONS *
1192              **************************/
1193
1194             r30              = _mm256_mul_ps(rsq30,rinv30);
1195
1196             /* Compute parameters for interactions between i and j atoms */
1197             qq30             = _mm256_mul_ps(iq3,jq0);
1198
1199             /* EWALD ELECTROSTATICS */
1200             
1201             /* Analytical PME correction */
1202             zeta2            = _mm256_mul_ps(beta2,rsq30);
1203             rinv3            = _mm256_mul_ps(rinvsq30,rinv30);
1204             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
1205             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1206             felec            = _mm256_mul_ps(qq30,felec);
1207             
1208             fscal            = felec;
1209
1210             /* Calculate temporary vectorial force */
1211             tx               = _mm256_mul_ps(fscal,dx30);
1212             ty               = _mm256_mul_ps(fscal,dy30);
1213             tz               = _mm256_mul_ps(fscal,dz30);
1214
1215             /* Update vectorial force */
1216             fix3             = _mm256_add_ps(fix3,tx);
1217             fiy3             = _mm256_add_ps(fiy3,ty);
1218             fiz3             = _mm256_add_ps(fiz3,tz);
1219
1220             fjx0             = _mm256_add_ps(fjx0,tx);
1221             fjy0             = _mm256_add_ps(fjy0,ty);
1222             fjz0             = _mm256_add_ps(fjz0,tz);
1223
1224             fjptrA             = f+j_coord_offsetA;
1225             fjptrB             = f+j_coord_offsetB;
1226             fjptrC             = f+j_coord_offsetC;
1227             fjptrD             = f+j_coord_offsetD;
1228             fjptrE             = f+j_coord_offsetE;
1229             fjptrF             = f+j_coord_offsetF;
1230             fjptrG             = f+j_coord_offsetG;
1231             fjptrH             = f+j_coord_offsetH;
1232
1233             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
1234
1235             /* Inner loop uses 217 flops */
1236         }
1237
1238         if(jidx<j_index_end)
1239         {
1240
1241             /* Get j neighbor index, and coordinate index */
1242             jnrlistA         = jjnr[jidx];
1243             jnrlistB         = jjnr[jidx+1];
1244             jnrlistC         = jjnr[jidx+2];
1245             jnrlistD         = jjnr[jidx+3];
1246             jnrlistE         = jjnr[jidx+4];
1247             jnrlistF         = jjnr[jidx+5];
1248             jnrlistG         = jjnr[jidx+6];
1249             jnrlistH         = jjnr[jidx+7];
1250             /* Sign of each element will be negative for non-real atoms.
1251              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1252              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1253              */
1254             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
1255                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
1256                                             
1257             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
1258             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
1259             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
1260             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
1261             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
1262             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
1263             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
1264             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
1265             j_coord_offsetA  = DIM*jnrA;
1266             j_coord_offsetB  = DIM*jnrB;
1267             j_coord_offsetC  = DIM*jnrC;
1268             j_coord_offsetD  = DIM*jnrD;
1269             j_coord_offsetE  = DIM*jnrE;
1270             j_coord_offsetF  = DIM*jnrF;
1271             j_coord_offsetG  = DIM*jnrG;
1272             j_coord_offsetH  = DIM*jnrH;
1273
1274             /* load j atom coordinates */
1275             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1276                                                  x+j_coord_offsetC,x+j_coord_offsetD,
1277                                                  x+j_coord_offsetE,x+j_coord_offsetF,
1278                                                  x+j_coord_offsetG,x+j_coord_offsetH,
1279                                                  &jx0,&jy0,&jz0);
1280
1281             /* Calculate displacement vector */
1282             dx00             = _mm256_sub_ps(ix0,jx0);
1283             dy00             = _mm256_sub_ps(iy0,jy0);
1284             dz00             = _mm256_sub_ps(iz0,jz0);
1285             dx10             = _mm256_sub_ps(ix1,jx0);
1286             dy10             = _mm256_sub_ps(iy1,jy0);
1287             dz10             = _mm256_sub_ps(iz1,jz0);
1288             dx20             = _mm256_sub_ps(ix2,jx0);
1289             dy20             = _mm256_sub_ps(iy2,jy0);
1290             dz20             = _mm256_sub_ps(iz2,jz0);
1291             dx30             = _mm256_sub_ps(ix3,jx0);
1292             dy30             = _mm256_sub_ps(iy3,jy0);
1293             dz30             = _mm256_sub_ps(iz3,jz0);
1294
1295             /* Calculate squared distance and things based on it */
1296             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
1297             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
1298             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
1299             rsq30            = gmx_mm256_calc_rsq_ps(dx30,dy30,dz30);
1300
1301             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
1302             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
1303             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
1304             rinv30           = gmx_mm256_invsqrt_ps(rsq30);
1305
1306             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
1307             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
1308             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
1309             rinvsq30         = _mm256_mul_ps(rinv30,rinv30);
1310
1311             /* Load parameters for j particles */
1312             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1313                                                                  charge+jnrC+0,charge+jnrD+0,
1314                                                                  charge+jnrE+0,charge+jnrF+0,
1315                                                                  charge+jnrG+0,charge+jnrH+0);
1316             vdwjidx0A        = 2*vdwtype[jnrA+0];
1317             vdwjidx0B        = 2*vdwtype[jnrB+0];
1318             vdwjidx0C        = 2*vdwtype[jnrC+0];
1319             vdwjidx0D        = 2*vdwtype[jnrD+0];
1320             vdwjidx0E        = 2*vdwtype[jnrE+0];
1321             vdwjidx0F        = 2*vdwtype[jnrF+0];
1322             vdwjidx0G        = 2*vdwtype[jnrG+0];
1323             vdwjidx0H        = 2*vdwtype[jnrH+0];
1324
1325             fjx0             = _mm256_setzero_ps();
1326             fjy0             = _mm256_setzero_ps();
1327             fjz0             = _mm256_setzero_ps();
1328
1329             /**************************
1330              * CALCULATE INTERACTIONS *
1331              **************************/
1332
1333             r00              = _mm256_mul_ps(rsq00,rinv00);
1334             r00              = _mm256_andnot_ps(dummy_mask,r00);
1335
1336             /* Compute parameters for interactions between i and j atoms */
1337             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
1338                                             vdwioffsetptr0+vdwjidx0B,
1339                                             vdwioffsetptr0+vdwjidx0C,
1340                                             vdwioffsetptr0+vdwjidx0D,
1341                                             vdwioffsetptr0+vdwjidx0E,
1342                                             vdwioffsetptr0+vdwjidx0F,
1343                                             vdwioffsetptr0+vdwjidx0G,
1344                                             vdwioffsetptr0+vdwjidx0H,
1345                                             &c6_00,&c12_00);
1346
1347             c6grid_00       = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
1348                                                                   vdwgridioffsetptr0+vdwjidx0B,
1349                                                                   vdwgridioffsetptr0+vdwjidx0C,
1350                                                                   vdwgridioffsetptr0+vdwjidx0D,
1351                                                                   vdwgridioffsetptr0+vdwjidx0E,
1352                                                                   vdwgridioffsetptr0+vdwjidx0F,
1353                                                                   vdwgridioffsetptr0+vdwjidx0G,
1354                                                                   vdwgridioffsetptr0+vdwjidx0H);
1355
1356             /* Analytical LJ-PME */
1357             rinvsix          = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
1358             ewcljrsq         = _mm256_mul_ps(ewclj2,rsq00);
1359             ewclj6           = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
1360             exponent         = gmx_simd_exp_r(ewcljrsq);
1361             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
1362             poly             = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
1363             /* f6A = 6 * C6grid * (1 - poly) */
1364             f6A              = _mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly));
1365             /* f6B = C6grid * exponent * beta^6 */
1366             f6B              = _mm256_mul_ps(_mm256_mul_ps(c6grid_00,one_sixth),_mm256_mul_ps(exponent,ewclj6));
1367             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
1368             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);
1369
1370             fscal            = fvdw;
1371
1372             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1373
1374             /* Calculate temporary vectorial force */
1375             tx               = _mm256_mul_ps(fscal,dx00);
1376             ty               = _mm256_mul_ps(fscal,dy00);
1377             tz               = _mm256_mul_ps(fscal,dz00);
1378
1379             /* Update vectorial force */
1380             fix0             = _mm256_add_ps(fix0,tx);
1381             fiy0             = _mm256_add_ps(fiy0,ty);
1382             fiz0             = _mm256_add_ps(fiz0,tz);
1383
1384             fjx0             = _mm256_add_ps(fjx0,tx);
1385             fjy0             = _mm256_add_ps(fjy0,ty);
1386             fjz0             = _mm256_add_ps(fjz0,tz);
1387
1388             /**************************
1389              * CALCULATE INTERACTIONS *
1390              **************************/
1391
1392             r10              = _mm256_mul_ps(rsq10,rinv10);
1393             r10              = _mm256_andnot_ps(dummy_mask,r10);
1394
1395             /* Compute parameters for interactions between i and j atoms */
1396             qq10             = _mm256_mul_ps(iq1,jq0);
1397
1398             /* EWALD ELECTROSTATICS */
1399             
1400             /* Analytical PME correction */
1401             zeta2            = _mm256_mul_ps(beta2,rsq10);
1402             rinv3            = _mm256_mul_ps(rinvsq10,rinv10);
1403             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
1404             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1405             felec            = _mm256_mul_ps(qq10,felec);
1406             
1407             fscal            = felec;
1408
1409             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1410
1411             /* Calculate temporary vectorial force */
1412             tx               = _mm256_mul_ps(fscal,dx10);
1413             ty               = _mm256_mul_ps(fscal,dy10);
1414             tz               = _mm256_mul_ps(fscal,dz10);
1415
1416             /* Update vectorial force */
1417             fix1             = _mm256_add_ps(fix1,tx);
1418             fiy1             = _mm256_add_ps(fiy1,ty);
1419             fiz1             = _mm256_add_ps(fiz1,tz);
1420
1421             fjx0             = _mm256_add_ps(fjx0,tx);
1422             fjy0             = _mm256_add_ps(fjy0,ty);
1423             fjz0             = _mm256_add_ps(fjz0,tz);
1424
1425             /**************************
1426              * CALCULATE INTERACTIONS *
1427              **************************/
1428
1429             r20              = _mm256_mul_ps(rsq20,rinv20);
1430             r20              = _mm256_andnot_ps(dummy_mask,r20);
1431
1432             /* Compute parameters for interactions between i and j atoms */
1433             qq20             = _mm256_mul_ps(iq2,jq0);
1434
1435             /* EWALD ELECTROSTATICS */
1436             
1437             /* Analytical PME correction */
1438             zeta2            = _mm256_mul_ps(beta2,rsq20);
1439             rinv3            = _mm256_mul_ps(rinvsq20,rinv20);
1440             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
1441             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1442             felec            = _mm256_mul_ps(qq20,felec);
1443             
1444             fscal            = felec;
1445
1446             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1447
1448             /* Calculate temporary vectorial force */
1449             tx               = _mm256_mul_ps(fscal,dx20);
1450             ty               = _mm256_mul_ps(fscal,dy20);
1451             tz               = _mm256_mul_ps(fscal,dz20);
1452
1453             /* Update vectorial force */
1454             fix2             = _mm256_add_ps(fix2,tx);
1455             fiy2             = _mm256_add_ps(fiy2,ty);
1456             fiz2             = _mm256_add_ps(fiz2,tz);
1457
1458             fjx0             = _mm256_add_ps(fjx0,tx);
1459             fjy0             = _mm256_add_ps(fjy0,ty);
1460             fjz0             = _mm256_add_ps(fjz0,tz);
1461
1462             /**************************
1463              * CALCULATE INTERACTIONS *
1464              **************************/
1465
1466             r30              = _mm256_mul_ps(rsq30,rinv30);
1467             r30              = _mm256_andnot_ps(dummy_mask,r30);
1468
1469             /* Compute parameters for interactions between i and j atoms */
1470             qq30             = _mm256_mul_ps(iq3,jq0);
1471
1472             /* EWALD ELECTROSTATICS */
1473             
1474             /* Analytical PME correction */
1475             zeta2            = _mm256_mul_ps(beta2,rsq30);
1476             rinv3            = _mm256_mul_ps(rinvsq30,rinv30);
1477             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
1478             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1479             felec            = _mm256_mul_ps(qq30,felec);
1480             
1481             fscal            = felec;
1482
1483             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1484
1485             /* Calculate temporary vectorial force */
1486             tx               = _mm256_mul_ps(fscal,dx30);
1487             ty               = _mm256_mul_ps(fscal,dy30);
1488             tz               = _mm256_mul_ps(fscal,dz30);
1489
1490             /* Update vectorial force */
1491             fix3             = _mm256_add_ps(fix3,tx);
1492             fiy3             = _mm256_add_ps(fiy3,ty);
1493             fiz3             = _mm256_add_ps(fiz3,tz);
1494
1495             fjx0             = _mm256_add_ps(fjx0,tx);
1496             fjy0             = _mm256_add_ps(fjy0,ty);
1497             fjz0             = _mm256_add_ps(fjz0,tz);
1498
1499             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1500             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1501             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1502             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1503             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1504             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1505             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1506             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1507
1508             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
1509
1510             /* Inner loop uses 221 flops */
1511         }
1512
1513         /* End of innermost loop */
1514
1515         gmx_mm256_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1516                                                  f+i_coord_offset,fshift+i_shift_offset);
1517
1518         /* Increment number of inner iterations */
1519         inneriter                  += j_index_end - j_index_start;
1520
1521         /* Outer loop uses 24 flops */
1522     }
1523
1524     /* Increment number of outer iterations */
1525     outeriter        += nri;
1526
1527     /* Update outer/inner flops */
1528
1529     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*221);
1530 }