Remove all unnecessary HAVE_CONFIG_H
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sse2_single / nb_kernel_ElecRFCut_VdwNone_GeomW4P1_sse2_single.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS sse2_single kernel generator.
37  */
38 #include "config.h"
39
40 #include <math.h>
41
42 #include "../nb_kernel.h"
43 #include "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_ElecRFCut_VdwNone_GeomW4P1_VF_sse2_single
52  * Electrostatics interaction: ReactionField
53  * VdW interaction:            None
54  * Geometry:                   Water4-Particle
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecRFCut_VdwNone_GeomW4P1_VF_sse2_single
59                     (t_nblist                    * gmx_restrict       nlist,
60                      rvec                        * gmx_restrict          xx,
61                      rvec                        * gmx_restrict          ff,
62                      t_forcerec                  * gmx_restrict          fr,
63                      t_mdatoms                   * gmx_restrict     mdatoms,
64                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
65                      t_nrnb                      * gmx_restrict        nrnb)
66 {
67     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
68      * just 0 for non-waters.
69      * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
70      * jnr indices corresponding to data put in the four positions in the SIMD register.
71      */
72     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
73     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
74     int              jnrA,jnrB,jnrC,jnrD;
75     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
76     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
77     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
78     real             rcutoff_scalar;
79     real             *shiftvec,*fshift,*x,*f;
80     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
81     real             scratch[4*DIM];
82     __m128           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
83     int              vdwioffset1;
84     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
85     int              vdwioffset2;
86     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
87     int              vdwioffset3;
88     __m128           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
89     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
90     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
91     __m128           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
92     __m128           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
93     __m128           dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
94     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
95     real             *charge;
96     __m128           dummy_mask,cutoff_mask;
97     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
98     __m128           one     = _mm_set1_ps(1.0);
99     __m128           two     = _mm_set1_ps(2.0);
100     x                = xx[0];
101     f                = ff[0];
102
103     nri              = nlist->nri;
104     iinr             = nlist->iinr;
105     jindex           = nlist->jindex;
106     jjnr             = nlist->jjnr;
107     shiftidx         = nlist->shift;
108     gid              = nlist->gid;
109     shiftvec         = fr->shift_vec[0];
110     fshift           = fr->fshift[0];
111     facel            = _mm_set1_ps(fr->epsfac);
112     charge           = mdatoms->chargeA;
113     krf              = _mm_set1_ps(fr->ic->k_rf);
114     krf2             = _mm_set1_ps(fr->ic->k_rf*2.0);
115     crf              = _mm_set1_ps(fr->ic->c_rf);
116
117     /* Setup water-specific parameters */
118     inr              = nlist->iinr[0];
119     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
120     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
121     iq3              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
122
123     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
124     rcutoff_scalar   = fr->rcoulomb;
125     rcutoff          = _mm_set1_ps(rcutoff_scalar);
126     rcutoff2         = _mm_mul_ps(rcutoff,rcutoff);
127
128     /* Avoid stupid compiler warnings */
129     jnrA = jnrB = jnrC = jnrD = 0;
130     j_coord_offsetA = 0;
131     j_coord_offsetB = 0;
132     j_coord_offsetC = 0;
133     j_coord_offsetD = 0;
134
135     outeriter        = 0;
136     inneriter        = 0;
137
138     for(iidx=0;iidx<4*DIM;iidx++)
139     {
140         scratch[iidx] = 0.0;
141     }  
142
143     /* Start outer loop over neighborlists */
144     for(iidx=0; iidx<nri; iidx++)
145     {
146         /* Load shift vector for this list */
147         i_shift_offset   = DIM*shiftidx[iidx];
148
149         /* Load limits for loop over neighbors */
150         j_index_start    = jindex[iidx];
151         j_index_end      = jindex[iidx+1];
152
153         /* Get outer coordinate index */
154         inr              = iinr[iidx];
155         i_coord_offset   = DIM*inr;
156
157         /* Load i particle coords and add shift vector */
158         gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
159                                                  &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
160         
161         fix1             = _mm_setzero_ps();
162         fiy1             = _mm_setzero_ps();
163         fiz1             = _mm_setzero_ps();
164         fix2             = _mm_setzero_ps();
165         fiy2             = _mm_setzero_ps();
166         fiz2             = _mm_setzero_ps();
167         fix3             = _mm_setzero_ps();
168         fiy3             = _mm_setzero_ps();
169         fiz3             = _mm_setzero_ps();
170
171         /* Reset potential sums */
172         velecsum         = _mm_setzero_ps();
173
174         /* Start inner kernel loop */
175         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
176         {
177
178             /* Get j neighbor index, and coordinate index */
179             jnrA             = jjnr[jidx];
180             jnrB             = jjnr[jidx+1];
181             jnrC             = jjnr[jidx+2];
182             jnrD             = jjnr[jidx+3];
183             j_coord_offsetA  = DIM*jnrA;
184             j_coord_offsetB  = DIM*jnrB;
185             j_coord_offsetC  = DIM*jnrC;
186             j_coord_offsetD  = DIM*jnrD;
187
188             /* load j atom coordinates */
189             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
190                                               x+j_coord_offsetC,x+j_coord_offsetD,
191                                               &jx0,&jy0,&jz0);
192
193             /* Calculate displacement vector */
194             dx10             = _mm_sub_ps(ix1,jx0);
195             dy10             = _mm_sub_ps(iy1,jy0);
196             dz10             = _mm_sub_ps(iz1,jz0);
197             dx20             = _mm_sub_ps(ix2,jx0);
198             dy20             = _mm_sub_ps(iy2,jy0);
199             dz20             = _mm_sub_ps(iz2,jz0);
200             dx30             = _mm_sub_ps(ix3,jx0);
201             dy30             = _mm_sub_ps(iy3,jy0);
202             dz30             = _mm_sub_ps(iz3,jz0);
203
204             /* Calculate squared distance and things based on it */
205             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
206             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
207             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
208
209             rinv10           = gmx_mm_invsqrt_ps(rsq10);
210             rinv20           = gmx_mm_invsqrt_ps(rsq20);
211             rinv30           = gmx_mm_invsqrt_ps(rsq30);
212
213             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
214             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
215             rinvsq30         = _mm_mul_ps(rinv30,rinv30);
216
217             /* Load parameters for j particles */
218             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
219                                                               charge+jnrC+0,charge+jnrD+0);
220
221             fjx0             = _mm_setzero_ps();
222             fjy0             = _mm_setzero_ps();
223             fjz0             = _mm_setzero_ps();
224
225             /**************************
226              * CALCULATE INTERACTIONS *
227              **************************/
228
229             if (gmx_mm_any_lt(rsq10,rcutoff2))
230             {
231
232             /* Compute parameters for interactions between i and j atoms */
233             qq10             = _mm_mul_ps(iq1,jq0);
234
235             /* REACTION-FIELD ELECTROSTATICS */
236             velec            = _mm_mul_ps(qq10,_mm_sub_ps(_mm_add_ps(rinv10,_mm_mul_ps(krf,rsq10)),crf));
237             felec            = _mm_mul_ps(qq10,_mm_sub_ps(_mm_mul_ps(rinv10,rinvsq10),krf2));
238
239             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
240
241             /* Update potential sum for this i atom from the interaction with this j atom. */
242             velec            = _mm_and_ps(velec,cutoff_mask);
243             velecsum         = _mm_add_ps(velecsum,velec);
244
245             fscal            = felec;
246
247             fscal            = _mm_and_ps(fscal,cutoff_mask);
248
249             /* Calculate temporary vectorial force */
250             tx               = _mm_mul_ps(fscal,dx10);
251             ty               = _mm_mul_ps(fscal,dy10);
252             tz               = _mm_mul_ps(fscal,dz10);
253
254             /* Update vectorial force */
255             fix1             = _mm_add_ps(fix1,tx);
256             fiy1             = _mm_add_ps(fiy1,ty);
257             fiz1             = _mm_add_ps(fiz1,tz);
258
259             fjx0             = _mm_add_ps(fjx0,tx);
260             fjy0             = _mm_add_ps(fjy0,ty);
261             fjz0             = _mm_add_ps(fjz0,tz);
262             
263             }
264
265             /**************************
266              * CALCULATE INTERACTIONS *
267              **************************/
268
269             if (gmx_mm_any_lt(rsq20,rcutoff2))
270             {
271
272             /* Compute parameters for interactions between i and j atoms */
273             qq20             = _mm_mul_ps(iq2,jq0);
274
275             /* REACTION-FIELD ELECTROSTATICS */
276             velec            = _mm_mul_ps(qq20,_mm_sub_ps(_mm_add_ps(rinv20,_mm_mul_ps(krf,rsq20)),crf));
277             felec            = _mm_mul_ps(qq20,_mm_sub_ps(_mm_mul_ps(rinv20,rinvsq20),krf2));
278
279             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
280
281             /* Update potential sum for this i atom from the interaction with this j atom. */
282             velec            = _mm_and_ps(velec,cutoff_mask);
283             velecsum         = _mm_add_ps(velecsum,velec);
284
285             fscal            = felec;
286
287             fscal            = _mm_and_ps(fscal,cutoff_mask);
288
289             /* Calculate temporary vectorial force */
290             tx               = _mm_mul_ps(fscal,dx20);
291             ty               = _mm_mul_ps(fscal,dy20);
292             tz               = _mm_mul_ps(fscal,dz20);
293
294             /* Update vectorial force */
295             fix2             = _mm_add_ps(fix2,tx);
296             fiy2             = _mm_add_ps(fiy2,ty);
297             fiz2             = _mm_add_ps(fiz2,tz);
298
299             fjx0             = _mm_add_ps(fjx0,tx);
300             fjy0             = _mm_add_ps(fjy0,ty);
301             fjz0             = _mm_add_ps(fjz0,tz);
302             
303             }
304
305             /**************************
306              * CALCULATE INTERACTIONS *
307              **************************/
308
309             if (gmx_mm_any_lt(rsq30,rcutoff2))
310             {
311
312             /* Compute parameters for interactions between i and j atoms */
313             qq30             = _mm_mul_ps(iq3,jq0);
314
315             /* REACTION-FIELD ELECTROSTATICS */
316             velec            = _mm_mul_ps(qq30,_mm_sub_ps(_mm_add_ps(rinv30,_mm_mul_ps(krf,rsq30)),crf));
317             felec            = _mm_mul_ps(qq30,_mm_sub_ps(_mm_mul_ps(rinv30,rinvsq30),krf2));
318
319             cutoff_mask      = _mm_cmplt_ps(rsq30,rcutoff2);
320
321             /* Update potential sum for this i atom from the interaction with this j atom. */
322             velec            = _mm_and_ps(velec,cutoff_mask);
323             velecsum         = _mm_add_ps(velecsum,velec);
324
325             fscal            = felec;
326
327             fscal            = _mm_and_ps(fscal,cutoff_mask);
328
329             /* Calculate temporary vectorial force */
330             tx               = _mm_mul_ps(fscal,dx30);
331             ty               = _mm_mul_ps(fscal,dy30);
332             tz               = _mm_mul_ps(fscal,dz30);
333
334             /* Update vectorial force */
335             fix3             = _mm_add_ps(fix3,tx);
336             fiy3             = _mm_add_ps(fiy3,ty);
337             fiz3             = _mm_add_ps(fiz3,tz);
338
339             fjx0             = _mm_add_ps(fjx0,tx);
340             fjy0             = _mm_add_ps(fjy0,ty);
341             fjz0             = _mm_add_ps(fjz0,tz);
342             
343             }
344
345             fjptrA             = f+j_coord_offsetA;
346             fjptrB             = f+j_coord_offsetB;
347             fjptrC             = f+j_coord_offsetC;
348             fjptrD             = f+j_coord_offsetD;
349
350             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
351
352             /* Inner loop uses 108 flops */
353         }
354
355         if(jidx<j_index_end)
356         {
357
358             /* Get j neighbor index, and coordinate index */
359             jnrlistA         = jjnr[jidx];
360             jnrlistB         = jjnr[jidx+1];
361             jnrlistC         = jjnr[jidx+2];
362             jnrlistD         = jjnr[jidx+3];
363             /* Sign of each element will be negative for non-real atoms.
364              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
365              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
366              */
367             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
368             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
369             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
370             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
371             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
372             j_coord_offsetA  = DIM*jnrA;
373             j_coord_offsetB  = DIM*jnrB;
374             j_coord_offsetC  = DIM*jnrC;
375             j_coord_offsetD  = DIM*jnrD;
376
377             /* load j atom coordinates */
378             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
379                                               x+j_coord_offsetC,x+j_coord_offsetD,
380                                               &jx0,&jy0,&jz0);
381
382             /* Calculate displacement vector */
383             dx10             = _mm_sub_ps(ix1,jx0);
384             dy10             = _mm_sub_ps(iy1,jy0);
385             dz10             = _mm_sub_ps(iz1,jz0);
386             dx20             = _mm_sub_ps(ix2,jx0);
387             dy20             = _mm_sub_ps(iy2,jy0);
388             dz20             = _mm_sub_ps(iz2,jz0);
389             dx30             = _mm_sub_ps(ix3,jx0);
390             dy30             = _mm_sub_ps(iy3,jy0);
391             dz30             = _mm_sub_ps(iz3,jz0);
392
393             /* Calculate squared distance and things based on it */
394             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
395             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
396             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
397
398             rinv10           = gmx_mm_invsqrt_ps(rsq10);
399             rinv20           = gmx_mm_invsqrt_ps(rsq20);
400             rinv30           = gmx_mm_invsqrt_ps(rsq30);
401
402             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
403             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
404             rinvsq30         = _mm_mul_ps(rinv30,rinv30);
405
406             /* Load parameters for j particles */
407             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
408                                                               charge+jnrC+0,charge+jnrD+0);
409
410             fjx0             = _mm_setzero_ps();
411             fjy0             = _mm_setzero_ps();
412             fjz0             = _mm_setzero_ps();
413
414             /**************************
415              * CALCULATE INTERACTIONS *
416              **************************/
417
418             if (gmx_mm_any_lt(rsq10,rcutoff2))
419             {
420
421             /* Compute parameters for interactions between i and j atoms */
422             qq10             = _mm_mul_ps(iq1,jq0);
423
424             /* REACTION-FIELD ELECTROSTATICS */
425             velec            = _mm_mul_ps(qq10,_mm_sub_ps(_mm_add_ps(rinv10,_mm_mul_ps(krf,rsq10)),crf));
426             felec            = _mm_mul_ps(qq10,_mm_sub_ps(_mm_mul_ps(rinv10,rinvsq10),krf2));
427
428             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
429
430             /* Update potential sum for this i atom from the interaction with this j atom. */
431             velec            = _mm_and_ps(velec,cutoff_mask);
432             velec            = _mm_andnot_ps(dummy_mask,velec);
433             velecsum         = _mm_add_ps(velecsum,velec);
434
435             fscal            = felec;
436
437             fscal            = _mm_and_ps(fscal,cutoff_mask);
438
439             fscal            = _mm_andnot_ps(dummy_mask,fscal);
440
441             /* Calculate temporary vectorial force */
442             tx               = _mm_mul_ps(fscal,dx10);
443             ty               = _mm_mul_ps(fscal,dy10);
444             tz               = _mm_mul_ps(fscal,dz10);
445
446             /* Update vectorial force */
447             fix1             = _mm_add_ps(fix1,tx);
448             fiy1             = _mm_add_ps(fiy1,ty);
449             fiz1             = _mm_add_ps(fiz1,tz);
450
451             fjx0             = _mm_add_ps(fjx0,tx);
452             fjy0             = _mm_add_ps(fjy0,ty);
453             fjz0             = _mm_add_ps(fjz0,tz);
454             
455             }
456
457             /**************************
458              * CALCULATE INTERACTIONS *
459              **************************/
460
461             if (gmx_mm_any_lt(rsq20,rcutoff2))
462             {
463
464             /* Compute parameters for interactions between i and j atoms */
465             qq20             = _mm_mul_ps(iq2,jq0);
466
467             /* REACTION-FIELD ELECTROSTATICS */
468             velec            = _mm_mul_ps(qq20,_mm_sub_ps(_mm_add_ps(rinv20,_mm_mul_ps(krf,rsq20)),crf));
469             felec            = _mm_mul_ps(qq20,_mm_sub_ps(_mm_mul_ps(rinv20,rinvsq20),krf2));
470
471             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
472
473             /* Update potential sum for this i atom from the interaction with this j atom. */
474             velec            = _mm_and_ps(velec,cutoff_mask);
475             velec            = _mm_andnot_ps(dummy_mask,velec);
476             velecsum         = _mm_add_ps(velecsum,velec);
477
478             fscal            = felec;
479
480             fscal            = _mm_and_ps(fscal,cutoff_mask);
481
482             fscal            = _mm_andnot_ps(dummy_mask,fscal);
483
484             /* Calculate temporary vectorial force */
485             tx               = _mm_mul_ps(fscal,dx20);
486             ty               = _mm_mul_ps(fscal,dy20);
487             tz               = _mm_mul_ps(fscal,dz20);
488
489             /* Update vectorial force */
490             fix2             = _mm_add_ps(fix2,tx);
491             fiy2             = _mm_add_ps(fiy2,ty);
492             fiz2             = _mm_add_ps(fiz2,tz);
493
494             fjx0             = _mm_add_ps(fjx0,tx);
495             fjy0             = _mm_add_ps(fjy0,ty);
496             fjz0             = _mm_add_ps(fjz0,tz);
497             
498             }
499
500             /**************************
501              * CALCULATE INTERACTIONS *
502              **************************/
503
504             if (gmx_mm_any_lt(rsq30,rcutoff2))
505             {
506
507             /* Compute parameters for interactions between i and j atoms */
508             qq30             = _mm_mul_ps(iq3,jq0);
509
510             /* REACTION-FIELD ELECTROSTATICS */
511             velec            = _mm_mul_ps(qq30,_mm_sub_ps(_mm_add_ps(rinv30,_mm_mul_ps(krf,rsq30)),crf));
512             felec            = _mm_mul_ps(qq30,_mm_sub_ps(_mm_mul_ps(rinv30,rinvsq30),krf2));
513
514             cutoff_mask      = _mm_cmplt_ps(rsq30,rcutoff2);
515
516             /* Update potential sum for this i atom from the interaction with this j atom. */
517             velec            = _mm_and_ps(velec,cutoff_mask);
518             velec            = _mm_andnot_ps(dummy_mask,velec);
519             velecsum         = _mm_add_ps(velecsum,velec);
520
521             fscal            = felec;
522
523             fscal            = _mm_and_ps(fscal,cutoff_mask);
524
525             fscal            = _mm_andnot_ps(dummy_mask,fscal);
526
527             /* Calculate temporary vectorial force */
528             tx               = _mm_mul_ps(fscal,dx30);
529             ty               = _mm_mul_ps(fscal,dy30);
530             tz               = _mm_mul_ps(fscal,dz30);
531
532             /* Update vectorial force */
533             fix3             = _mm_add_ps(fix3,tx);
534             fiy3             = _mm_add_ps(fiy3,ty);
535             fiz3             = _mm_add_ps(fiz3,tz);
536
537             fjx0             = _mm_add_ps(fjx0,tx);
538             fjy0             = _mm_add_ps(fjy0,ty);
539             fjz0             = _mm_add_ps(fjz0,tz);
540             
541             }
542
543             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
544             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
545             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
546             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
547
548             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
549
550             /* Inner loop uses 108 flops */
551         }
552
553         /* End of innermost loop */
554
555         gmx_mm_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
556                                               f+i_coord_offset+DIM,fshift+i_shift_offset);
557
558         ggid                        = gid[iidx];
559         /* Update potential energies */
560         gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
561
562         /* Increment number of inner iterations */
563         inneriter                  += j_index_end - j_index_start;
564
565         /* Outer loop uses 19 flops */
566     }
567
568     /* Increment number of outer iterations */
569     outeriter        += nri;
570
571     /* Update outer/inner flops */
572
573     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_VF,outeriter*19 + inneriter*108);
574 }
575 /*
576  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwNone_GeomW4P1_F_sse2_single
577  * Electrostatics interaction: ReactionField
578  * VdW interaction:            None
579  * Geometry:                   Water4-Particle
580  * Calculate force/pot:        Force
581  */
582 void
583 nb_kernel_ElecRFCut_VdwNone_GeomW4P1_F_sse2_single
584                     (t_nblist                    * gmx_restrict       nlist,
585                      rvec                        * gmx_restrict          xx,
586                      rvec                        * gmx_restrict          ff,
587                      t_forcerec                  * gmx_restrict          fr,
588                      t_mdatoms                   * gmx_restrict     mdatoms,
589                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
590                      t_nrnb                      * gmx_restrict        nrnb)
591 {
592     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
593      * just 0 for non-waters.
594      * Suffixes A,B,C,D refer to j loop unrolling done with SSE, e.g. for the four different
595      * jnr indices corresponding to data put in the four positions in the SIMD register.
596      */
597     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
598     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
599     int              jnrA,jnrB,jnrC,jnrD;
600     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
601     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
602     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
603     real             rcutoff_scalar;
604     real             *shiftvec,*fshift,*x,*f;
605     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
606     real             scratch[4*DIM];
607     __m128           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
608     int              vdwioffset1;
609     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
610     int              vdwioffset2;
611     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
612     int              vdwioffset3;
613     __m128           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
614     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
615     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
616     __m128           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
617     __m128           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
618     __m128           dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
619     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
620     real             *charge;
621     __m128           dummy_mask,cutoff_mask;
622     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
623     __m128           one     = _mm_set1_ps(1.0);
624     __m128           two     = _mm_set1_ps(2.0);
625     x                = xx[0];
626     f                = ff[0];
627
628     nri              = nlist->nri;
629     iinr             = nlist->iinr;
630     jindex           = nlist->jindex;
631     jjnr             = nlist->jjnr;
632     shiftidx         = nlist->shift;
633     gid              = nlist->gid;
634     shiftvec         = fr->shift_vec[0];
635     fshift           = fr->fshift[0];
636     facel            = _mm_set1_ps(fr->epsfac);
637     charge           = mdatoms->chargeA;
638     krf              = _mm_set1_ps(fr->ic->k_rf);
639     krf2             = _mm_set1_ps(fr->ic->k_rf*2.0);
640     crf              = _mm_set1_ps(fr->ic->c_rf);
641
642     /* Setup water-specific parameters */
643     inr              = nlist->iinr[0];
644     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
645     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
646     iq3              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
647
648     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
649     rcutoff_scalar   = fr->rcoulomb;
650     rcutoff          = _mm_set1_ps(rcutoff_scalar);
651     rcutoff2         = _mm_mul_ps(rcutoff,rcutoff);
652
653     /* Avoid stupid compiler warnings */
654     jnrA = jnrB = jnrC = jnrD = 0;
655     j_coord_offsetA = 0;
656     j_coord_offsetB = 0;
657     j_coord_offsetC = 0;
658     j_coord_offsetD = 0;
659
660     outeriter        = 0;
661     inneriter        = 0;
662
663     for(iidx=0;iidx<4*DIM;iidx++)
664     {
665         scratch[iidx] = 0.0;
666     }  
667
668     /* Start outer loop over neighborlists */
669     for(iidx=0; iidx<nri; iidx++)
670     {
671         /* Load shift vector for this list */
672         i_shift_offset   = DIM*shiftidx[iidx];
673
674         /* Load limits for loop over neighbors */
675         j_index_start    = jindex[iidx];
676         j_index_end      = jindex[iidx+1];
677
678         /* Get outer coordinate index */
679         inr              = iinr[iidx];
680         i_coord_offset   = DIM*inr;
681
682         /* Load i particle coords and add shift vector */
683         gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
684                                                  &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
685         
686         fix1             = _mm_setzero_ps();
687         fiy1             = _mm_setzero_ps();
688         fiz1             = _mm_setzero_ps();
689         fix2             = _mm_setzero_ps();
690         fiy2             = _mm_setzero_ps();
691         fiz2             = _mm_setzero_ps();
692         fix3             = _mm_setzero_ps();
693         fiy3             = _mm_setzero_ps();
694         fiz3             = _mm_setzero_ps();
695
696         /* Start inner kernel loop */
697         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
698         {
699
700             /* Get j neighbor index, and coordinate index */
701             jnrA             = jjnr[jidx];
702             jnrB             = jjnr[jidx+1];
703             jnrC             = jjnr[jidx+2];
704             jnrD             = jjnr[jidx+3];
705             j_coord_offsetA  = DIM*jnrA;
706             j_coord_offsetB  = DIM*jnrB;
707             j_coord_offsetC  = DIM*jnrC;
708             j_coord_offsetD  = DIM*jnrD;
709
710             /* load j atom coordinates */
711             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
712                                               x+j_coord_offsetC,x+j_coord_offsetD,
713                                               &jx0,&jy0,&jz0);
714
715             /* Calculate displacement vector */
716             dx10             = _mm_sub_ps(ix1,jx0);
717             dy10             = _mm_sub_ps(iy1,jy0);
718             dz10             = _mm_sub_ps(iz1,jz0);
719             dx20             = _mm_sub_ps(ix2,jx0);
720             dy20             = _mm_sub_ps(iy2,jy0);
721             dz20             = _mm_sub_ps(iz2,jz0);
722             dx30             = _mm_sub_ps(ix3,jx0);
723             dy30             = _mm_sub_ps(iy3,jy0);
724             dz30             = _mm_sub_ps(iz3,jz0);
725
726             /* Calculate squared distance and things based on it */
727             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
728             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
729             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
730
731             rinv10           = gmx_mm_invsqrt_ps(rsq10);
732             rinv20           = gmx_mm_invsqrt_ps(rsq20);
733             rinv30           = gmx_mm_invsqrt_ps(rsq30);
734
735             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
736             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
737             rinvsq30         = _mm_mul_ps(rinv30,rinv30);
738
739             /* Load parameters for j particles */
740             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
741                                                               charge+jnrC+0,charge+jnrD+0);
742
743             fjx0             = _mm_setzero_ps();
744             fjy0             = _mm_setzero_ps();
745             fjz0             = _mm_setzero_ps();
746
747             /**************************
748              * CALCULATE INTERACTIONS *
749              **************************/
750
751             if (gmx_mm_any_lt(rsq10,rcutoff2))
752             {
753
754             /* Compute parameters for interactions between i and j atoms */
755             qq10             = _mm_mul_ps(iq1,jq0);
756
757             /* REACTION-FIELD ELECTROSTATICS */
758             felec            = _mm_mul_ps(qq10,_mm_sub_ps(_mm_mul_ps(rinv10,rinvsq10),krf2));
759
760             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
761
762             fscal            = felec;
763
764             fscal            = _mm_and_ps(fscal,cutoff_mask);
765
766             /* Calculate temporary vectorial force */
767             tx               = _mm_mul_ps(fscal,dx10);
768             ty               = _mm_mul_ps(fscal,dy10);
769             tz               = _mm_mul_ps(fscal,dz10);
770
771             /* Update vectorial force */
772             fix1             = _mm_add_ps(fix1,tx);
773             fiy1             = _mm_add_ps(fiy1,ty);
774             fiz1             = _mm_add_ps(fiz1,tz);
775
776             fjx0             = _mm_add_ps(fjx0,tx);
777             fjy0             = _mm_add_ps(fjy0,ty);
778             fjz0             = _mm_add_ps(fjz0,tz);
779             
780             }
781
782             /**************************
783              * CALCULATE INTERACTIONS *
784              **************************/
785
786             if (gmx_mm_any_lt(rsq20,rcutoff2))
787             {
788
789             /* Compute parameters for interactions between i and j atoms */
790             qq20             = _mm_mul_ps(iq2,jq0);
791
792             /* REACTION-FIELD ELECTROSTATICS */
793             felec            = _mm_mul_ps(qq20,_mm_sub_ps(_mm_mul_ps(rinv20,rinvsq20),krf2));
794
795             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
796
797             fscal            = felec;
798
799             fscal            = _mm_and_ps(fscal,cutoff_mask);
800
801             /* Calculate temporary vectorial force */
802             tx               = _mm_mul_ps(fscal,dx20);
803             ty               = _mm_mul_ps(fscal,dy20);
804             tz               = _mm_mul_ps(fscal,dz20);
805
806             /* Update vectorial force */
807             fix2             = _mm_add_ps(fix2,tx);
808             fiy2             = _mm_add_ps(fiy2,ty);
809             fiz2             = _mm_add_ps(fiz2,tz);
810
811             fjx0             = _mm_add_ps(fjx0,tx);
812             fjy0             = _mm_add_ps(fjy0,ty);
813             fjz0             = _mm_add_ps(fjz0,tz);
814             
815             }
816
817             /**************************
818              * CALCULATE INTERACTIONS *
819              **************************/
820
821             if (gmx_mm_any_lt(rsq30,rcutoff2))
822             {
823
824             /* Compute parameters for interactions between i and j atoms */
825             qq30             = _mm_mul_ps(iq3,jq0);
826
827             /* REACTION-FIELD ELECTROSTATICS */
828             felec            = _mm_mul_ps(qq30,_mm_sub_ps(_mm_mul_ps(rinv30,rinvsq30),krf2));
829
830             cutoff_mask      = _mm_cmplt_ps(rsq30,rcutoff2);
831
832             fscal            = felec;
833
834             fscal            = _mm_and_ps(fscal,cutoff_mask);
835
836             /* Calculate temporary vectorial force */
837             tx               = _mm_mul_ps(fscal,dx30);
838             ty               = _mm_mul_ps(fscal,dy30);
839             tz               = _mm_mul_ps(fscal,dz30);
840
841             /* Update vectorial force */
842             fix3             = _mm_add_ps(fix3,tx);
843             fiy3             = _mm_add_ps(fiy3,ty);
844             fiz3             = _mm_add_ps(fiz3,tz);
845
846             fjx0             = _mm_add_ps(fjx0,tx);
847             fjy0             = _mm_add_ps(fjy0,ty);
848             fjz0             = _mm_add_ps(fjz0,tz);
849             
850             }
851
852             fjptrA             = f+j_coord_offsetA;
853             fjptrB             = f+j_coord_offsetB;
854             fjptrC             = f+j_coord_offsetC;
855             fjptrD             = f+j_coord_offsetD;
856
857             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
858
859             /* Inner loop uses 90 flops */
860         }
861
862         if(jidx<j_index_end)
863         {
864
865             /* Get j neighbor index, and coordinate index */
866             jnrlistA         = jjnr[jidx];
867             jnrlistB         = jjnr[jidx+1];
868             jnrlistC         = jjnr[jidx+2];
869             jnrlistD         = jjnr[jidx+3];
870             /* Sign of each element will be negative for non-real atoms.
871              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
872              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
873              */
874             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
875             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
876             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
877             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
878             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
879             j_coord_offsetA  = DIM*jnrA;
880             j_coord_offsetB  = DIM*jnrB;
881             j_coord_offsetC  = DIM*jnrC;
882             j_coord_offsetD  = DIM*jnrD;
883
884             /* load j atom coordinates */
885             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
886                                               x+j_coord_offsetC,x+j_coord_offsetD,
887                                               &jx0,&jy0,&jz0);
888
889             /* Calculate displacement vector */
890             dx10             = _mm_sub_ps(ix1,jx0);
891             dy10             = _mm_sub_ps(iy1,jy0);
892             dz10             = _mm_sub_ps(iz1,jz0);
893             dx20             = _mm_sub_ps(ix2,jx0);
894             dy20             = _mm_sub_ps(iy2,jy0);
895             dz20             = _mm_sub_ps(iz2,jz0);
896             dx30             = _mm_sub_ps(ix3,jx0);
897             dy30             = _mm_sub_ps(iy3,jy0);
898             dz30             = _mm_sub_ps(iz3,jz0);
899
900             /* Calculate squared distance and things based on it */
901             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
902             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
903             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
904
905             rinv10           = gmx_mm_invsqrt_ps(rsq10);
906             rinv20           = gmx_mm_invsqrt_ps(rsq20);
907             rinv30           = gmx_mm_invsqrt_ps(rsq30);
908
909             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
910             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
911             rinvsq30         = _mm_mul_ps(rinv30,rinv30);
912
913             /* Load parameters for j particles */
914             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
915                                                               charge+jnrC+0,charge+jnrD+0);
916
917             fjx0             = _mm_setzero_ps();
918             fjy0             = _mm_setzero_ps();
919             fjz0             = _mm_setzero_ps();
920
921             /**************************
922              * CALCULATE INTERACTIONS *
923              **************************/
924
925             if (gmx_mm_any_lt(rsq10,rcutoff2))
926             {
927
928             /* Compute parameters for interactions between i and j atoms */
929             qq10             = _mm_mul_ps(iq1,jq0);
930
931             /* REACTION-FIELD ELECTROSTATICS */
932             felec            = _mm_mul_ps(qq10,_mm_sub_ps(_mm_mul_ps(rinv10,rinvsq10),krf2));
933
934             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
935
936             fscal            = felec;
937
938             fscal            = _mm_and_ps(fscal,cutoff_mask);
939
940             fscal            = _mm_andnot_ps(dummy_mask,fscal);
941
942             /* Calculate temporary vectorial force */
943             tx               = _mm_mul_ps(fscal,dx10);
944             ty               = _mm_mul_ps(fscal,dy10);
945             tz               = _mm_mul_ps(fscal,dz10);
946
947             /* Update vectorial force */
948             fix1             = _mm_add_ps(fix1,tx);
949             fiy1             = _mm_add_ps(fiy1,ty);
950             fiz1             = _mm_add_ps(fiz1,tz);
951
952             fjx0             = _mm_add_ps(fjx0,tx);
953             fjy0             = _mm_add_ps(fjy0,ty);
954             fjz0             = _mm_add_ps(fjz0,tz);
955             
956             }
957
958             /**************************
959              * CALCULATE INTERACTIONS *
960              **************************/
961
962             if (gmx_mm_any_lt(rsq20,rcutoff2))
963             {
964
965             /* Compute parameters for interactions between i and j atoms */
966             qq20             = _mm_mul_ps(iq2,jq0);
967
968             /* REACTION-FIELD ELECTROSTATICS */
969             felec            = _mm_mul_ps(qq20,_mm_sub_ps(_mm_mul_ps(rinv20,rinvsq20),krf2));
970
971             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
972
973             fscal            = felec;
974
975             fscal            = _mm_and_ps(fscal,cutoff_mask);
976
977             fscal            = _mm_andnot_ps(dummy_mask,fscal);
978
979             /* Calculate temporary vectorial force */
980             tx               = _mm_mul_ps(fscal,dx20);
981             ty               = _mm_mul_ps(fscal,dy20);
982             tz               = _mm_mul_ps(fscal,dz20);
983
984             /* Update vectorial force */
985             fix2             = _mm_add_ps(fix2,tx);
986             fiy2             = _mm_add_ps(fiy2,ty);
987             fiz2             = _mm_add_ps(fiz2,tz);
988
989             fjx0             = _mm_add_ps(fjx0,tx);
990             fjy0             = _mm_add_ps(fjy0,ty);
991             fjz0             = _mm_add_ps(fjz0,tz);
992             
993             }
994
995             /**************************
996              * CALCULATE INTERACTIONS *
997              **************************/
998
999             if (gmx_mm_any_lt(rsq30,rcutoff2))
1000             {
1001
1002             /* Compute parameters for interactions between i and j atoms */
1003             qq30             = _mm_mul_ps(iq3,jq0);
1004
1005             /* REACTION-FIELD ELECTROSTATICS */
1006             felec            = _mm_mul_ps(qq30,_mm_sub_ps(_mm_mul_ps(rinv30,rinvsq30),krf2));
1007
1008             cutoff_mask      = _mm_cmplt_ps(rsq30,rcutoff2);
1009
1010             fscal            = felec;
1011
1012             fscal            = _mm_and_ps(fscal,cutoff_mask);
1013
1014             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1015
1016             /* Calculate temporary vectorial force */
1017             tx               = _mm_mul_ps(fscal,dx30);
1018             ty               = _mm_mul_ps(fscal,dy30);
1019             tz               = _mm_mul_ps(fscal,dz30);
1020
1021             /* Update vectorial force */
1022             fix3             = _mm_add_ps(fix3,tx);
1023             fiy3             = _mm_add_ps(fiy3,ty);
1024             fiz3             = _mm_add_ps(fiz3,tz);
1025
1026             fjx0             = _mm_add_ps(fjx0,tx);
1027             fjy0             = _mm_add_ps(fjy0,ty);
1028             fjz0             = _mm_add_ps(fjz0,tz);
1029             
1030             }
1031
1032             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1033             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1034             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1035             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1036
1037             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1038
1039             /* Inner loop uses 90 flops */
1040         }
1041
1042         /* End of innermost loop */
1043
1044         gmx_mm_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1045                                               f+i_coord_offset+DIM,fshift+i_shift_offset);
1046
1047         /* Increment number of inner iterations */
1048         inneriter                  += j_index_end - j_index_start;
1049
1050         /* Outer loop uses 18 flops */
1051     }
1052
1053     /* Increment number of outer iterations */
1054     outeriter        += nri;
1055
1056     /* Update outer/inner flops */
1057
1058     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_F,outeriter*18 + inneriter*90);
1059 }