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