ba86d93aa2f9b5c026b60884dfc9ccff074ff16a
[alexxy/gromacs.git] / src / gmxlib / nonbonded / nb_kernel_sse2_single / nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_sse2_single.c
1 /*
2  * Note: this file was generated by the Gromacs sse2_single kernel generator.
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  * Copyright (c) 2001-2012, The GROMACS Development Team
9  *
10  * Gromacs is a library for molecular simulation and trajectory analysis,
11  * written by Erik Lindahl, David van der Spoel, Berk Hess, and others - for
12  * a full list of developers and information, check out http://www.gromacs.org
13  *
14  * This program is free software; you can redistribute it and/or modify it under
15  * the terms of the GNU Lesser General Public License as published by the Free
16  * Software Foundation; either version 2 of the License, or (at your option) any
17  * later version.
18  *
19  * To help fund GROMACS development, we humbly ask that you cite
20  * the papers people have written on it - you can find them on the website.
21  */
22 #ifdef HAVE_CONFIG_H
23 #include <config.h>
24 #endif
25
26 #include <math.h>
27
28 #include "../nb_kernel.h"
29 #include "types/simple.h"
30 #include "vec.h"
31 #include "nrnb.h"
32
33 #include "gmx_math_x86_sse2_single.h"
34 #include "kernelutil_x86_sse2_single.h"
35
36 /*
37  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_VF_sse2_single
38  * Electrostatics interaction: ReactionField
39  * VdW interaction:            CubicSplineTable
40  * Geometry:                   Water3-Particle
41  * Calculate force/pot:        PotentialAndForce
42  */
43 void
44 nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_VF_sse2_single
45                     (t_nblist * gmx_restrict                nlist,
46                      rvec * gmx_restrict                    xx,
47                      rvec * gmx_restrict                    ff,
48                      t_forcerec * gmx_restrict              fr,
49                      t_mdatoms * gmx_restrict               mdatoms,
50                      nb_kernel_data_t * gmx_restrict        kernel_data,
51                      t_nrnb * gmx_restrict                  nrnb)
52 {
53     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
54      * just 0 for non-waters.
55      * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
56      * jnr indices corresponding to data put in the four positions in the SIMD register.
57      */
58     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
59     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
60     int              jnrA,jnrB,jnrC,jnrD;
61     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
62     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
63     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
64     real             rcutoff_scalar;
65     real             *shiftvec,*fshift,*x,*f;
66     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
67     real             scratch[4*DIM];
68     __m128           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
69     int              vdwioffset0;
70     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
71     int              vdwioffset1;
72     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
73     int              vdwioffset2;
74     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
75     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
76     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
77     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
78     __m128           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
79     __m128           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
80     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
81     real             *charge;
82     int              nvdwtype;
83     __m128           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
84     int              *vdwtype;
85     real             *vdwparam;
86     __m128           one_sixth   = _mm_set1_ps(1.0/6.0);
87     __m128           one_twelfth = _mm_set1_ps(1.0/12.0);
88     __m128i          vfitab;
89     __m128i          ifour       = _mm_set1_epi32(4);
90     __m128           rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
91     real             *vftab;
92     __m128           dummy_mask,cutoff_mask;
93     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
94     __m128           one     = _mm_set1_ps(1.0);
95     __m128           two     = _mm_set1_ps(2.0);
96     x                = xx[0];
97     f                = ff[0];
98
99     nri              = nlist->nri;
100     iinr             = nlist->iinr;
101     jindex           = nlist->jindex;
102     jjnr             = nlist->jjnr;
103     shiftidx         = nlist->shift;
104     gid              = nlist->gid;
105     shiftvec         = fr->shift_vec[0];
106     fshift           = fr->fshift[0];
107     facel            = _mm_set1_ps(fr->epsfac);
108     charge           = mdatoms->chargeA;
109     krf              = _mm_set1_ps(fr->ic->k_rf);
110     krf2             = _mm_set1_ps(fr->ic->k_rf*2.0);
111     crf              = _mm_set1_ps(fr->ic->c_rf);
112     nvdwtype         = fr->ntype;
113     vdwparam         = fr->nbfp;
114     vdwtype          = mdatoms->typeA;
115
116     vftab            = kernel_data->table_vdw->data;
117     vftabscale       = _mm_set1_ps(kernel_data->table_vdw->scale);
118
119     /* Setup water-specific parameters */
120     inr              = nlist->iinr[0];
121     iq0              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+0]));
122     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
123     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
124     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
125
126     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
127     rcutoff_scalar   = fr->rcoulomb;
128     rcutoff          = _mm_set1_ps(rcutoff_scalar);
129     rcutoff2         = _mm_mul_ps(rcutoff,rcutoff);
130
131     /* Avoid stupid compiler warnings */
132     jnrA = jnrB = jnrC = jnrD = 0;
133     j_coord_offsetA = 0;
134     j_coord_offsetB = 0;
135     j_coord_offsetC = 0;
136     j_coord_offsetD = 0;
137
138     outeriter        = 0;
139     inneriter        = 0;
140
141     for(iidx=0;iidx<4*DIM;iidx++)
142     {
143         scratch[iidx] = 0.0;
144     }  
145
146     /* Start outer loop over neighborlists */
147     for(iidx=0; iidx<nri; iidx++)
148     {
149         /* Load shift vector for this list */
150         i_shift_offset   = DIM*shiftidx[iidx];
151
152         /* Load limits for loop over neighbors */
153         j_index_start    = jindex[iidx];
154         j_index_end      = jindex[iidx+1];
155
156         /* Get outer coordinate index */
157         inr              = iinr[iidx];
158         i_coord_offset   = DIM*inr;
159
160         /* Load i particle coords and add shift vector */
161         gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
162                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
163         
164         fix0             = _mm_setzero_ps();
165         fiy0             = _mm_setzero_ps();
166         fiz0             = _mm_setzero_ps();
167         fix1             = _mm_setzero_ps();
168         fiy1             = _mm_setzero_ps();
169         fiz1             = _mm_setzero_ps();
170         fix2             = _mm_setzero_ps();
171         fiy2             = _mm_setzero_ps();
172         fiz2             = _mm_setzero_ps();
173
174         /* Reset potential sums */
175         velecsum         = _mm_setzero_ps();
176         vvdwsum          = _mm_setzero_ps();
177
178         /* Start inner kernel loop */
179         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
180         {
181
182             /* Get j neighbor index, and coordinate index */
183             jnrA             = jjnr[jidx];
184             jnrB             = jjnr[jidx+1];
185             jnrC             = jjnr[jidx+2];
186             jnrD             = jjnr[jidx+3];
187             j_coord_offsetA  = DIM*jnrA;
188             j_coord_offsetB  = DIM*jnrB;
189             j_coord_offsetC  = DIM*jnrC;
190             j_coord_offsetD  = DIM*jnrD;
191
192             /* load j atom coordinates */
193             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
194                                               x+j_coord_offsetC,x+j_coord_offsetD,
195                                               &jx0,&jy0,&jz0);
196
197             /* Calculate displacement vector */
198             dx00             = _mm_sub_ps(ix0,jx0);
199             dy00             = _mm_sub_ps(iy0,jy0);
200             dz00             = _mm_sub_ps(iz0,jz0);
201             dx10             = _mm_sub_ps(ix1,jx0);
202             dy10             = _mm_sub_ps(iy1,jy0);
203             dz10             = _mm_sub_ps(iz1,jz0);
204             dx20             = _mm_sub_ps(ix2,jx0);
205             dy20             = _mm_sub_ps(iy2,jy0);
206             dz20             = _mm_sub_ps(iz2,jz0);
207
208             /* Calculate squared distance and things based on it */
209             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
210             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
211             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
212
213             rinv00           = gmx_mm_invsqrt_ps(rsq00);
214             rinv10           = gmx_mm_invsqrt_ps(rsq10);
215             rinv20           = gmx_mm_invsqrt_ps(rsq20);
216
217             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
218             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
219             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
220
221             /* Load parameters for j particles */
222             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
223                                                               charge+jnrC+0,charge+jnrD+0);
224             vdwjidx0A        = 2*vdwtype[jnrA+0];
225             vdwjidx0B        = 2*vdwtype[jnrB+0];
226             vdwjidx0C        = 2*vdwtype[jnrC+0];
227             vdwjidx0D        = 2*vdwtype[jnrD+0];
228
229             fjx0             = _mm_setzero_ps();
230             fjy0             = _mm_setzero_ps();
231             fjz0             = _mm_setzero_ps();
232
233             /**************************
234              * CALCULATE INTERACTIONS *
235              **************************/
236
237             if (gmx_mm_any_lt(rsq00,rcutoff2))
238             {
239
240             r00              = _mm_mul_ps(rsq00,rinv00);
241
242             /* Compute parameters for interactions between i and j atoms */
243             qq00             = _mm_mul_ps(iq0,jq0);
244             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
245                                          vdwparam+vdwioffset0+vdwjidx0B,
246                                          vdwparam+vdwioffset0+vdwjidx0C,
247                                          vdwparam+vdwioffset0+vdwjidx0D,
248                                          &c6_00,&c12_00);
249
250             /* Calculate table index by multiplying r with table scale and truncate to integer */
251             rt               = _mm_mul_ps(r00,vftabscale);
252             vfitab           = _mm_cvttps_epi32(rt);
253             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
254             vfitab           = _mm_slli_epi32(vfitab,3);
255
256             /* REACTION-FIELD ELECTROSTATICS */
257             velec            = _mm_mul_ps(qq00,_mm_sub_ps(_mm_add_ps(rinv00,_mm_mul_ps(krf,rsq00)),crf));
258             felec            = _mm_mul_ps(qq00,_mm_sub_ps(_mm_mul_ps(rinv00,rinvsq00),krf2));
259
260             /* CUBIC SPLINE TABLE DISPERSION */
261             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
262             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
263             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
264             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
265             _MM_TRANSPOSE4_PS(Y,F,G,H);
266             Heps             = _mm_mul_ps(vfeps,H);
267             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
268             VV               = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
269             vvdw6            = _mm_mul_ps(c6_00,VV);
270             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
271             fvdw6            = _mm_mul_ps(c6_00,FF);
272
273             /* CUBIC SPLINE TABLE REPULSION */
274             vfitab           = _mm_add_epi32(vfitab,ifour);
275             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
276             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
277             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
278             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
279             _MM_TRANSPOSE4_PS(Y,F,G,H);
280             Heps             = _mm_mul_ps(vfeps,H);
281             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
282             VV               = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
283             vvdw12           = _mm_mul_ps(c12_00,VV);
284             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
285             fvdw12           = _mm_mul_ps(c12_00,FF);
286             vvdw             = _mm_add_ps(vvdw12,vvdw6);
287             fvdw             = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
288
289             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
290
291             /* Update potential sum for this i atom from the interaction with this j atom. */
292             velec            = _mm_and_ps(velec,cutoff_mask);
293             velecsum         = _mm_add_ps(velecsum,velec);
294             vvdw             = _mm_and_ps(vvdw,cutoff_mask);
295             vvdwsum          = _mm_add_ps(vvdwsum,vvdw);
296
297             fscal            = _mm_add_ps(felec,fvdw);
298
299             fscal            = _mm_and_ps(fscal,cutoff_mask);
300
301             /* Calculate temporary vectorial force */
302             tx               = _mm_mul_ps(fscal,dx00);
303             ty               = _mm_mul_ps(fscal,dy00);
304             tz               = _mm_mul_ps(fscal,dz00);
305
306             /* Update vectorial force */
307             fix0             = _mm_add_ps(fix0,tx);
308             fiy0             = _mm_add_ps(fiy0,ty);
309             fiz0             = _mm_add_ps(fiz0,tz);
310
311             fjx0             = _mm_add_ps(fjx0,tx);
312             fjy0             = _mm_add_ps(fjy0,ty);
313             fjz0             = _mm_add_ps(fjz0,tz);
314             
315             }
316
317             /**************************
318              * CALCULATE INTERACTIONS *
319              **************************/
320
321             if (gmx_mm_any_lt(rsq10,rcutoff2))
322             {
323
324             /* Compute parameters for interactions between i and j atoms */
325             qq10             = _mm_mul_ps(iq1,jq0);
326
327             /* REACTION-FIELD ELECTROSTATICS */
328             velec            = _mm_mul_ps(qq10,_mm_sub_ps(_mm_add_ps(rinv10,_mm_mul_ps(krf,rsq10)),crf));
329             felec            = _mm_mul_ps(qq10,_mm_sub_ps(_mm_mul_ps(rinv10,rinvsq10),krf2));
330
331             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
332
333             /* Update potential sum for this i atom from the interaction with this j atom. */
334             velec            = _mm_and_ps(velec,cutoff_mask);
335             velecsum         = _mm_add_ps(velecsum,velec);
336
337             fscal            = felec;
338
339             fscal            = _mm_and_ps(fscal,cutoff_mask);
340
341             /* Calculate temporary vectorial force */
342             tx               = _mm_mul_ps(fscal,dx10);
343             ty               = _mm_mul_ps(fscal,dy10);
344             tz               = _mm_mul_ps(fscal,dz10);
345
346             /* Update vectorial force */
347             fix1             = _mm_add_ps(fix1,tx);
348             fiy1             = _mm_add_ps(fiy1,ty);
349             fiz1             = _mm_add_ps(fiz1,tz);
350
351             fjx0             = _mm_add_ps(fjx0,tx);
352             fjy0             = _mm_add_ps(fjy0,ty);
353             fjz0             = _mm_add_ps(fjz0,tz);
354             
355             }
356
357             /**************************
358              * CALCULATE INTERACTIONS *
359              **************************/
360
361             if (gmx_mm_any_lt(rsq20,rcutoff2))
362             {
363
364             /* Compute parameters for interactions between i and j atoms */
365             qq20             = _mm_mul_ps(iq2,jq0);
366
367             /* REACTION-FIELD ELECTROSTATICS */
368             velec            = _mm_mul_ps(qq20,_mm_sub_ps(_mm_add_ps(rinv20,_mm_mul_ps(krf,rsq20)),crf));
369             felec            = _mm_mul_ps(qq20,_mm_sub_ps(_mm_mul_ps(rinv20,rinvsq20),krf2));
370
371             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
372
373             /* Update potential sum for this i atom from the interaction with this j atom. */
374             velec            = _mm_and_ps(velec,cutoff_mask);
375             velecsum         = _mm_add_ps(velecsum,velec);
376
377             fscal            = felec;
378
379             fscal            = _mm_and_ps(fscal,cutoff_mask);
380
381             /* Calculate temporary vectorial force */
382             tx               = _mm_mul_ps(fscal,dx20);
383             ty               = _mm_mul_ps(fscal,dy20);
384             tz               = _mm_mul_ps(fscal,dz20);
385
386             /* Update vectorial force */
387             fix2             = _mm_add_ps(fix2,tx);
388             fiy2             = _mm_add_ps(fiy2,ty);
389             fiz2             = _mm_add_ps(fiz2,tz);
390
391             fjx0             = _mm_add_ps(fjx0,tx);
392             fjy0             = _mm_add_ps(fjy0,ty);
393             fjz0             = _mm_add_ps(fjz0,tz);
394             
395             }
396
397             fjptrA             = f+j_coord_offsetA;
398             fjptrB             = f+j_coord_offsetB;
399             fjptrC             = f+j_coord_offsetC;
400             fjptrD             = f+j_coord_offsetD;
401
402             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
403
404             /* Inner loop uses 144 flops */
405         }
406
407         if(jidx<j_index_end)
408         {
409
410             /* Get j neighbor index, and coordinate index */
411             jnrlistA         = jjnr[jidx];
412             jnrlistB         = jjnr[jidx+1];
413             jnrlistC         = jjnr[jidx+2];
414             jnrlistD         = jjnr[jidx+3];
415             /* Sign of each element will be negative for non-real atoms.
416              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
417              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
418              */
419             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
420             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
421             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
422             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
423             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
424             j_coord_offsetA  = DIM*jnrA;
425             j_coord_offsetB  = DIM*jnrB;
426             j_coord_offsetC  = DIM*jnrC;
427             j_coord_offsetD  = DIM*jnrD;
428
429             /* load j atom coordinates */
430             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
431                                               x+j_coord_offsetC,x+j_coord_offsetD,
432                                               &jx0,&jy0,&jz0);
433
434             /* Calculate displacement vector */
435             dx00             = _mm_sub_ps(ix0,jx0);
436             dy00             = _mm_sub_ps(iy0,jy0);
437             dz00             = _mm_sub_ps(iz0,jz0);
438             dx10             = _mm_sub_ps(ix1,jx0);
439             dy10             = _mm_sub_ps(iy1,jy0);
440             dz10             = _mm_sub_ps(iz1,jz0);
441             dx20             = _mm_sub_ps(ix2,jx0);
442             dy20             = _mm_sub_ps(iy2,jy0);
443             dz20             = _mm_sub_ps(iz2,jz0);
444
445             /* Calculate squared distance and things based on it */
446             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
447             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
448             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
449
450             rinv00           = gmx_mm_invsqrt_ps(rsq00);
451             rinv10           = gmx_mm_invsqrt_ps(rsq10);
452             rinv20           = gmx_mm_invsqrt_ps(rsq20);
453
454             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
455             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
456             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
457
458             /* Load parameters for j particles */
459             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
460                                                               charge+jnrC+0,charge+jnrD+0);
461             vdwjidx0A        = 2*vdwtype[jnrA+0];
462             vdwjidx0B        = 2*vdwtype[jnrB+0];
463             vdwjidx0C        = 2*vdwtype[jnrC+0];
464             vdwjidx0D        = 2*vdwtype[jnrD+0];
465
466             fjx0             = _mm_setzero_ps();
467             fjy0             = _mm_setzero_ps();
468             fjz0             = _mm_setzero_ps();
469
470             /**************************
471              * CALCULATE INTERACTIONS *
472              **************************/
473
474             if (gmx_mm_any_lt(rsq00,rcutoff2))
475             {
476
477             r00              = _mm_mul_ps(rsq00,rinv00);
478             r00              = _mm_andnot_ps(dummy_mask,r00);
479
480             /* Compute parameters for interactions between i and j atoms */
481             qq00             = _mm_mul_ps(iq0,jq0);
482             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
483                                          vdwparam+vdwioffset0+vdwjidx0B,
484                                          vdwparam+vdwioffset0+vdwjidx0C,
485                                          vdwparam+vdwioffset0+vdwjidx0D,
486                                          &c6_00,&c12_00);
487
488             /* Calculate table index by multiplying r with table scale and truncate to integer */
489             rt               = _mm_mul_ps(r00,vftabscale);
490             vfitab           = _mm_cvttps_epi32(rt);
491             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
492             vfitab           = _mm_slli_epi32(vfitab,3);
493
494             /* REACTION-FIELD ELECTROSTATICS */
495             velec            = _mm_mul_ps(qq00,_mm_sub_ps(_mm_add_ps(rinv00,_mm_mul_ps(krf,rsq00)),crf));
496             felec            = _mm_mul_ps(qq00,_mm_sub_ps(_mm_mul_ps(rinv00,rinvsq00),krf2));
497
498             /* CUBIC SPLINE TABLE DISPERSION */
499             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
500             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
501             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
502             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
503             _MM_TRANSPOSE4_PS(Y,F,G,H);
504             Heps             = _mm_mul_ps(vfeps,H);
505             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
506             VV               = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
507             vvdw6            = _mm_mul_ps(c6_00,VV);
508             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
509             fvdw6            = _mm_mul_ps(c6_00,FF);
510
511             /* CUBIC SPLINE TABLE REPULSION */
512             vfitab           = _mm_add_epi32(vfitab,ifour);
513             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
514             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
515             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
516             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
517             _MM_TRANSPOSE4_PS(Y,F,G,H);
518             Heps             = _mm_mul_ps(vfeps,H);
519             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
520             VV               = _mm_add_ps(Y,_mm_mul_ps(vfeps,Fp));
521             vvdw12           = _mm_mul_ps(c12_00,VV);
522             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
523             fvdw12           = _mm_mul_ps(c12_00,FF);
524             vvdw             = _mm_add_ps(vvdw12,vvdw6);
525             fvdw             = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
526
527             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
528
529             /* Update potential sum for this i atom from the interaction with this j atom. */
530             velec            = _mm_and_ps(velec,cutoff_mask);
531             velec            = _mm_andnot_ps(dummy_mask,velec);
532             velecsum         = _mm_add_ps(velecsum,velec);
533             vvdw             = _mm_and_ps(vvdw,cutoff_mask);
534             vvdw             = _mm_andnot_ps(dummy_mask,vvdw);
535             vvdwsum          = _mm_add_ps(vvdwsum,vvdw);
536
537             fscal            = _mm_add_ps(felec,fvdw);
538
539             fscal            = _mm_and_ps(fscal,cutoff_mask);
540
541             fscal            = _mm_andnot_ps(dummy_mask,fscal);
542
543             /* Calculate temporary vectorial force */
544             tx               = _mm_mul_ps(fscal,dx00);
545             ty               = _mm_mul_ps(fscal,dy00);
546             tz               = _mm_mul_ps(fscal,dz00);
547
548             /* Update vectorial force */
549             fix0             = _mm_add_ps(fix0,tx);
550             fiy0             = _mm_add_ps(fiy0,ty);
551             fiz0             = _mm_add_ps(fiz0,tz);
552
553             fjx0             = _mm_add_ps(fjx0,tx);
554             fjy0             = _mm_add_ps(fjy0,ty);
555             fjz0             = _mm_add_ps(fjz0,tz);
556             
557             }
558
559             /**************************
560              * CALCULATE INTERACTIONS *
561              **************************/
562
563             if (gmx_mm_any_lt(rsq10,rcutoff2))
564             {
565
566             /* Compute parameters for interactions between i and j atoms */
567             qq10             = _mm_mul_ps(iq1,jq0);
568
569             /* REACTION-FIELD ELECTROSTATICS */
570             velec            = _mm_mul_ps(qq10,_mm_sub_ps(_mm_add_ps(rinv10,_mm_mul_ps(krf,rsq10)),crf));
571             felec            = _mm_mul_ps(qq10,_mm_sub_ps(_mm_mul_ps(rinv10,rinvsq10),krf2));
572
573             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
574
575             /* Update potential sum for this i atom from the interaction with this j atom. */
576             velec            = _mm_and_ps(velec,cutoff_mask);
577             velec            = _mm_andnot_ps(dummy_mask,velec);
578             velecsum         = _mm_add_ps(velecsum,velec);
579
580             fscal            = felec;
581
582             fscal            = _mm_and_ps(fscal,cutoff_mask);
583
584             fscal            = _mm_andnot_ps(dummy_mask,fscal);
585
586             /* Calculate temporary vectorial force */
587             tx               = _mm_mul_ps(fscal,dx10);
588             ty               = _mm_mul_ps(fscal,dy10);
589             tz               = _mm_mul_ps(fscal,dz10);
590
591             /* Update vectorial force */
592             fix1             = _mm_add_ps(fix1,tx);
593             fiy1             = _mm_add_ps(fiy1,ty);
594             fiz1             = _mm_add_ps(fiz1,tz);
595
596             fjx0             = _mm_add_ps(fjx0,tx);
597             fjy0             = _mm_add_ps(fjy0,ty);
598             fjz0             = _mm_add_ps(fjz0,tz);
599             
600             }
601
602             /**************************
603              * CALCULATE INTERACTIONS *
604              **************************/
605
606             if (gmx_mm_any_lt(rsq20,rcutoff2))
607             {
608
609             /* Compute parameters for interactions between i and j atoms */
610             qq20             = _mm_mul_ps(iq2,jq0);
611
612             /* REACTION-FIELD ELECTROSTATICS */
613             velec            = _mm_mul_ps(qq20,_mm_sub_ps(_mm_add_ps(rinv20,_mm_mul_ps(krf,rsq20)),crf));
614             felec            = _mm_mul_ps(qq20,_mm_sub_ps(_mm_mul_ps(rinv20,rinvsq20),krf2));
615
616             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
617
618             /* Update potential sum for this i atom from the interaction with this j atom. */
619             velec            = _mm_and_ps(velec,cutoff_mask);
620             velec            = _mm_andnot_ps(dummy_mask,velec);
621             velecsum         = _mm_add_ps(velecsum,velec);
622
623             fscal            = felec;
624
625             fscal            = _mm_and_ps(fscal,cutoff_mask);
626
627             fscal            = _mm_andnot_ps(dummy_mask,fscal);
628
629             /* Calculate temporary vectorial force */
630             tx               = _mm_mul_ps(fscal,dx20);
631             ty               = _mm_mul_ps(fscal,dy20);
632             tz               = _mm_mul_ps(fscal,dz20);
633
634             /* Update vectorial force */
635             fix2             = _mm_add_ps(fix2,tx);
636             fiy2             = _mm_add_ps(fiy2,ty);
637             fiz2             = _mm_add_ps(fiz2,tz);
638
639             fjx0             = _mm_add_ps(fjx0,tx);
640             fjy0             = _mm_add_ps(fjy0,ty);
641             fjz0             = _mm_add_ps(fjz0,tz);
642             
643             }
644
645             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
646             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
647             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
648             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
649
650             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
651
652             /* Inner loop uses 145 flops */
653         }
654
655         /* End of innermost loop */
656
657         gmx_mm_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
658                                               f+i_coord_offset,fshift+i_shift_offset);
659
660         ggid                        = gid[iidx];
661         /* Update potential energies */
662         gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
663         gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
664
665         /* Increment number of inner iterations */
666         inneriter                  += j_index_end - j_index_start;
667
668         /* Outer loop uses 20 flops */
669     }
670
671     /* Increment number of outer iterations */
672     outeriter        += nri;
673
674     /* Update outer/inner flops */
675
676     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_VF,outeriter*20 + inneriter*145);
677 }
678 /*
679  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_F_sse2_single
680  * Electrostatics interaction: ReactionField
681  * VdW interaction:            CubicSplineTable
682  * Geometry:                   Water3-Particle
683  * Calculate force/pot:        Force
684  */
685 void
686 nb_kernel_ElecRFCut_VdwCSTab_GeomW3P1_F_sse2_single
687                     (t_nblist * gmx_restrict                nlist,
688                      rvec * gmx_restrict                    xx,
689                      rvec * gmx_restrict                    ff,
690                      t_forcerec * gmx_restrict              fr,
691                      t_mdatoms * gmx_restrict               mdatoms,
692                      nb_kernel_data_t * gmx_restrict        kernel_data,
693                      t_nrnb * gmx_restrict                  nrnb)
694 {
695     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
696      * just 0 for non-waters.
697      * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
698      * jnr indices corresponding to data put in the four positions in the SIMD register.
699      */
700     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
701     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
702     int              jnrA,jnrB,jnrC,jnrD;
703     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
704     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
705     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
706     real             rcutoff_scalar;
707     real             *shiftvec,*fshift,*x,*f;
708     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
709     real             scratch[4*DIM];
710     __m128           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
711     int              vdwioffset0;
712     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
713     int              vdwioffset1;
714     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
715     int              vdwioffset2;
716     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
717     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
718     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
719     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
720     __m128           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
721     __m128           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
722     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
723     real             *charge;
724     int              nvdwtype;
725     __m128           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
726     int              *vdwtype;
727     real             *vdwparam;
728     __m128           one_sixth   = _mm_set1_ps(1.0/6.0);
729     __m128           one_twelfth = _mm_set1_ps(1.0/12.0);
730     __m128i          vfitab;
731     __m128i          ifour       = _mm_set1_epi32(4);
732     __m128           rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
733     real             *vftab;
734     __m128           dummy_mask,cutoff_mask;
735     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
736     __m128           one     = _mm_set1_ps(1.0);
737     __m128           two     = _mm_set1_ps(2.0);
738     x                = xx[0];
739     f                = ff[0];
740
741     nri              = nlist->nri;
742     iinr             = nlist->iinr;
743     jindex           = nlist->jindex;
744     jjnr             = nlist->jjnr;
745     shiftidx         = nlist->shift;
746     gid              = nlist->gid;
747     shiftvec         = fr->shift_vec[0];
748     fshift           = fr->fshift[0];
749     facel            = _mm_set1_ps(fr->epsfac);
750     charge           = mdatoms->chargeA;
751     krf              = _mm_set1_ps(fr->ic->k_rf);
752     krf2             = _mm_set1_ps(fr->ic->k_rf*2.0);
753     crf              = _mm_set1_ps(fr->ic->c_rf);
754     nvdwtype         = fr->ntype;
755     vdwparam         = fr->nbfp;
756     vdwtype          = mdatoms->typeA;
757
758     vftab            = kernel_data->table_vdw->data;
759     vftabscale       = _mm_set1_ps(kernel_data->table_vdw->scale);
760
761     /* Setup water-specific parameters */
762     inr              = nlist->iinr[0];
763     iq0              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+0]));
764     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
765     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
766     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
767
768     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
769     rcutoff_scalar   = fr->rcoulomb;
770     rcutoff          = _mm_set1_ps(rcutoff_scalar);
771     rcutoff2         = _mm_mul_ps(rcutoff,rcutoff);
772
773     /* Avoid stupid compiler warnings */
774     jnrA = jnrB = jnrC = jnrD = 0;
775     j_coord_offsetA = 0;
776     j_coord_offsetB = 0;
777     j_coord_offsetC = 0;
778     j_coord_offsetD = 0;
779
780     outeriter        = 0;
781     inneriter        = 0;
782
783     for(iidx=0;iidx<4*DIM;iidx++)
784     {
785         scratch[iidx] = 0.0;
786     }  
787
788     /* Start outer loop over neighborlists */
789     for(iidx=0; iidx<nri; iidx++)
790     {
791         /* Load shift vector for this list */
792         i_shift_offset   = DIM*shiftidx[iidx];
793
794         /* Load limits for loop over neighbors */
795         j_index_start    = jindex[iidx];
796         j_index_end      = jindex[iidx+1];
797
798         /* Get outer coordinate index */
799         inr              = iinr[iidx];
800         i_coord_offset   = DIM*inr;
801
802         /* Load i particle coords and add shift vector */
803         gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
804                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
805         
806         fix0             = _mm_setzero_ps();
807         fiy0             = _mm_setzero_ps();
808         fiz0             = _mm_setzero_ps();
809         fix1             = _mm_setzero_ps();
810         fiy1             = _mm_setzero_ps();
811         fiz1             = _mm_setzero_ps();
812         fix2             = _mm_setzero_ps();
813         fiy2             = _mm_setzero_ps();
814         fiz2             = _mm_setzero_ps();
815
816         /* Start inner kernel loop */
817         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
818         {
819
820             /* Get j neighbor index, and coordinate index */
821             jnrA             = jjnr[jidx];
822             jnrB             = jjnr[jidx+1];
823             jnrC             = jjnr[jidx+2];
824             jnrD             = jjnr[jidx+3];
825             j_coord_offsetA  = DIM*jnrA;
826             j_coord_offsetB  = DIM*jnrB;
827             j_coord_offsetC  = DIM*jnrC;
828             j_coord_offsetD  = DIM*jnrD;
829
830             /* load j atom coordinates */
831             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
832                                               x+j_coord_offsetC,x+j_coord_offsetD,
833                                               &jx0,&jy0,&jz0);
834
835             /* Calculate displacement vector */
836             dx00             = _mm_sub_ps(ix0,jx0);
837             dy00             = _mm_sub_ps(iy0,jy0);
838             dz00             = _mm_sub_ps(iz0,jz0);
839             dx10             = _mm_sub_ps(ix1,jx0);
840             dy10             = _mm_sub_ps(iy1,jy0);
841             dz10             = _mm_sub_ps(iz1,jz0);
842             dx20             = _mm_sub_ps(ix2,jx0);
843             dy20             = _mm_sub_ps(iy2,jy0);
844             dz20             = _mm_sub_ps(iz2,jz0);
845
846             /* Calculate squared distance and things based on it */
847             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
848             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
849             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
850
851             rinv00           = gmx_mm_invsqrt_ps(rsq00);
852             rinv10           = gmx_mm_invsqrt_ps(rsq10);
853             rinv20           = gmx_mm_invsqrt_ps(rsq20);
854
855             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
856             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
857             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
858
859             /* Load parameters for j particles */
860             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
861                                                               charge+jnrC+0,charge+jnrD+0);
862             vdwjidx0A        = 2*vdwtype[jnrA+0];
863             vdwjidx0B        = 2*vdwtype[jnrB+0];
864             vdwjidx0C        = 2*vdwtype[jnrC+0];
865             vdwjidx0D        = 2*vdwtype[jnrD+0];
866
867             fjx0             = _mm_setzero_ps();
868             fjy0             = _mm_setzero_ps();
869             fjz0             = _mm_setzero_ps();
870
871             /**************************
872              * CALCULATE INTERACTIONS *
873              **************************/
874
875             if (gmx_mm_any_lt(rsq00,rcutoff2))
876             {
877
878             r00              = _mm_mul_ps(rsq00,rinv00);
879
880             /* Compute parameters for interactions between i and j atoms */
881             qq00             = _mm_mul_ps(iq0,jq0);
882             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
883                                          vdwparam+vdwioffset0+vdwjidx0B,
884                                          vdwparam+vdwioffset0+vdwjidx0C,
885                                          vdwparam+vdwioffset0+vdwjidx0D,
886                                          &c6_00,&c12_00);
887
888             /* Calculate table index by multiplying r with table scale and truncate to integer */
889             rt               = _mm_mul_ps(r00,vftabscale);
890             vfitab           = _mm_cvttps_epi32(rt);
891             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
892             vfitab           = _mm_slli_epi32(vfitab,3);
893
894             /* REACTION-FIELD ELECTROSTATICS */
895             felec            = _mm_mul_ps(qq00,_mm_sub_ps(_mm_mul_ps(rinv00,rinvsq00),krf2));
896
897             /* CUBIC SPLINE TABLE DISPERSION */
898             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
899             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
900             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
901             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
902             _MM_TRANSPOSE4_PS(Y,F,G,H);
903             Heps             = _mm_mul_ps(vfeps,H);
904             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
905             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
906             fvdw6            = _mm_mul_ps(c6_00,FF);
907
908             /* CUBIC SPLINE TABLE REPULSION */
909             vfitab           = _mm_add_epi32(vfitab,ifour);
910             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
911             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
912             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
913             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
914             _MM_TRANSPOSE4_PS(Y,F,G,H);
915             Heps             = _mm_mul_ps(vfeps,H);
916             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
917             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
918             fvdw12           = _mm_mul_ps(c12_00,FF);
919             fvdw             = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
920
921             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
922
923             fscal            = _mm_add_ps(felec,fvdw);
924
925             fscal            = _mm_and_ps(fscal,cutoff_mask);
926
927             /* Calculate temporary vectorial force */
928             tx               = _mm_mul_ps(fscal,dx00);
929             ty               = _mm_mul_ps(fscal,dy00);
930             tz               = _mm_mul_ps(fscal,dz00);
931
932             /* Update vectorial force */
933             fix0             = _mm_add_ps(fix0,tx);
934             fiy0             = _mm_add_ps(fiy0,ty);
935             fiz0             = _mm_add_ps(fiz0,tz);
936
937             fjx0             = _mm_add_ps(fjx0,tx);
938             fjy0             = _mm_add_ps(fjy0,ty);
939             fjz0             = _mm_add_ps(fjz0,tz);
940             
941             }
942
943             /**************************
944              * CALCULATE INTERACTIONS *
945              **************************/
946
947             if (gmx_mm_any_lt(rsq10,rcutoff2))
948             {
949
950             /* Compute parameters for interactions between i and j atoms */
951             qq10             = _mm_mul_ps(iq1,jq0);
952
953             /* REACTION-FIELD ELECTROSTATICS */
954             felec            = _mm_mul_ps(qq10,_mm_sub_ps(_mm_mul_ps(rinv10,rinvsq10),krf2));
955
956             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
957
958             fscal            = felec;
959
960             fscal            = _mm_and_ps(fscal,cutoff_mask);
961
962             /* Calculate temporary vectorial force */
963             tx               = _mm_mul_ps(fscal,dx10);
964             ty               = _mm_mul_ps(fscal,dy10);
965             tz               = _mm_mul_ps(fscal,dz10);
966
967             /* Update vectorial force */
968             fix1             = _mm_add_ps(fix1,tx);
969             fiy1             = _mm_add_ps(fiy1,ty);
970             fiz1             = _mm_add_ps(fiz1,tz);
971
972             fjx0             = _mm_add_ps(fjx0,tx);
973             fjy0             = _mm_add_ps(fjy0,ty);
974             fjz0             = _mm_add_ps(fjz0,tz);
975             
976             }
977
978             /**************************
979              * CALCULATE INTERACTIONS *
980              **************************/
981
982             if (gmx_mm_any_lt(rsq20,rcutoff2))
983             {
984
985             /* Compute parameters for interactions between i and j atoms */
986             qq20             = _mm_mul_ps(iq2,jq0);
987
988             /* REACTION-FIELD ELECTROSTATICS */
989             felec            = _mm_mul_ps(qq20,_mm_sub_ps(_mm_mul_ps(rinv20,rinvsq20),krf2));
990
991             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
992
993             fscal            = felec;
994
995             fscal            = _mm_and_ps(fscal,cutoff_mask);
996
997             /* Calculate temporary vectorial force */
998             tx               = _mm_mul_ps(fscal,dx20);
999             ty               = _mm_mul_ps(fscal,dy20);
1000             tz               = _mm_mul_ps(fscal,dz20);
1001
1002             /* Update vectorial force */
1003             fix2             = _mm_add_ps(fix2,tx);
1004             fiy2             = _mm_add_ps(fiy2,ty);
1005             fiz2             = _mm_add_ps(fiz2,tz);
1006
1007             fjx0             = _mm_add_ps(fjx0,tx);
1008             fjy0             = _mm_add_ps(fjy0,ty);
1009             fjz0             = _mm_add_ps(fjz0,tz);
1010             
1011             }
1012
1013             fjptrA             = f+j_coord_offsetA;
1014             fjptrB             = f+j_coord_offsetB;
1015             fjptrC             = f+j_coord_offsetC;
1016             fjptrD             = f+j_coord_offsetD;
1017
1018             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1019
1020             /* Inner loop uses 117 flops */
1021         }
1022
1023         if(jidx<j_index_end)
1024         {
1025
1026             /* Get j neighbor index, and coordinate index */
1027             jnrlistA         = jjnr[jidx];
1028             jnrlistB         = jjnr[jidx+1];
1029             jnrlistC         = jjnr[jidx+2];
1030             jnrlistD         = jjnr[jidx+3];
1031             /* Sign of each element will be negative for non-real atoms.
1032              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1033              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1034              */
1035             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
1036             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
1037             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
1038             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
1039             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
1040             j_coord_offsetA  = DIM*jnrA;
1041             j_coord_offsetB  = DIM*jnrB;
1042             j_coord_offsetC  = DIM*jnrC;
1043             j_coord_offsetD  = DIM*jnrD;
1044
1045             /* load j atom coordinates */
1046             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1047                                               x+j_coord_offsetC,x+j_coord_offsetD,
1048                                               &jx0,&jy0,&jz0);
1049
1050             /* Calculate displacement vector */
1051             dx00             = _mm_sub_ps(ix0,jx0);
1052             dy00             = _mm_sub_ps(iy0,jy0);
1053             dz00             = _mm_sub_ps(iz0,jz0);
1054             dx10             = _mm_sub_ps(ix1,jx0);
1055             dy10             = _mm_sub_ps(iy1,jy0);
1056             dz10             = _mm_sub_ps(iz1,jz0);
1057             dx20             = _mm_sub_ps(ix2,jx0);
1058             dy20             = _mm_sub_ps(iy2,jy0);
1059             dz20             = _mm_sub_ps(iz2,jz0);
1060
1061             /* Calculate squared distance and things based on it */
1062             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
1063             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
1064             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
1065
1066             rinv00           = gmx_mm_invsqrt_ps(rsq00);
1067             rinv10           = gmx_mm_invsqrt_ps(rsq10);
1068             rinv20           = gmx_mm_invsqrt_ps(rsq20);
1069
1070             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
1071             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
1072             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
1073
1074             /* Load parameters for j particles */
1075             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1076                                                               charge+jnrC+0,charge+jnrD+0);
1077             vdwjidx0A        = 2*vdwtype[jnrA+0];
1078             vdwjidx0B        = 2*vdwtype[jnrB+0];
1079             vdwjidx0C        = 2*vdwtype[jnrC+0];
1080             vdwjidx0D        = 2*vdwtype[jnrD+0];
1081
1082             fjx0             = _mm_setzero_ps();
1083             fjy0             = _mm_setzero_ps();
1084             fjz0             = _mm_setzero_ps();
1085
1086             /**************************
1087              * CALCULATE INTERACTIONS *
1088              **************************/
1089
1090             if (gmx_mm_any_lt(rsq00,rcutoff2))
1091             {
1092
1093             r00              = _mm_mul_ps(rsq00,rinv00);
1094             r00              = _mm_andnot_ps(dummy_mask,r00);
1095
1096             /* Compute parameters for interactions between i and j atoms */
1097             qq00             = _mm_mul_ps(iq0,jq0);
1098             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
1099                                          vdwparam+vdwioffset0+vdwjidx0B,
1100                                          vdwparam+vdwioffset0+vdwjidx0C,
1101                                          vdwparam+vdwioffset0+vdwjidx0D,
1102                                          &c6_00,&c12_00);
1103
1104             /* Calculate table index by multiplying r with table scale and truncate to integer */
1105             rt               = _mm_mul_ps(r00,vftabscale);
1106             vfitab           = _mm_cvttps_epi32(rt);
1107             vfeps            = _mm_sub_ps(rt,_mm_cvtepi32_ps(vfitab));
1108             vfitab           = _mm_slli_epi32(vfitab,3);
1109
1110             /* REACTION-FIELD ELECTROSTATICS */
1111             felec            = _mm_mul_ps(qq00,_mm_sub_ps(_mm_mul_ps(rinv00,rinvsq00),krf2));
1112
1113             /* CUBIC SPLINE TABLE DISPERSION */
1114             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
1115             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
1116             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
1117             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
1118             _MM_TRANSPOSE4_PS(Y,F,G,H);
1119             Heps             = _mm_mul_ps(vfeps,H);
1120             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
1121             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
1122             fvdw6            = _mm_mul_ps(c6_00,FF);
1123
1124             /* CUBIC SPLINE TABLE REPULSION */
1125             vfitab           = _mm_add_epi32(vfitab,ifour);
1126             Y                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,0) );
1127             F                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,1) );
1128             G                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,2) );
1129             H                = _mm_load_ps( vftab + gmx_mm_extract_epi32(vfitab,3) );
1130             _MM_TRANSPOSE4_PS(Y,F,G,H);
1131             Heps             = _mm_mul_ps(vfeps,H);
1132             Fp               = _mm_add_ps(F,_mm_mul_ps(vfeps,_mm_add_ps(G,Heps)));
1133             FF               = _mm_add_ps(Fp,_mm_mul_ps(vfeps,_mm_add_ps(G,_mm_add_ps(Heps,Heps))));
1134             fvdw12           = _mm_mul_ps(c12_00,FF);
1135             fvdw             = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv00)));
1136
1137             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
1138
1139             fscal            = _mm_add_ps(felec,fvdw);
1140
1141             fscal            = _mm_and_ps(fscal,cutoff_mask);
1142
1143             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1144
1145             /* Calculate temporary vectorial force */
1146             tx               = _mm_mul_ps(fscal,dx00);
1147             ty               = _mm_mul_ps(fscal,dy00);
1148             tz               = _mm_mul_ps(fscal,dz00);
1149
1150             /* Update vectorial force */
1151             fix0             = _mm_add_ps(fix0,tx);
1152             fiy0             = _mm_add_ps(fiy0,ty);
1153             fiz0             = _mm_add_ps(fiz0,tz);
1154
1155             fjx0             = _mm_add_ps(fjx0,tx);
1156             fjy0             = _mm_add_ps(fjy0,ty);
1157             fjz0             = _mm_add_ps(fjz0,tz);
1158             
1159             }
1160
1161             /**************************
1162              * CALCULATE INTERACTIONS *
1163              **************************/
1164
1165             if (gmx_mm_any_lt(rsq10,rcutoff2))
1166             {
1167
1168             /* Compute parameters for interactions between i and j atoms */
1169             qq10             = _mm_mul_ps(iq1,jq0);
1170
1171             /* REACTION-FIELD ELECTROSTATICS */
1172             felec            = _mm_mul_ps(qq10,_mm_sub_ps(_mm_mul_ps(rinv10,rinvsq10),krf2));
1173
1174             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
1175
1176             fscal            = felec;
1177
1178             fscal            = _mm_and_ps(fscal,cutoff_mask);
1179
1180             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1181
1182             /* Calculate temporary vectorial force */
1183             tx               = _mm_mul_ps(fscal,dx10);
1184             ty               = _mm_mul_ps(fscal,dy10);
1185             tz               = _mm_mul_ps(fscal,dz10);
1186
1187             /* Update vectorial force */
1188             fix1             = _mm_add_ps(fix1,tx);
1189             fiy1             = _mm_add_ps(fiy1,ty);
1190             fiz1             = _mm_add_ps(fiz1,tz);
1191
1192             fjx0             = _mm_add_ps(fjx0,tx);
1193             fjy0             = _mm_add_ps(fjy0,ty);
1194             fjz0             = _mm_add_ps(fjz0,tz);
1195             
1196             }
1197
1198             /**************************
1199              * CALCULATE INTERACTIONS *
1200              **************************/
1201
1202             if (gmx_mm_any_lt(rsq20,rcutoff2))
1203             {
1204
1205             /* Compute parameters for interactions between i and j atoms */
1206             qq20             = _mm_mul_ps(iq2,jq0);
1207
1208             /* REACTION-FIELD ELECTROSTATICS */
1209             felec            = _mm_mul_ps(qq20,_mm_sub_ps(_mm_mul_ps(rinv20,rinvsq20),krf2));
1210
1211             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
1212
1213             fscal            = felec;
1214
1215             fscal            = _mm_and_ps(fscal,cutoff_mask);
1216
1217             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1218
1219             /* Calculate temporary vectorial force */
1220             tx               = _mm_mul_ps(fscal,dx20);
1221             ty               = _mm_mul_ps(fscal,dy20);
1222             tz               = _mm_mul_ps(fscal,dz20);
1223
1224             /* Update vectorial force */
1225             fix2             = _mm_add_ps(fix2,tx);
1226             fiy2             = _mm_add_ps(fiy2,ty);
1227             fiz2             = _mm_add_ps(fiz2,tz);
1228
1229             fjx0             = _mm_add_ps(fjx0,tx);
1230             fjy0             = _mm_add_ps(fjy0,ty);
1231             fjz0             = _mm_add_ps(fjz0,tz);
1232             
1233             }
1234
1235             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1236             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1237             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1238             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1239
1240             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1241
1242             /* Inner loop uses 118 flops */
1243         }
1244
1245         /* End of innermost loop */
1246
1247         gmx_mm_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1248                                               f+i_coord_offset,fshift+i_shift_offset);
1249
1250         /* Increment number of inner iterations */
1251         inneriter                  += j_index_end - j_index_start;
1252
1253         /* Outer loop uses 18 flops */
1254     }
1255
1256     /* Increment number of outer iterations */
1257     outeriter        += nri;
1258
1259     /* Update outer/inner flops */
1260
1261     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3_F,outeriter*18 + inneriter*118);
1262 }