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