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