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