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