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