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