e78c3a328796c1e24dbec782cdf0a7d706975b9a
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sse2_single / nb_kernel_ElecCSTab_VdwCSTab_GeomW4P1_sse2_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 sse2_single kernel generator.
37  */
38 #include "config.h"
39
40 #include <math.h>
41
42 #include "../nb_kernel.h"
43 #include "gromacs/legacyheaders/types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/legacyheaders/nrnb.h"
46
47 #include "gromacs/simd/math_x86_sse2_single.h"
48 #include "kernelutil_x86_sse2_single.h"
49
50 /*
51  * Gromacs nonbonded kernel:   nb_kernel_ElecCSTab_VdwCSTab_GeomW4P1_VF_sse2_single
52  * Electrostatics interaction: CubicSplineTable
53  * VdW interaction:            CubicSplineTable
54  * Geometry:                   Water4-Particle
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecCSTab_VdwCSTab_GeomW4P1_VF_sse2_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              vdwioffset1;
86     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
87     int              vdwioffset2;
88     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
89     int              vdwioffset3;
90     __m128           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
91     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
92     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
93     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
94     __m128           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
95     __m128           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
96     __m128           dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
97     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
98     real             *charge;
99     int              nvdwtype;
100     __m128           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
101     int              *vdwtype;
102     real             *vdwparam;
103     __m128           one_sixth   = _mm_set1_ps(1.0/6.0);
104     __m128           one_twelfth = _mm_set1_ps(1.0/12.0);
105     __m128i          vfitab;
106     __m128i          ifour       = _mm_set1_epi32(4);
107     __m128           rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
108     real             *vftab;
109     __m128           dummy_mask,cutoff_mask;
110     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
111     __m128           one     = _mm_set1_ps(1.0);
112     __m128           two     = _mm_set1_ps(2.0);
113     x                = xx[0];
114     f                = ff[0];
115
116     nri              = nlist->nri;
117     iinr             = nlist->iinr;
118     jindex           = nlist->jindex;
119     jjnr             = nlist->jjnr;
120     shiftidx         = nlist->shift;
121     gid              = nlist->gid;
122     shiftvec         = fr->shift_vec[0];
123     fshift           = fr->fshift[0];
124     facel            = _mm_set1_ps(fr->epsfac);
125     charge           = mdatoms->chargeA;
126     nvdwtype         = fr->ntype;
127     vdwparam         = fr->nbfp;
128     vdwtype          = mdatoms->typeA;
129
130     vftab            = kernel_data->table_elec_vdw->data;
131     vftabscale       = _mm_set1_ps(kernel_data->table_elec_vdw->scale);
132
133     /* Setup water-specific parameters */
134     inr              = nlist->iinr[0];
135     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
136     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
137     iq3              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
138     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
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_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
171                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
172         
173         fix0             = _mm_setzero_ps();
174         fiy0             = _mm_setzero_ps();
175         fiz0             = _mm_setzero_ps();
176         fix1             = _mm_setzero_ps();
177         fiy1             = _mm_setzero_ps();
178         fiz1             = _mm_setzero_ps();
179         fix2             = _mm_setzero_ps();
180         fiy2             = _mm_setzero_ps();
181         fiz2             = _mm_setzero_ps();
182         fix3             = _mm_setzero_ps();
183         fiy3             = _mm_setzero_ps();
184         fiz3             = _mm_setzero_ps();
185
186         /* Reset potential sums */
187         velecsum         = _mm_setzero_ps();
188         vvdwsum          = _mm_setzero_ps();
189
190         /* Start inner kernel loop */
191         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
192         {
193
194             /* Get j neighbor index, and coordinate index */
195             jnrA             = jjnr[jidx];
196             jnrB             = jjnr[jidx+1];
197             jnrC             = jjnr[jidx+2];
198             jnrD             = jjnr[jidx+3];
199             j_coord_offsetA  = DIM*jnrA;
200             j_coord_offsetB  = DIM*jnrB;
201             j_coord_offsetC  = DIM*jnrC;
202             j_coord_offsetD  = DIM*jnrD;
203
204             /* load j atom coordinates */
205             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
206                                               x+j_coord_offsetC,x+j_coord_offsetD,
207                                               &jx0,&jy0,&jz0);
208
209             /* Calculate displacement vector */
210             dx00             = _mm_sub_ps(ix0,jx0);
211             dy00             = _mm_sub_ps(iy0,jy0);
212             dz00             = _mm_sub_ps(iz0,jz0);
213             dx10             = _mm_sub_ps(ix1,jx0);
214             dy10             = _mm_sub_ps(iy1,jy0);
215             dz10             = _mm_sub_ps(iz1,jz0);
216             dx20             = _mm_sub_ps(ix2,jx0);
217             dy20             = _mm_sub_ps(iy2,jy0);
218             dz20             = _mm_sub_ps(iz2,jz0);
219             dx30             = _mm_sub_ps(ix3,jx0);
220             dy30             = _mm_sub_ps(iy3,jy0);
221             dz30             = _mm_sub_ps(iz3,jz0);
222
223             /* Calculate squared distance and things based on it */
224             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
225             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
226             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
227             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
228
229             rinv00           = gmx_mm_invsqrt_ps(rsq00);
230             rinv10           = gmx_mm_invsqrt_ps(rsq10);
231             rinv20           = gmx_mm_invsqrt_ps(rsq20);
232             rinv30           = gmx_mm_invsqrt_ps(rsq30);
233
234             /* Load parameters for j particles */
235             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
236                                                               charge+jnrC+0,charge+jnrD+0);
237             vdwjidx0A        = 2*vdwtype[jnrA+0];
238             vdwjidx0B        = 2*vdwtype[jnrB+0];
239             vdwjidx0C        = 2*vdwtype[jnrC+0];
240             vdwjidx0D        = 2*vdwtype[jnrD+0];
241
242             fjx0             = _mm_setzero_ps();
243             fjy0             = _mm_setzero_ps();
244             fjz0             = _mm_setzero_ps();
245
246             /**************************
247              * CALCULATE INTERACTIONS *
248              **************************/
249
250             r00              = _mm_mul_ps(rsq00,rinv00);
251
252             /* Compute parameters for interactions between i and j atoms */
253             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
254                                          vdwparam+vdwioffset0+vdwjidx0B,
255                                          vdwparam+vdwioffset0+vdwjidx0C,
256                                          vdwparam+vdwioffset0+vdwjidx0D,
257                                          &c6_00,&c12_00);
258
259             /* Calculate table index by multiplying r with table scale and truncate to integer */
260             rt               = _mm_mul_ps(r00,vftabscale);
261             vfitab           = _mm_cvttps_epi32(rt);
262             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
263             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
264
265             /* CUBIC SPLINE TABLE DISPERSION */
266             vfitab           = _mm_add_epi32(vfitab,ifour);
267             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
268             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
269             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
270             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
271             _MM_TRANSPOSE4_PS(Y,F,G,H);
272             Heps             = _mm_mul_ps(vfeps,H);
273             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
274             VV               = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
275             vvdw6            = _mm_mul_ps(c6_00,VV);
276             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
277             fvdw6            = _mm_mul_ps(c6_00,FF);
278
279             /* CUBIC SPLINE TABLE REPULSION */
280             vfitab           = _mm_add_epi32(vfitab,ifour);
281             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
282             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
283             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
284             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
285             _MM_TRANSPOSE4_PS(Y,F,G,H);
286             Heps             = _mm_mul_ps(vfeps,H);
287             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
288             VV               = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
289             vvdw12           = _mm_mul_ps(c12_00,VV);
290             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
291             fvdw12           = _mm_mul_ps(c12_00,FF);
292             vvdw             = _mm_add_ps(vvdw12,vvdw6);
293             fvdw             = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
294
295             /* Update potential sum for this i atom from the interaction with this j atom. */
296             vvdwsum          = _mm_add_ps(vvdwsum,vvdw);
297
298             fscal            = fvdw;
299
300             /* Calculate temporary vectorial force */
301             tx               = _mm_mul_ps(fscal,dx00);
302             ty               = _mm_mul_ps(fscal,dy00);
303             tz               = _mm_mul_ps(fscal,dz00);
304
305             /* Update vectorial force */
306             fix0             = _mm_add_ps(fix0,tx);
307             fiy0             = _mm_add_ps(fiy0,ty);
308             fiz0             = _mm_add_ps(fiz0,tz);
309
310             fjx0             = _mm_add_ps(fjx0,tx);
311             fjy0             = _mm_add_ps(fjy0,ty);
312             fjz0             = _mm_add_ps(fjz0,tz);
313             
314             /**************************
315              * CALCULATE INTERACTIONS *
316              **************************/
317
318             r10              = _mm_mul_ps(rsq10,rinv10);
319
320             /* Compute parameters for interactions between i and j atoms */
321             qq10             = _mm_mul_ps(iq1,jq0);
322
323             /* Calculate table index by multiplying r with table scale and truncate to integer */
324             rt               = _mm_mul_ps(r10,vftabscale);
325             vfitab           = _mm_cvttps_epi32(rt);
326             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
327             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
328
329             /* CUBIC SPLINE TABLE ELECTROSTATICS */
330             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
331             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
332             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
333             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
334             _MM_TRANSPOSE4_PS(Y,F,G,H);
335             Heps             = _mm_mul_ps(vfeps,H);
336             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
337             VV               = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
338             velec            = _mm_mul_ps(qq10,VV);
339             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
340             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq10,FF),_mm_mul_ps(vftabscale,rinv10)));
341
342             /* Update potential sum for this i atom from the interaction with this j atom. */
343             velecsum         = _mm_add_ps(velecsum,velec);
344
345             fscal            = felec;
346
347             /* Calculate temporary vectorial force */
348             tx               = _mm_mul_ps(fscal,dx10);
349             ty               = _mm_mul_ps(fscal,dy10);
350             tz               = _mm_mul_ps(fscal,dz10);
351
352             /* Update vectorial force */
353             fix1             = _mm_add_ps(fix1,tx);
354             fiy1             = _mm_add_ps(fiy1,ty);
355             fiz1             = _mm_add_ps(fiz1,tz);
356
357             fjx0             = _mm_add_ps(fjx0,tx);
358             fjy0             = _mm_add_ps(fjy0,ty);
359             fjz0             = _mm_add_ps(fjz0,tz);
360             
361             /**************************
362              * CALCULATE INTERACTIONS *
363              **************************/
364
365             r20              = _mm_mul_ps(rsq20,rinv20);
366
367             /* Compute parameters for interactions between i and j atoms */
368             qq20             = _mm_mul_ps(iq2,jq0);
369
370             /* Calculate table index by multiplying r with table scale and truncate to integer */
371             rt               = _mm_mul_ps(r20,vftabscale);
372             vfitab           = _mm_cvttps_epi32(rt);
373             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
374             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
375
376             /* CUBIC SPLINE TABLE ELECTROSTATICS */
377             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
378             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
379             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
380             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
381             _MM_TRANSPOSE4_PS(Y,F,G,H);
382             Heps             = _mm_mul_ps(vfeps,H);
383             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
384             VV               = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
385             velec            = _mm_mul_ps(qq20,VV);
386             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
387             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq20,FF),_mm_mul_ps(vftabscale,rinv20)));
388
389             /* Update potential sum for this i atom from the interaction with this j atom. */
390             velecsum         = _mm_add_ps(velecsum,velec);
391
392             fscal            = felec;
393
394             /* Calculate temporary vectorial force */
395             tx               = _mm_mul_ps(fscal,dx20);
396             ty               = _mm_mul_ps(fscal,dy20);
397             tz               = _mm_mul_ps(fscal,dz20);
398
399             /* Update vectorial force */
400             fix2             = _mm_add_ps(fix2,tx);
401             fiy2             = _mm_add_ps(fiy2,ty);
402             fiz2             = _mm_add_ps(fiz2,tz);
403
404             fjx0             = _mm_add_ps(fjx0,tx);
405             fjy0             = _mm_add_ps(fjy0,ty);
406             fjz0             = _mm_add_ps(fjz0,tz);
407             
408             /**************************
409              * CALCULATE INTERACTIONS *
410              **************************/
411
412             r30              = _mm_mul_ps(rsq30,rinv30);
413
414             /* Compute parameters for interactions between i and j atoms */
415             qq30             = _mm_mul_ps(iq3,jq0);
416
417             /* Calculate table index by multiplying r with table scale and truncate to integer */
418             rt               = _mm_mul_ps(r30,vftabscale);
419             vfitab           = _mm_cvttps_epi32(rt);
420             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
421             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
422
423             /* CUBIC SPLINE TABLE ELECTROSTATICS */
424             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
425             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
426             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
427             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
428             _MM_TRANSPOSE4_PS(Y,F,G,H);
429             Heps             = _mm_mul_ps(vfeps,H);
430             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
431             VV               = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
432             velec            = _mm_mul_ps(qq30,VV);
433             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
434             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq30,FF),_mm_mul_ps(vftabscale,rinv30)));
435
436             /* Update potential sum for this i atom from the interaction with this j atom. */
437             velecsum         = _mm_add_ps(velecsum,velec);
438
439             fscal            = felec;
440
441             /* Calculate temporary vectorial force */
442             tx               = _mm_mul_ps(fscal,dx30);
443             ty               = _mm_mul_ps(fscal,dy30);
444             tz               = _mm_mul_ps(fscal,dz30);
445
446             /* Update vectorial force */
447             fix3             = _mm_add_ps(fix3,tx);
448             fiy3             = _mm_add_ps(fiy3,ty);
449             fiz3             = _mm_add_ps(fiz3,tz);
450
451             fjx0             = _mm_add_ps(fjx0,tx);
452             fjy0             = _mm_add_ps(fjy0,ty);
453             fjz0             = _mm_add_ps(fjz0,tz);
454             
455             fjptrA             = f+j_coord_offsetA;
456             fjptrB             = f+j_coord_offsetB;
457             fjptrC             = f+j_coord_offsetC;
458             fjptrD             = f+j_coord_offsetD;
459
460             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
461
462             /* Inner loop uses 185 flops */
463         }
464
465         if(jidx<j_index_end)
466         {
467
468             /* Get j neighbor index, and coordinate index */
469             jnrlistA         = jjnr[jidx];
470             jnrlistB         = jjnr[jidx+1];
471             jnrlistC         = jjnr[jidx+2];
472             jnrlistD         = jjnr[jidx+3];
473             /* Sign of each element will be negative for non-real atoms.
474              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
475              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
476              */
477             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
478             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
479             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
480             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
481             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
482             j_coord_offsetA  = DIM*jnrA;
483             j_coord_offsetB  = DIM*jnrB;
484             j_coord_offsetC  = DIM*jnrC;
485             j_coord_offsetD  = DIM*jnrD;
486
487             /* load j atom coordinates */
488             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
489                                               x+j_coord_offsetC,x+j_coord_offsetD,
490                                               &jx0,&jy0,&jz0);
491
492             /* Calculate displacement vector */
493             dx00             = _mm_sub_ps(ix0,jx0);
494             dy00             = _mm_sub_ps(iy0,jy0);
495             dz00             = _mm_sub_ps(iz0,jz0);
496             dx10             = _mm_sub_ps(ix1,jx0);
497             dy10             = _mm_sub_ps(iy1,jy0);
498             dz10             = _mm_sub_ps(iz1,jz0);
499             dx20             = _mm_sub_ps(ix2,jx0);
500             dy20             = _mm_sub_ps(iy2,jy0);
501             dz20             = _mm_sub_ps(iz2,jz0);
502             dx30             = _mm_sub_ps(ix3,jx0);
503             dy30             = _mm_sub_ps(iy3,jy0);
504             dz30             = _mm_sub_ps(iz3,jz0);
505
506             /* Calculate squared distance and things based on it */
507             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
508             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
509             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
510             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
511
512             rinv00           = gmx_mm_invsqrt_ps(rsq00);
513             rinv10           = gmx_mm_invsqrt_ps(rsq10);
514             rinv20           = gmx_mm_invsqrt_ps(rsq20);
515             rinv30           = gmx_mm_invsqrt_ps(rsq30);
516
517             /* Load parameters for j particles */
518             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
519                                                               charge+jnrC+0,charge+jnrD+0);
520             vdwjidx0A        = 2*vdwtype[jnrA+0];
521             vdwjidx0B        = 2*vdwtype[jnrB+0];
522             vdwjidx0C        = 2*vdwtype[jnrC+0];
523             vdwjidx0D        = 2*vdwtype[jnrD+0];
524
525             fjx0             = _mm_setzero_ps();
526             fjy0             = _mm_setzero_ps();
527             fjz0             = _mm_setzero_ps();
528
529             /**************************
530              * CALCULATE INTERACTIONS *
531              **************************/
532
533             r00              = _mm_mul_ps(rsq00,rinv00);
534             r00              = _mm_andnot_ps(dummy_mask,r00);
535
536             /* Compute parameters for interactions between i and j atoms */
537             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
538                                          vdwparam+vdwioffset0+vdwjidx0B,
539                                          vdwparam+vdwioffset0+vdwjidx0C,
540                                          vdwparam+vdwioffset0+vdwjidx0D,
541                                          &c6_00,&c12_00);
542
543             /* Calculate table index by multiplying r with table scale and truncate to integer */
544             rt               = _mm_mul_ps(r00,vftabscale);
545             vfitab           = _mm_cvttps_epi32(rt);
546             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
547             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
548
549             /* CUBIC SPLINE TABLE DISPERSION */
550             vfitab           = _mm_add_epi32(vfitab,ifour);
551             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
552             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
553             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
554             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
555             _MM_TRANSPOSE4_PS(Y,F,G,H);
556             Heps             = _mm_mul_ps(vfeps,H);
557             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
558             VV               = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
559             vvdw6            = _mm_mul_ps(c6_00,VV);
560             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
561             fvdw6            = _mm_mul_ps(c6_00,FF);
562
563             /* CUBIC SPLINE TABLE REPULSION */
564             vfitab           = _mm_add_epi32(vfitab,ifour);
565             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
566             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
567             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
568             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
569             _MM_TRANSPOSE4_PS(Y,F,G,H);
570             Heps             = _mm_mul_ps(vfeps,H);
571             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
572             VV               = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
573             vvdw12           = _mm_mul_ps(c12_00,VV);
574             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
575             fvdw12           = _mm_mul_ps(c12_00,FF);
576             vvdw             = _mm_add_ps(vvdw12,vvdw6);
577             fvdw             = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
578
579             /* Update potential sum for this i atom from the interaction with this j atom. */
580             vvdw             = _mm_andnot_ps(dummy_mask,vvdw);
581             vvdwsum          = _mm_add_ps(vvdwsum,vvdw);
582
583             fscal            = fvdw;
584
585             fscal            = _mm_andnot_ps(dummy_mask,fscal);
586
587             /* Calculate temporary vectorial force */
588             tx               = _mm_mul_ps(fscal,dx00);
589             ty               = _mm_mul_ps(fscal,dy00);
590             tz               = _mm_mul_ps(fscal,dz00);
591
592             /* Update vectorial force */
593             fix0             = _mm_add_ps(fix0,tx);
594             fiy0             = _mm_add_ps(fiy0,ty);
595             fiz0             = _mm_add_ps(fiz0,tz);
596
597             fjx0             = _mm_add_ps(fjx0,tx);
598             fjy0             = _mm_add_ps(fjy0,ty);
599             fjz0             = _mm_add_ps(fjz0,tz);
600             
601             /**************************
602              * CALCULATE INTERACTIONS *
603              **************************/
604
605             r10              = _mm_mul_ps(rsq10,rinv10);
606             r10              = _mm_andnot_ps(dummy_mask,r10);
607
608             /* Compute parameters for interactions between i and j atoms */
609             qq10             = _mm_mul_ps(iq1,jq0);
610
611             /* Calculate table index by multiplying r with table scale and truncate to integer */
612             rt               = _mm_mul_ps(r10,vftabscale);
613             vfitab           = _mm_cvttps_epi32(rt);
614             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
615             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
616
617             /* CUBIC SPLINE TABLE ELECTROSTATICS */
618             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
619             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
620             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
621             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
622             _MM_TRANSPOSE4_PS(Y,F,G,H);
623             Heps             = _mm_mul_ps(vfeps,H);
624             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
625             VV               = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
626             velec            = _mm_mul_ps(qq10,VV);
627             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
628             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq10,FF),_mm_mul_ps(vftabscale,rinv10)));
629
630             /* Update potential sum for this i atom from the interaction with this j atom. */
631             velec            = _mm_andnot_ps(dummy_mask,velec);
632             velecsum         = _mm_add_ps(velecsum,velec);
633
634             fscal            = felec;
635
636             fscal            = _mm_andnot_ps(dummy_mask,fscal);
637
638             /* Calculate temporary vectorial force */
639             tx               = _mm_mul_ps(fscal,dx10);
640             ty               = _mm_mul_ps(fscal,dy10);
641             tz               = _mm_mul_ps(fscal,dz10);
642
643             /* Update vectorial force */
644             fix1             = _mm_add_ps(fix1,tx);
645             fiy1             = _mm_add_ps(fiy1,ty);
646             fiz1             = _mm_add_ps(fiz1,tz);
647
648             fjx0             = _mm_add_ps(fjx0,tx);
649             fjy0             = _mm_add_ps(fjy0,ty);
650             fjz0             = _mm_add_ps(fjz0,tz);
651             
652             /**************************
653              * CALCULATE INTERACTIONS *
654              **************************/
655
656             r20              = _mm_mul_ps(rsq20,rinv20);
657             r20              = _mm_andnot_ps(dummy_mask,r20);
658
659             /* Compute parameters for interactions between i and j atoms */
660             qq20             = _mm_mul_ps(iq2,jq0);
661
662             /* Calculate table index by multiplying r with table scale and truncate to integer */
663             rt               = _mm_mul_ps(r20,vftabscale);
664             vfitab           = _mm_cvttps_epi32(rt);
665             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
666             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
667
668             /* CUBIC SPLINE TABLE ELECTROSTATICS */
669             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
670             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
671             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
672             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
673             _MM_TRANSPOSE4_PS(Y,F,G,H);
674             Heps             = _mm_mul_ps(vfeps,H);
675             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
676             VV               = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
677             velec            = _mm_mul_ps(qq20,VV);
678             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
679             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq20,FF),_mm_mul_ps(vftabscale,rinv20)));
680
681             /* Update potential sum for this i atom from the interaction with this j atom. */
682             velec            = _mm_andnot_ps(dummy_mask,velec);
683             velecsum         = _mm_add_ps(velecsum,velec);
684
685             fscal            = felec;
686
687             fscal            = _mm_andnot_ps(dummy_mask,fscal);
688
689             /* Calculate temporary vectorial force */
690             tx               = _mm_mul_ps(fscal,dx20);
691             ty               = _mm_mul_ps(fscal,dy20);
692             tz               = _mm_mul_ps(fscal,dz20);
693
694             /* Update vectorial force */
695             fix2             = _mm_add_ps(fix2,tx);
696             fiy2             = _mm_add_ps(fiy2,ty);
697             fiz2             = _mm_add_ps(fiz2,tz);
698
699             fjx0             = _mm_add_ps(fjx0,tx);
700             fjy0             = _mm_add_ps(fjy0,ty);
701             fjz0             = _mm_add_ps(fjz0,tz);
702             
703             /**************************
704              * CALCULATE INTERACTIONS *
705              **************************/
706
707             r30              = _mm_mul_ps(rsq30,rinv30);
708             r30              = _mm_andnot_ps(dummy_mask,r30);
709
710             /* Compute parameters for interactions between i and j atoms */
711             qq30             = _mm_mul_ps(iq3,jq0);
712
713             /* Calculate table index by multiplying r with table scale and truncate to integer */
714             rt               = _mm_mul_ps(r30,vftabscale);
715             vfitab           = _mm_cvttps_epi32(rt);
716             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
717             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
718
719             /* CUBIC SPLINE TABLE ELECTROSTATICS */
720             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
721             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
722             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
723             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
724             _MM_TRANSPOSE4_PS(Y,F,G,H);
725             Heps             = _mm_mul_ps(vfeps,H);
726             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
727             VV               = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
728             velec            = _mm_mul_ps(qq30,VV);
729             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
730             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq30,FF),_mm_mul_ps(vftabscale,rinv30)));
731
732             /* Update potential sum for this i atom from the interaction with this j atom. */
733             velec            = _mm_andnot_ps(dummy_mask,velec);
734             velecsum         = _mm_add_ps(velecsum,velec);
735
736             fscal            = felec;
737
738             fscal            = _mm_andnot_ps(dummy_mask,fscal);
739
740             /* Calculate temporary vectorial force */
741             tx               = _mm_mul_ps(fscal,dx30);
742             ty               = _mm_mul_ps(fscal,dy30);
743             tz               = _mm_mul_ps(fscal,dz30);
744
745             /* Update vectorial force */
746             fix3             = _mm_add_ps(fix3,tx);
747             fiy3             = _mm_add_ps(fiy3,ty);
748             fiz3             = _mm_add_ps(fiz3,tz);
749
750             fjx0             = _mm_add_ps(fjx0,tx);
751             fjy0             = _mm_add_ps(fjy0,ty);
752             fjz0             = _mm_add_ps(fjz0,tz);
753             
754             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
755             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
756             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
757             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
758
759             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
760
761             /* Inner loop uses 189 flops */
762         }
763
764         /* End of innermost loop */
765
766         gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
767                                               f+i_coord_offset,fshift+i_shift_offset);
768
769         ggid                        = gid[iidx];
770         /* Update potential energies */
771         gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
772         gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
773
774         /* Increment number of inner iterations */
775         inneriter                  += j_index_end - j_index_start;
776
777         /* Outer loop uses 26 flops */
778     }
779
780     /* Increment number of outer iterations */
781     outeriter        += nri;
782
783     /* Update outer/inner flops */
784
785     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*189);
786 }
787 /*
788  * Gromacs nonbonded kernel:   nb_kernel_ElecCSTab_VdwCSTab_GeomW4P1_F_sse2_single
789  * Electrostatics interaction: CubicSplineTable
790  * VdW interaction:            CubicSplineTable
791  * Geometry:                   Water4-Particle
792  * Calculate force/pot:        Force
793  */
794 void
795 nb_kernel_ElecCSTab_VdwCSTab_GeomW4P1_F_sse2_single
796                     (t_nblist                    * gmx_restrict       nlist,
797                      rvec                        * gmx_restrict          xx,
798                      rvec                        * gmx_restrict          ff,
799                      t_forcerec                  * gmx_restrict          fr,
800                      t_mdatoms                   * gmx_restrict     mdatoms,
801                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
802                      t_nrnb                      * gmx_restrict        nrnb)
803 {
804     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
805      * just 0 for non-waters.
806      * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
807      * jnr indices corresponding to data put in the four positions in the SIMD register.
808      */
809     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
810     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
811     int              jnrA,jnrB,jnrC,jnrD;
812     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
813     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
814     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
815     real             rcutoff_scalar;
816     real             *shiftvec,*fshift,*x,*f;
817     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
818     real             scratch[4*DIM];
819     __m128           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
820     int              vdwioffset0;
821     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
822     int              vdwioffset1;
823     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
824     int              vdwioffset2;
825     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
826     int              vdwioffset3;
827     __m128           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
828     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
829     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
830     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
831     __m128           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
832     __m128           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
833     __m128           dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
834     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
835     real             *charge;
836     int              nvdwtype;
837     __m128           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
838     int              *vdwtype;
839     real             *vdwparam;
840     __m128           one_sixth   = _mm_set1_ps(1.0/6.0);
841     __m128           one_twelfth = _mm_set1_ps(1.0/12.0);
842     __m128i          vfitab;
843     __m128i          ifour       = _mm_set1_epi32(4);
844     __m128           rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
845     real             *vftab;
846     __m128           dummy_mask,cutoff_mask;
847     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
848     __m128           one     = _mm_set1_ps(1.0);
849     __m128           two     = _mm_set1_ps(2.0);
850     x                = xx[0];
851     f                = ff[0];
852
853     nri              = nlist->nri;
854     iinr             = nlist->iinr;
855     jindex           = nlist->jindex;
856     jjnr             = nlist->jjnr;
857     shiftidx         = nlist->shift;
858     gid              = nlist->gid;
859     shiftvec         = fr->shift_vec[0];
860     fshift           = fr->fshift[0];
861     facel            = _mm_set1_ps(fr->epsfac);
862     charge           = mdatoms->chargeA;
863     nvdwtype         = fr->ntype;
864     vdwparam         = fr->nbfp;
865     vdwtype          = mdatoms->typeA;
866
867     vftab            = kernel_data->table_elec_vdw->data;
868     vftabscale       = _mm_set1_ps(kernel_data->table_elec_vdw->scale);
869
870     /* Setup water-specific parameters */
871     inr              = nlist->iinr[0];
872     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
873     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
874     iq3              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
875     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
876
877     /* Avoid stupid compiler warnings */
878     jnrA = jnrB = jnrC = jnrD = 0;
879     j_coord_offsetA = 0;
880     j_coord_offsetB = 0;
881     j_coord_offsetC = 0;
882     j_coord_offsetD = 0;
883
884     outeriter        = 0;
885     inneriter        = 0;
886
887     for(iidx=0;iidx<4*DIM;iidx++)
888     {
889         scratch[iidx] = 0.0;
890     }  
891
892     /* Start outer loop over neighborlists */
893     for(iidx=0; iidx<nri; iidx++)
894     {
895         /* Load shift vector for this list */
896         i_shift_offset   = DIM*shiftidx[iidx];
897
898         /* Load limits for loop over neighbors */
899         j_index_start    = jindex[iidx];
900         j_index_end      = jindex[iidx+1];
901
902         /* Get outer coordinate index */
903         inr              = iinr[iidx];
904         i_coord_offset   = DIM*inr;
905
906         /* Load i particle coords and add shift vector */
907         gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
908                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
909         
910         fix0             = _mm_setzero_ps();
911         fiy0             = _mm_setzero_ps();
912         fiz0             = _mm_setzero_ps();
913         fix1             = _mm_setzero_ps();
914         fiy1             = _mm_setzero_ps();
915         fiz1             = _mm_setzero_ps();
916         fix2             = _mm_setzero_ps();
917         fiy2             = _mm_setzero_ps();
918         fiz2             = _mm_setzero_ps();
919         fix3             = _mm_setzero_ps();
920         fiy3             = _mm_setzero_ps();
921         fiz3             = _mm_setzero_ps();
922
923         /* Start inner kernel loop */
924         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
925         {
926
927             /* Get j neighbor index, and coordinate index */
928             jnrA             = jjnr[jidx];
929             jnrB             = jjnr[jidx+1];
930             jnrC             = jjnr[jidx+2];
931             jnrD             = jjnr[jidx+3];
932             j_coord_offsetA  = DIM*jnrA;
933             j_coord_offsetB  = DIM*jnrB;
934             j_coord_offsetC  = DIM*jnrC;
935             j_coord_offsetD  = DIM*jnrD;
936
937             /* load j atom coordinates */
938             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
939                                               x+j_coord_offsetC,x+j_coord_offsetD,
940                                               &jx0,&jy0,&jz0);
941
942             /* Calculate displacement vector */
943             dx00             = _mm_sub_ps(ix0,jx0);
944             dy00             = _mm_sub_ps(iy0,jy0);
945             dz00             = _mm_sub_ps(iz0,jz0);
946             dx10             = _mm_sub_ps(ix1,jx0);
947             dy10             = _mm_sub_ps(iy1,jy0);
948             dz10             = _mm_sub_ps(iz1,jz0);
949             dx20             = _mm_sub_ps(ix2,jx0);
950             dy20             = _mm_sub_ps(iy2,jy0);
951             dz20             = _mm_sub_ps(iz2,jz0);
952             dx30             = _mm_sub_ps(ix3,jx0);
953             dy30             = _mm_sub_ps(iy3,jy0);
954             dz30             = _mm_sub_ps(iz3,jz0);
955
956             /* Calculate squared distance and things based on it */
957             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
958             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
959             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
960             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
961
962             rinv00           = gmx_mm_invsqrt_ps(rsq00);
963             rinv10           = gmx_mm_invsqrt_ps(rsq10);
964             rinv20           = gmx_mm_invsqrt_ps(rsq20);
965             rinv30           = gmx_mm_invsqrt_ps(rsq30);
966
967             /* Load parameters for j particles */
968             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
969                                                               charge+jnrC+0,charge+jnrD+0);
970             vdwjidx0A        = 2*vdwtype[jnrA+0];
971             vdwjidx0B        = 2*vdwtype[jnrB+0];
972             vdwjidx0C        = 2*vdwtype[jnrC+0];
973             vdwjidx0D        = 2*vdwtype[jnrD+0];
974
975             fjx0             = _mm_setzero_ps();
976             fjy0             = _mm_setzero_ps();
977             fjz0             = _mm_setzero_ps();
978
979             /**************************
980              * CALCULATE INTERACTIONS *
981              **************************/
982
983             r00              = _mm_mul_ps(rsq00,rinv00);
984
985             /* Compute parameters for interactions between i and j atoms */
986             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
987                                          vdwparam+vdwioffset0+vdwjidx0B,
988                                          vdwparam+vdwioffset0+vdwjidx0C,
989                                          vdwparam+vdwioffset0+vdwjidx0D,
990                                          &c6_00,&c12_00);
991
992             /* Calculate table index by multiplying r with table scale and truncate to integer */
993             rt               = _mm_mul_ps(r00,vftabscale);
994             vfitab           = _mm_cvttps_epi32(rt);
995             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
996             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
997
998             /* CUBIC SPLINE TABLE DISPERSION */
999             vfitab           = _mm_add_epi32(vfitab,ifour);
1000             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
1001             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
1002             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
1003             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
1004             _MM_TRANSPOSE4_PS(Y,F,G,H);
1005             Heps             = _mm_mul_ps(vfeps,H);
1006             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
1007             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
1008             fvdw6            = _mm_mul_ps(c6_00,FF);
1009
1010             /* CUBIC SPLINE TABLE REPULSION */
1011             vfitab           = _mm_add_epi32(vfitab,ifour);
1012             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
1013             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
1014             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
1015             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
1016             _MM_TRANSPOSE4_PS(Y,F,G,H);
1017             Heps             = _mm_mul_ps(vfeps,H);
1018             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
1019             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
1020             fvdw12           = _mm_mul_ps(c12_00,FF);
1021             fvdw             = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
1022
1023             fscal            = fvdw;
1024
1025             /* Calculate temporary vectorial force */
1026             tx               = _mm_mul_ps(fscal,dx00);
1027             ty               = _mm_mul_ps(fscal,dy00);
1028             tz               = _mm_mul_ps(fscal,dz00);
1029
1030             /* Update vectorial force */
1031             fix0             = _mm_add_ps(fix0,tx);
1032             fiy0             = _mm_add_ps(fiy0,ty);
1033             fiz0             = _mm_add_ps(fiz0,tz);
1034
1035             fjx0             = _mm_add_ps(fjx0,tx);
1036             fjy0             = _mm_add_ps(fjy0,ty);
1037             fjz0             = _mm_add_ps(fjz0,tz);
1038             
1039             /**************************
1040              * CALCULATE INTERACTIONS *
1041              **************************/
1042
1043             r10              = _mm_mul_ps(rsq10,rinv10);
1044
1045             /* Compute parameters for interactions between i and j atoms */
1046             qq10             = _mm_mul_ps(iq1,jq0);
1047
1048             /* Calculate table index by multiplying r with table scale and truncate to integer */
1049             rt               = _mm_mul_ps(r10,vftabscale);
1050             vfitab           = _mm_cvttps_epi32(rt);
1051             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
1052             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
1053
1054             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1055             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
1056             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
1057             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
1058             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
1059             _MM_TRANSPOSE4_PS(Y,F,G,H);
1060             Heps             = _mm_mul_ps(vfeps,H);
1061             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
1062             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
1063             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq10,FF),_mm_mul_ps(vftabscale,rinv10)));
1064
1065             fscal            = felec;
1066
1067             /* Calculate temporary vectorial force */
1068             tx               = _mm_mul_ps(fscal,dx10);
1069             ty               = _mm_mul_ps(fscal,dy10);
1070             tz               = _mm_mul_ps(fscal,dz10);
1071
1072             /* Update vectorial force */
1073             fix1             = _mm_add_ps(fix1,tx);
1074             fiy1             = _mm_add_ps(fiy1,ty);
1075             fiz1             = _mm_add_ps(fiz1,tz);
1076
1077             fjx0             = _mm_add_ps(fjx0,tx);
1078             fjy0             = _mm_add_ps(fjy0,ty);
1079             fjz0             = _mm_add_ps(fjz0,tz);
1080             
1081             /**************************
1082              * CALCULATE INTERACTIONS *
1083              **************************/
1084
1085             r20              = _mm_mul_ps(rsq20,rinv20);
1086
1087             /* Compute parameters for interactions between i and j atoms */
1088             qq20             = _mm_mul_ps(iq2,jq0);
1089
1090             /* Calculate table index by multiplying r with table scale and truncate to integer */
1091             rt               = _mm_mul_ps(r20,vftabscale);
1092             vfitab           = _mm_cvttps_epi32(rt);
1093             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
1094             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
1095
1096             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1097             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
1098             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
1099             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
1100             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
1101             _MM_TRANSPOSE4_PS(Y,F,G,H);
1102             Heps             = _mm_mul_ps(vfeps,H);
1103             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
1104             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
1105             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq20,FF),_mm_mul_ps(vftabscale,rinv20)));
1106
1107             fscal            = felec;
1108
1109             /* Calculate temporary vectorial force */
1110             tx               = _mm_mul_ps(fscal,dx20);
1111             ty               = _mm_mul_ps(fscal,dy20);
1112             tz               = _mm_mul_ps(fscal,dz20);
1113
1114             /* Update vectorial force */
1115             fix2             = _mm_add_ps(fix2,tx);
1116             fiy2             = _mm_add_ps(fiy2,ty);
1117             fiz2             = _mm_add_ps(fiz2,tz);
1118
1119             fjx0             = _mm_add_ps(fjx0,tx);
1120             fjy0             = _mm_add_ps(fjy0,ty);
1121             fjz0             = _mm_add_ps(fjz0,tz);
1122             
1123             /**************************
1124              * CALCULATE INTERACTIONS *
1125              **************************/
1126
1127             r30              = _mm_mul_ps(rsq30,rinv30);
1128
1129             /* Compute parameters for interactions between i and j atoms */
1130             qq30             = _mm_mul_ps(iq3,jq0);
1131
1132             /* Calculate table index by multiplying r with table scale and truncate to integer */
1133             rt               = _mm_mul_ps(r30,vftabscale);
1134             vfitab           = _mm_cvttps_epi32(rt);
1135             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
1136             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
1137
1138             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1139             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
1140             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
1141             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
1142             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
1143             _MM_TRANSPOSE4_PS(Y,F,G,H);
1144             Heps             = _mm_mul_ps(vfeps,H);
1145             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
1146             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
1147             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq30,FF),_mm_mul_ps(vftabscale,rinv30)));
1148
1149             fscal            = felec;
1150
1151             /* Calculate temporary vectorial force */
1152             tx               = _mm_mul_ps(fscal,dx30);
1153             ty               = _mm_mul_ps(fscal,dy30);
1154             tz               = _mm_mul_ps(fscal,dz30);
1155
1156             /* Update vectorial force */
1157             fix3             = _mm_add_ps(fix3,tx);
1158             fiy3             = _mm_add_ps(fiy3,ty);
1159             fiz3             = _mm_add_ps(fiz3,tz);
1160
1161             fjx0             = _mm_add_ps(fjx0,tx);
1162             fjy0             = _mm_add_ps(fjy0,ty);
1163             fjz0             = _mm_add_ps(fjz0,tz);
1164             
1165             fjptrA             = f+j_coord_offsetA;
1166             fjptrB             = f+j_coord_offsetB;
1167             fjptrC             = f+j_coord_offsetC;
1168             fjptrD             = f+j_coord_offsetD;
1169
1170             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1171
1172             /* Inner loop uses 165 flops */
1173         }
1174
1175         if(jidx<j_index_end)
1176         {
1177
1178             /* Get j neighbor index, and coordinate index */
1179             jnrlistA         = jjnr[jidx];
1180             jnrlistB         = jjnr[jidx+1];
1181             jnrlistC         = jjnr[jidx+2];
1182             jnrlistD         = jjnr[jidx+3];
1183             /* Sign of each element will be negative for non-real atoms.
1184              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1185              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1186              */
1187             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
1188             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
1189             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
1190             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
1191             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
1192             j_coord_offsetA  = DIM*jnrA;
1193             j_coord_offsetB  = DIM*jnrB;
1194             j_coord_offsetC  = DIM*jnrC;
1195             j_coord_offsetD  = DIM*jnrD;
1196
1197             /* load j atom coordinates */
1198             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1199                                               x+j_coord_offsetC,x+j_coord_offsetD,
1200                                               &jx0,&jy0,&jz0);
1201
1202             /* Calculate displacement vector */
1203             dx00             = _mm_sub_ps(ix0,jx0);
1204             dy00             = _mm_sub_ps(iy0,jy0);
1205             dz00             = _mm_sub_ps(iz0,jz0);
1206             dx10             = _mm_sub_ps(ix1,jx0);
1207             dy10             = _mm_sub_ps(iy1,jy0);
1208             dz10             = _mm_sub_ps(iz1,jz0);
1209             dx20             = _mm_sub_ps(ix2,jx0);
1210             dy20             = _mm_sub_ps(iy2,jy0);
1211             dz20             = _mm_sub_ps(iz2,jz0);
1212             dx30             = _mm_sub_ps(ix3,jx0);
1213             dy30             = _mm_sub_ps(iy3,jy0);
1214             dz30             = _mm_sub_ps(iz3,jz0);
1215
1216             /* Calculate squared distance and things based on it */
1217             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
1218             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
1219             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
1220             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
1221
1222             rinv00           = gmx_mm_invsqrt_ps(rsq00);
1223             rinv10           = gmx_mm_invsqrt_ps(rsq10);
1224             rinv20           = gmx_mm_invsqrt_ps(rsq20);
1225             rinv30           = gmx_mm_invsqrt_ps(rsq30);
1226
1227             /* Load parameters for j particles */
1228             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1229                                                               charge+jnrC+0,charge+jnrD+0);
1230             vdwjidx0A        = 2*vdwtype[jnrA+0];
1231             vdwjidx0B        = 2*vdwtype[jnrB+0];
1232             vdwjidx0C        = 2*vdwtype[jnrC+0];
1233             vdwjidx0D        = 2*vdwtype[jnrD+0];
1234
1235             fjx0             = _mm_setzero_ps();
1236             fjy0             = _mm_setzero_ps();
1237             fjz0             = _mm_setzero_ps();
1238
1239             /**************************
1240              * CALCULATE INTERACTIONS *
1241              **************************/
1242
1243             r00              = _mm_mul_ps(rsq00,rinv00);
1244             r00              = _mm_andnot_ps(dummy_mask,r00);
1245
1246             /* Compute parameters for interactions between i and j atoms */
1247             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
1248                                          vdwparam+vdwioffset0+vdwjidx0B,
1249                                          vdwparam+vdwioffset0+vdwjidx0C,
1250                                          vdwparam+vdwioffset0+vdwjidx0D,
1251                                          &c6_00,&c12_00);
1252
1253             /* Calculate table index by multiplying r with table scale and truncate to integer */
1254             rt               = _mm_mul_ps(r00,vftabscale);
1255             vfitab           = _mm_cvttps_epi32(rt);
1256             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
1257             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
1258
1259             /* CUBIC SPLINE TABLE DISPERSION */
1260             vfitab           = _mm_add_epi32(vfitab,ifour);
1261             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
1262             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
1263             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
1264             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
1265             _MM_TRANSPOSE4_PS(Y,F,G,H);
1266             Heps             = _mm_mul_ps(vfeps,H);
1267             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
1268             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
1269             fvdw6            = _mm_mul_ps(c6_00,FF);
1270
1271             /* CUBIC SPLINE TABLE REPULSION */
1272             vfitab           = _mm_add_epi32(vfitab,ifour);
1273             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
1274             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
1275             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
1276             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
1277             _MM_TRANSPOSE4_PS(Y,F,G,H);
1278             Heps             = _mm_mul_ps(vfeps,H);
1279             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
1280             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
1281             fvdw12           = _mm_mul_ps(c12_00,FF);
1282             fvdw             = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
1283
1284             fscal            = fvdw;
1285
1286             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1287
1288             /* Calculate temporary vectorial force */
1289             tx               = _mm_mul_ps(fscal,dx00);
1290             ty               = _mm_mul_ps(fscal,dy00);
1291             tz               = _mm_mul_ps(fscal,dz00);
1292
1293             /* Update vectorial force */
1294             fix0             = _mm_add_ps(fix0,tx);
1295             fiy0             = _mm_add_ps(fiy0,ty);
1296             fiz0             = _mm_add_ps(fiz0,tz);
1297
1298             fjx0             = _mm_add_ps(fjx0,tx);
1299             fjy0             = _mm_add_ps(fjy0,ty);
1300             fjz0             = _mm_add_ps(fjz0,tz);
1301             
1302             /**************************
1303              * CALCULATE INTERACTIONS *
1304              **************************/
1305
1306             r10              = _mm_mul_ps(rsq10,rinv10);
1307             r10              = _mm_andnot_ps(dummy_mask,r10);
1308
1309             /* Compute parameters for interactions between i and j atoms */
1310             qq10             = _mm_mul_ps(iq1,jq0);
1311
1312             /* Calculate table index by multiplying r with table scale and truncate to integer */
1313             rt               = _mm_mul_ps(r10,vftabscale);
1314             vfitab           = _mm_cvttps_epi32(rt);
1315             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
1316             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
1317
1318             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1319             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
1320             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
1321             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
1322             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
1323             _MM_TRANSPOSE4_PS(Y,F,G,H);
1324             Heps             = _mm_mul_ps(vfeps,H);
1325             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
1326             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
1327             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq10,FF),_mm_mul_ps(vftabscale,rinv10)));
1328
1329             fscal            = felec;
1330
1331             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1332
1333             /* Calculate temporary vectorial force */
1334             tx               = _mm_mul_ps(fscal,dx10);
1335             ty               = _mm_mul_ps(fscal,dy10);
1336             tz               = _mm_mul_ps(fscal,dz10);
1337
1338             /* Update vectorial force */
1339             fix1             = _mm_add_ps(fix1,tx);
1340             fiy1             = _mm_add_ps(fiy1,ty);
1341             fiz1             = _mm_add_ps(fiz1,tz);
1342
1343             fjx0             = _mm_add_ps(fjx0,tx);
1344             fjy0             = _mm_add_ps(fjy0,ty);
1345             fjz0             = _mm_add_ps(fjz0,tz);
1346             
1347             /**************************
1348              * CALCULATE INTERACTIONS *
1349              **************************/
1350
1351             r20              = _mm_mul_ps(rsq20,rinv20);
1352             r20              = _mm_andnot_ps(dummy_mask,r20);
1353
1354             /* Compute parameters for interactions between i and j atoms */
1355             qq20             = _mm_mul_ps(iq2,jq0);
1356
1357             /* Calculate table index by multiplying r with table scale and truncate to integer */
1358             rt               = _mm_mul_ps(r20,vftabscale);
1359             vfitab           = _mm_cvttps_epi32(rt);
1360             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
1361             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
1362
1363             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1364             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
1365             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
1366             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
1367             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
1368             _MM_TRANSPOSE4_PS(Y,F,G,H);
1369             Heps             = _mm_mul_ps(vfeps,H);
1370             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
1371             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
1372             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq20,FF),_mm_mul_ps(vftabscale,rinv20)));
1373
1374             fscal            = felec;
1375
1376             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1377
1378             /* Calculate temporary vectorial force */
1379             tx               = _mm_mul_ps(fscal,dx20);
1380             ty               = _mm_mul_ps(fscal,dy20);
1381             tz               = _mm_mul_ps(fscal,dz20);
1382
1383             /* Update vectorial force */
1384             fix2             = _mm_add_ps(fix2,tx);
1385             fiy2             = _mm_add_ps(fiy2,ty);
1386             fiz2             = _mm_add_ps(fiz2,tz);
1387
1388             fjx0             = _mm_add_ps(fjx0,tx);
1389             fjy0             = _mm_add_ps(fjy0,ty);
1390             fjz0             = _mm_add_ps(fjz0,tz);
1391             
1392             /**************************
1393              * CALCULATE INTERACTIONS *
1394              **************************/
1395
1396             r30              = _mm_mul_ps(rsq30,rinv30);
1397             r30              = _mm_andnot_ps(dummy_mask,r30);
1398
1399             /* Compute parameters for interactions between i and j atoms */
1400             qq30             = _mm_mul_ps(iq3,jq0);
1401
1402             /* Calculate table index by multiplying r with table scale and truncate to integer */
1403             rt               = _mm_mul_ps(r30,vftabscale);
1404             vfitab           = _mm_cvttps_epi32(rt);
1405             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
1406             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
1407
1408             /* CUBIC SPLINE TABLE ELECTROSTATICS */
1409             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
1410             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
1411             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
1412             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
1413             _MM_TRANSPOSE4_PS(Y,F,G,H);
1414             Heps             = _mm_mul_ps(vfeps,H);
1415             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
1416             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
1417             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq30,FF),_mm_mul_ps(vftabscale,rinv30)));
1418
1419             fscal            = felec;
1420
1421             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1422
1423             /* Calculate temporary vectorial force */
1424             tx               = _mm_mul_ps(fscal,dx30);
1425             ty               = _mm_mul_ps(fscal,dy30);
1426             tz               = _mm_mul_ps(fscal,dz30);
1427
1428             /* Update vectorial force */
1429             fix3             = _mm_add_ps(fix3,tx);
1430             fiy3             = _mm_add_ps(fiy3,ty);
1431             fiz3             = _mm_add_ps(fiz3,tz);
1432
1433             fjx0             = _mm_add_ps(fjx0,tx);
1434             fjy0             = _mm_add_ps(fjy0,ty);
1435             fjz0             = _mm_add_ps(fjz0,tz);
1436             
1437             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1438             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1439             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1440             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1441
1442             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1443
1444             /* Inner loop uses 169 flops */
1445         }
1446
1447         /* End of innermost loop */
1448
1449         gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1450                                               f+i_coord_offset,fshift+i_shift_offset);
1451
1452         /* Increment number of inner iterations */
1453         inneriter                  += j_index_end - j_index_start;
1454
1455         /* Outer loop uses 24 flops */
1456     }
1457
1458     /* Increment number of outer iterations */
1459     outeriter        += nri;
1460
1461     /* Update outer/inner flops */
1462
1463     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*169);
1464 }