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