Introduce gmxpre.h for truly global definitions
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_single / nb_kernel_ElecEw_VdwLJEw_GeomP1P1_avx_128_fma_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_128_fma_single kernel generator.
37  */
38 #include "gmxpre.h"
39
40 #include "config.h"
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "gromacs/legacyheaders/types/simple.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/legacyheaders/nrnb.h"
48
49 #include "gromacs/simd/math_x86_avx_128_fma_single.h"
50 #include "kernelutil_x86_avx_128_fma_single.h"
51
52 /*
53  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwLJEw_GeomP1P1_VF_avx_128_fma_single
54  * Electrostatics interaction: Ewald
55  * VdW interaction:            LJEwald
56  * Geometry:                   Particle-Particle
57  * Calculate force/pot:        PotentialAndForce
58  */
59 void
60 nb_kernel_ElecEw_VdwLJEw_GeomP1P1_VF_avx_128_fma_single
61                     (t_nblist                    * gmx_restrict       nlist,
62                      rvec                        * gmx_restrict          xx,
63                      rvec                        * gmx_restrict          ff,
64                      t_forcerec                  * gmx_restrict          fr,
65                      t_mdatoms                   * gmx_restrict     mdatoms,
66                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67                      t_nrnb                      * gmx_restrict        nrnb)
68 {
69     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
70      * just 0 for non-waters.
71      * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
72      * jnr indices corresponding to data put in the four positions in the SIMD register.
73      */
74     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
75     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
76     int              jnrA,jnrB,jnrC,jnrD;
77     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
78     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
79     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
80     real             rcutoff_scalar;
81     real             *shiftvec,*fshift,*x,*f;
82     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
83     real             scratch[4*DIM];
84     __m128           fscal,rcutoff,rcutoff2,jidxall;
85     int              vdwioffset0;
86     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
87     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
88     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
89     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
90     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
91     real             *charge;
92     int              nvdwtype;
93     __m128           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
94     int              *vdwtype;
95     real             *vdwparam;
96     __m128           one_sixth   = _mm_set1_ps(1.0/6.0);
97     __m128           one_twelfth = _mm_set1_ps(1.0/12.0);
98     __m128           c6grid_00;
99     real             *vdwgridparam;
100     __m128           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
101     __m128           one_half = _mm_set1_ps(0.5);
102     __m128           minus_one = _mm_set1_ps(-1.0);
103     __m128i          ewitab;
104     __m128           ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
105     __m128           beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
106     real             *ewtab;
107     __m128           dummy_mask,cutoff_mask;
108     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
109     __m128           one     = _mm_set1_ps(1.0);
110     __m128           two     = _mm_set1_ps(2.0);
111     x                = xx[0];
112     f                = ff[0];
113
114     nri              = nlist->nri;
115     iinr             = nlist->iinr;
116     jindex           = nlist->jindex;
117     jjnr             = nlist->jjnr;
118     shiftidx         = nlist->shift;
119     gid              = nlist->gid;
120     shiftvec         = fr->shift_vec[0];
121     fshift           = fr->fshift[0];
122     facel            = _mm_set1_ps(fr->epsfac);
123     charge           = mdatoms->chargeA;
124     nvdwtype         = fr->ntype;
125     vdwparam         = fr->nbfp;
126     vdwtype          = mdatoms->typeA;
127     vdwgridparam     = fr->ljpme_c6grid;
128     sh_lj_ewald      = _mm_set1_ps(fr->ic->sh_lj_ewald);
129     ewclj            = _mm_set1_ps(fr->ewaldcoeff_lj);
130     ewclj2           = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
131
132     sh_ewald         = _mm_set1_ps(fr->ic->sh_ewald);
133     beta             = _mm_set1_ps(fr->ic->ewaldcoeff_q);
134     beta2            = _mm_mul_ps(beta,beta);
135     beta3            = _mm_mul_ps(beta,beta2);
136     ewtab            = fr->ic->tabq_coul_FDV0;
137     ewtabscale       = _mm_set1_ps(fr->ic->tabq_scale);
138     ewtabhalfspace   = _mm_set1_ps(0.5/fr->ic->tabq_scale);
139
140     /* Avoid stupid compiler warnings */
141     jnrA = jnrB = jnrC = jnrD = 0;
142     j_coord_offsetA = 0;
143     j_coord_offsetB = 0;
144     j_coord_offsetC = 0;
145     j_coord_offsetD = 0;
146
147     outeriter        = 0;
148     inneriter        = 0;
149
150     for(iidx=0;iidx<4*DIM;iidx++)
151     {
152         scratch[iidx] = 0.0;
153     }
154
155     /* Start outer loop over neighborlists */
156     for(iidx=0; iidx<nri; iidx++)
157     {
158         /* Load shift vector for this list */
159         i_shift_offset   = DIM*shiftidx[iidx];
160
161         /* Load limits for loop over neighbors */
162         j_index_start    = jindex[iidx];
163         j_index_end      = jindex[iidx+1];
164
165         /* Get outer coordinate index */
166         inr              = iinr[iidx];
167         i_coord_offset   = DIM*inr;
168
169         /* Load i particle coords and add shift vector */
170         gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
171
172         fix0             = _mm_setzero_ps();
173         fiy0             = _mm_setzero_ps();
174         fiz0             = _mm_setzero_ps();
175
176         /* Load parameters for i particles */
177         iq0              = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
178         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
179
180         /* Reset potential sums */
181         velecsum         = _mm_setzero_ps();
182         vvdwsum          = _mm_setzero_ps();
183
184         /* Start inner kernel loop */
185         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
186         {
187
188             /* Get j neighbor index, and coordinate index */
189             jnrA             = jjnr[jidx];
190             jnrB             = jjnr[jidx+1];
191             jnrC             = jjnr[jidx+2];
192             jnrD             = jjnr[jidx+3];
193             j_coord_offsetA  = DIM*jnrA;
194             j_coord_offsetB  = DIM*jnrB;
195             j_coord_offsetC  = DIM*jnrC;
196             j_coord_offsetD  = DIM*jnrD;
197
198             /* load j atom coordinates */
199             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
200                                               x+j_coord_offsetC,x+j_coord_offsetD,
201                                               &jx0,&jy0,&jz0);
202
203             /* Calculate displacement vector */
204             dx00             = _mm_sub_ps(ix0,jx0);
205             dy00             = _mm_sub_ps(iy0,jy0);
206             dz00             = _mm_sub_ps(iz0,jz0);
207
208             /* Calculate squared distance and things based on it */
209             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
210
211             rinv00           = gmx_mm_invsqrt_ps(rsq00);
212
213             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
214
215             /* Load parameters for j particles */
216             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
217                                                               charge+jnrC+0,charge+jnrD+0);
218             vdwjidx0A        = 2*vdwtype[jnrA+0];
219             vdwjidx0B        = 2*vdwtype[jnrB+0];
220             vdwjidx0C        = 2*vdwtype[jnrC+0];
221             vdwjidx0D        = 2*vdwtype[jnrD+0];
222
223             /**************************
224              * CALCULATE INTERACTIONS *
225              **************************/
226
227             r00              = _mm_mul_ps(rsq00,rinv00);
228
229             /* Compute parameters for interactions between i and j atoms */
230             qq00             = _mm_mul_ps(iq0,jq0);
231             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
232                                          vdwparam+vdwioffset0+vdwjidx0B,
233                                          vdwparam+vdwioffset0+vdwjidx0C,
234                                          vdwparam+vdwioffset0+vdwjidx0D,
235                                          &c6_00,&c12_00);
236
237             c6grid_00       = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
238                                                                vdwgridparam+vdwioffset0+vdwjidx0B,
239                                                                vdwgridparam+vdwioffset0+vdwjidx0C,
240                                                                vdwgridparam+vdwioffset0+vdwjidx0D);
241
242             /* EWALD ELECTROSTATICS */
243
244             /* Analytical PME correction */
245             zeta2            = _mm_mul_ps(beta2,rsq00);
246             rinv3            = _mm_mul_ps(rinvsq00,rinv00);
247             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
248             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
249             felec            = _mm_mul_ps(qq00,felec);
250             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
251             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv00);
252             velec            = _mm_mul_ps(qq00,velec);
253
254             /* Analytical LJ-PME */
255             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
256             ewcljrsq         = _mm_mul_ps(ewclj2,rsq00);
257             ewclj6           = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
258             exponent         = gmx_simd_exp_r(ewcljrsq);
259             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
260             poly             = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
261             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
262             vvdw6            = _mm_mul_ps(_mm_macc_ps(-c6grid_00,_mm_sub_ps(one,poly),c6_00),rinvsix);
263             vvdw12           = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
264             vvdw             = _mm_msub_ps(vvdw12,one_twelfth,_mm_mul_ps(vvdw6,one_sixth));
265             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
266             fvdw             = _mm_mul_ps(_mm_add_ps(vvdw12,_mm_msub_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6),vvdw6)),rinvsq00);
267
268             /* Update potential sum for this i atom from the interaction with this j atom. */
269             velecsum         = _mm_add_ps(velecsum,velec);
270             vvdwsum          = _mm_add_ps(vvdwsum,vvdw);
271
272             fscal            = _mm_add_ps(felec,fvdw);
273
274              /* Update vectorial force */
275             fix0             = _mm_macc_ps(dx00,fscal,fix0);
276             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
277             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
278
279             fjptrA             = f+j_coord_offsetA;
280             fjptrB             = f+j_coord_offsetB;
281             fjptrC             = f+j_coord_offsetC;
282             fjptrD             = f+j_coord_offsetD;
283             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
284                                                    _mm_mul_ps(dx00,fscal),
285                                                    _mm_mul_ps(dy00,fscal),
286                                                    _mm_mul_ps(dz00,fscal));
287
288             /* Inner loop uses 53 flops */
289         }
290
291         if(jidx<j_index_end)
292         {
293
294             /* Get j neighbor index, and coordinate index */
295             jnrlistA         = jjnr[jidx];
296             jnrlistB         = jjnr[jidx+1];
297             jnrlistC         = jjnr[jidx+2];
298             jnrlistD         = jjnr[jidx+3];
299             /* Sign of each element will be negative for non-real atoms.
300              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
301              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
302              */
303             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
304             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
305             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
306             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
307             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
308             j_coord_offsetA  = DIM*jnrA;
309             j_coord_offsetB  = DIM*jnrB;
310             j_coord_offsetC  = DIM*jnrC;
311             j_coord_offsetD  = DIM*jnrD;
312
313             /* load j atom coordinates */
314             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
315                                               x+j_coord_offsetC,x+j_coord_offsetD,
316                                               &jx0,&jy0,&jz0);
317
318             /* Calculate displacement vector */
319             dx00             = _mm_sub_ps(ix0,jx0);
320             dy00             = _mm_sub_ps(iy0,jy0);
321             dz00             = _mm_sub_ps(iz0,jz0);
322
323             /* Calculate squared distance and things based on it */
324             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
325
326             rinv00           = gmx_mm_invsqrt_ps(rsq00);
327
328             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
329
330             /* Load parameters for j particles */
331             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
332                                                               charge+jnrC+0,charge+jnrD+0);
333             vdwjidx0A        = 2*vdwtype[jnrA+0];
334             vdwjidx0B        = 2*vdwtype[jnrB+0];
335             vdwjidx0C        = 2*vdwtype[jnrC+0];
336             vdwjidx0D        = 2*vdwtype[jnrD+0];
337
338             /**************************
339              * CALCULATE INTERACTIONS *
340              **************************/
341
342             r00              = _mm_mul_ps(rsq00,rinv00);
343             r00              = _mm_andnot_ps(dummy_mask,r00);
344
345             /* Compute parameters for interactions between i and j atoms */
346             qq00             = _mm_mul_ps(iq0,jq0);
347             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
348                                          vdwparam+vdwioffset0+vdwjidx0B,
349                                          vdwparam+vdwioffset0+vdwjidx0C,
350                                          vdwparam+vdwioffset0+vdwjidx0D,
351                                          &c6_00,&c12_00);
352
353             c6grid_00       = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
354                                                                vdwgridparam+vdwioffset0+vdwjidx0B,
355                                                                vdwgridparam+vdwioffset0+vdwjidx0C,
356                                                                vdwgridparam+vdwioffset0+vdwjidx0D);
357
358             /* EWALD ELECTROSTATICS */
359
360             /* Analytical PME correction */
361             zeta2            = _mm_mul_ps(beta2,rsq00);
362             rinv3            = _mm_mul_ps(rinvsq00,rinv00);
363             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
364             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
365             felec            = _mm_mul_ps(qq00,felec);
366             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
367             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv00);
368             velec            = _mm_mul_ps(qq00,velec);
369
370             /* Analytical LJ-PME */
371             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
372             ewcljrsq         = _mm_mul_ps(ewclj2,rsq00);
373             ewclj6           = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
374             exponent         = gmx_simd_exp_r(ewcljrsq);
375             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
376             poly             = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
377             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
378             vvdw6            = _mm_mul_ps(_mm_macc_ps(-c6grid_00,_mm_sub_ps(one,poly),c6_00),rinvsix);
379             vvdw12           = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
380             vvdw             = _mm_msub_ps(vvdw12,one_twelfth,_mm_mul_ps(vvdw6,one_sixth));
381             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
382             fvdw             = _mm_mul_ps(_mm_add_ps(vvdw12,_mm_msub_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6),vvdw6)),rinvsq00);
383
384             /* Update potential sum for this i atom from the interaction with this j atom. */
385             velec            = _mm_andnot_ps(dummy_mask,velec);
386             velecsum         = _mm_add_ps(velecsum,velec);
387             vvdw             = _mm_andnot_ps(dummy_mask,vvdw);
388             vvdwsum          = _mm_add_ps(vvdwsum,vvdw);
389
390             fscal            = _mm_add_ps(felec,fvdw);
391
392             fscal            = _mm_andnot_ps(dummy_mask,fscal);
393
394              /* Update vectorial force */
395             fix0             = _mm_macc_ps(dx00,fscal,fix0);
396             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
397             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
398
399             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
400             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
401             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
402             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
403             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
404                                                    _mm_mul_ps(dx00,fscal),
405                                                    _mm_mul_ps(dy00,fscal),
406                                                    _mm_mul_ps(dz00,fscal));
407
408             /* Inner loop uses 54 flops */
409         }
410
411         /* End of innermost loop */
412
413         gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
414                                               f+i_coord_offset,fshift+i_shift_offset);
415
416         ggid                        = gid[iidx];
417         /* Update potential energies */
418         gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
419         gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
420
421         /* Increment number of inner iterations */
422         inneriter                  += j_index_end - j_index_start;
423
424         /* Outer loop uses 9 flops */
425     }
426
427     /* Increment number of outer iterations */
428     outeriter        += nri;
429
430     /* Update outer/inner flops */
431
432     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*54);
433 }
434 /*
435  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwLJEw_GeomP1P1_F_avx_128_fma_single
436  * Electrostatics interaction: Ewald
437  * VdW interaction:            LJEwald
438  * Geometry:                   Particle-Particle
439  * Calculate force/pot:        Force
440  */
441 void
442 nb_kernel_ElecEw_VdwLJEw_GeomP1P1_F_avx_128_fma_single
443                     (t_nblist                    * gmx_restrict       nlist,
444                      rvec                        * gmx_restrict          xx,
445                      rvec                        * gmx_restrict          ff,
446                      t_forcerec                  * gmx_restrict          fr,
447                      t_mdatoms                   * gmx_restrict     mdatoms,
448                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
449                      t_nrnb                      * gmx_restrict        nrnb)
450 {
451     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
452      * just 0 for non-waters.
453      * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
454      * jnr indices corresponding to data put in the four positions in the SIMD register.
455      */
456     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
457     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
458     int              jnrA,jnrB,jnrC,jnrD;
459     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
460     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
461     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
462     real             rcutoff_scalar;
463     real             *shiftvec,*fshift,*x,*f;
464     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
465     real             scratch[4*DIM];
466     __m128           fscal,rcutoff,rcutoff2,jidxall;
467     int              vdwioffset0;
468     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
469     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
470     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
471     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
472     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
473     real             *charge;
474     int              nvdwtype;
475     __m128           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
476     int              *vdwtype;
477     real             *vdwparam;
478     __m128           one_sixth   = _mm_set1_ps(1.0/6.0);
479     __m128           one_twelfth = _mm_set1_ps(1.0/12.0);
480     __m128           c6grid_00;
481     real             *vdwgridparam;
482     __m128           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
483     __m128           one_half = _mm_set1_ps(0.5);
484     __m128           minus_one = _mm_set1_ps(-1.0);
485     __m128i          ewitab;
486     __m128           ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
487     __m128           beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
488     real             *ewtab;
489     __m128           dummy_mask,cutoff_mask;
490     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
491     __m128           one     = _mm_set1_ps(1.0);
492     __m128           two     = _mm_set1_ps(2.0);
493     x                = xx[0];
494     f                = ff[0];
495
496     nri              = nlist->nri;
497     iinr             = nlist->iinr;
498     jindex           = nlist->jindex;
499     jjnr             = nlist->jjnr;
500     shiftidx         = nlist->shift;
501     gid              = nlist->gid;
502     shiftvec         = fr->shift_vec[0];
503     fshift           = fr->fshift[0];
504     facel            = _mm_set1_ps(fr->epsfac);
505     charge           = mdatoms->chargeA;
506     nvdwtype         = fr->ntype;
507     vdwparam         = fr->nbfp;
508     vdwtype          = mdatoms->typeA;
509     vdwgridparam     = fr->ljpme_c6grid;
510     sh_lj_ewald      = _mm_set1_ps(fr->ic->sh_lj_ewald);
511     ewclj            = _mm_set1_ps(fr->ewaldcoeff_lj);
512     ewclj2           = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
513
514     sh_ewald         = _mm_set1_ps(fr->ic->sh_ewald);
515     beta             = _mm_set1_ps(fr->ic->ewaldcoeff_q);
516     beta2            = _mm_mul_ps(beta,beta);
517     beta3            = _mm_mul_ps(beta,beta2);
518     ewtab            = fr->ic->tabq_coul_F;
519     ewtabscale       = _mm_set1_ps(fr->ic->tabq_scale);
520     ewtabhalfspace   = _mm_set1_ps(0.5/fr->ic->tabq_scale);
521
522     /* Avoid stupid compiler warnings */
523     jnrA = jnrB = jnrC = jnrD = 0;
524     j_coord_offsetA = 0;
525     j_coord_offsetB = 0;
526     j_coord_offsetC = 0;
527     j_coord_offsetD = 0;
528
529     outeriter        = 0;
530     inneriter        = 0;
531
532     for(iidx=0;iidx<4*DIM;iidx++)
533     {
534         scratch[iidx] = 0.0;
535     }
536
537     /* Start outer loop over neighborlists */
538     for(iidx=0; iidx<nri; iidx++)
539     {
540         /* Load shift vector for this list */
541         i_shift_offset   = DIM*shiftidx[iidx];
542
543         /* Load limits for loop over neighbors */
544         j_index_start    = jindex[iidx];
545         j_index_end      = jindex[iidx+1];
546
547         /* Get outer coordinate index */
548         inr              = iinr[iidx];
549         i_coord_offset   = DIM*inr;
550
551         /* Load i particle coords and add shift vector */
552         gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
553
554         fix0             = _mm_setzero_ps();
555         fiy0             = _mm_setzero_ps();
556         fiz0             = _mm_setzero_ps();
557
558         /* Load parameters for i particles */
559         iq0              = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+0));
560         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
561
562         /* Start inner kernel loop */
563         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
564         {
565
566             /* Get j neighbor index, and coordinate index */
567             jnrA             = jjnr[jidx];
568             jnrB             = jjnr[jidx+1];
569             jnrC             = jjnr[jidx+2];
570             jnrD             = jjnr[jidx+3];
571             j_coord_offsetA  = DIM*jnrA;
572             j_coord_offsetB  = DIM*jnrB;
573             j_coord_offsetC  = DIM*jnrC;
574             j_coord_offsetD  = DIM*jnrD;
575
576             /* load j atom coordinates */
577             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
578                                               x+j_coord_offsetC,x+j_coord_offsetD,
579                                               &jx0,&jy0,&jz0);
580
581             /* Calculate displacement vector */
582             dx00             = _mm_sub_ps(ix0,jx0);
583             dy00             = _mm_sub_ps(iy0,jy0);
584             dz00             = _mm_sub_ps(iz0,jz0);
585
586             /* Calculate squared distance and things based on it */
587             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
588
589             rinv00           = gmx_mm_invsqrt_ps(rsq00);
590
591             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
592
593             /* Load parameters for j particles */
594             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
595                                                               charge+jnrC+0,charge+jnrD+0);
596             vdwjidx0A        = 2*vdwtype[jnrA+0];
597             vdwjidx0B        = 2*vdwtype[jnrB+0];
598             vdwjidx0C        = 2*vdwtype[jnrC+0];
599             vdwjidx0D        = 2*vdwtype[jnrD+0];
600
601             /**************************
602              * CALCULATE INTERACTIONS *
603              **************************/
604
605             r00              = _mm_mul_ps(rsq00,rinv00);
606
607             /* Compute parameters for interactions between i and j atoms */
608             qq00             = _mm_mul_ps(iq0,jq0);
609             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
610                                          vdwparam+vdwioffset0+vdwjidx0B,
611                                          vdwparam+vdwioffset0+vdwjidx0C,
612                                          vdwparam+vdwioffset0+vdwjidx0D,
613                                          &c6_00,&c12_00);
614
615             c6grid_00       = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
616                                                                vdwgridparam+vdwioffset0+vdwjidx0B,
617                                                                vdwgridparam+vdwioffset0+vdwjidx0C,
618                                                                vdwgridparam+vdwioffset0+vdwjidx0D);
619
620             /* EWALD ELECTROSTATICS */
621
622             /* Analytical PME correction */
623             zeta2            = _mm_mul_ps(beta2,rsq00);
624             rinv3            = _mm_mul_ps(rinvsq00,rinv00);
625             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
626             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
627             felec            = _mm_mul_ps(qq00,felec);
628
629             /* Analytical LJ-PME */
630             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
631             ewcljrsq         = _mm_mul_ps(ewclj2,rsq00);
632             ewclj6           = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
633             exponent         = gmx_simd_exp_r(ewcljrsq);
634             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
635             poly             = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
636             /* f6A = 6 * C6grid * (1 - poly) */
637             f6A              = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
638             /* f6B = C6grid * exponent * beta^6 */
639             f6B              = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
640             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
641             fvdw              = _mm_mul_ps(_mm_macc_ps(_mm_msub_ps(c12_00,rinvsix,_mm_sub_ps(c6_00,f6A)),rinvsix,f6B),rinvsq00);
642
643             fscal            = _mm_add_ps(felec,fvdw);
644
645              /* Update vectorial force */
646             fix0             = _mm_macc_ps(dx00,fscal,fix0);
647             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
648             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
649
650             fjptrA             = f+j_coord_offsetA;
651             fjptrB             = f+j_coord_offsetB;
652             fjptrC             = f+j_coord_offsetC;
653             fjptrD             = f+j_coord_offsetD;
654             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
655                                                    _mm_mul_ps(dx00,fscal),
656                                                    _mm_mul_ps(dy00,fscal),
657                                                    _mm_mul_ps(dz00,fscal));
658
659             /* Inner loop uses 49 flops */
660         }
661
662         if(jidx<j_index_end)
663         {
664
665             /* Get j neighbor index, and coordinate index */
666             jnrlistA         = jjnr[jidx];
667             jnrlistB         = jjnr[jidx+1];
668             jnrlistC         = jjnr[jidx+2];
669             jnrlistD         = jjnr[jidx+3];
670             /* Sign of each element will be negative for non-real atoms.
671              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
672              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
673              */
674             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
675             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
676             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
677             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
678             jnrD       = (jnrlistD>=0) ? jnrlistD : 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
684             /* load j atom coordinates */
685             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
686                                               x+j_coord_offsetC,x+j_coord_offsetD,
687                                               &jx0,&jy0,&jz0);
688
689             /* Calculate displacement vector */
690             dx00             = _mm_sub_ps(ix0,jx0);
691             dy00             = _mm_sub_ps(iy0,jy0);
692             dz00             = _mm_sub_ps(iz0,jz0);
693
694             /* Calculate squared distance and things based on it */
695             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
696
697             rinv00           = gmx_mm_invsqrt_ps(rsq00);
698
699             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
700
701             /* Load parameters for j particles */
702             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
703                                                               charge+jnrC+0,charge+jnrD+0);
704             vdwjidx0A        = 2*vdwtype[jnrA+0];
705             vdwjidx0B        = 2*vdwtype[jnrB+0];
706             vdwjidx0C        = 2*vdwtype[jnrC+0];
707             vdwjidx0D        = 2*vdwtype[jnrD+0];
708
709             /**************************
710              * CALCULATE INTERACTIONS *
711              **************************/
712
713             r00              = _mm_mul_ps(rsq00,rinv00);
714             r00              = _mm_andnot_ps(dummy_mask,r00);
715
716             /* Compute parameters for interactions between i and j atoms */
717             qq00             = _mm_mul_ps(iq0,jq0);
718             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
719                                          vdwparam+vdwioffset0+vdwjidx0B,
720                                          vdwparam+vdwioffset0+vdwjidx0C,
721                                          vdwparam+vdwioffset0+vdwjidx0D,
722                                          &c6_00,&c12_00);
723
724             c6grid_00       = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset0+vdwjidx0A,
725                                                                vdwgridparam+vdwioffset0+vdwjidx0B,
726                                                                vdwgridparam+vdwioffset0+vdwjidx0C,
727                                                                vdwgridparam+vdwioffset0+vdwjidx0D);
728
729             /* EWALD ELECTROSTATICS */
730
731             /* Analytical PME correction */
732             zeta2            = _mm_mul_ps(beta2,rsq00);
733             rinv3            = _mm_mul_ps(rinvsq00,rinv00);
734             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
735             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
736             felec            = _mm_mul_ps(qq00,felec);
737
738             /* Analytical LJ-PME */
739             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
740             ewcljrsq         = _mm_mul_ps(ewclj2,rsq00);
741             ewclj6           = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
742             exponent         = gmx_simd_exp_r(ewcljrsq);
743             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
744             poly             = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
745             /* f6A = 6 * C6grid * (1 - poly) */
746             f6A              = _mm_mul_ps(c6grid_00,_mm_sub_ps(one,poly));
747             /* f6B = C6grid * exponent * beta^6 */
748             f6B              = _mm_mul_ps(_mm_mul_ps(c6grid_00,one_sixth),_mm_mul_ps(exponent,ewclj6));
749             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
750             fvdw              = _mm_mul_ps(_mm_macc_ps(_mm_msub_ps(c12_00,rinvsix,_mm_sub_ps(c6_00,f6A)),rinvsix,f6B),rinvsq00);
751
752             fscal            = _mm_add_ps(felec,fvdw);
753
754             fscal            = _mm_andnot_ps(dummy_mask,fscal);
755
756              /* Update vectorial force */
757             fix0             = _mm_macc_ps(dx00,fscal,fix0);
758             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
759             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
760
761             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
762             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
763             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
764             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
765             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
766                                                    _mm_mul_ps(dx00,fscal),
767                                                    _mm_mul_ps(dy00,fscal),
768                                                    _mm_mul_ps(dz00,fscal));
769
770             /* Inner loop uses 50 flops */
771         }
772
773         /* End of innermost loop */
774
775         gmx_mm_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*50);
790 }