Introduce gmxpre.h for truly global definitions
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_double / nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_avx_256_double.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_double 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_double.h"
50 #include "kernelutil_x86_avx_256_double.h"
51
52 /*
53  * Gromacs nonbonded kernel:   nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_avx_256_double
54  * Electrostatics interaction: None
55  * VdW interaction:            LJEwald
56  * Geometry:                   Particle-Particle
57  * Calculate force/pot:        PotentialAndForce
58  */
59 void
60 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_VF_avx_256_double
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 refer to j loop unrolling done with AVX, e.g. for the four 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              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
78     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
79     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
80     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
81     real             rcutoff_scalar;
82     real             *shiftvec,*fshift,*x,*f;
83     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
84     real             scratch[4*DIM];
85     __m256d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
86     real *           vdwioffsetptr0;
87     real *           vdwgridioffsetptr0;
88     __m256d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
89     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
90     __m256d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
91     __m256d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
92     int              nvdwtype;
93     __m256d          rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
94     int              *vdwtype;
95     real             *vdwparam;
96     __m256d          one_sixth   = _mm256_set1_pd(1.0/6.0);
97     __m256d          one_twelfth = _mm256_set1_pd(1.0/12.0);
98     __m256d           c6grid_00;
99     real             *vdwgridparam;
100     __m256d           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
101     __m256d           one_half  = _mm256_set1_pd(0.5);
102     __m256d           minus_one = _mm256_set1_pd(-1.0);
103     __m256d          dummy_mask,cutoff_mask;
104     __m128           tmpmask0,tmpmask1;
105     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
106     __m256d          one     = _mm256_set1_pd(1.0);
107     __m256d          two     = _mm256_set1_pd(2.0);
108     x                = xx[0];
109     f                = ff[0];
110
111     nri              = nlist->nri;
112     iinr             = nlist->iinr;
113     jindex           = nlist->jindex;
114     jjnr             = nlist->jjnr;
115     shiftidx         = nlist->shift;
116     gid              = nlist->gid;
117     shiftvec         = fr->shift_vec[0];
118     fshift           = fr->fshift[0];
119     nvdwtype         = fr->ntype;
120     vdwparam         = fr->nbfp;
121     vdwtype          = mdatoms->typeA;
122     vdwgridparam     = fr->ljpme_c6grid;
123     sh_lj_ewald      = _mm256_set1_pd(fr->ic->sh_lj_ewald);
124     ewclj            = _mm256_set1_pd(fr->ewaldcoeff_lj);
125     ewclj2           = _mm256_mul_pd(minus_one,_mm256_mul_pd(ewclj,ewclj));
126
127     rcutoff_scalar   = fr->rvdw;
128     rcutoff          = _mm256_set1_pd(rcutoff_scalar);
129     rcutoff2         = _mm256_mul_pd(rcutoff,rcutoff);
130
131     sh_vdw_invrcut6  = _mm256_set1_pd(fr->ic->sh_invrc6);
132     rvdw             = _mm256_set1_pd(fr->rvdw);
133
134     /* Avoid stupid compiler warnings */
135     jnrA = jnrB = jnrC = jnrD = 0;
136     j_coord_offsetA = 0;
137     j_coord_offsetB = 0;
138     j_coord_offsetC = 0;
139     j_coord_offsetD = 0;
140
141     outeriter        = 0;
142     inneriter        = 0;
143
144     for(iidx=0;iidx<4*DIM;iidx++)
145     {
146         scratch[iidx] = 0.0;
147     }
148
149     /* Start outer loop over neighborlists */
150     for(iidx=0; iidx<nri; iidx++)
151     {
152         /* Load shift vector for this list */
153         i_shift_offset   = DIM*shiftidx[iidx];
154
155         /* Load limits for loop over neighbors */
156         j_index_start    = jindex[iidx];
157         j_index_end      = jindex[iidx+1];
158
159         /* Get outer coordinate index */
160         inr              = iinr[iidx];
161         i_coord_offset   = DIM*inr;
162
163         /* Load i particle coords and add shift vector */
164         gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
165
166         fix0             = _mm256_setzero_pd();
167         fiy0             = _mm256_setzero_pd();
168         fiz0             = _mm256_setzero_pd();
169
170         /* Load parameters for i particles */
171         vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
172         vdwgridioffsetptr0 = vdwgridparam+2*nvdwtype*vdwtype[inr+0];
173
174         /* Reset potential sums */
175         vvdwsum          = _mm256_setzero_pd();
176
177         /* Start inner kernel loop */
178         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
179         {
180
181             /* Get j neighbor index, and coordinate index */
182             jnrA             = jjnr[jidx];
183             jnrB             = jjnr[jidx+1];
184             jnrC             = jjnr[jidx+2];
185             jnrD             = jjnr[jidx+3];
186             j_coord_offsetA  = DIM*jnrA;
187             j_coord_offsetB  = DIM*jnrB;
188             j_coord_offsetC  = DIM*jnrC;
189             j_coord_offsetD  = DIM*jnrD;
190
191             /* load j atom coordinates */
192             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
193                                                  x+j_coord_offsetC,x+j_coord_offsetD,
194                                                  &jx0,&jy0,&jz0);
195
196             /* Calculate displacement vector */
197             dx00             = _mm256_sub_pd(ix0,jx0);
198             dy00             = _mm256_sub_pd(iy0,jy0);
199             dz00             = _mm256_sub_pd(iz0,jz0);
200
201             /* Calculate squared distance and things based on it */
202             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
203
204             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
205
206             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
207
208             /* Load parameters for j particles */
209             vdwjidx0A        = 2*vdwtype[jnrA+0];
210             vdwjidx0B        = 2*vdwtype[jnrB+0];
211             vdwjidx0C        = 2*vdwtype[jnrC+0];
212             vdwjidx0D        = 2*vdwtype[jnrD+0];
213
214             /**************************
215              * CALCULATE INTERACTIONS *
216              **************************/
217
218             if (gmx_mm256_any_lt(rsq00,rcutoff2))
219             {
220
221             r00              = _mm256_mul_pd(rsq00,rinv00);
222
223             /* Compute parameters for interactions between i and j atoms */
224             gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
225                                             vdwioffsetptr0+vdwjidx0B,
226                                             vdwioffsetptr0+vdwjidx0C,
227                                             vdwioffsetptr0+vdwjidx0D,
228                                             &c6_00,&c12_00);
229
230             c6grid_00       = gmx_mm256_load_4real_swizzle_pd(vdwgridioffsetptr0+vdwjidx0A,
231                                                                   vdwgridioffsetptr0+vdwjidx0B,
232                                                                   vdwgridioffsetptr0+vdwjidx0C,
233                                                                   vdwgridioffsetptr0+vdwjidx0D);
234
235             /* Analytical LJ-PME */
236             rinvsix          = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
237             ewcljrsq         = _mm256_mul_pd(ewclj2,rsq00);
238             ewclj6           = _mm256_mul_pd(ewclj2,_mm256_mul_pd(ewclj2,ewclj2));
239             exponent         = gmx_simd_exp_d(ewcljrsq);
240             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
241             poly             = _mm256_mul_pd(exponent,_mm256_add_pd(_mm256_sub_pd(one,ewcljrsq),_mm256_mul_pd(_mm256_mul_pd(ewcljrsq,ewcljrsq),one_half)));
242             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
243             vvdw6            = _mm256_mul_pd(_mm256_sub_pd(c6_00,_mm256_mul_pd(c6grid_00,_mm256_sub_pd(one,poly))),rinvsix);
244             vvdw12           = _mm256_mul_pd(c12_00,_mm256_mul_pd(rinvsix,rinvsix));
245             vvdw             = _mm256_sub_pd(_mm256_mul_pd( _mm256_sub_pd(vvdw12 , _mm256_mul_pd(c12_00,_mm256_mul_pd(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
246                                           _mm256_mul_pd( _mm256_sub_pd(vvdw6,_mm256_add_pd(_mm256_mul_pd(c6_00,sh_vdw_invrcut6),_mm256_mul_pd(c6grid_00,sh_lj_ewald))),one_sixth));
247             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
248             fvdw             = _mm256_mul_pd(_mm256_sub_pd(vvdw12,_mm256_sub_pd(vvdw6,_mm256_mul_pd(_mm256_mul_pd(c6grid_00,one_sixth),_mm256_mul_pd(exponent,ewclj6)))),rinvsq00);
249
250             cutoff_mask      = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
251
252             /* Update potential sum for this i atom from the interaction with this j atom. */
253             vvdw             = _mm256_and_pd(vvdw,cutoff_mask);
254             vvdwsum          = _mm256_add_pd(vvdwsum,vvdw);
255
256             fscal            = fvdw;
257
258             fscal            = _mm256_and_pd(fscal,cutoff_mask);
259
260             /* Calculate temporary vectorial force */
261             tx               = _mm256_mul_pd(fscal,dx00);
262             ty               = _mm256_mul_pd(fscal,dy00);
263             tz               = _mm256_mul_pd(fscal,dz00);
264
265             /* Update vectorial force */
266             fix0             = _mm256_add_pd(fix0,tx);
267             fiy0             = _mm256_add_pd(fiy0,ty);
268             fiz0             = _mm256_add_pd(fiz0,tz);
269
270             fjptrA             = f+j_coord_offsetA;
271             fjptrB             = f+j_coord_offsetB;
272             fjptrC             = f+j_coord_offsetC;
273             fjptrD             = f+j_coord_offsetD;
274             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
275
276             }
277
278             /* Inner loop uses 62 flops */
279         }
280
281         if(jidx<j_index_end)
282         {
283
284             /* Get j neighbor index, and coordinate index */
285             jnrlistA         = jjnr[jidx];
286             jnrlistB         = jjnr[jidx+1];
287             jnrlistC         = jjnr[jidx+2];
288             jnrlistD         = jjnr[jidx+3];
289             /* Sign of each element will be negative for non-real atoms.
290              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
291              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
292              */
293             tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
294
295             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
296             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
297             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
298
299             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
300             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
301             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
302             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
303             j_coord_offsetA  = DIM*jnrA;
304             j_coord_offsetB  = DIM*jnrB;
305             j_coord_offsetC  = DIM*jnrC;
306             j_coord_offsetD  = DIM*jnrD;
307
308             /* load j atom coordinates */
309             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
310                                                  x+j_coord_offsetC,x+j_coord_offsetD,
311                                                  &jx0,&jy0,&jz0);
312
313             /* Calculate displacement vector */
314             dx00             = _mm256_sub_pd(ix0,jx0);
315             dy00             = _mm256_sub_pd(iy0,jy0);
316             dz00             = _mm256_sub_pd(iz0,jz0);
317
318             /* Calculate squared distance and things based on it */
319             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
320
321             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
322
323             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
324
325             /* Load parameters for j particles */
326             vdwjidx0A        = 2*vdwtype[jnrA+0];
327             vdwjidx0B        = 2*vdwtype[jnrB+0];
328             vdwjidx0C        = 2*vdwtype[jnrC+0];
329             vdwjidx0D        = 2*vdwtype[jnrD+0];
330
331             /**************************
332              * CALCULATE INTERACTIONS *
333              **************************/
334
335             if (gmx_mm256_any_lt(rsq00,rcutoff2))
336             {
337
338             r00              = _mm256_mul_pd(rsq00,rinv00);
339             r00              = _mm256_andnot_pd(dummy_mask,r00);
340
341             /* Compute parameters for interactions between i and j atoms */
342             gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
343                                             vdwioffsetptr0+vdwjidx0B,
344                                             vdwioffsetptr0+vdwjidx0C,
345                                             vdwioffsetptr0+vdwjidx0D,
346                                             &c6_00,&c12_00);
347
348             c6grid_00       = gmx_mm256_load_4real_swizzle_pd(vdwgridioffsetptr0+vdwjidx0A,
349                                                                   vdwgridioffsetptr0+vdwjidx0B,
350                                                                   vdwgridioffsetptr0+vdwjidx0C,
351                                                                   vdwgridioffsetptr0+vdwjidx0D);
352
353             /* Analytical LJ-PME */
354             rinvsix          = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
355             ewcljrsq         = _mm256_mul_pd(ewclj2,rsq00);
356             ewclj6           = _mm256_mul_pd(ewclj2,_mm256_mul_pd(ewclj2,ewclj2));
357             exponent         = gmx_simd_exp_d(ewcljrsq);
358             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
359             poly             = _mm256_mul_pd(exponent,_mm256_add_pd(_mm256_sub_pd(one,ewcljrsq),_mm256_mul_pd(_mm256_mul_pd(ewcljrsq,ewcljrsq),one_half)));
360             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
361             vvdw6            = _mm256_mul_pd(_mm256_sub_pd(c6_00,_mm256_mul_pd(c6grid_00,_mm256_sub_pd(one,poly))),rinvsix);
362             vvdw12           = _mm256_mul_pd(c12_00,_mm256_mul_pd(rinvsix,rinvsix));
363             vvdw             = _mm256_sub_pd(_mm256_mul_pd( _mm256_sub_pd(vvdw12 , _mm256_mul_pd(c12_00,_mm256_mul_pd(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
364                                           _mm256_mul_pd( _mm256_sub_pd(vvdw6,_mm256_add_pd(_mm256_mul_pd(c6_00,sh_vdw_invrcut6),_mm256_mul_pd(c6grid_00,sh_lj_ewald))),one_sixth));
365             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
366             fvdw             = _mm256_mul_pd(_mm256_sub_pd(vvdw12,_mm256_sub_pd(vvdw6,_mm256_mul_pd(_mm256_mul_pd(c6grid_00,one_sixth),_mm256_mul_pd(exponent,ewclj6)))),rinvsq00);
367
368             cutoff_mask      = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
369
370             /* Update potential sum for this i atom from the interaction with this j atom. */
371             vvdw             = _mm256_and_pd(vvdw,cutoff_mask);
372             vvdw             = _mm256_andnot_pd(dummy_mask,vvdw);
373             vvdwsum          = _mm256_add_pd(vvdwsum,vvdw);
374
375             fscal            = fvdw;
376
377             fscal            = _mm256_and_pd(fscal,cutoff_mask);
378
379             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
380
381             /* Calculate temporary vectorial force */
382             tx               = _mm256_mul_pd(fscal,dx00);
383             ty               = _mm256_mul_pd(fscal,dy00);
384             tz               = _mm256_mul_pd(fscal,dz00);
385
386             /* Update vectorial force */
387             fix0             = _mm256_add_pd(fix0,tx);
388             fiy0             = _mm256_add_pd(fiy0,ty);
389             fiz0             = _mm256_add_pd(fiz0,tz);
390
391             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
392             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
393             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
394             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
395             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
396
397             }
398
399             /* Inner loop uses 63 flops */
400         }
401
402         /* End of innermost loop */
403
404         gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
405                                                  f+i_coord_offset,fshift+i_shift_offset);
406
407         ggid                        = gid[iidx];
408         /* Update potential energies */
409         gmx_mm256_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
410
411         /* Increment number of inner iterations */
412         inneriter                  += j_index_end - j_index_start;
413
414         /* Outer loop uses 7 flops */
415     }
416
417     /* Increment number of outer iterations */
418     outeriter        += nri;
419
420     /* Update outer/inner flops */
421
422     inc_nrnb(nrnb,eNR_NBKERNEL_VDW_VF,outeriter*7 + inneriter*63);
423 }
424 /*
425  * Gromacs nonbonded kernel:   nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_avx_256_double
426  * Electrostatics interaction: None
427  * VdW interaction:            LJEwald
428  * Geometry:                   Particle-Particle
429  * Calculate force/pot:        Force
430  */
431 void
432 nb_kernel_ElecNone_VdwLJEwSh_GeomP1P1_F_avx_256_double
433                     (t_nblist                    * gmx_restrict       nlist,
434                      rvec                        * gmx_restrict          xx,
435                      rvec                        * gmx_restrict          ff,
436                      t_forcerec                  * gmx_restrict          fr,
437                      t_mdatoms                   * gmx_restrict     mdatoms,
438                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
439                      t_nrnb                      * gmx_restrict        nrnb)
440 {
441     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
442      * just 0 for non-waters.
443      * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
444      * jnr indices corresponding to data put in the four positions in the SIMD register.
445      */
446     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
447     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
448     int              jnrA,jnrB,jnrC,jnrD;
449     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
450     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
451     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
452     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
453     real             rcutoff_scalar;
454     real             *shiftvec,*fshift,*x,*f;
455     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
456     real             scratch[4*DIM];
457     __m256d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
458     real *           vdwioffsetptr0;
459     real *           vdwgridioffsetptr0;
460     __m256d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
461     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
462     __m256d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
463     __m256d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
464     int              nvdwtype;
465     __m256d          rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
466     int              *vdwtype;
467     real             *vdwparam;
468     __m256d          one_sixth   = _mm256_set1_pd(1.0/6.0);
469     __m256d          one_twelfth = _mm256_set1_pd(1.0/12.0);
470     __m256d           c6grid_00;
471     real             *vdwgridparam;
472     __m256d           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
473     __m256d           one_half  = _mm256_set1_pd(0.5);
474     __m256d           minus_one = _mm256_set1_pd(-1.0);
475     __m256d          dummy_mask,cutoff_mask;
476     __m128           tmpmask0,tmpmask1;
477     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
478     __m256d          one     = _mm256_set1_pd(1.0);
479     __m256d          two     = _mm256_set1_pd(2.0);
480     x                = xx[0];
481     f                = ff[0];
482
483     nri              = nlist->nri;
484     iinr             = nlist->iinr;
485     jindex           = nlist->jindex;
486     jjnr             = nlist->jjnr;
487     shiftidx         = nlist->shift;
488     gid              = nlist->gid;
489     shiftvec         = fr->shift_vec[0];
490     fshift           = fr->fshift[0];
491     nvdwtype         = fr->ntype;
492     vdwparam         = fr->nbfp;
493     vdwtype          = mdatoms->typeA;
494     vdwgridparam     = fr->ljpme_c6grid;
495     sh_lj_ewald      = _mm256_set1_pd(fr->ic->sh_lj_ewald);
496     ewclj            = _mm256_set1_pd(fr->ewaldcoeff_lj);
497     ewclj2           = _mm256_mul_pd(minus_one,_mm256_mul_pd(ewclj,ewclj));
498
499     rcutoff_scalar   = fr->rvdw;
500     rcutoff          = _mm256_set1_pd(rcutoff_scalar);
501     rcutoff2         = _mm256_mul_pd(rcutoff,rcutoff);
502
503     sh_vdw_invrcut6  = _mm256_set1_pd(fr->ic->sh_invrc6);
504     rvdw             = _mm256_set1_pd(fr->rvdw);
505
506     /* Avoid stupid compiler warnings */
507     jnrA = jnrB = jnrC = jnrD = 0;
508     j_coord_offsetA = 0;
509     j_coord_offsetB = 0;
510     j_coord_offsetC = 0;
511     j_coord_offsetD = 0;
512
513     outeriter        = 0;
514     inneriter        = 0;
515
516     for(iidx=0;iidx<4*DIM;iidx++)
517     {
518         scratch[iidx] = 0.0;
519     }
520
521     /* Start outer loop over neighborlists */
522     for(iidx=0; iidx<nri; iidx++)
523     {
524         /* Load shift vector for this list */
525         i_shift_offset   = DIM*shiftidx[iidx];
526
527         /* Load limits for loop over neighbors */
528         j_index_start    = jindex[iidx];
529         j_index_end      = jindex[iidx+1];
530
531         /* Get outer coordinate index */
532         inr              = iinr[iidx];
533         i_coord_offset   = DIM*inr;
534
535         /* Load i particle coords and add shift vector */
536         gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
537
538         fix0             = _mm256_setzero_pd();
539         fiy0             = _mm256_setzero_pd();
540         fiz0             = _mm256_setzero_pd();
541
542         /* Load parameters for i particles */
543         vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
544         vdwgridioffsetptr0 = vdwgridparam+2*nvdwtype*vdwtype[inr+0];
545
546         /* Start inner kernel loop */
547         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
548         {
549
550             /* Get j neighbor index, and coordinate index */
551             jnrA             = jjnr[jidx];
552             jnrB             = jjnr[jidx+1];
553             jnrC             = jjnr[jidx+2];
554             jnrD             = jjnr[jidx+3];
555             j_coord_offsetA  = DIM*jnrA;
556             j_coord_offsetB  = DIM*jnrB;
557             j_coord_offsetC  = DIM*jnrC;
558             j_coord_offsetD  = DIM*jnrD;
559
560             /* load j atom coordinates */
561             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
562                                                  x+j_coord_offsetC,x+j_coord_offsetD,
563                                                  &jx0,&jy0,&jz0);
564
565             /* Calculate displacement vector */
566             dx00             = _mm256_sub_pd(ix0,jx0);
567             dy00             = _mm256_sub_pd(iy0,jy0);
568             dz00             = _mm256_sub_pd(iz0,jz0);
569
570             /* Calculate squared distance and things based on it */
571             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
572
573             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
574
575             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
576
577             /* Load parameters for j particles */
578             vdwjidx0A        = 2*vdwtype[jnrA+0];
579             vdwjidx0B        = 2*vdwtype[jnrB+0];
580             vdwjidx0C        = 2*vdwtype[jnrC+0];
581             vdwjidx0D        = 2*vdwtype[jnrD+0];
582
583             /**************************
584              * CALCULATE INTERACTIONS *
585              **************************/
586
587             if (gmx_mm256_any_lt(rsq00,rcutoff2))
588             {
589
590             r00              = _mm256_mul_pd(rsq00,rinv00);
591
592             /* Compute parameters for interactions between i and j atoms */
593             gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
594                                             vdwioffsetptr0+vdwjidx0B,
595                                             vdwioffsetptr0+vdwjidx0C,
596                                             vdwioffsetptr0+vdwjidx0D,
597                                             &c6_00,&c12_00);
598
599             c6grid_00       = gmx_mm256_load_4real_swizzle_pd(vdwgridioffsetptr0+vdwjidx0A,
600                                                                   vdwgridioffsetptr0+vdwjidx0B,
601                                                                   vdwgridioffsetptr0+vdwjidx0C,
602                                                                   vdwgridioffsetptr0+vdwjidx0D);
603
604             /* Analytical LJ-PME */
605             rinvsix          = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
606             ewcljrsq         = _mm256_mul_pd(ewclj2,rsq00);
607             ewclj6           = _mm256_mul_pd(ewclj2,_mm256_mul_pd(ewclj2,ewclj2));
608             exponent         = gmx_simd_exp_d(ewcljrsq);
609             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
610             poly             = _mm256_mul_pd(exponent,_mm256_add_pd(_mm256_sub_pd(one,ewcljrsq),_mm256_mul_pd(_mm256_mul_pd(ewcljrsq,ewcljrsq),one_half)));
611             /* f6A = 6 * C6grid * (1 - poly) */
612             f6A              = _mm256_mul_pd(c6grid_00,_mm256_sub_pd(one,poly));
613             /* f6B = C6grid * exponent * beta^6 */
614             f6B              = _mm256_mul_pd(_mm256_mul_pd(c6grid_00,one_sixth),_mm256_mul_pd(exponent,ewclj6));
615             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
616             fvdw              = _mm256_mul_pd(_mm256_add_pd(_mm256_mul_pd(_mm256_sub_pd(_mm256_mul_pd(c12_00,rinvsix),_mm256_sub_pd(c6_00,f6A)),rinvsix),f6B),rinvsq00);
617
618             cutoff_mask      = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
619
620             fscal            = fvdw;
621
622             fscal            = _mm256_and_pd(fscal,cutoff_mask);
623
624             /* Calculate temporary vectorial force */
625             tx               = _mm256_mul_pd(fscal,dx00);
626             ty               = _mm256_mul_pd(fscal,dy00);
627             tz               = _mm256_mul_pd(fscal,dz00);
628
629             /* Update vectorial force */
630             fix0             = _mm256_add_pd(fix0,tx);
631             fiy0             = _mm256_add_pd(fiy0,ty);
632             fiz0             = _mm256_add_pd(fiz0,tz);
633
634             fjptrA             = f+j_coord_offsetA;
635             fjptrB             = f+j_coord_offsetB;
636             fjptrC             = f+j_coord_offsetC;
637             fjptrD             = f+j_coord_offsetD;
638             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
639
640             }
641
642             /* Inner loop uses 49 flops */
643         }
644
645         if(jidx<j_index_end)
646         {
647
648             /* Get j neighbor index, and coordinate index */
649             jnrlistA         = jjnr[jidx];
650             jnrlistB         = jjnr[jidx+1];
651             jnrlistC         = jjnr[jidx+2];
652             jnrlistD         = jjnr[jidx+3];
653             /* Sign of each element will be negative for non-real atoms.
654              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
655              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
656              */
657             tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
658
659             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
660             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
661             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
662
663             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
664             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
665             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
666             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
667             j_coord_offsetA  = DIM*jnrA;
668             j_coord_offsetB  = DIM*jnrB;
669             j_coord_offsetC  = DIM*jnrC;
670             j_coord_offsetD  = DIM*jnrD;
671
672             /* load j atom coordinates */
673             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
674                                                  x+j_coord_offsetC,x+j_coord_offsetD,
675                                                  &jx0,&jy0,&jz0);
676
677             /* Calculate displacement vector */
678             dx00             = _mm256_sub_pd(ix0,jx0);
679             dy00             = _mm256_sub_pd(iy0,jy0);
680             dz00             = _mm256_sub_pd(iz0,jz0);
681
682             /* Calculate squared distance and things based on it */
683             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
684
685             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
686
687             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
688
689             /* Load parameters for j particles */
690             vdwjidx0A        = 2*vdwtype[jnrA+0];
691             vdwjidx0B        = 2*vdwtype[jnrB+0];
692             vdwjidx0C        = 2*vdwtype[jnrC+0];
693             vdwjidx0D        = 2*vdwtype[jnrD+0];
694
695             /**************************
696              * CALCULATE INTERACTIONS *
697              **************************/
698
699             if (gmx_mm256_any_lt(rsq00,rcutoff2))
700             {
701
702             r00              = _mm256_mul_pd(rsq00,rinv00);
703             r00              = _mm256_andnot_pd(dummy_mask,r00);
704
705             /* Compute parameters for interactions between i and j atoms */
706             gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr0+vdwjidx0A,
707                                             vdwioffsetptr0+vdwjidx0B,
708                                             vdwioffsetptr0+vdwjidx0C,
709                                             vdwioffsetptr0+vdwjidx0D,
710                                             &c6_00,&c12_00);
711
712             c6grid_00       = gmx_mm256_load_4real_swizzle_pd(vdwgridioffsetptr0+vdwjidx0A,
713                                                                   vdwgridioffsetptr0+vdwjidx0B,
714                                                                   vdwgridioffsetptr0+vdwjidx0C,
715                                                                   vdwgridioffsetptr0+vdwjidx0D);
716
717             /* Analytical LJ-PME */
718             rinvsix          = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
719             ewcljrsq         = _mm256_mul_pd(ewclj2,rsq00);
720             ewclj6           = _mm256_mul_pd(ewclj2,_mm256_mul_pd(ewclj2,ewclj2));
721             exponent         = gmx_simd_exp_d(ewcljrsq);
722             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
723             poly             = _mm256_mul_pd(exponent,_mm256_add_pd(_mm256_sub_pd(one,ewcljrsq),_mm256_mul_pd(_mm256_mul_pd(ewcljrsq,ewcljrsq),one_half)));
724             /* f6A = 6 * C6grid * (1 - poly) */
725             f6A              = _mm256_mul_pd(c6grid_00,_mm256_sub_pd(one,poly));
726             /* f6B = C6grid * exponent * beta^6 */
727             f6B              = _mm256_mul_pd(_mm256_mul_pd(c6grid_00,one_sixth),_mm256_mul_pd(exponent,ewclj6));
728             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
729             fvdw              = _mm256_mul_pd(_mm256_add_pd(_mm256_mul_pd(_mm256_sub_pd(_mm256_mul_pd(c12_00,rinvsix),_mm256_sub_pd(c6_00,f6A)),rinvsix),f6B),rinvsq00);
730
731             cutoff_mask      = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
732
733             fscal            = fvdw;
734
735             fscal            = _mm256_and_pd(fscal,cutoff_mask);
736
737             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
738
739             /* Calculate temporary vectorial force */
740             tx               = _mm256_mul_pd(fscal,dx00);
741             ty               = _mm256_mul_pd(fscal,dy00);
742             tz               = _mm256_mul_pd(fscal,dz00);
743
744             /* Update vectorial force */
745             fix0             = _mm256_add_pd(fix0,tx);
746             fiy0             = _mm256_add_pd(fiy0,ty);
747             fiz0             = _mm256_add_pd(fiz0,tz);
748
749             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
750             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
751             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
752             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
753             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
754
755             }
756
757             /* Inner loop uses 50 flops */
758         }
759
760         /* End of innermost loop */
761
762         gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
763                                                  f+i_coord_offset,fshift+i_shift_offset);
764
765         /* Increment number of inner iterations */
766         inneriter                  += j_index_end - j_index_start;
767
768         /* Outer loop uses 6 flops */
769     }
770
771     /* Increment number of outer iterations */
772     outeriter        += nri;
773
774     /* Update outer/inner flops */
775
776     inc_nrnb(nrnb,eNR_NBKERNEL_VDW_F,outeriter*6 + inneriter*50);
777 }