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