Remove no-inline-max-size and suppress remark
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_single / nb_kernel_ElecEwSh_VdwLJEwSh_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 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "types/simple.h"
46 #include "vec.h"
47 #include "nrnb.h"
48
49 #include "gromacs/simd/math_x86_avx_256_single.h"
50 #include "kernelutil_x86_avx_256_single.h"
51
52 /*
53  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwLJEwSh_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_ElecEwSh_VdwLJEwSh_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     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
170     rcutoff_scalar   = fr->rcoulomb;
171     rcutoff          = _mm256_set1_ps(rcutoff_scalar);
172     rcutoff2         = _mm256_mul_ps(rcutoff,rcutoff);
173
174     sh_vdw_invrcut6  = _mm256_set1_ps(fr->ic->sh_invrc6);
175     rvdw             = _mm256_set1_ps(fr->rvdw);
176
177     /* Avoid stupid compiler warnings */
178     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
179     j_coord_offsetA = 0;
180     j_coord_offsetB = 0;
181     j_coord_offsetC = 0;
182     j_coord_offsetD = 0;
183     j_coord_offsetE = 0;
184     j_coord_offsetF = 0;
185     j_coord_offsetG = 0;
186     j_coord_offsetH = 0;
187
188     outeriter        = 0;
189     inneriter        = 0;
190
191     for(iidx=0;iidx<4*DIM;iidx++)
192     {
193         scratch[iidx] = 0.0;
194     }
195
196     /* Start outer loop over neighborlists */
197     for(iidx=0; iidx<nri; iidx++)
198     {
199         /* Load shift vector for this list */
200         i_shift_offset   = DIM*shiftidx[iidx];
201
202         /* Load limits for loop over neighbors */
203         j_index_start    = jindex[iidx];
204         j_index_end      = jindex[iidx+1];
205
206         /* Get outer coordinate index */
207         inr              = iinr[iidx];
208         i_coord_offset   = DIM*inr;
209
210         /* Load i particle coords and add shift vector */
211         gmx_mm256_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
212                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
213
214         fix0             = _mm256_setzero_ps();
215         fiy0             = _mm256_setzero_ps();
216         fiz0             = _mm256_setzero_ps();
217         fix1             = _mm256_setzero_ps();
218         fiy1             = _mm256_setzero_ps();
219         fiz1             = _mm256_setzero_ps();
220         fix2             = _mm256_setzero_ps();
221         fiy2             = _mm256_setzero_ps();
222         fiz2             = _mm256_setzero_ps();
223         fix3             = _mm256_setzero_ps();
224         fiy3             = _mm256_setzero_ps();
225         fiz3             = _mm256_setzero_ps();
226
227         /* Reset potential sums */
228         velecsum         = _mm256_setzero_ps();
229         vvdwsum          = _mm256_setzero_ps();
230
231         /* Start inner kernel loop */
232         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
233         {
234
235             /* Get j neighbor index, and coordinate index */
236             jnrA             = jjnr[jidx];
237             jnrB             = jjnr[jidx+1];
238             jnrC             = jjnr[jidx+2];
239             jnrD             = jjnr[jidx+3];
240             jnrE             = jjnr[jidx+4];
241             jnrF             = jjnr[jidx+5];
242             jnrG             = jjnr[jidx+6];
243             jnrH             = jjnr[jidx+7];
244             j_coord_offsetA  = DIM*jnrA;
245             j_coord_offsetB  = DIM*jnrB;
246             j_coord_offsetC  = DIM*jnrC;
247             j_coord_offsetD  = DIM*jnrD;
248             j_coord_offsetE  = DIM*jnrE;
249             j_coord_offsetF  = DIM*jnrF;
250             j_coord_offsetG  = DIM*jnrG;
251             j_coord_offsetH  = DIM*jnrH;
252
253             /* load j atom coordinates */
254             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
255                                                  x+j_coord_offsetC,x+j_coord_offsetD,
256                                                  x+j_coord_offsetE,x+j_coord_offsetF,
257                                                  x+j_coord_offsetG,x+j_coord_offsetH,
258                                                  &jx0,&jy0,&jz0);
259
260             /* Calculate displacement vector */
261             dx00             = _mm256_sub_ps(ix0,jx0);
262             dy00             = _mm256_sub_ps(iy0,jy0);
263             dz00             = _mm256_sub_ps(iz0,jz0);
264             dx10             = _mm256_sub_ps(ix1,jx0);
265             dy10             = _mm256_sub_ps(iy1,jy0);
266             dz10             = _mm256_sub_ps(iz1,jz0);
267             dx20             = _mm256_sub_ps(ix2,jx0);
268             dy20             = _mm256_sub_ps(iy2,jy0);
269             dz20             = _mm256_sub_ps(iz2,jz0);
270             dx30             = _mm256_sub_ps(ix3,jx0);
271             dy30             = _mm256_sub_ps(iy3,jy0);
272             dz30             = _mm256_sub_ps(iz3,jz0);
273
274             /* Calculate squared distance and things based on it */
275             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
276             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
277             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
278             rsq30            = gmx_mm256_calc_rsq_ps(dx30,dy30,dz30);
279
280             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
281             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
282             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
283             rinv30           = gmx_mm256_invsqrt_ps(rsq30);
284
285             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
286             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
287             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
288             rinvsq30         = _mm256_mul_ps(rinv30,rinv30);
289
290             /* Load parameters for j particles */
291             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
292                                                                  charge+jnrC+0,charge+jnrD+0,
293                                                                  charge+jnrE+0,charge+jnrF+0,
294                                                                  charge+jnrG+0,charge+jnrH+0);
295             vdwjidx0A        = 2*vdwtype[jnrA+0];
296             vdwjidx0B        = 2*vdwtype[jnrB+0];
297             vdwjidx0C        = 2*vdwtype[jnrC+0];
298             vdwjidx0D        = 2*vdwtype[jnrD+0];
299             vdwjidx0E        = 2*vdwtype[jnrE+0];
300             vdwjidx0F        = 2*vdwtype[jnrF+0];
301             vdwjidx0G        = 2*vdwtype[jnrG+0];
302             vdwjidx0H        = 2*vdwtype[jnrH+0];
303
304             fjx0             = _mm256_setzero_ps();
305             fjy0             = _mm256_setzero_ps();
306             fjz0             = _mm256_setzero_ps();
307
308             /**************************
309              * CALCULATE INTERACTIONS *
310              **************************/
311
312             if (gmx_mm256_any_lt(rsq00,rcutoff2))
313             {
314
315             r00              = _mm256_mul_ps(rsq00,rinv00);
316
317             /* Compute parameters for interactions between i and j atoms */
318             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
319                                             vdwioffsetptr0+vdwjidx0B,
320                                             vdwioffsetptr0+vdwjidx0C,
321                                             vdwioffsetptr0+vdwjidx0D,
322                                             vdwioffsetptr0+vdwjidx0E,
323                                             vdwioffsetptr0+vdwjidx0F,
324                                             vdwioffsetptr0+vdwjidx0G,
325                                             vdwioffsetptr0+vdwjidx0H,
326                                             &c6_00,&c12_00);
327
328             c6grid_00       = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
329                                                                   vdwgridioffsetptr0+vdwjidx0B,
330                                                                   vdwgridioffsetptr0+vdwjidx0C,
331                                                                   vdwgridioffsetptr0+vdwjidx0D,
332                                                                   vdwgridioffsetptr0+vdwjidx0E,
333                                                                   vdwgridioffsetptr0+vdwjidx0F,
334                                                                   vdwgridioffsetptr0+vdwjidx0G,
335                                                                   vdwgridioffsetptr0+vdwjidx0H);
336
337             /* Analytical LJ-PME */
338             rinvsix          = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
339             ewcljrsq         = _mm256_mul_ps(ewclj2,rsq00);
340             ewclj6           = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
341             exponent         = gmx_simd_exp_r(ewcljrsq);
342             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
343             poly             = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
344             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
345             vvdw6            = _mm256_mul_ps(_mm256_sub_ps(c6_00,_mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly))),rinvsix);
346             vvdw12           = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
347             vvdw             = _mm256_sub_ps(_mm256_mul_ps( _mm256_sub_ps(vvdw12 , _mm256_mul_ps(c12_00,_mm256_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
348                                           _mm256_mul_ps( _mm256_sub_ps(vvdw6,_mm256_add_ps(_mm256_mul_ps(c6_00,sh_vdw_invrcut6),_mm256_mul_ps(c6grid_00,sh_lj_ewald))),one_sixth));
349             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
350             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);
351
352             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
353
354             /* Update potential sum for this i atom from the interaction with this j atom. */
355             vvdw             = _mm256_and_ps(vvdw,cutoff_mask);
356             vvdwsum          = _mm256_add_ps(vvdwsum,vvdw);
357
358             fscal            = fvdw;
359
360             fscal            = _mm256_and_ps(fscal,cutoff_mask);
361
362             /* Calculate temporary vectorial force */
363             tx               = _mm256_mul_ps(fscal,dx00);
364             ty               = _mm256_mul_ps(fscal,dy00);
365             tz               = _mm256_mul_ps(fscal,dz00);
366
367             /* Update vectorial force */
368             fix0             = _mm256_add_ps(fix0,tx);
369             fiy0             = _mm256_add_ps(fiy0,ty);
370             fiz0             = _mm256_add_ps(fiz0,tz);
371
372             fjx0             = _mm256_add_ps(fjx0,tx);
373             fjy0             = _mm256_add_ps(fjy0,ty);
374             fjz0             = _mm256_add_ps(fjz0,tz);
375
376             }
377
378             /**************************
379              * CALCULATE INTERACTIONS *
380              **************************/
381
382             if (gmx_mm256_any_lt(rsq10,rcutoff2))
383             {
384
385             r10              = _mm256_mul_ps(rsq10,rinv10);
386
387             /* Compute parameters for interactions between i and j atoms */
388             qq10             = _mm256_mul_ps(iq1,jq0);
389
390             /* EWALD ELECTROSTATICS */
391             
392             /* Analytical PME correction */
393             zeta2            = _mm256_mul_ps(beta2,rsq10);
394             rinv3            = _mm256_mul_ps(rinvsq10,rinv10);
395             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
396             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
397             felec            = _mm256_mul_ps(qq10,felec);
398             pmecorrV         = gmx_mm256_pmecorrV_ps(zeta2);
399             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
400             velec            = _mm256_sub_ps(_mm256_sub_ps(rinv10,sh_ewald),pmecorrV);
401             velec            = _mm256_mul_ps(qq10,velec);
402             
403             cutoff_mask      = _mm256_cmp_ps(rsq10,rcutoff2,_CMP_LT_OQ);
404
405             /* Update potential sum for this i atom from the interaction with this j atom. */
406             velec            = _mm256_and_ps(velec,cutoff_mask);
407             velecsum         = _mm256_add_ps(velecsum,velec);
408
409             fscal            = felec;
410
411             fscal            = _mm256_and_ps(fscal,cutoff_mask);
412
413             /* Calculate temporary vectorial force */
414             tx               = _mm256_mul_ps(fscal,dx10);
415             ty               = _mm256_mul_ps(fscal,dy10);
416             tz               = _mm256_mul_ps(fscal,dz10);
417
418             /* Update vectorial force */
419             fix1             = _mm256_add_ps(fix1,tx);
420             fiy1             = _mm256_add_ps(fiy1,ty);
421             fiz1             = _mm256_add_ps(fiz1,tz);
422
423             fjx0             = _mm256_add_ps(fjx0,tx);
424             fjy0             = _mm256_add_ps(fjy0,ty);
425             fjz0             = _mm256_add_ps(fjz0,tz);
426
427             }
428
429             /**************************
430              * CALCULATE INTERACTIONS *
431              **************************/
432
433             if (gmx_mm256_any_lt(rsq20,rcutoff2))
434             {
435
436             r20              = _mm256_mul_ps(rsq20,rinv20);
437
438             /* Compute parameters for interactions between i and j atoms */
439             qq20             = _mm256_mul_ps(iq2,jq0);
440
441             /* EWALD ELECTROSTATICS */
442             
443             /* Analytical PME correction */
444             zeta2            = _mm256_mul_ps(beta2,rsq20);
445             rinv3            = _mm256_mul_ps(rinvsq20,rinv20);
446             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
447             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
448             felec            = _mm256_mul_ps(qq20,felec);
449             pmecorrV         = gmx_mm256_pmecorrV_ps(zeta2);
450             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
451             velec            = _mm256_sub_ps(_mm256_sub_ps(rinv20,sh_ewald),pmecorrV);
452             velec            = _mm256_mul_ps(qq20,velec);
453             
454             cutoff_mask      = _mm256_cmp_ps(rsq20,rcutoff2,_CMP_LT_OQ);
455
456             /* Update potential sum for this i atom from the interaction with this j atom. */
457             velec            = _mm256_and_ps(velec,cutoff_mask);
458             velecsum         = _mm256_add_ps(velecsum,velec);
459
460             fscal            = felec;
461
462             fscal            = _mm256_and_ps(fscal,cutoff_mask);
463
464             /* Calculate temporary vectorial force */
465             tx               = _mm256_mul_ps(fscal,dx20);
466             ty               = _mm256_mul_ps(fscal,dy20);
467             tz               = _mm256_mul_ps(fscal,dz20);
468
469             /* Update vectorial force */
470             fix2             = _mm256_add_ps(fix2,tx);
471             fiy2             = _mm256_add_ps(fiy2,ty);
472             fiz2             = _mm256_add_ps(fiz2,tz);
473
474             fjx0             = _mm256_add_ps(fjx0,tx);
475             fjy0             = _mm256_add_ps(fjy0,ty);
476             fjz0             = _mm256_add_ps(fjz0,tz);
477
478             }
479
480             /**************************
481              * CALCULATE INTERACTIONS *
482              **************************/
483
484             if (gmx_mm256_any_lt(rsq30,rcutoff2))
485             {
486
487             r30              = _mm256_mul_ps(rsq30,rinv30);
488
489             /* Compute parameters for interactions between i and j atoms */
490             qq30             = _mm256_mul_ps(iq3,jq0);
491
492             /* EWALD ELECTROSTATICS */
493             
494             /* Analytical PME correction */
495             zeta2            = _mm256_mul_ps(beta2,rsq30);
496             rinv3            = _mm256_mul_ps(rinvsq30,rinv30);
497             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
498             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
499             felec            = _mm256_mul_ps(qq30,felec);
500             pmecorrV         = gmx_mm256_pmecorrV_ps(zeta2);
501             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
502             velec            = _mm256_sub_ps(_mm256_sub_ps(rinv30,sh_ewald),pmecorrV);
503             velec            = _mm256_mul_ps(qq30,velec);
504             
505             cutoff_mask      = _mm256_cmp_ps(rsq30,rcutoff2,_CMP_LT_OQ);
506
507             /* Update potential sum for this i atom from the interaction with this j atom. */
508             velec            = _mm256_and_ps(velec,cutoff_mask);
509             velecsum         = _mm256_add_ps(velecsum,velec);
510
511             fscal            = felec;
512
513             fscal            = _mm256_and_ps(fscal,cutoff_mask);
514
515             /* Calculate temporary vectorial force */
516             tx               = _mm256_mul_ps(fscal,dx30);
517             ty               = _mm256_mul_ps(fscal,dy30);
518             tz               = _mm256_mul_ps(fscal,dz30);
519
520             /* Update vectorial force */
521             fix3             = _mm256_add_ps(fix3,tx);
522             fiy3             = _mm256_add_ps(fiy3,ty);
523             fiz3             = _mm256_add_ps(fiz3,tz);
524
525             fjx0             = _mm256_add_ps(fjx0,tx);
526             fjy0             = _mm256_add_ps(fjy0,ty);
527             fjz0             = _mm256_add_ps(fjz0,tz);
528
529             }
530
531             fjptrA             = f+j_coord_offsetA;
532             fjptrB             = f+j_coord_offsetB;
533             fjptrC             = f+j_coord_offsetC;
534             fjptrD             = f+j_coord_offsetD;
535             fjptrE             = f+j_coord_offsetE;
536             fjptrF             = f+j_coord_offsetF;
537             fjptrG             = f+j_coord_offsetG;
538             fjptrH             = f+j_coord_offsetH;
539
540             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
541
542             /* Inner loop uses 392 flops */
543         }
544
545         if(jidx<j_index_end)
546         {
547
548             /* Get j neighbor index, and coordinate index */
549             jnrlistA         = jjnr[jidx];
550             jnrlistB         = jjnr[jidx+1];
551             jnrlistC         = jjnr[jidx+2];
552             jnrlistD         = jjnr[jidx+3];
553             jnrlistE         = jjnr[jidx+4];
554             jnrlistF         = jjnr[jidx+5];
555             jnrlistG         = jjnr[jidx+6];
556             jnrlistH         = jjnr[jidx+7];
557             /* Sign of each element will be negative for non-real atoms.
558              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
559              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
560              */
561             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
562                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
563                                             
564             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
565             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
566             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
567             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
568             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
569             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
570             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
571             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
572             j_coord_offsetA  = DIM*jnrA;
573             j_coord_offsetB  = DIM*jnrB;
574             j_coord_offsetC  = DIM*jnrC;
575             j_coord_offsetD  = DIM*jnrD;
576             j_coord_offsetE  = DIM*jnrE;
577             j_coord_offsetF  = DIM*jnrF;
578             j_coord_offsetG  = DIM*jnrG;
579             j_coord_offsetH  = DIM*jnrH;
580
581             /* load j atom coordinates */
582             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
583                                                  x+j_coord_offsetC,x+j_coord_offsetD,
584                                                  x+j_coord_offsetE,x+j_coord_offsetF,
585                                                  x+j_coord_offsetG,x+j_coord_offsetH,
586                                                  &jx0,&jy0,&jz0);
587
588             /* Calculate displacement vector */
589             dx00             = _mm256_sub_ps(ix0,jx0);
590             dy00             = _mm256_sub_ps(iy0,jy0);
591             dz00             = _mm256_sub_ps(iz0,jz0);
592             dx10             = _mm256_sub_ps(ix1,jx0);
593             dy10             = _mm256_sub_ps(iy1,jy0);
594             dz10             = _mm256_sub_ps(iz1,jz0);
595             dx20             = _mm256_sub_ps(ix2,jx0);
596             dy20             = _mm256_sub_ps(iy2,jy0);
597             dz20             = _mm256_sub_ps(iz2,jz0);
598             dx30             = _mm256_sub_ps(ix3,jx0);
599             dy30             = _mm256_sub_ps(iy3,jy0);
600             dz30             = _mm256_sub_ps(iz3,jz0);
601
602             /* Calculate squared distance and things based on it */
603             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
604             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
605             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
606             rsq30            = gmx_mm256_calc_rsq_ps(dx30,dy30,dz30);
607
608             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
609             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
610             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
611             rinv30           = gmx_mm256_invsqrt_ps(rsq30);
612
613             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
614             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
615             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
616             rinvsq30         = _mm256_mul_ps(rinv30,rinv30);
617
618             /* Load parameters for j particles */
619             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
620                                                                  charge+jnrC+0,charge+jnrD+0,
621                                                                  charge+jnrE+0,charge+jnrF+0,
622                                                                  charge+jnrG+0,charge+jnrH+0);
623             vdwjidx0A        = 2*vdwtype[jnrA+0];
624             vdwjidx0B        = 2*vdwtype[jnrB+0];
625             vdwjidx0C        = 2*vdwtype[jnrC+0];
626             vdwjidx0D        = 2*vdwtype[jnrD+0];
627             vdwjidx0E        = 2*vdwtype[jnrE+0];
628             vdwjidx0F        = 2*vdwtype[jnrF+0];
629             vdwjidx0G        = 2*vdwtype[jnrG+0];
630             vdwjidx0H        = 2*vdwtype[jnrH+0];
631
632             fjx0             = _mm256_setzero_ps();
633             fjy0             = _mm256_setzero_ps();
634             fjz0             = _mm256_setzero_ps();
635
636             /**************************
637              * CALCULATE INTERACTIONS *
638              **************************/
639
640             if (gmx_mm256_any_lt(rsq00,rcutoff2))
641             {
642
643             r00              = _mm256_mul_ps(rsq00,rinv00);
644             r00              = _mm256_andnot_ps(dummy_mask,r00);
645
646             /* Compute parameters for interactions between i and j atoms */
647             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
648                                             vdwioffsetptr0+vdwjidx0B,
649                                             vdwioffsetptr0+vdwjidx0C,
650                                             vdwioffsetptr0+vdwjidx0D,
651                                             vdwioffsetptr0+vdwjidx0E,
652                                             vdwioffsetptr0+vdwjidx0F,
653                                             vdwioffsetptr0+vdwjidx0G,
654                                             vdwioffsetptr0+vdwjidx0H,
655                                             &c6_00,&c12_00);
656
657             c6grid_00       = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
658                                                                   vdwgridioffsetptr0+vdwjidx0B,
659                                                                   vdwgridioffsetptr0+vdwjidx0C,
660                                                                   vdwgridioffsetptr0+vdwjidx0D,
661                                                                   vdwgridioffsetptr0+vdwjidx0E,
662                                                                   vdwgridioffsetptr0+vdwjidx0F,
663                                                                   vdwgridioffsetptr0+vdwjidx0G,
664                                                                   vdwgridioffsetptr0+vdwjidx0H);
665
666             /* Analytical LJ-PME */
667             rinvsix          = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
668             ewcljrsq         = _mm256_mul_ps(ewclj2,rsq00);
669             ewclj6           = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
670             exponent         = gmx_simd_exp_r(ewcljrsq);
671             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
672             poly             = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
673             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
674             vvdw6            = _mm256_mul_ps(_mm256_sub_ps(c6_00,_mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly))),rinvsix);
675             vvdw12           = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
676             vvdw             = _mm256_sub_ps(_mm256_mul_ps( _mm256_sub_ps(vvdw12 , _mm256_mul_ps(c12_00,_mm256_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
677                                           _mm256_mul_ps( _mm256_sub_ps(vvdw6,_mm256_add_ps(_mm256_mul_ps(c6_00,sh_vdw_invrcut6),_mm256_mul_ps(c6grid_00,sh_lj_ewald))),one_sixth));
678             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
679             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);
680
681             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
682
683             /* Update potential sum for this i atom from the interaction with this j atom. */
684             vvdw             = _mm256_and_ps(vvdw,cutoff_mask);
685             vvdw             = _mm256_andnot_ps(dummy_mask,vvdw);
686             vvdwsum          = _mm256_add_ps(vvdwsum,vvdw);
687
688             fscal            = fvdw;
689
690             fscal            = _mm256_and_ps(fscal,cutoff_mask);
691
692             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
693
694             /* Calculate temporary vectorial force */
695             tx               = _mm256_mul_ps(fscal,dx00);
696             ty               = _mm256_mul_ps(fscal,dy00);
697             tz               = _mm256_mul_ps(fscal,dz00);
698
699             /* Update vectorial force */
700             fix0             = _mm256_add_ps(fix0,tx);
701             fiy0             = _mm256_add_ps(fiy0,ty);
702             fiz0             = _mm256_add_ps(fiz0,tz);
703
704             fjx0             = _mm256_add_ps(fjx0,tx);
705             fjy0             = _mm256_add_ps(fjy0,ty);
706             fjz0             = _mm256_add_ps(fjz0,tz);
707
708             }
709
710             /**************************
711              * CALCULATE INTERACTIONS *
712              **************************/
713
714             if (gmx_mm256_any_lt(rsq10,rcutoff2))
715             {
716
717             r10              = _mm256_mul_ps(rsq10,rinv10);
718             r10              = _mm256_andnot_ps(dummy_mask,r10);
719
720             /* Compute parameters for interactions between i and j atoms */
721             qq10             = _mm256_mul_ps(iq1,jq0);
722
723             /* EWALD ELECTROSTATICS */
724             
725             /* Analytical PME correction */
726             zeta2            = _mm256_mul_ps(beta2,rsq10);
727             rinv3            = _mm256_mul_ps(rinvsq10,rinv10);
728             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
729             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
730             felec            = _mm256_mul_ps(qq10,felec);
731             pmecorrV         = gmx_mm256_pmecorrV_ps(zeta2);
732             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
733             velec            = _mm256_sub_ps(_mm256_sub_ps(rinv10,sh_ewald),pmecorrV);
734             velec            = _mm256_mul_ps(qq10,velec);
735             
736             cutoff_mask      = _mm256_cmp_ps(rsq10,rcutoff2,_CMP_LT_OQ);
737
738             /* Update potential sum for this i atom from the interaction with this j atom. */
739             velec            = _mm256_and_ps(velec,cutoff_mask);
740             velec            = _mm256_andnot_ps(dummy_mask,velec);
741             velecsum         = _mm256_add_ps(velecsum,velec);
742
743             fscal            = felec;
744
745             fscal            = _mm256_and_ps(fscal,cutoff_mask);
746
747             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
748
749             /* Calculate temporary vectorial force */
750             tx               = _mm256_mul_ps(fscal,dx10);
751             ty               = _mm256_mul_ps(fscal,dy10);
752             tz               = _mm256_mul_ps(fscal,dz10);
753
754             /* Update vectorial force */
755             fix1             = _mm256_add_ps(fix1,tx);
756             fiy1             = _mm256_add_ps(fiy1,ty);
757             fiz1             = _mm256_add_ps(fiz1,tz);
758
759             fjx0             = _mm256_add_ps(fjx0,tx);
760             fjy0             = _mm256_add_ps(fjy0,ty);
761             fjz0             = _mm256_add_ps(fjz0,tz);
762
763             }
764
765             /**************************
766              * CALCULATE INTERACTIONS *
767              **************************/
768
769             if (gmx_mm256_any_lt(rsq20,rcutoff2))
770             {
771
772             r20              = _mm256_mul_ps(rsq20,rinv20);
773             r20              = _mm256_andnot_ps(dummy_mask,r20);
774
775             /* Compute parameters for interactions between i and j atoms */
776             qq20             = _mm256_mul_ps(iq2,jq0);
777
778             /* EWALD ELECTROSTATICS */
779             
780             /* Analytical PME correction */
781             zeta2            = _mm256_mul_ps(beta2,rsq20);
782             rinv3            = _mm256_mul_ps(rinvsq20,rinv20);
783             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
784             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
785             felec            = _mm256_mul_ps(qq20,felec);
786             pmecorrV         = gmx_mm256_pmecorrV_ps(zeta2);
787             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
788             velec            = _mm256_sub_ps(_mm256_sub_ps(rinv20,sh_ewald),pmecorrV);
789             velec            = _mm256_mul_ps(qq20,velec);
790             
791             cutoff_mask      = _mm256_cmp_ps(rsq20,rcutoff2,_CMP_LT_OQ);
792
793             /* Update potential sum for this i atom from the interaction with this j atom. */
794             velec            = _mm256_and_ps(velec,cutoff_mask);
795             velec            = _mm256_andnot_ps(dummy_mask,velec);
796             velecsum         = _mm256_add_ps(velecsum,velec);
797
798             fscal            = felec;
799
800             fscal            = _mm256_and_ps(fscal,cutoff_mask);
801
802             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
803
804             /* Calculate temporary vectorial force */
805             tx               = _mm256_mul_ps(fscal,dx20);
806             ty               = _mm256_mul_ps(fscal,dy20);
807             tz               = _mm256_mul_ps(fscal,dz20);
808
809             /* Update vectorial force */
810             fix2             = _mm256_add_ps(fix2,tx);
811             fiy2             = _mm256_add_ps(fiy2,ty);
812             fiz2             = _mm256_add_ps(fiz2,tz);
813
814             fjx0             = _mm256_add_ps(fjx0,tx);
815             fjy0             = _mm256_add_ps(fjy0,ty);
816             fjz0             = _mm256_add_ps(fjz0,tz);
817
818             }
819
820             /**************************
821              * CALCULATE INTERACTIONS *
822              **************************/
823
824             if (gmx_mm256_any_lt(rsq30,rcutoff2))
825             {
826
827             r30              = _mm256_mul_ps(rsq30,rinv30);
828             r30              = _mm256_andnot_ps(dummy_mask,r30);
829
830             /* Compute parameters for interactions between i and j atoms */
831             qq30             = _mm256_mul_ps(iq3,jq0);
832
833             /* EWALD ELECTROSTATICS */
834             
835             /* Analytical PME correction */
836             zeta2            = _mm256_mul_ps(beta2,rsq30);
837             rinv3            = _mm256_mul_ps(rinvsq30,rinv30);
838             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
839             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
840             felec            = _mm256_mul_ps(qq30,felec);
841             pmecorrV         = gmx_mm256_pmecorrV_ps(zeta2);
842             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
843             velec            = _mm256_sub_ps(_mm256_sub_ps(rinv30,sh_ewald),pmecorrV);
844             velec            = _mm256_mul_ps(qq30,velec);
845             
846             cutoff_mask      = _mm256_cmp_ps(rsq30,rcutoff2,_CMP_LT_OQ);
847
848             /* Update potential sum for this i atom from the interaction with this j atom. */
849             velec            = _mm256_and_ps(velec,cutoff_mask);
850             velec            = _mm256_andnot_ps(dummy_mask,velec);
851             velecsum         = _mm256_add_ps(velecsum,velec);
852
853             fscal            = felec;
854
855             fscal            = _mm256_and_ps(fscal,cutoff_mask);
856
857             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
858
859             /* Calculate temporary vectorial force */
860             tx               = _mm256_mul_ps(fscal,dx30);
861             ty               = _mm256_mul_ps(fscal,dy30);
862             tz               = _mm256_mul_ps(fscal,dz30);
863
864             /* Update vectorial force */
865             fix3             = _mm256_add_ps(fix3,tx);
866             fiy3             = _mm256_add_ps(fiy3,ty);
867             fiz3             = _mm256_add_ps(fiz3,tz);
868
869             fjx0             = _mm256_add_ps(fjx0,tx);
870             fjy0             = _mm256_add_ps(fjy0,ty);
871             fjz0             = _mm256_add_ps(fjz0,tz);
872
873             }
874
875             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
876             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
877             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
878             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
879             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
880             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
881             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
882             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
883
884             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
885
886             /* Inner loop uses 396 flops */
887         }
888
889         /* End of innermost loop */
890
891         gmx_mm256_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
892                                                  f+i_coord_offset,fshift+i_shift_offset);
893
894         ggid                        = gid[iidx];
895         /* Update potential energies */
896         gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
897         gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
898
899         /* Increment number of inner iterations */
900         inneriter                  += j_index_end - j_index_start;
901
902         /* Outer loop uses 26 flops */
903     }
904
905     /* Increment number of outer iterations */
906     outeriter        += nri;
907
908     /* Update outer/inner flops */
909
910     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*396);
911 }
912 /*
913  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwLJEwSh_GeomW4P1_F_avx_256_single
914  * Electrostatics interaction: Ewald
915  * VdW interaction:            LJEwald
916  * Geometry:                   Water4-Particle
917  * Calculate force/pot:        Force
918  */
919 void
920 nb_kernel_ElecEwSh_VdwLJEwSh_GeomW4P1_F_avx_256_single
921                     (t_nblist                    * gmx_restrict       nlist,
922                      rvec                        * gmx_restrict          xx,
923                      rvec                        * gmx_restrict          ff,
924                      t_forcerec                  * gmx_restrict          fr,
925                      t_mdatoms                   * gmx_restrict     mdatoms,
926                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
927                      t_nrnb                      * gmx_restrict        nrnb)
928 {
929     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
930      * just 0 for non-waters.
931      * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
932      * jnr indices corresponding to data put in the four positions in the SIMD register.
933      */
934     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
935     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
936     int              jnrA,jnrB,jnrC,jnrD;
937     int              jnrE,jnrF,jnrG,jnrH;
938     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
939     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
940     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
941     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
942     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
943     real             rcutoff_scalar;
944     real             *shiftvec,*fshift,*x,*f;
945     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
946     real             scratch[4*DIM];
947     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
948     real *           vdwioffsetptr0;
949     real *           vdwgridioffsetptr0;
950     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
951     real *           vdwioffsetptr1;
952     real *           vdwgridioffsetptr1;
953     __m256           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
954     real *           vdwioffsetptr2;
955     real *           vdwgridioffsetptr2;
956     __m256           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
957     real *           vdwioffsetptr3;
958     real *           vdwgridioffsetptr3;
959     __m256           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
960     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
961     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
962     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
963     __m256           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
964     __m256           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
965     __m256           dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
966     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
967     real             *charge;
968     int              nvdwtype;
969     __m256           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
970     int              *vdwtype;
971     real             *vdwparam;
972     __m256           one_sixth   = _mm256_set1_ps(1.0/6.0);
973     __m256           one_twelfth = _mm256_set1_ps(1.0/12.0);
974     __m256           c6grid_00;
975     __m256           c6grid_10;
976     __m256           c6grid_20;
977     __m256           c6grid_30;
978     real             *vdwgridparam;
979     __m256           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
980     __m256           one_half  = _mm256_set1_ps(0.5);
981     __m256           minus_one = _mm256_set1_ps(-1.0);
982     __m256i          ewitab;
983     __m128i          ewitab_lo,ewitab_hi;
984     __m256           ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
985     __m256           beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
986     real             *ewtab;
987     __m256           dummy_mask,cutoff_mask;
988     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
989     __m256           one     = _mm256_set1_ps(1.0);
990     __m256           two     = _mm256_set1_ps(2.0);
991     x                = xx[0];
992     f                = ff[0];
993
994     nri              = nlist->nri;
995     iinr             = nlist->iinr;
996     jindex           = nlist->jindex;
997     jjnr             = nlist->jjnr;
998     shiftidx         = nlist->shift;
999     gid              = nlist->gid;
1000     shiftvec         = fr->shift_vec[0];
1001     fshift           = fr->fshift[0];
1002     facel            = _mm256_set1_ps(fr->epsfac);
1003     charge           = mdatoms->chargeA;
1004     nvdwtype         = fr->ntype;
1005     vdwparam         = fr->nbfp;
1006     vdwtype          = mdatoms->typeA;
1007     vdwgridparam     = fr->ljpme_c6grid;
1008     sh_lj_ewald      = _mm256_set1_ps(fr->ic->sh_lj_ewald);
1009     ewclj            = _mm256_set1_ps(fr->ewaldcoeff_lj);
1010     ewclj2           = _mm256_mul_ps(minus_one,_mm256_mul_ps(ewclj,ewclj));
1011
1012     sh_ewald         = _mm256_set1_ps(fr->ic->sh_ewald);
1013     beta             = _mm256_set1_ps(fr->ic->ewaldcoeff_q);
1014     beta2            = _mm256_mul_ps(beta,beta);
1015     beta3            = _mm256_mul_ps(beta,beta2);
1016
1017     ewtab            = fr->ic->tabq_coul_F;
1018     ewtabscale       = _mm256_set1_ps(fr->ic->tabq_scale);
1019     ewtabhalfspace   = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
1020
1021     /* Setup water-specific parameters */
1022     inr              = nlist->iinr[0];
1023     iq1              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
1024     iq2              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
1025     iq3              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+3]));
1026     vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
1027     vdwgridioffsetptr0 = vdwgridparam+2*nvdwtype*vdwtype[inr+0];
1028
1029     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
1030     rcutoff_scalar   = fr->rcoulomb;
1031     rcutoff          = _mm256_set1_ps(rcutoff_scalar);
1032     rcutoff2         = _mm256_mul_ps(rcutoff,rcutoff);
1033
1034     sh_vdw_invrcut6  = _mm256_set1_ps(fr->ic->sh_invrc6);
1035     rvdw             = _mm256_set1_ps(fr->rvdw);
1036
1037     /* Avoid stupid compiler warnings */
1038     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
1039     j_coord_offsetA = 0;
1040     j_coord_offsetB = 0;
1041     j_coord_offsetC = 0;
1042     j_coord_offsetD = 0;
1043     j_coord_offsetE = 0;
1044     j_coord_offsetF = 0;
1045     j_coord_offsetG = 0;
1046     j_coord_offsetH = 0;
1047
1048     outeriter        = 0;
1049     inneriter        = 0;
1050
1051     for(iidx=0;iidx<4*DIM;iidx++)
1052     {
1053         scratch[iidx] = 0.0;
1054     }
1055
1056     /* Start outer loop over neighborlists */
1057     for(iidx=0; iidx<nri; iidx++)
1058     {
1059         /* Load shift vector for this list */
1060         i_shift_offset   = DIM*shiftidx[iidx];
1061
1062         /* Load limits for loop over neighbors */
1063         j_index_start    = jindex[iidx];
1064         j_index_end      = jindex[iidx+1];
1065
1066         /* Get outer coordinate index */
1067         inr              = iinr[iidx];
1068         i_coord_offset   = DIM*inr;
1069
1070         /* Load i particle coords and add shift vector */
1071         gmx_mm256_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
1072                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
1073
1074         fix0             = _mm256_setzero_ps();
1075         fiy0             = _mm256_setzero_ps();
1076         fiz0             = _mm256_setzero_ps();
1077         fix1             = _mm256_setzero_ps();
1078         fiy1             = _mm256_setzero_ps();
1079         fiz1             = _mm256_setzero_ps();
1080         fix2             = _mm256_setzero_ps();
1081         fiy2             = _mm256_setzero_ps();
1082         fiz2             = _mm256_setzero_ps();
1083         fix3             = _mm256_setzero_ps();
1084         fiy3             = _mm256_setzero_ps();
1085         fiz3             = _mm256_setzero_ps();
1086
1087         /* Start inner kernel loop */
1088         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
1089         {
1090
1091             /* Get j neighbor index, and coordinate index */
1092             jnrA             = jjnr[jidx];
1093             jnrB             = jjnr[jidx+1];
1094             jnrC             = jjnr[jidx+2];
1095             jnrD             = jjnr[jidx+3];
1096             jnrE             = jjnr[jidx+4];
1097             jnrF             = jjnr[jidx+5];
1098             jnrG             = jjnr[jidx+6];
1099             jnrH             = jjnr[jidx+7];
1100             j_coord_offsetA  = DIM*jnrA;
1101             j_coord_offsetB  = DIM*jnrB;
1102             j_coord_offsetC  = DIM*jnrC;
1103             j_coord_offsetD  = DIM*jnrD;
1104             j_coord_offsetE  = DIM*jnrE;
1105             j_coord_offsetF  = DIM*jnrF;
1106             j_coord_offsetG  = DIM*jnrG;
1107             j_coord_offsetH  = DIM*jnrH;
1108
1109             /* load j atom coordinates */
1110             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1111                                                  x+j_coord_offsetC,x+j_coord_offsetD,
1112                                                  x+j_coord_offsetE,x+j_coord_offsetF,
1113                                                  x+j_coord_offsetG,x+j_coord_offsetH,
1114                                                  &jx0,&jy0,&jz0);
1115
1116             /* Calculate displacement vector */
1117             dx00             = _mm256_sub_ps(ix0,jx0);
1118             dy00             = _mm256_sub_ps(iy0,jy0);
1119             dz00             = _mm256_sub_ps(iz0,jz0);
1120             dx10             = _mm256_sub_ps(ix1,jx0);
1121             dy10             = _mm256_sub_ps(iy1,jy0);
1122             dz10             = _mm256_sub_ps(iz1,jz0);
1123             dx20             = _mm256_sub_ps(ix2,jx0);
1124             dy20             = _mm256_sub_ps(iy2,jy0);
1125             dz20             = _mm256_sub_ps(iz2,jz0);
1126             dx30             = _mm256_sub_ps(ix3,jx0);
1127             dy30             = _mm256_sub_ps(iy3,jy0);
1128             dz30             = _mm256_sub_ps(iz3,jz0);
1129
1130             /* Calculate squared distance and things based on it */
1131             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
1132             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
1133             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
1134             rsq30            = gmx_mm256_calc_rsq_ps(dx30,dy30,dz30);
1135
1136             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
1137             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
1138             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
1139             rinv30           = gmx_mm256_invsqrt_ps(rsq30);
1140
1141             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
1142             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
1143             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
1144             rinvsq30         = _mm256_mul_ps(rinv30,rinv30);
1145
1146             /* Load parameters for j particles */
1147             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1148                                                                  charge+jnrC+0,charge+jnrD+0,
1149                                                                  charge+jnrE+0,charge+jnrF+0,
1150                                                                  charge+jnrG+0,charge+jnrH+0);
1151             vdwjidx0A        = 2*vdwtype[jnrA+0];
1152             vdwjidx0B        = 2*vdwtype[jnrB+0];
1153             vdwjidx0C        = 2*vdwtype[jnrC+0];
1154             vdwjidx0D        = 2*vdwtype[jnrD+0];
1155             vdwjidx0E        = 2*vdwtype[jnrE+0];
1156             vdwjidx0F        = 2*vdwtype[jnrF+0];
1157             vdwjidx0G        = 2*vdwtype[jnrG+0];
1158             vdwjidx0H        = 2*vdwtype[jnrH+0];
1159
1160             fjx0             = _mm256_setzero_ps();
1161             fjy0             = _mm256_setzero_ps();
1162             fjz0             = _mm256_setzero_ps();
1163
1164             /**************************
1165              * CALCULATE INTERACTIONS *
1166              **************************/
1167
1168             if (gmx_mm256_any_lt(rsq00,rcutoff2))
1169             {
1170
1171             r00              = _mm256_mul_ps(rsq00,rinv00);
1172
1173             /* Compute parameters for interactions between i and j atoms */
1174             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
1175                                             vdwioffsetptr0+vdwjidx0B,
1176                                             vdwioffsetptr0+vdwjidx0C,
1177                                             vdwioffsetptr0+vdwjidx0D,
1178                                             vdwioffsetptr0+vdwjidx0E,
1179                                             vdwioffsetptr0+vdwjidx0F,
1180                                             vdwioffsetptr0+vdwjidx0G,
1181                                             vdwioffsetptr0+vdwjidx0H,
1182                                             &c6_00,&c12_00);
1183
1184             c6grid_00       = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
1185                                                                   vdwgridioffsetptr0+vdwjidx0B,
1186                                                                   vdwgridioffsetptr0+vdwjidx0C,
1187                                                                   vdwgridioffsetptr0+vdwjidx0D,
1188                                                                   vdwgridioffsetptr0+vdwjidx0E,
1189                                                                   vdwgridioffsetptr0+vdwjidx0F,
1190                                                                   vdwgridioffsetptr0+vdwjidx0G,
1191                                                                   vdwgridioffsetptr0+vdwjidx0H);
1192
1193             /* Analytical LJ-PME */
1194             rinvsix          = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
1195             ewcljrsq         = _mm256_mul_ps(ewclj2,rsq00);
1196             ewclj6           = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
1197             exponent         = gmx_simd_exp_r(ewcljrsq);
1198             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
1199             poly             = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
1200             /* f6A = 6 * C6grid * (1 - poly) */
1201             f6A              = _mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly));
1202             /* f6B = C6grid * exponent * beta^6 */
1203             f6B              = _mm256_mul_ps(_mm256_mul_ps(c6grid_00,one_sixth),_mm256_mul_ps(exponent,ewclj6));
1204             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
1205             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);
1206
1207             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
1208
1209             fscal            = fvdw;
1210
1211             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1212
1213             /* Calculate temporary vectorial force */
1214             tx               = _mm256_mul_ps(fscal,dx00);
1215             ty               = _mm256_mul_ps(fscal,dy00);
1216             tz               = _mm256_mul_ps(fscal,dz00);
1217
1218             /* Update vectorial force */
1219             fix0             = _mm256_add_ps(fix0,tx);
1220             fiy0             = _mm256_add_ps(fiy0,ty);
1221             fiz0             = _mm256_add_ps(fiz0,tz);
1222
1223             fjx0             = _mm256_add_ps(fjx0,tx);
1224             fjy0             = _mm256_add_ps(fjy0,ty);
1225             fjz0             = _mm256_add_ps(fjz0,tz);
1226
1227             }
1228
1229             /**************************
1230              * CALCULATE INTERACTIONS *
1231              **************************/
1232
1233             if (gmx_mm256_any_lt(rsq10,rcutoff2))
1234             {
1235
1236             r10              = _mm256_mul_ps(rsq10,rinv10);
1237
1238             /* Compute parameters for interactions between i and j atoms */
1239             qq10             = _mm256_mul_ps(iq1,jq0);
1240
1241             /* EWALD ELECTROSTATICS */
1242             
1243             /* Analytical PME correction */
1244             zeta2            = _mm256_mul_ps(beta2,rsq10);
1245             rinv3            = _mm256_mul_ps(rinvsq10,rinv10);
1246             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
1247             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1248             felec            = _mm256_mul_ps(qq10,felec);
1249             
1250             cutoff_mask      = _mm256_cmp_ps(rsq10,rcutoff2,_CMP_LT_OQ);
1251
1252             fscal            = felec;
1253
1254             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1255
1256             /* Calculate temporary vectorial force */
1257             tx               = _mm256_mul_ps(fscal,dx10);
1258             ty               = _mm256_mul_ps(fscal,dy10);
1259             tz               = _mm256_mul_ps(fscal,dz10);
1260
1261             /* Update vectorial force */
1262             fix1             = _mm256_add_ps(fix1,tx);
1263             fiy1             = _mm256_add_ps(fiy1,ty);
1264             fiz1             = _mm256_add_ps(fiz1,tz);
1265
1266             fjx0             = _mm256_add_ps(fjx0,tx);
1267             fjy0             = _mm256_add_ps(fjy0,ty);
1268             fjz0             = _mm256_add_ps(fjz0,tz);
1269
1270             }
1271
1272             /**************************
1273              * CALCULATE INTERACTIONS *
1274              **************************/
1275
1276             if (gmx_mm256_any_lt(rsq20,rcutoff2))
1277             {
1278
1279             r20              = _mm256_mul_ps(rsq20,rinv20);
1280
1281             /* Compute parameters for interactions between i and j atoms */
1282             qq20             = _mm256_mul_ps(iq2,jq0);
1283
1284             /* EWALD ELECTROSTATICS */
1285             
1286             /* Analytical PME correction */
1287             zeta2            = _mm256_mul_ps(beta2,rsq20);
1288             rinv3            = _mm256_mul_ps(rinvsq20,rinv20);
1289             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
1290             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1291             felec            = _mm256_mul_ps(qq20,felec);
1292             
1293             cutoff_mask      = _mm256_cmp_ps(rsq20,rcutoff2,_CMP_LT_OQ);
1294
1295             fscal            = felec;
1296
1297             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1298
1299             /* Calculate temporary vectorial force */
1300             tx               = _mm256_mul_ps(fscal,dx20);
1301             ty               = _mm256_mul_ps(fscal,dy20);
1302             tz               = _mm256_mul_ps(fscal,dz20);
1303
1304             /* Update vectorial force */
1305             fix2             = _mm256_add_ps(fix2,tx);
1306             fiy2             = _mm256_add_ps(fiy2,ty);
1307             fiz2             = _mm256_add_ps(fiz2,tz);
1308
1309             fjx0             = _mm256_add_ps(fjx0,tx);
1310             fjy0             = _mm256_add_ps(fjy0,ty);
1311             fjz0             = _mm256_add_ps(fjz0,tz);
1312
1313             }
1314
1315             /**************************
1316              * CALCULATE INTERACTIONS *
1317              **************************/
1318
1319             if (gmx_mm256_any_lt(rsq30,rcutoff2))
1320             {
1321
1322             r30              = _mm256_mul_ps(rsq30,rinv30);
1323
1324             /* Compute parameters for interactions between i and j atoms */
1325             qq30             = _mm256_mul_ps(iq3,jq0);
1326
1327             /* EWALD ELECTROSTATICS */
1328             
1329             /* Analytical PME correction */
1330             zeta2            = _mm256_mul_ps(beta2,rsq30);
1331             rinv3            = _mm256_mul_ps(rinvsq30,rinv30);
1332             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
1333             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1334             felec            = _mm256_mul_ps(qq30,felec);
1335             
1336             cutoff_mask      = _mm256_cmp_ps(rsq30,rcutoff2,_CMP_LT_OQ);
1337
1338             fscal            = felec;
1339
1340             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1341
1342             /* Calculate temporary vectorial force */
1343             tx               = _mm256_mul_ps(fscal,dx30);
1344             ty               = _mm256_mul_ps(fscal,dy30);
1345             tz               = _mm256_mul_ps(fscal,dz30);
1346
1347             /* Update vectorial force */
1348             fix3             = _mm256_add_ps(fix3,tx);
1349             fiy3             = _mm256_add_ps(fiy3,ty);
1350             fiz3             = _mm256_add_ps(fiz3,tz);
1351
1352             fjx0             = _mm256_add_ps(fjx0,tx);
1353             fjy0             = _mm256_add_ps(fjy0,ty);
1354             fjz0             = _mm256_add_ps(fjz0,tz);
1355
1356             }
1357
1358             fjptrA             = f+j_coord_offsetA;
1359             fjptrB             = f+j_coord_offsetB;
1360             fjptrC             = f+j_coord_offsetC;
1361             fjptrD             = f+j_coord_offsetD;
1362             fjptrE             = f+j_coord_offsetE;
1363             fjptrF             = f+j_coord_offsetF;
1364             fjptrG             = f+j_coord_offsetG;
1365             fjptrH             = f+j_coord_offsetH;
1366
1367             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
1368
1369             /* Inner loop uses 229 flops */
1370         }
1371
1372         if(jidx<j_index_end)
1373         {
1374
1375             /* Get j neighbor index, and coordinate index */
1376             jnrlistA         = jjnr[jidx];
1377             jnrlistB         = jjnr[jidx+1];
1378             jnrlistC         = jjnr[jidx+2];
1379             jnrlistD         = jjnr[jidx+3];
1380             jnrlistE         = jjnr[jidx+4];
1381             jnrlistF         = jjnr[jidx+5];
1382             jnrlistG         = jjnr[jidx+6];
1383             jnrlistH         = jjnr[jidx+7];
1384             /* Sign of each element will be negative for non-real atoms.
1385              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1386              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1387              */
1388             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
1389                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
1390                                             
1391             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
1392             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
1393             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
1394             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
1395             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
1396             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
1397             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
1398             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
1399             j_coord_offsetA  = DIM*jnrA;
1400             j_coord_offsetB  = DIM*jnrB;
1401             j_coord_offsetC  = DIM*jnrC;
1402             j_coord_offsetD  = DIM*jnrD;
1403             j_coord_offsetE  = DIM*jnrE;
1404             j_coord_offsetF  = DIM*jnrF;
1405             j_coord_offsetG  = DIM*jnrG;
1406             j_coord_offsetH  = DIM*jnrH;
1407
1408             /* load j atom coordinates */
1409             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1410                                                  x+j_coord_offsetC,x+j_coord_offsetD,
1411                                                  x+j_coord_offsetE,x+j_coord_offsetF,
1412                                                  x+j_coord_offsetG,x+j_coord_offsetH,
1413                                                  &jx0,&jy0,&jz0);
1414
1415             /* Calculate displacement vector */
1416             dx00             = _mm256_sub_ps(ix0,jx0);
1417             dy00             = _mm256_sub_ps(iy0,jy0);
1418             dz00             = _mm256_sub_ps(iz0,jz0);
1419             dx10             = _mm256_sub_ps(ix1,jx0);
1420             dy10             = _mm256_sub_ps(iy1,jy0);
1421             dz10             = _mm256_sub_ps(iz1,jz0);
1422             dx20             = _mm256_sub_ps(ix2,jx0);
1423             dy20             = _mm256_sub_ps(iy2,jy0);
1424             dz20             = _mm256_sub_ps(iz2,jz0);
1425             dx30             = _mm256_sub_ps(ix3,jx0);
1426             dy30             = _mm256_sub_ps(iy3,jy0);
1427             dz30             = _mm256_sub_ps(iz3,jz0);
1428
1429             /* Calculate squared distance and things based on it */
1430             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
1431             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
1432             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
1433             rsq30            = gmx_mm256_calc_rsq_ps(dx30,dy30,dz30);
1434
1435             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
1436             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
1437             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
1438             rinv30           = gmx_mm256_invsqrt_ps(rsq30);
1439
1440             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
1441             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
1442             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
1443             rinvsq30         = _mm256_mul_ps(rinv30,rinv30);
1444
1445             /* Load parameters for j particles */
1446             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1447                                                                  charge+jnrC+0,charge+jnrD+0,
1448                                                                  charge+jnrE+0,charge+jnrF+0,
1449                                                                  charge+jnrG+0,charge+jnrH+0);
1450             vdwjidx0A        = 2*vdwtype[jnrA+0];
1451             vdwjidx0B        = 2*vdwtype[jnrB+0];
1452             vdwjidx0C        = 2*vdwtype[jnrC+0];
1453             vdwjidx0D        = 2*vdwtype[jnrD+0];
1454             vdwjidx0E        = 2*vdwtype[jnrE+0];
1455             vdwjidx0F        = 2*vdwtype[jnrF+0];
1456             vdwjidx0G        = 2*vdwtype[jnrG+0];
1457             vdwjidx0H        = 2*vdwtype[jnrH+0];
1458
1459             fjx0             = _mm256_setzero_ps();
1460             fjy0             = _mm256_setzero_ps();
1461             fjz0             = _mm256_setzero_ps();
1462
1463             /**************************
1464              * CALCULATE INTERACTIONS *
1465              **************************/
1466
1467             if (gmx_mm256_any_lt(rsq00,rcutoff2))
1468             {
1469
1470             r00              = _mm256_mul_ps(rsq00,rinv00);
1471             r00              = _mm256_andnot_ps(dummy_mask,r00);
1472
1473             /* Compute parameters for interactions between i and j atoms */
1474             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr0+vdwjidx0A,
1475                                             vdwioffsetptr0+vdwjidx0B,
1476                                             vdwioffsetptr0+vdwjidx0C,
1477                                             vdwioffsetptr0+vdwjidx0D,
1478                                             vdwioffsetptr0+vdwjidx0E,
1479                                             vdwioffsetptr0+vdwjidx0F,
1480                                             vdwioffsetptr0+vdwjidx0G,
1481                                             vdwioffsetptr0+vdwjidx0H,
1482                                             &c6_00,&c12_00);
1483
1484             c6grid_00       = gmx_mm256_load_8real_swizzle_ps(vdwgridioffsetptr0+vdwjidx0A,
1485                                                                   vdwgridioffsetptr0+vdwjidx0B,
1486                                                                   vdwgridioffsetptr0+vdwjidx0C,
1487                                                                   vdwgridioffsetptr0+vdwjidx0D,
1488                                                                   vdwgridioffsetptr0+vdwjidx0E,
1489                                                                   vdwgridioffsetptr0+vdwjidx0F,
1490                                                                   vdwgridioffsetptr0+vdwjidx0G,
1491                                                                   vdwgridioffsetptr0+vdwjidx0H);
1492
1493             /* Analytical LJ-PME */
1494             rinvsix          = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
1495             ewcljrsq         = _mm256_mul_ps(ewclj2,rsq00);
1496             ewclj6           = _mm256_mul_ps(ewclj2,_mm256_mul_ps(ewclj2,ewclj2));
1497             exponent         = gmx_simd_exp_r(ewcljrsq);
1498             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
1499             poly             = _mm256_mul_ps(exponent,_mm256_add_ps(_mm256_sub_ps(one,ewcljrsq),_mm256_mul_ps(_mm256_mul_ps(ewcljrsq,ewcljrsq),one_half)));
1500             /* f6A = 6 * C6grid * (1 - poly) */
1501             f6A              = _mm256_mul_ps(c6grid_00,_mm256_sub_ps(one,poly));
1502             /* f6B = C6grid * exponent * beta^6 */
1503             f6B              = _mm256_mul_ps(_mm256_mul_ps(c6grid_00,one_sixth),_mm256_mul_ps(exponent,ewclj6));
1504             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
1505             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);
1506
1507             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
1508
1509             fscal            = fvdw;
1510
1511             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1512
1513             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1514
1515             /* Calculate temporary vectorial force */
1516             tx               = _mm256_mul_ps(fscal,dx00);
1517             ty               = _mm256_mul_ps(fscal,dy00);
1518             tz               = _mm256_mul_ps(fscal,dz00);
1519
1520             /* Update vectorial force */
1521             fix0             = _mm256_add_ps(fix0,tx);
1522             fiy0             = _mm256_add_ps(fiy0,ty);
1523             fiz0             = _mm256_add_ps(fiz0,tz);
1524
1525             fjx0             = _mm256_add_ps(fjx0,tx);
1526             fjy0             = _mm256_add_ps(fjy0,ty);
1527             fjz0             = _mm256_add_ps(fjz0,tz);
1528
1529             }
1530
1531             /**************************
1532              * CALCULATE INTERACTIONS *
1533              **************************/
1534
1535             if (gmx_mm256_any_lt(rsq10,rcutoff2))
1536             {
1537
1538             r10              = _mm256_mul_ps(rsq10,rinv10);
1539             r10              = _mm256_andnot_ps(dummy_mask,r10);
1540
1541             /* Compute parameters for interactions between i and j atoms */
1542             qq10             = _mm256_mul_ps(iq1,jq0);
1543
1544             /* EWALD ELECTROSTATICS */
1545             
1546             /* Analytical PME correction */
1547             zeta2            = _mm256_mul_ps(beta2,rsq10);
1548             rinv3            = _mm256_mul_ps(rinvsq10,rinv10);
1549             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
1550             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1551             felec            = _mm256_mul_ps(qq10,felec);
1552             
1553             cutoff_mask      = _mm256_cmp_ps(rsq10,rcutoff2,_CMP_LT_OQ);
1554
1555             fscal            = felec;
1556
1557             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1558
1559             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1560
1561             /* Calculate temporary vectorial force */
1562             tx               = _mm256_mul_ps(fscal,dx10);
1563             ty               = _mm256_mul_ps(fscal,dy10);
1564             tz               = _mm256_mul_ps(fscal,dz10);
1565
1566             /* Update vectorial force */
1567             fix1             = _mm256_add_ps(fix1,tx);
1568             fiy1             = _mm256_add_ps(fiy1,ty);
1569             fiz1             = _mm256_add_ps(fiz1,tz);
1570
1571             fjx0             = _mm256_add_ps(fjx0,tx);
1572             fjy0             = _mm256_add_ps(fjy0,ty);
1573             fjz0             = _mm256_add_ps(fjz0,tz);
1574
1575             }
1576
1577             /**************************
1578              * CALCULATE INTERACTIONS *
1579              **************************/
1580
1581             if (gmx_mm256_any_lt(rsq20,rcutoff2))
1582             {
1583
1584             r20              = _mm256_mul_ps(rsq20,rinv20);
1585             r20              = _mm256_andnot_ps(dummy_mask,r20);
1586
1587             /* Compute parameters for interactions between i and j atoms */
1588             qq20             = _mm256_mul_ps(iq2,jq0);
1589
1590             /* EWALD ELECTROSTATICS */
1591             
1592             /* Analytical PME correction */
1593             zeta2            = _mm256_mul_ps(beta2,rsq20);
1594             rinv3            = _mm256_mul_ps(rinvsq20,rinv20);
1595             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
1596             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1597             felec            = _mm256_mul_ps(qq20,felec);
1598             
1599             cutoff_mask      = _mm256_cmp_ps(rsq20,rcutoff2,_CMP_LT_OQ);
1600
1601             fscal            = felec;
1602
1603             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1604
1605             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1606
1607             /* Calculate temporary vectorial force */
1608             tx               = _mm256_mul_ps(fscal,dx20);
1609             ty               = _mm256_mul_ps(fscal,dy20);
1610             tz               = _mm256_mul_ps(fscal,dz20);
1611
1612             /* Update vectorial force */
1613             fix2             = _mm256_add_ps(fix2,tx);
1614             fiy2             = _mm256_add_ps(fiy2,ty);
1615             fiz2             = _mm256_add_ps(fiz2,tz);
1616
1617             fjx0             = _mm256_add_ps(fjx0,tx);
1618             fjy0             = _mm256_add_ps(fjy0,ty);
1619             fjz0             = _mm256_add_ps(fjz0,tz);
1620
1621             }
1622
1623             /**************************
1624              * CALCULATE INTERACTIONS *
1625              **************************/
1626
1627             if (gmx_mm256_any_lt(rsq30,rcutoff2))
1628             {
1629
1630             r30              = _mm256_mul_ps(rsq30,rinv30);
1631             r30              = _mm256_andnot_ps(dummy_mask,r30);
1632
1633             /* Compute parameters for interactions between i and j atoms */
1634             qq30             = _mm256_mul_ps(iq3,jq0);
1635
1636             /* EWALD ELECTROSTATICS */
1637             
1638             /* Analytical PME correction */
1639             zeta2            = _mm256_mul_ps(beta2,rsq30);
1640             rinv3            = _mm256_mul_ps(rinvsq30,rinv30);
1641             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
1642             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
1643             felec            = _mm256_mul_ps(qq30,felec);
1644             
1645             cutoff_mask      = _mm256_cmp_ps(rsq30,rcutoff2,_CMP_LT_OQ);
1646
1647             fscal            = felec;
1648
1649             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1650
1651             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1652
1653             /* Calculate temporary vectorial force */
1654             tx               = _mm256_mul_ps(fscal,dx30);
1655             ty               = _mm256_mul_ps(fscal,dy30);
1656             tz               = _mm256_mul_ps(fscal,dz30);
1657
1658             /* Update vectorial force */
1659             fix3             = _mm256_add_ps(fix3,tx);
1660             fiy3             = _mm256_add_ps(fiy3,ty);
1661             fiz3             = _mm256_add_ps(fiz3,tz);
1662
1663             fjx0             = _mm256_add_ps(fjx0,tx);
1664             fjy0             = _mm256_add_ps(fjy0,ty);
1665             fjz0             = _mm256_add_ps(fjz0,tz);
1666
1667             }
1668
1669             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1670             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1671             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1672             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1673             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1674             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1675             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1676             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1677
1678             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
1679
1680             /* Inner loop uses 233 flops */
1681         }
1682
1683         /* End of innermost loop */
1684
1685         gmx_mm256_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1686                                                  f+i_coord_offset,fshift+i_shift_offset);
1687
1688         /* Increment number of inner iterations */
1689         inneriter                  += j_index_end - j_index_start;
1690
1691         /* Outer loop uses 24 flops */
1692     }
1693
1694     /* Increment number of outer iterations */
1695     outeriter        += nri;
1696
1697     /* Update outer/inner flops */
1698
1699     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*233);
1700 }