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