cf7ab9299c806439d36e1705daffc29d794bfcaa
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_single / nb_kernel_ElecRF_VdwNone_GeomW4W4_avx_128_fma_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 avx_128_fma_single kernel generator.
37  */
38 #include "config.h"
39
40 #include <math.h>
41
42 #include "../nb_kernel.h"
43 #include "gromacs/legacyheaders/types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "gromacs/legacyheaders/nrnb.h"
46
47 #include "gromacs/simd/math_x86_avx_128_fma_single.h"
48 #include "kernelutil_x86_avx_128_fma_single.h"
49
50 /*
51  * Gromacs nonbonded kernel:   nb_kernel_ElecRF_VdwNone_GeomW4W4_VF_avx_128_fma_single
52  * Electrostatics interaction: ReactionField
53  * VdW interaction:            None
54  * Geometry:                   Water4-Water4
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecRF_VdwNone_GeomW4W4_VF_avx_128_fma_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 AVX_128, 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           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              vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D;
90     __m128           jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
91     int              vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D;
92     __m128           jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
93     int              vdwjidx3A,vdwjidx3B,vdwjidx3C,vdwjidx3D;
94     __m128           jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
95     __m128           dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
96     __m128           dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
97     __m128           dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
98     __m128           dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
99     __m128           dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
100     __m128           dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
101     __m128           dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
102     __m128           dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
103     __m128           dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
104     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
105     real             *charge;
106     __m128           dummy_mask,cutoff_mask;
107     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
108     __m128           one     = _mm_set1_ps(1.0);
109     __m128           two     = _mm_set1_ps(2.0);
110     x                = xx[0];
111     f                = ff[0];
112
113     nri              = nlist->nri;
114     iinr             = nlist->iinr;
115     jindex           = nlist->jindex;
116     jjnr             = nlist->jjnr;
117     shiftidx         = nlist->shift;
118     gid              = nlist->gid;
119     shiftvec         = fr->shift_vec[0];
120     fshift           = fr->fshift[0];
121     facel            = _mm_set1_ps(fr->epsfac);
122     charge           = mdatoms->chargeA;
123     krf              = _mm_set1_ps(fr->ic->k_rf);
124     krf2             = _mm_set1_ps(fr->ic->k_rf*2.0);
125     crf              = _mm_set1_ps(fr->ic->c_rf);
126
127     /* Setup water-specific parameters */
128     inr              = nlist->iinr[0];
129     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
130     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
131     iq3              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
132
133     jq1              = _mm_set1_ps(charge[inr+1]);
134     jq2              = _mm_set1_ps(charge[inr+2]);
135     jq3              = _mm_set1_ps(charge[inr+3]);
136     qq11             = _mm_mul_ps(iq1,jq1);
137     qq12             = _mm_mul_ps(iq1,jq2);
138     qq13             = _mm_mul_ps(iq1,jq3);
139     qq21             = _mm_mul_ps(iq2,jq1);
140     qq22             = _mm_mul_ps(iq2,jq2);
141     qq23             = _mm_mul_ps(iq2,jq3);
142     qq31             = _mm_mul_ps(iq3,jq1);
143     qq32             = _mm_mul_ps(iq3,jq2);
144     qq33             = _mm_mul_ps(iq3,jq3);
145
146     /* Avoid stupid compiler warnings */
147     jnrA = jnrB = jnrC = jnrD = 0;
148     j_coord_offsetA = 0;
149     j_coord_offsetB = 0;
150     j_coord_offsetC = 0;
151     j_coord_offsetD = 0;
152
153     outeriter        = 0;
154     inneriter        = 0;
155
156     for(iidx=0;iidx<4*DIM;iidx++)
157     {
158         scratch[iidx] = 0.0;
159     }
160
161     /* Start outer loop over neighborlists */
162     for(iidx=0; iidx<nri; iidx++)
163     {
164         /* Load shift vector for this list */
165         i_shift_offset   = DIM*shiftidx[iidx];
166
167         /* Load limits for loop over neighbors */
168         j_index_start    = jindex[iidx];
169         j_index_end      = jindex[iidx+1];
170
171         /* Get outer coordinate index */
172         inr              = iinr[iidx];
173         i_coord_offset   = DIM*inr;
174
175         /* Load i particle coords and add shift vector */
176         gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
177                                                  &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
178
179         fix1             = _mm_setzero_ps();
180         fiy1             = _mm_setzero_ps();
181         fiz1             = _mm_setzero_ps();
182         fix2             = _mm_setzero_ps();
183         fiy2             = _mm_setzero_ps();
184         fiz2             = _mm_setzero_ps();
185         fix3             = _mm_setzero_ps();
186         fiy3             = _mm_setzero_ps();
187         fiz3             = _mm_setzero_ps();
188
189         /* Reset potential sums */
190         velecsum         = _mm_setzero_ps();
191
192         /* Start inner kernel loop */
193         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
194         {
195
196             /* Get j neighbor index, and coordinate index */
197             jnrA             = jjnr[jidx];
198             jnrB             = jjnr[jidx+1];
199             jnrC             = jjnr[jidx+2];
200             jnrD             = jjnr[jidx+3];
201             j_coord_offsetA  = DIM*jnrA;
202             j_coord_offsetB  = DIM*jnrB;
203             j_coord_offsetC  = DIM*jnrC;
204             j_coord_offsetD  = DIM*jnrD;
205
206             /* load j atom coordinates */
207             gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
208                                               x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
209                                               &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
210
211             /* Calculate displacement vector */
212             dx11             = _mm_sub_ps(ix1,jx1);
213             dy11             = _mm_sub_ps(iy1,jy1);
214             dz11             = _mm_sub_ps(iz1,jz1);
215             dx12             = _mm_sub_ps(ix1,jx2);
216             dy12             = _mm_sub_ps(iy1,jy2);
217             dz12             = _mm_sub_ps(iz1,jz2);
218             dx13             = _mm_sub_ps(ix1,jx3);
219             dy13             = _mm_sub_ps(iy1,jy3);
220             dz13             = _mm_sub_ps(iz1,jz3);
221             dx21             = _mm_sub_ps(ix2,jx1);
222             dy21             = _mm_sub_ps(iy2,jy1);
223             dz21             = _mm_sub_ps(iz2,jz1);
224             dx22             = _mm_sub_ps(ix2,jx2);
225             dy22             = _mm_sub_ps(iy2,jy2);
226             dz22             = _mm_sub_ps(iz2,jz2);
227             dx23             = _mm_sub_ps(ix2,jx3);
228             dy23             = _mm_sub_ps(iy2,jy3);
229             dz23             = _mm_sub_ps(iz2,jz3);
230             dx31             = _mm_sub_ps(ix3,jx1);
231             dy31             = _mm_sub_ps(iy3,jy1);
232             dz31             = _mm_sub_ps(iz3,jz1);
233             dx32             = _mm_sub_ps(ix3,jx2);
234             dy32             = _mm_sub_ps(iy3,jy2);
235             dz32             = _mm_sub_ps(iz3,jz2);
236             dx33             = _mm_sub_ps(ix3,jx3);
237             dy33             = _mm_sub_ps(iy3,jy3);
238             dz33             = _mm_sub_ps(iz3,jz3);
239
240             /* Calculate squared distance and things based on it */
241             rsq11            = gmx_mm_calc_rsq_ps(dx11,dy11,dz11);
242             rsq12            = gmx_mm_calc_rsq_ps(dx12,dy12,dz12);
243             rsq13            = gmx_mm_calc_rsq_ps(dx13,dy13,dz13);
244             rsq21            = gmx_mm_calc_rsq_ps(dx21,dy21,dz21);
245             rsq22            = gmx_mm_calc_rsq_ps(dx22,dy22,dz22);
246             rsq23            = gmx_mm_calc_rsq_ps(dx23,dy23,dz23);
247             rsq31            = gmx_mm_calc_rsq_ps(dx31,dy31,dz31);
248             rsq32            = gmx_mm_calc_rsq_ps(dx32,dy32,dz32);
249             rsq33            = gmx_mm_calc_rsq_ps(dx33,dy33,dz33);
250
251             rinv11           = gmx_mm_invsqrt_ps(rsq11);
252             rinv12           = gmx_mm_invsqrt_ps(rsq12);
253             rinv13           = gmx_mm_invsqrt_ps(rsq13);
254             rinv21           = gmx_mm_invsqrt_ps(rsq21);
255             rinv22           = gmx_mm_invsqrt_ps(rsq22);
256             rinv23           = gmx_mm_invsqrt_ps(rsq23);
257             rinv31           = gmx_mm_invsqrt_ps(rsq31);
258             rinv32           = gmx_mm_invsqrt_ps(rsq32);
259             rinv33           = gmx_mm_invsqrt_ps(rsq33);
260
261             rinvsq11         = _mm_mul_ps(rinv11,rinv11);
262             rinvsq12         = _mm_mul_ps(rinv12,rinv12);
263             rinvsq13         = _mm_mul_ps(rinv13,rinv13);
264             rinvsq21         = _mm_mul_ps(rinv21,rinv21);
265             rinvsq22         = _mm_mul_ps(rinv22,rinv22);
266             rinvsq23         = _mm_mul_ps(rinv23,rinv23);
267             rinvsq31         = _mm_mul_ps(rinv31,rinv31);
268             rinvsq32         = _mm_mul_ps(rinv32,rinv32);
269             rinvsq33         = _mm_mul_ps(rinv33,rinv33);
270
271             fjx1             = _mm_setzero_ps();
272             fjy1             = _mm_setzero_ps();
273             fjz1             = _mm_setzero_ps();
274             fjx2             = _mm_setzero_ps();
275             fjy2             = _mm_setzero_ps();
276             fjz2             = _mm_setzero_ps();
277             fjx3             = _mm_setzero_ps();
278             fjy3             = _mm_setzero_ps();
279             fjz3             = _mm_setzero_ps();
280
281             /**************************
282              * CALCULATE INTERACTIONS *
283              **************************/
284
285             /* REACTION-FIELD ELECTROSTATICS */
286             velec            = _mm_mul_ps(qq11,_mm_sub_ps(_mm_macc_ps(krf,rsq11,rinv11),crf));
287             felec            = _mm_mul_ps(qq11,_mm_msub_ps(rinv11,rinvsq11,krf2));
288
289             /* Update potential sum for this i atom from the interaction with this j atom. */
290             velecsum         = _mm_add_ps(velecsum,velec);
291
292             fscal            = felec;
293
294              /* Update vectorial force */
295             fix1             = _mm_macc_ps(dx11,fscal,fix1);
296             fiy1             = _mm_macc_ps(dy11,fscal,fiy1);
297             fiz1             = _mm_macc_ps(dz11,fscal,fiz1);
298
299             fjx1             = _mm_macc_ps(dx11,fscal,fjx1);
300             fjy1             = _mm_macc_ps(dy11,fscal,fjy1);
301             fjz1             = _mm_macc_ps(dz11,fscal,fjz1);
302
303             /**************************
304              * CALCULATE INTERACTIONS *
305              **************************/
306
307             /* REACTION-FIELD ELECTROSTATICS */
308             velec            = _mm_mul_ps(qq12,_mm_sub_ps(_mm_macc_ps(krf,rsq12,rinv12),crf));
309             felec            = _mm_mul_ps(qq12,_mm_msub_ps(rinv12,rinvsq12,krf2));
310
311             /* Update potential sum for this i atom from the interaction with this j atom. */
312             velecsum         = _mm_add_ps(velecsum,velec);
313
314             fscal            = felec;
315
316              /* Update vectorial force */
317             fix1             = _mm_macc_ps(dx12,fscal,fix1);
318             fiy1             = _mm_macc_ps(dy12,fscal,fiy1);
319             fiz1             = _mm_macc_ps(dz12,fscal,fiz1);
320
321             fjx2             = _mm_macc_ps(dx12,fscal,fjx2);
322             fjy2             = _mm_macc_ps(dy12,fscal,fjy2);
323             fjz2             = _mm_macc_ps(dz12,fscal,fjz2);
324
325             /**************************
326              * CALCULATE INTERACTIONS *
327              **************************/
328
329             /* REACTION-FIELD ELECTROSTATICS */
330             velec            = _mm_mul_ps(qq13,_mm_sub_ps(_mm_macc_ps(krf,rsq13,rinv13),crf));
331             felec            = _mm_mul_ps(qq13,_mm_msub_ps(rinv13,rinvsq13,krf2));
332
333             /* Update potential sum for this i atom from the interaction with this j atom. */
334             velecsum         = _mm_add_ps(velecsum,velec);
335
336             fscal            = felec;
337
338              /* Update vectorial force */
339             fix1             = _mm_macc_ps(dx13,fscal,fix1);
340             fiy1             = _mm_macc_ps(dy13,fscal,fiy1);
341             fiz1             = _mm_macc_ps(dz13,fscal,fiz1);
342
343             fjx3             = _mm_macc_ps(dx13,fscal,fjx3);
344             fjy3             = _mm_macc_ps(dy13,fscal,fjy3);
345             fjz3             = _mm_macc_ps(dz13,fscal,fjz3);
346
347             /**************************
348              * CALCULATE INTERACTIONS *
349              **************************/
350
351             /* REACTION-FIELD ELECTROSTATICS */
352             velec            = _mm_mul_ps(qq21,_mm_sub_ps(_mm_macc_ps(krf,rsq21,rinv21),crf));
353             felec            = _mm_mul_ps(qq21,_mm_msub_ps(rinv21,rinvsq21,krf2));
354
355             /* Update potential sum for this i atom from the interaction with this j atom. */
356             velecsum         = _mm_add_ps(velecsum,velec);
357
358             fscal            = felec;
359
360              /* Update vectorial force */
361             fix2             = _mm_macc_ps(dx21,fscal,fix2);
362             fiy2             = _mm_macc_ps(dy21,fscal,fiy2);
363             fiz2             = _mm_macc_ps(dz21,fscal,fiz2);
364
365             fjx1             = _mm_macc_ps(dx21,fscal,fjx1);
366             fjy1             = _mm_macc_ps(dy21,fscal,fjy1);
367             fjz1             = _mm_macc_ps(dz21,fscal,fjz1);
368
369             /**************************
370              * CALCULATE INTERACTIONS *
371              **************************/
372
373             /* REACTION-FIELD ELECTROSTATICS */
374             velec            = _mm_mul_ps(qq22,_mm_sub_ps(_mm_macc_ps(krf,rsq22,rinv22),crf));
375             felec            = _mm_mul_ps(qq22,_mm_msub_ps(rinv22,rinvsq22,krf2));
376
377             /* Update potential sum for this i atom from the interaction with this j atom. */
378             velecsum         = _mm_add_ps(velecsum,velec);
379
380             fscal            = felec;
381
382              /* Update vectorial force */
383             fix2             = _mm_macc_ps(dx22,fscal,fix2);
384             fiy2             = _mm_macc_ps(dy22,fscal,fiy2);
385             fiz2             = _mm_macc_ps(dz22,fscal,fiz2);
386
387             fjx2             = _mm_macc_ps(dx22,fscal,fjx2);
388             fjy2             = _mm_macc_ps(dy22,fscal,fjy2);
389             fjz2             = _mm_macc_ps(dz22,fscal,fjz2);
390
391             /**************************
392              * CALCULATE INTERACTIONS *
393              **************************/
394
395             /* REACTION-FIELD ELECTROSTATICS */
396             velec            = _mm_mul_ps(qq23,_mm_sub_ps(_mm_macc_ps(krf,rsq23,rinv23),crf));
397             felec            = _mm_mul_ps(qq23,_mm_msub_ps(rinv23,rinvsq23,krf2));
398
399             /* Update potential sum for this i atom from the interaction with this j atom. */
400             velecsum         = _mm_add_ps(velecsum,velec);
401
402             fscal            = felec;
403
404              /* Update vectorial force */
405             fix2             = _mm_macc_ps(dx23,fscal,fix2);
406             fiy2             = _mm_macc_ps(dy23,fscal,fiy2);
407             fiz2             = _mm_macc_ps(dz23,fscal,fiz2);
408
409             fjx3             = _mm_macc_ps(dx23,fscal,fjx3);
410             fjy3             = _mm_macc_ps(dy23,fscal,fjy3);
411             fjz3             = _mm_macc_ps(dz23,fscal,fjz3);
412
413             /**************************
414              * CALCULATE INTERACTIONS *
415              **************************/
416
417             /* REACTION-FIELD ELECTROSTATICS */
418             velec            = _mm_mul_ps(qq31,_mm_sub_ps(_mm_macc_ps(krf,rsq31,rinv31),crf));
419             felec            = _mm_mul_ps(qq31,_mm_msub_ps(rinv31,rinvsq31,krf2));
420
421             /* Update potential sum for this i atom from the interaction with this j atom. */
422             velecsum         = _mm_add_ps(velecsum,velec);
423
424             fscal            = felec;
425
426              /* Update vectorial force */
427             fix3             = _mm_macc_ps(dx31,fscal,fix3);
428             fiy3             = _mm_macc_ps(dy31,fscal,fiy3);
429             fiz3             = _mm_macc_ps(dz31,fscal,fiz3);
430
431             fjx1             = _mm_macc_ps(dx31,fscal,fjx1);
432             fjy1             = _mm_macc_ps(dy31,fscal,fjy1);
433             fjz1             = _mm_macc_ps(dz31,fscal,fjz1);
434
435             /**************************
436              * CALCULATE INTERACTIONS *
437              **************************/
438
439             /* REACTION-FIELD ELECTROSTATICS */
440             velec            = _mm_mul_ps(qq32,_mm_sub_ps(_mm_macc_ps(krf,rsq32,rinv32),crf));
441             felec            = _mm_mul_ps(qq32,_mm_msub_ps(rinv32,rinvsq32,krf2));
442
443             /* Update potential sum for this i atom from the interaction with this j atom. */
444             velecsum         = _mm_add_ps(velecsum,velec);
445
446             fscal            = felec;
447
448              /* Update vectorial force */
449             fix3             = _mm_macc_ps(dx32,fscal,fix3);
450             fiy3             = _mm_macc_ps(dy32,fscal,fiy3);
451             fiz3             = _mm_macc_ps(dz32,fscal,fiz3);
452
453             fjx2             = _mm_macc_ps(dx32,fscal,fjx2);
454             fjy2             = _mm_macc_ps(dy32,fscal,fjy2);
455             fjz2             = _mm_macc_ps(dz32,fscal,fjz2);
456
457             /**************************
458              * CALCULATE INTERACTIONS *
459              **************************/
460
461             /* REACTION-FIELD ELECTROSTATICS */
462             velec            = _mm_mul_ps(qq33,_mm_sub_ps(_mm_macc_ps(krf,rsq33,rinv33),crf));
463             felec            = _mm_mul_ps(qq33,_mm_msub_ps(rinv33,rinvsq33,krf2));
464
465             /* Update potential sum for this i atom from the interaction with this j atom. */
466             velecsum         = _mm_add_ps(velecsum,velec);
467
468             fscal            = felec;
469
470              /* Update vectorial force */
471             fix3             = _mm_macc_ps(dx33,fscal,fix3);
472             fiy3             = _mm_macc_ps(dy33,fscal,fiy3);
473             fiz3             = _mm_macc_ps(dz33,fscal,fiz3);
474
475             fjx3             = _mm_macc_ps(dx33,fscal,fjx3);
476             fjy3             = _mm_macc_ps(dy33,fscal,fjy3);
477             fjz3             = _mm_macc_ps(dz33,fscal,fjz3);
478
479             fjptrA             = f+j_coord_offsetA;
480             fjptrB             = f+j_coord_offsetB;
481             fjptrC             = f+j_coord_offsetC;
482             fjptrD             = f+j_coord_offsetD;
483
484             gmx_mm_decrement_3rvec_4ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
485                                                    fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
486
487             /* Inner loop uses 315 flops */
488         }
489
490         if(jidx<j_index_end)
491         {
492
493             /* Get j neighbor index, and coordinate index */
494             jnrlistA         = jjnr[jidx];
495             jnrlistB         = jjnr[jidx+1];
496             jnrlistC         = jjnr[jidx+2];
497             jnrlistD         = jjnr[jidx+3];
498             /* Sign of each element will be negative for non-real atoms.
499              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
500              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
501              */
502             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
503             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
504             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
505             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
506             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
507             j_coord_offsetA  = DIM*jnrA;
508             j_coord_offsetB  = DIM*jnrB;
509             j_coord_offsetC  = DIM*jnrC;
510             j_coord_offsetD  = DIM*jnrD;
511
512             /* load j atom coordinates */
513             gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
514                                               x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
515                                               &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
516
517             /* Calculate displacement vector */
518             dx11             = _mm_sub_ps(ix1,jx1);
519             dy11             = _mm_sub_ps(iy1,jy1);
520             dz11             = _mm_sub_ps(iz1,jz1);
521             dx12             = _mm_sub_ps(ix1,jx2);
522             dy12             = _mm_sub_ps(iy1,jy2);
523             dz12             = _mm_sub_ps(iz1,jz2);
524             dx13             = _mm_sub_ps(ix1,jx3);
525             dy13             = _mm_sub_ps(iy1,jy3);
526             dz13             = _mm_sub_ps(iz1,jz3);
527             dx21             = _mm_sub_ps(ix2,jx1);
528             dy21             = _mm_sub_ps(iy2,jy1);
529             dz21             = _mm_sub_ps(iz2,jz1);
530             dx22             = _mm_sub_ps(ix2,jx2);
531             dy22             = _mm_sub_ps(iy2,jy2);
532             dz22             = _mm_sub_ps(iz2,jz2);
533             dx23             = _mm_sub_ps(ix2,jx3);
534             dy23             = _mm_sub_ps(iy2,jy3);
535             dz23             = _mm_sub_ps(iz2,jz3);
536             dx31             = _mm_sub_ps(ix3,jx1);
537             dy31             = _mm_sub_ps(iy3,jy1);
538             dz31             = _mm_sub_ps(iz3,jz1);
539             dx32             = _mm_sub_ps(ix3,jx2);
540             dy32             = _mm_sub_ps(iy3,jy2);
541             dz32             = _mm_sub_ps(iz3,jz2);
542             dx33             = _mm_sub_ps(ix3,jx3);
543             dy33             = _mm_sub_ps(iy3,jy3);
544             dz33             = _mm_sub_ps(iz3,jz3);
545
546             /* Calculate squared distance and things based on it */
547             rsq11            = gmx_mm_calc_rsq_ps(dx11,dy11,dz11);
548             rsq12            = gmx_mm_calc_rsq_ps(dx12,dy12,dz12);
549             rsq13            = gmx_mm_calc_rsq_ps(dx13,dy13,dz13);
550             rsq21            = gmx_mm_calc_rsq_ps(dx21,dy21,dz21);
551             rsq22            = gmx_mm_calc_rsq_ps(dx22,dy22,dz22);
552             rsq23            = gmx_mm_calc_rsq_ps(dx23,dy23,dz23);
553             rsq31            = gmx_mm_calc_rsq_ps(dx31,dy31,dz31);
554             rsq32            = gmx_mm_calc_rsq_ps(dx32,dy32,dz32);
555             rsq33            = gmx_mm_calc_rsq_ps(dx33,dy33,dz33);
556
557             rinv11           = gmx_mm_invsqrt_ps(rsq11);
558             rinv12           = gmx_mm_invsqrt_ps(rsq12);
559             rinv13           = gmx_mm_invsqrt_ps(rsq13);
560             rinv21           = gmx_mm_invsqrt_ps(rsq21);
561             rinv22           = gmx_mm_invsqrt_ps(rsq22);
562             rinv23           = gmx_mm_invsqrt_ps(rsq23);
563             rinv31           = gmx_mm_invsqrt_ps(rsq31);
564             rinv32           = gmx_mm_invsqrt_ps(rsq32);
565             rinv33           = gmx_mm_invsqrt_ps(rsq33);
566
567             rinvsq11         = _mm_mul_ps(rinv11,rinv11);
568             rinvsq12         = _mm_mul_ps(rinv12,rinv12);
569             rinvsq13         = _mm_mul_ps(rinv13,rinv13);
570             rinvsq21         = _mm_mul_ps(rinv21,rinv21);
571             rinvsq22         = _mm_mul_ps(rinv22,rinv22);
572             rinvsq23         = _mm_mul_ps(rinv23,rinv23);
573             rinvsq31         = _mm_mul_ps(rinv31,rinv31);
574             rinvsq32         = _mm_mul_ps(rinv32,rinv32);
575             rinvsq33         = _mm_mul_ps(rinv33,rinv33);
576
577             fjx1             = _mm_setzero_ps();
578             fjy1             = _mm_setzero_ps();
579             fjz1             = _mm_setzero_ps();
580             fjx2             = _mm_setzero_ps();
581             fjy2             = _mm_setzero_ps();
582             fjz2             = _mm_setzero_ps();
583             fjx3             = _mm_setzero_ps();
584             fjy3             = _mm_setzero_ps();
585             fjz3             = _mm_setzero_ps();
586
587             /**************************
588              * CALCULATE INTERACTIONS *
589              **************************/
590
591             /* REACTION-FIELD ELECTROSTATICS */
592             velec            = _mm_mul_ps(qq11,_mm_sub_ps(_mm_macc_ps(krf,rsq11,rinv11),crf));
593             felec            = _mm_mul_ps(qq11,_mm_msub_ps(rinv11,rinvsq11,krf2));
594
595             /* Update potential sum for this i atom from the interaction with this j atom. */
596             velec            = _mm_andnot_ps(dummy_mask,velec);
597             velecsum         = _mm_add_ps(velecsum,velec);
598
599             fscal            = felec;
600
601             fscal            = _mm_andnot_ps(dummy_mask,fscal);
602
603              /* Update vectorial force */
604             fix1             = _mm_macc_ps(dx11,fscal,fix1);
605             fiy1             = _mm_macc_ps(dy11,fscal,fiy1);
606             fiz1             = _mm_macc_ps(dz11,fscal,fiz1);
607
608             fjx1             = _mm_macc_ps(dx11,fscal,fjx1);
609             fjy1             = _mm_macc_ps(dy11,fscal,fjy1);
610             fjz1             = _mm_macc_ps(dz11,fscal,fjz1);
611
612             /**************************
613              * CALCULATE INTERACTIONS *
614              **************************/
615
616             /* REACTION-FIELD ELECTROSTATICS */
617             velec            = _mm_mul_ps(qq12,_mm_sub_ps(_mm_macc_ps(krf,rsq12,rinv12),crf));
618             felec            = _mm_mul_ps(qq12,_mm_msub_ps(rinv12,rinvsq12,krf2));
619
620             /* Update potential sum for this i atom from the interaction with this j atom. */
621             velec            = _mm_andnot_ps(dummy_mask,velec);
622             velecsum         = _mm_add_ps(velecsum,velec);
623
624             fscal            = felec;
625
626             fscal            = _mm_andnot_ps(dummy_mask,fscal);
627
628              /* Update vectorial force */
629             fix1             = _mm_macc_ps(dx12,fscal,fix1);
630             fiy1             = _mm_macc_ps(dy12,fscal,fiy1);
631             fiz1             = _mm_macc_ps(dz12,fscal,fiz1);
632
633             fjx2             = _mm_macc_ps(dx12,fscal,fjx2);
634             fjy2             = _mm_macc_ps(dy12,fscal,fjy2);
635             fjz2             = _mm_macc_ps(dz12,fscal,fjz2);
636
637             /**************************
638              * CALCULATE INTERACTIONS *
639              **************************/
640
641             /* REACTION-FIELD ELECTROSTATICS */
642             velec            = _mm_mul_ps(qq13,_mm_sub_ps(_mm_macc_ps(krf,rsq13,rinv13),crf));
643             felec            = _mm_mul_ps(qq13,_mm_msub_ps(rinv13,rinvsq13,krf2));
644
645             /* Update potential sum for this i atom from the interaction with this j atom. */
646             velec            = _mm_andnot_ps(dummy_mask,velec);
647             velecsum         = _mm_add_ps(velecsum,velec);
648
649             fscal            = felec;
650
651             fscal            = _mm_andnot_ps(dummy_mask,fscal);
652
653              /* Update vectorial force */
654             fix1             = _mm_macc_ps(dx13,fscal,fix1);
655             fiy1             = _mm_macc_ps(dy13,fscal,fiy1);
656             fiz1             = _mm_macc_ps(dz13,fscal,fiz1);
657
658             fjx3             = _mm_macc_ps(dx13,fscal,fjx3);
659             fjy3             = _mm_macc_ps(dy13,fscal,fjy3);
660             fjz3             = _mm_macc_ps(dz13,fscal,fjz3);
661
662             /**************************
663              * CALCULATE INTERACTIONS *
664              **************************/
665
666             /* REACTION-FIELD ELECTROSTATICS */
667             velec            = _mm_mul_ps(qq21,_mm_sub_ps(_mm_macc_ps(krf,rsq21,rinv21),crf));
668             felec            = _mm_mul_ps(qq21,_mm_msub_ps(rinv21,rinvsq21,krf2));
669
670             /* Update potential sum for this i atom from the interaction with this j atom. */
671             velec            = _mm_andnot_ps(dummy_mask,velec);
672             velecsum         = _mm_add_ps(velecsum,velec);
673
674             fscal            = felec;
675
676             fscal            = _mm_andnot_ps(dummy_mask,fscal);
677
678              /* Update vectorial force */
679             fix2             = _mm_macc_ps(dx21,fscal,fix2);
680             fiy2             = _mm_macc_ps(dy21,fscal,fiy2);
681             fiz2             = _mm_macc_ps(dz21,fscal,fiz2);
682
683             fjx1             = _mm_macc_ps(dx21,fscal,fjx1);
684             fjy1             = _mm_macc_ps(dy21,fscal,fjy1);
685             fjz1             = _mm_macc_ps(dz21,fscal,fjz1);
686
687             /**************************
688              * CALCULATE INTERACTIONS *
689              **************************/
690
691             /* REACTION-FIELD ELECTROSTATICS */
692             velec            = _mm_mul_ps(qq22,_mm_sub_ps(_mm_macc_ps(krf,rsq22,rinv22),crf));
693             felec            = _mm_mul_ps(qq22,_mm_msub_ps(rinv22,rinvsq22,krf2));
694
695             /* Update potential sum for this i atom from the interaction with this j atom. */
696             velec            = _mm_andnot_ps(dummy_mask,velec);
697             velecsum         = _mm_add_ps(velecsum,velec);
698
699             fscal            = felec;
700
701             fscal            = _mm_andnot_ps(dummy_mask,fscal);
702
703              /* Update vectorial force */
704             fix2             = _mm_macc_ps(dx22,fscal,fix2);
705             fiy2             = _mm_macc_ps(dy22,fscal,fiy2);
706             fiz2             = _mm_macc_ps(dz22,fscal,fiz2);
707
708             fjx2             = _mm_macc_ps(dx22,fscal,fjx2);
709             fjy2             = _mm_macc_ps(dy22,fscal,fjy2);
710             fjz2             = _mm_macc_ps(dz22,fscal,fjz2);
711
712             /**************************
713              * CALCULATE INTERACTIONS *
714              **************************/
715
716             /* REACTION-FIELD ELECTROSTATICS */
717             velec            = _mm_mul_ps(qq23,_mm_sub_ps(_mm_macc_ps(krf,rsq23,rinv23),crf));
718             felec            = _mm_mul_ps(qq23,_mm_msub_ps(rinv23,rinvsq23,krf2));
719
720             /* Update potential sum for this i atom from the interaction with this j atom. */
721             velec            = _mm_andnot_ps(dummy_mask,velec);
722             velecsum         = _mm_add_ps(velecsum,velec);
723
724             fscal            = felec;
725
726             fscal            = _mm_andnot_ps(dummy_mask,fscal);
727
728              /* Update vectorial force */
729             fix2             = _mm_macc_ps(dx23,fscal,fix2);
730             fiy2             = _mm_macc_ps(dy23,fscal,fiy2);
731             fiz2             = _mm_macc_ps(dz23,fscal,fiz2);
732
733             fjx3             = _mm_macc_ps(dx23,fscal,fjx3);
734             fjy3             = _mm_macc_ps(dy23,fscal,fjy3);
735             fjz3             = _mm_macc_ps(dz23,fscal,fjz3);
736
737             /**************************
738              * CALCULATE INTERACTIONS *
739              **************************/
740
741             /* REACTION-FIELD ELECTROSTATICS */
742             velec            = _mm_mul_ps(qq31,_mm_sub_ps(_mm_macc_ps(krf,rsq31,rinv31),crf));
743             felec            = _mm_mul_ps(qq31,_mm_msub_ps(rinv31,rinvsq31,krf2));
744
745             /* Update potential sum for this i atom from the interaction with this j atom. */
746             velec            = _mm_andnot_ps(dummy_mask,velec);
747             velecsum         = _mm_add_ps(velecsum,velec);
748
749             fscal            = felec;
750
751             fscal            = _mm_andnot_ps(dummy_mask,fscal);
752
753              /* Update vectorial force */
754             fix3             = _mm_macc_ps(dx31,fscal,fix3);
755             fiy3             = _mm_macc_ps(dy31,fscal,fiy3);
756             fiz3             = _mm_macc_ps(dz31,fscal,fiz3);
757
758             fjx1             = _mm_macc_ps(dx31,fscal,fjx1);
759             fjy1             = _mm_macc_ps(dy31,fscal,fjy1);
760             fjz1             = _mm_macc_ps(dz31,fscal,fjz1);
761
762             /**************************
763              * CALCULATE INTERACTIONS *
764              **************************/
765
766             /* REACTION-FIELD ELECTROSTATICS */
767             velec            = _mm_mul_ps(qq32,_mm_sub_ps(_mm_macc_ps(krf,rsq32,rinv32),crf));
768             felec            = _mm_mul_ps(qq32,_mm_msub_ps(rinv32,rinvsq32,krf2));
769
770             /* Update potential sum for this i atom from the interaction with this j atom. */
771             velec            = _mm_andnot_ps(dummy_mask,velec);
772             velecsum         = _mm_add_ps(velecsum,velec);
773
774             fscal            = felec;
775
776             fscal            = _mm_andnot_ps(dummy_mask,fscal);
777
778              /* Update vectorial force */
779             fix3             = _mm_macc_ps(dx32,fscal,fix3);
780             fiy3             = _mm_macc_ps(dy32,fscal,fiy3);
781             fiz3             = _mm_macc_ps(dz32,fscal,fiz3);
782
783             fjx2             = _mm_macc_ps(dx32,fscal,fjx2);
784             fjy2             = _mm_macc_ps(dy32,fscal,fjy2);
785             fjz2             = _mm_macc_ps(dz32,fscal,fjz2);
786
787             /**************************
788              * CALCULATE INTERACTIONS *
789              **************************/
790
791             /* REACTION-FIELD ELECTROSTATICS */
792             velec            = _mm_mul_ps(qq33,_mm_sub_ps(_mm_macc_ps(krf,rsq33,rinv33),crf));
793             felec            = _mm_mul_ps(qq33,_mm_msub_ps(rinv33,rinvsq33,krf2));
794
795             /* Update potential sum for this i atom from the interaction with this j atom. */
796             velec            = _mm_andnot_ps(dummy_mask,velec);
797             velecsum         = _mm_add_ps(velecsum,velec);
798
799             fscal            = felec;
800
801             fscal            = _mm_andnot_ps(dummy_mask,fscal);
802
803              /* Update vectorial force */
804             fix3             = _mm_macc_ps(dx33,fscal,fix3);
805             fiy3             = _mm_macc_ps(dy33,fscal,fiy3);
806             fiz3             = _mm_macc_ps(dz33,fscal,fiz3);
807
808             fjx3             = _mm_macc_ps(dx33,fscal,fjx3);
809             fjy3             = _mm_macc_ps(dy33,fscal,fjy3);
810             fjz3             = _mm_macc_ps(dz33,fscal,fjz3);
811
812             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
813             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
814             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
815             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
816
817             gmx_mm_decrement_3rvec_4ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
818                                                    fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
819
820             /* Inner loop uses 315 flops */
821         }
822
823         /* End of innermost loop */
824
825         gmx_mm_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
826                                               f+i_coord_offset+DIM,fshift+i_shift_offset);
827
828         ggid                        = gid[iidx];
829         /* Update potential energies */
830         gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
831
832         /* Increment number of inner iterations */
833         inneriter                  += j_index_end - j_index_start;
834
835         /* Outer loop uses 19 flops */
836     }
837
838     /* Increment number of outer iterations */
839     outeriter        += nri;
840
841     /* Update outer/inner flops */
842
843     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4W4_VF,outeriter*19 + inneriter*315);
844 }
845 /*
846  * Gromacs nonbonded kernel:   nb_kernel_ElecRF_VdwNone_GeomW4W4_F_avx_128_fma_single
847  * Electrostatics interaction: ReactionField
848  * VdW interaction:            None
849  * Geometry:                   Water4-Water4
850  * Calculate force/pot:        Force
851  */
852 void
853 nb_kernel_ElecRF_VdwNone_GeomW4W4_F_avx_128_fma_single
854                     (t_nblist                    * gmx_restrict       nlist,
855                      rvec                        * gmx_restrict          xx,
856                      rvec                        * gmx_restrict          ff,
857                      t_forcerec                  * gmx_restrict          fr,
858                      t_mdatoms                   * gmx_restrict     mdatoms,
859                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
860                      t_nrnb                      * gmx_restrict        nrnb)
861 {
862     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
863      * just 0 for non-waters.
864      * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
865      * jnr indices corresponding to data put in the four positions in the SIMD register.
866      */
867     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
868     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
869     int              jnrA,jnrB,jnrC,jnrD;
870     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
871     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
872     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
873     real             rcutoff_scalar;
874     real             *shiftvec,*fshift,*x,*f;
875     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
876     real             scratch[4*DIM];
877     __m128           fscal,rcutoff,rcutoff2,jidxall;
878     int              vdwioffset1;
879     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
880     int              vdwioffset2;
881     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
882     int              vdwioffset3;
883     __m128           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
884     int              vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D;
885     __m128           jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
886     int              vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D;
887     __m128           jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
888     int              vdwjidx3A,vdwjidx3B,vdwjidx3C,vdwjidx3D;
889     __m128           jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
890     __m128           dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
891     __m128           dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
892     __m128           dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
893     __m128           dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
894     __m128           dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
895     __m128           dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
896     __m128           dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
897     __m128           dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
898     __m128           dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
899     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
900     real             *charge;
901     __m128           dummy_mask,cutoff_mask;
902     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
903     __m128           one     = _mm_set1_ps(1.0);
904     __m128           two     = _mm_set1_ps(2.0);
905     x                = xx[0];
906     f                = ff[0];
907
908     nri              = nlist->nri;
909     iinr             = nlist->iinr;
910     jindex           = nlist->jindex;
911     jjnr             = nlist->jjnr;
912     shiftidx         = nlist->shift;
913     gid              = nlist->gid;
914     shiftvec         = fr->shift_vec[0];
915     fshift           = fr->fshift[0];
916     facel            = _mm_set1_ps(fr->epsfac);
917     charge           = mdatoms->chargeA;
918     krf              = _mm_set1_ps(fr->ic->k_rf);
919     krf2             = _mm_set1_ps(fr->ic->k_rf*2.0);
920     crf              = _mm_set1_ps(fr->ic->c_rf);
921
922     /* Setup water-specific parameters */
923     inr              = nlist->iinr[0];
924     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
925     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
926     iq3              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
927
928     jq1              = _mm_set1_ps(charge[inr+1]);
929     jq2              = _mm_set1_ps(charge[inr+2]);
930     jq3              = _mm_set1_ps(charge[inr+3]);
931     qq11             = _mm_mul_ps(iq1,jq1);
932     qq12             = _mm_mul_ps(iq1,jq2);
933     qq13             = _mm_mul_ps(iq1,jq3);
934     qq21             = _mm_mul_ps(iq2,jq1);
935     qq22             = _mm_mul_ps(iq2,jq2);
936     qq23             = _mm_mul_ps(iq2,jq3);
937     qq31             = _mm_mul_ps(iq3,jq1);
938     qq32             = _mm_mul_ps(iq3,jq2);
939     qq33             = _mm_mul_ps(iq3,jq3);
940
941     /* Avoid stupid compiler warnings */
942     jnrA = jnrB = jnrC = jnrD = 0;
943     j_coord_offsetA = 0;
944     j_coord_offsetB = 0;
945     j_coord_offsetC = 0;
946     j_coord_offsetD = 0;
947
948     outeriter        = 0;
949     inneriter        = 0;
950
951     for(iidx=0;iidx<4*DIM;iidx++)
952     {
953         scratch[iidx] = 0.0;
954     }
955
956     /* Start outer loop over neighborlists */
957     for(iidx=0; iidx<nri; iidx++)
958     {
959         /* Load shift vector for this list */
960         i_shift_offset   = DIM*shiftidx[iidx];
961
962         /* Load limits for loop over neighbors */
963         j_index_start    = jindex[iidx];
964         j_index_end      = jindex[iidx+1];
965
966         /* Get outer coordinate index */
967         inr              = iinr[iidx];
968         i_coord_offset   = DIM*inr;
969
970         /* Load i particle coords and add shift vector */
971         gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
972                                                  &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
973
974         fix1             = _mm_setzero_ps();
975         fiy1             = _mm_setzero_ps();
976         fiz1             = _mm_setzero_ps();
977         fix2             = _mm_setzero_ps();
978         fiy2             = _mm_setzero_ps();
979         fiz2             = _mm_setzero_ps();
980         fix3             = _mm_setzero_ps();
981         fiy3             = _mm_setzero_ps();
982         fiz3             = _mm_setzero_ps();
983
984         /* Start inner kernel loop */
985         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
986         {
987
988             /* Get j neighbor index, and coordinate index */
989             jnrA             = jjnr[jidx];
990             jnrB             = jjnr[jidx+1];
991             jnrC             = jjnr[jidx+2];
992             jnrD             = jjnr[jidx+3];
993             j_coord_offsetA  = DIM*jnrA;
994             j_coord_offsetB  = DIM*jnrB;
995             j_coord_offsetC  = DIM*jnrC;
996             j_coord_offsetD  = DIM*jnrD;
997
998             /* load j atom coordinates */
999             gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
1000                                               x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
1001                                               &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
1002
1003             /* Calculate displacement vector */
1004             dx11             = _mm_sub_ps(ix1,jx1);
1005             dy11             = _mm_sub_ps(iy1,jy1);
1006             dz11             = _mm_sub_ps(iz1,jz1);
1007             dx12             = _mm_sub_ps(ix1,jx2);
1008             dy12             = _mm_sub_ps(iy1,jy2);
1009             dz12             = _mm_sub_ps(iz1,jz2);
1010             dx13             = _mm_sub_ps(ix1,jx3);
1011             dy13             = _mm_sub_ps(iy1,jy3);
1012             dz13             = _mm_sub_ps(iz1,jz3);
1013             dx21             = _mm_sub_ps(ix2,jx1);
1014             dy21             = _mm_sub_ps(iy2,jy1);
1015             dz21             = _mm_sub_ps(iz2,jz1);
1016             dx22             = _mm_sub_ps(ix2,jx2);
1017             dy22             = _mm_sub_ps(iy2,jy2);
1018             dz22             = _mm_sub_ps(iz2,jz2);
1019             dx23             = _mm_sub_ps(ix2,jx3);
1020             dy23             = _mm_sub_ps(iy2,jy3);
1021             dz23             = _mm_sub_ps(iz2,jz3);
1022             dx31             = _mm_sub_ps(ix3,jx1);
1023             dy31             = _mm_sub_ps(iy3,jy1);
1024             dz31             = _mm_sub_ps(iz3,jz1);
1025             dx32             = _mm_sub_ps(ix3,jx2);
1026             dy32             = _mm_sub_ps(iy3,jy2);
1027             dz32             = _mm_sub_ps(iz3,jz2);
1028             dx33             = _mm_sub_ps(ix3,jx3);
1029             dy33             = _mm_sub_ps(iy3,jy3);
1030             dz33             = _mm_sub_ps(iz3,jz3);
1031
1032             /* Calculate squared distance and things based on it */
1033             rsq11            = gmx_mm_calc_rsq_ps(dx11,dy11,dz11);
1034             rsq12            = gmx_mm_calc_rsq_ps(dx12,dy12,dz12);
1035             rsq13            = gmx_mm_calc_rsq_ps(dx13,dy13,dz13);
1036             rsq21            = gmx_mm_calc_rsq_ps(dx21,dy21,dz21);
1037             rsq22            = gmx_mm_calc_rsq_ps(dx22,dy22,dz22);
1038             rsq23            = gmx_mm_calc_rsq_ps(dx23,dy23,dz23);
1039             rsq31            = gmx_mm_calc_rsq_ps(dx31,dy31,dz31);
1040             rsq32            = gmx_mm_calc_rsq_ps(dx32,dy32,dz32);
1041             rsq33            = gmx_mm_calc_rsq_ps(dx33,dy33,dz33);
1042
1043             rinv11           = gmx_mm_invsqrt_ps(rsq11);
1044             rinv12           = gmx_mm_invsqrt_ps(rsq12);
1045             rinv13           = gmx_mm_invsqrt_ps(rsq13);
1046             rinv21           = gmx_mm_invsqrt_ps(rsq21);
1047             rinv22           = gmx_mm_invsqrt_ps(rsq22);
1048             rinv23           = gmx_mm_invsqrt_ps(rsq23);
1049             rinv31           = gmx_mm_invsqrt_ps(rsq31);
1050             rinv32           = gmx_mm_invsqrt_ps(rsq32);
1051             rinv33           = gmx_mm_invsqrt_ps(rsq33);
1052
1053             rinvsq11         = _mm_mul_ps(rinv11,rinv11);
1054             rinvsq12         = _mm_mul_ps(rinv12,rinv12);
1055             rinvsq13         = _mm_mul_ps(rinv13,rinv13);
1056             rinvsq21         = _mm_mul_ps(rinv21,rinv21);
1057             rinvsq22         = _mm_mul_ps(rinv22,rinv22);
1058             rinvsq23         = _mm_mul_ps(rinv23,rinv23);
1059             rinvsq31         = _mm_mul_ps(rinv31,rinv31);
1060             rinvsq32         = _mm_mul_ps(rinv32,rinv32);
1061             rinvsq33         = _mm_mul_ps(rinv33,rinv33);
1062
1063             fjx1             = _mm_setzero_ps();
1064             fjy1             = _mm_setzero_ps();
1065             fjz1             = _mm_setzero_ps();
1066             fjx2             = _mm_setzero_ps();
1067             fjy2             = _mm_setzero_ps();
1068             fjz2             = _mm_setzero_ps();
1069             fjx3             = _mm_setzero_ps();
1070             fjy3             = _mm_setzero_ps();
1071             fjz3             = _mm_setzero_ps();
1072
1073             /**************************
1074              * CALCULATE INTERACTIONS *
1075              **************************/
1076
1077             /* REACTION-FIELD ELECTROSTATICS */
1078             felec            = _mm_mul_ps(qq11,_mm_msub_ps(rinv11,rinvsq11,krf2));
1079
1080             fscal            = felec;
1081
1082              /* Update vectorial force */
1083             fix1             = _mm_macc_ps(dx11,fscal,fix1);
1084             fiy1             = _mm_macc_ps(dy11,fscal,fiy1);
1085             fiz1             = _mm_macc_ps(dz11,fscal,fiz1);
1086
1087             fjx1             = _mm_macc_ps(dx11,fscal,fjx1);
1088             fjy1             = _mm_macc_ps(dy11,fscal,fjy1);
1089             fjz1             = _mm_macc_ps(dz11,fscal,fjz1);
1090
1091             /**************************
1092              * CALCULATE INTERACTIONS *
1093              **************************/
1094
1095             /* REACTION-FIELD ELECTROSTATICS */
1096             felec            = _mm_mul_ps(qq12,_mm_msub_ps(rinv12,rinvsq12,krf2));
1097
1098             fscal            = felec;
1099
1100              /* Update vectorial force */
1101             fix1             = _mm_macc_ps(dx12,fscal,fix1);
1102             fiy1             = _mm_macc_ps(dy12,fscal,fiy1);
1103             fiz1             = _mm_macc_ps(dz12,fscal,fiz1);
1104
1105             fjx2             = _mm_macc_ps(dx12,fscal,fjx2);
1106             fjy2             = _mm_macc_ps(dy12,fscal,fjy2);
1107             fjz2             = _mm_macc_ps(dz12,fscal,fjz2);
1108
1109             /**************************
1110              * CALCULATE INTERACTIONS *
1111              **************************/
1112
1113             /* REACTION-FIELD ELECTROSTATICS */
1114             felec            = _mm_mul_ps(qq13,_mm_msub_ps(rinv13,rinvsq13,krf2));
1115
1116             fscal            = felec;
1117
1118              /* Update vectorial force */
1119             fix1             = _mm_macc_ps(dx13,fscal,fix1);
1120             fiy1             = _mm_macc_ps(dy13,fscal,fiy1);
1121             fiz1             = _mm_macc_ps(dz13,fscal,fiz1);
1122
1123             fjx3             = _mm_macc_ps(dx13,fscal,fjx3);
1124             fjy3             = _mm_macc_ps(dy13,fscal,fjy3);
1125             fjz3             = _mm_macc_ps(dz13,fscal,fjz3);
1126
1127             /**************************
1128              * CALCULATE INTERACTIONS *
1129              **************************/
1130
1131             /* REACTION-FIELD ELECTROSTATICS */
1132             felec            = _mm_mul_ps(qq21,_mm_msub_ps(rinv21,rinvsq21,krf2));
1133
1134             fscal            = felec;
1135
1136              /* Update vectorial force */
1137             fix2             = _mm_macc_ps(dx21,fscal,fix2);
1138             fiy2             = _mm_macc_ps(dy21,fscal,fiy2);
1139             fiz2             = _mm_macc_ps(dz21,fscal,fiz2);
1140
1141             fjx1             = _mm_macc_ps(dx21,fscal,fjx1);
1142             fjy1             = _mm_macc_ps(dy21,fscal,fjy1);
1143             fjz1             = _mm_macc_ps(dz21,fscal,fjz1);
1144
1145             /**************************
1146              * CALCULATE INTERACTIONS *
1147              **************************/
1148
1149             /* REACTION-FIELD ELECTROSTATICS */
1150             felec            = _mm_mul_ps(qq22,_mm_msub_ps(rinv22,rinvsq22,krf2));
1151
1152             fscal            = felec;
1153
1154              /* Update vectorial force */
1155             fix2             = _mm_macc_ps(dx22,fscal,fix2);
1156             fiy2             = _mm_macc_ps(dy22,fscal,fiy2);
1157             fiz2             = _mm_macc_ps(dz22,fscal,fiz2);
1158
1159             fjx2             = _mm_macc_ps(dx22,fscal,fjx2);
1160             fjy2             = _mm_macc_ps(dy22,fscal,fjy2);
1161             fjz2             = _mm_macc_ps(dz22,fscal,fjz2);
1162
1163             /**************************
1164              * CALCULATE INTERACTIONS *
1165              **************************/
1166
1167             /* REACTION-FIELD ELECTROSTATICS */
1168             felec            = _mm_mul_ps(qq23,_mm_msub_ps(rinv23,rinvsq23,krf2));
1169
1170             fscal            = felec;
1171
1172              /* Update vectorial force */
1173             fix2             = _mm_macc_ps(dx23,fscal,fix2);
1174             fiy2             = _mm_macc_ps(dy23,fscal,fiy2);
1175             fiz2             = _mm_macc_ps(dz23,fscal,fiz2);
1176
1177             fjx3             = _mm_macc_ps(dx23,fscal,fjx3);
1178             fjy3             = _mm_macc_ps(dy23,fscal,fjy3);
1179             fjz3             = _mm_macc_ps(dz23,fscal,fjz3);
1180
1181             /**************************
1182              * CALCULATE INTERACTIONS *
1183              **************************/
1184
1185             /* REACTION-FIELD ELECTROSTATICS */
1186             felec            = _mm_mul_ps(qq31,_mm_msub_ps(rinv31,rinvsq31,krf2));
1187
1188             fscal            = felec;
1189
1190              /* Update vectorial force */
1191             fix3             = _mm_macc_ps(dx31,fscal,fix3);
1192             fiy3             = _mm_macc_ps(dy31,fscal,fiy3);
1193             fiz3             = _mm_macc_ps(dz31,fscal,fiz3);
1194
1195             fjx1             = _mm_macc_ps(dx31,fscal,fjx1);
1196             fjy1             = _mm_macc_ps(dy31,fscal,fjy1);
1197             fjz1             = _mm_macc_ps(dz31,fscal,fjz1);
1198
1199             /**************************
1200              * CALCULATE INTERACTIONS *
1201              **************************/
1202
1203             /* REACTION-FIELD ELECTROSTATICS */
1204             felec            = _mm_mul_ps(qq32,_mm_msub_ps(rinv32,rinvsq32,krf2));
1205
1206             fscal            = felec;
1207
1208              /* Update vectorial force */
1209             fix3             = _mm_macc_ps(dx32,fscal,fix3);
1210             fiy3             = _mm_macc_ps(dy32,fscal,fiy3);
1211             fiz3             = _mm_macc_ps(dz32,fscal,fiz3);
1212
1213             fjx2             = _mm_macc_ps(dx32,fscal,fjx2);
1214             fjy2             = _mm_macc_ps(dy32,fscal,fjy2);
1215             fjz2             = _mm_macc_ps(dz32,fscal,fjz2);
1216
1217             /**************************
1218              * CALCULATE INTERACTIONS *
1219              **************************/
1220
1221             /* REACTION-FIELD ELECTROSTATICS */
1222             felec            = _mm_mul_ps(qq33,_mm_msub_ps(rinv33,rinvsq33,krf2));
1223
1224             fscal            = felec;
1225
1226              /* Update vectorial force */
1227             fix3             = _mm_macc_ps(dx33,fscal,fix3);
1228             fiy3             = _mm_macc_ps(dy33,fscal,fiy3);
1229             fiz3             = _mm_macc_ps(dz33,fscal,fiz3);
1230
1231             fjx3             = _mm_macc_ps(dx33,fscal,fjx3);
1232             fjy3             = _mm_macc_ps(dy33,fscal,fjy3);
1233             fjz3             = _mm_macc_ps(dz33,fscal,fjz3);
1234
1235             fjptrA             = f+j_coord_offsetA;
1236             fjptrB             = f+j_coord_offsetB;
1237             fjptrC             = f+j_coord_offsetC;
1238             fjptrD             = f+j_coord_offsetD;
1239
1240             gmx_mm_decrement_3rvec_4ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
1241                                                    fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1242
1243             /* Inner loop uses 270 flops */
1244         }
1245
1246         if(jidx<j_index_end)
1247         {
1248
1249             /* Get j neighbor index, and coordinate index */
1250             jnrlistA         = jjnr[jidx];
1251             jnrlistB         = jjnr[jidx+1];
1252             jnrlistC         = jjnr[jidx+2];
1253             jnrlistD         = jjnr[jidx+3];
1254             /* Sign of each element will be negative for non-real atoms.
1255              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1256              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1257              */
1258             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
1259             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
1260             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
1261             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
1262             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
1263             j_coord_offsetA  = DIM*jnrA;
1264             j_coord_offsetB  = DIM*jnrB;
1265             j_coord_offsetC  = DIM*jnrC;
1266             j_coord_offsetD  = DIM*jnrD;
1267
1268             /* load j atom coordinates */
1269             gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
1270                                               x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
1271                                               &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
1272
1273             /* Calculate displacement vector */
1274             dx11             = _mm_sub_ps(ix1,jx1);
1275             dy11             = _mm_sub_ps(iy1,jy1);
1276             dz11             = _mm_sub_ps(iz1,jz1);
1277             dx12             = _mm_sub_ps(ix1,jx2);
1278             dy12             = _mm_sub_ps(iy1,jy2);
1279             dz12             = _mm_sub_ps(iz1,jz2);
1280             dx13             = _mm_sub_ps(ix1,jx3);
1281             dy13             = _mm_sub_ps(iy1,jy3);
1282             dz13             = _mm_sub_ps(iz1,jz3);
1283             dx21             = _mm_sub_ps(ix2,jx1);
1284             dy21             = _mm_sub_ps(iy2,jy1);
1285             dz21             = _mm_sub_ps(iz2,jz1);
1286             dx22             = _mm_sub_ps(ix2,jx2);
1287             dy22             = _mm_sub_ps(iy2,jy2);
1288             dz22             = _mm_sub_ps(iz2,jz2);
1289             dx23             = _mm_sub_ps(ix2,jx3);
1290             dy23             = _mm_sub_ps(iy2,jy3);
1291             dz23             = _mm_sub_ps(iz2,jz3);
1292             dx31             = _mm_sub_ps(ix3,jx1);
1293             dy31             = _mm_sub_ps(iy3,jy1);
1294             dz31             = _mm_sub_ps(iz3,jz1);
1295             dx32             = _mm_sub_ps(ix3,jx2);
1296             dy32             = _mm_sub_ps(iy3,jy2);
1297             dz32             = _mm_sub_ps(iz3,jz2);
1298             dx33             = _mm_sub_ps(ix3,jx3);
1299             dy33             = _mm_sub_ps(iy3,jy3);
1300             dz33             = _mm_sub_ps(iz3,jz3);
1301
1302             /* Calculate squared distance and things based on it */
1303             rsq11            = gmx_mm_calc_rsq_ps(dx11,dy11,dz11);
1304             rsq12            = gmx_mm_calc_rsq_ps(dx12,dy12,dz12);
1305             rsq13            = gmx_mm_calc_rsq_ps(dx13,dy13,dz13);
1306             rsq21            = gmx_mm_calc_rsq_ps(dx21,dy21,dz21);
1307             rsq22            = gmx_mm_calc_rsq_ps(dx22,dy22,dz22);
1308             rsq23            = gmx_mm_calc_rsq_ps(dx23,dy23,dz23);
1309             rsq31            = gmx_mm_calc_rsq_ps(dx31,dy31,dz31);
1310             rsq32            = gmx_mm_calc_rsq_ps(dx32,dy32,dz32);
1311             rsq33            = gmx_mm_calc_rsq_ps(dx33,dy33,dz33);
1312
1313             rinv11           = gmx_mm_invsqrt_ps(rsq11);
1314             rinv12           = gmx_mm_invsqrt_ps(rsq12);
1315             rinv13           = gmx_mm_invsqrt_ps(rsq13);
1316             rinv21           = gmx_mm_invsqrt_ps(rsq21);
1317             rinv22           = gmx_mm_invsqrt_ps(rsq22);
1318             rinv23           = gmx_mm_invsqrt_ps(rsq23);
1319             rinv31           = gmx_mm_invsqrt_ps(rsq31);
1320             rinv32           = gmx_mm_invsqrt_ps(rsq32);
1321             rinv33           = gmx_mm_invsqrt_ps(rsq33);
1322
1323             rinvsq11         = _mm_mul_ps(rinv11,rinv11);
1324             rinvsq12         = _mm_mul_ps(rinv12,rinv12);
1325             rinvsq13         = _mm_mul_ps(rinv13,rinv13);
1326             rinvsq21         = _mm_mul_ps(rinv21,rinv21);
1327             rinvsq22         = _mm_mul_ps(rinv22,rinv22);
1328             rinvsq23         = _mm_mul_ps(rinv23,rinv23);
1329             rinvsq31         = _mm_mul_ps(rinv31,rinv31);
1330             rinvsq32         = _mm_mul_ps(rinv32,rinv32);
1331             rinvsq33         = _mm_mul_ps(rinv33,rinv33);
1332
1333             fjx1             = _mm_setzero_ps();
1334             fjy1             = _mm_setzero_ps();
1335             fjz1             = _mm_setzero_ps();
1336             fjx2             = _mm_setzero_ps();
1337             fjy2             = _mm_setzero_ps();
1338             fjz2             = _mm_setzero_ps();
1339             fjx3             = _mm_setzero_ps();
1340             fjy3             = _mm_setzero_ps();
1341             fjz3             = _mm_setzero_ps();
1342
1343             /**************************
1344              * CALCULATE INTERACTIONS *
1345              **************************/
1346
1347             /* REACTION-FIELD ELECTROSTATICS */
1348             felec            = _mm_mul_ps(qq11,_mm_msub_ps(rinv11,rinvsq11,krf2));
1349
1350             fscal            = felec;
1351
1352             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1353
1354              /* Update vectorial force */
1355             fix1             = _mm_macc_ps(dx11,fscal,fix1);
1356             fiy1             = _mm_macc_ps(dy11,fscal,fiy1);
1357             fiz1             = _mm_macc_ps(dz11,fscal,fiz1);
1358
1359             fjx1             = _mm_macc_ps(dx11,fscal,fjx1);
1360             fjy1             = _mm_macc_ps(dy11,fscal,fjy1);
1361             fjz1             = _mm_macc_ps(dz11,fscal,fjz1);
1362
1363             /**************************
1364              * CALCULATE INTERACTIONS *
1365              **************************/
1366
1367             /* REACTION-FIELD ELECTROSTATICS */
1368             felec            = _mm_mul_ps(qq12,_mm_msub_ps(rinv12,rinvsq12,krf2));
1369
1370             fscal            = felec;
1371
1372             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1373
1374              /* Update vectorial force */
1375             fix1             = _mm_macc_ps(dx12,fscal,fix1);
1376             fiy1             = _mm_macc_ps(dy12,fscal,fiy1);
1377             fiz1             = _mm_macc_ps(dz12,fscal,fiz1);
1378
1379             fjx2             = _mm_macc_ps(dx12,fscal,fjx2);
1380             fjy2             = _mm_macc_ps(dy12,fscal,fjy2);
1381             fjz2             = _mm_macc_ps(dz12,fscal,fjz2);
1382
1383             /**************************
1384              * CALCULATE INTERACTIONS *
1385              **************************/
1386
1387             /* REACTION-FIELD ELECTROSTATICS */
1388             felec            = _mm_mul_ps(qq13,_mm_msub_ps(rinv13,rinvsq13,krf2));
1389
1390             fscal            = felec;
1391
1392             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1393
1394              /* Update vectorial force */
1395             fix1             = _mm_macc_ps(dx13,fscal,fix1);
1396             fiy1             = _mm_macc_ps(dy13,fscal,fiy1);
1397             fiz1             = _mm_macc_ps(dz13,fscal,fiz1);
1398
1399             fjx3             = _mm_macc_ps(dx13,fscal,fjx3);
1400             fjy3             = _mm_macc_ps(dy13,fscal,fjy3);
1401             fjz3             = _mm_macc_ps(dz13,fscal,fjz3);
1402
1403             /**************************
1404              * CALCULATE INTERACTIONS *
1405              **************************/
1406
1407             /* REACTION-FIELD ELECTROSTATICS */
1408             felec            = _mm_mul_ps(qq21,_mm_msub_ps(rinv21,rinvsq21,krf2));
1409
1410             fscal            = felec;
1411
1412             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1413
1414              /* Update vectorial force */
1415             fix2             = _mm_macc_ps(dx21,fscal,fix2);
1416             fiy2             = _mm_macc_ps(dy21,fscal,fiy2);
1417             fiz2             = _mm_macc_ps(dz21,fscal,fiz2);
1418
1419             fjx1             = _mm_macc_ps(dx21,fscal,fjx1);
1420             fjy1             = _mm_macc_ps(dy21,fscal,fjy1);
1421             fjz1             = _mm_macc_ps(dz21,fscal,fjz1);
1422
1423             /**************************
1424              * CALCULATE INTERACTIONS *
1425              **************************/
1426
1427             /* REACTION-FIELD ELECTROSTATICS */
1428             felec            = _mm_mul_ps(qq22,_mm_msub_ps(rinv22,rinvsq22,krf2));
1429
1430             fscal            = felec;
1431
1432             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1433
1434              /* Update vectorial force */
1435             fix2             = _mm_macc_ps(dx22,fscal,fix2);
1436             fiy2             = _mm_macc_ps(dy22,fscal,fiy2);
1437             fiz2             = _mm_macc_ps(dz22,fscal,fiz2);
1438
1439             fjx2             = _mm_macc_ps(dx22,fscal,fjx2);
1440             fjy2             = _mm_macc_ps(dy22,fscal,fjy2);
1441             fjz2             = _mm_macc_ps(dz22,fscal,fjz2);
1442
1443             /**************************
1444              * CALCULATE INTERACTIONS *
1445              **************************/
1446
1447             /* REACTION-FIELD ELECTROSTATICS */
1448             felec            = _mm_mul_ps(qq23,_mm_msub_ps(rinv23,rinvsq23,krf2));
1449
1450             fscal            = felec;
1451
1452             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1453
1454              /* Update vectorial force */
1455             fix2             = _mm_macc_ps(dx23,fscal,fix2);
1456             fiy2             = _mm_macc_ps(dy23,fscal,fiy2);
1457             fiz2             = _mm_macc_ps(dz23,fscal,fiz2);
1458
1459             fjx3             = _mm_macc_ps(dx23,fscal,fjx3);
1460             fjy3             = _mm_macc_ps(dy23,fscal,fjy3);
1461             fjz3             = _mm_macc_ps(dz23,fscal,fjz3);
1462
1463             /**************************
1464              * CALCULATE INTERACTIONS *
1465              **************************/
1466
1467             /* REACTION-FIELD ELECTROSTATICS */
1468             felec            = _mm_mul_ps(qq31,_mm_msub_ps(rinv31,rinvsq31,krf2));
1469
1470             fscal            = felec;
1471
1472             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1473
1474              /* Update vectorial force */
1475             fix3             = _mm_macc_ps(dx31,fscal,fix3);
1476             fiy3             = _mm_macc_ps(dy31,fscal,fiy3);
1477             fiz3             = _mm_macc_ps(dz31,fscal,fiz3);
1478
1479             fjx1             = _mm_macc_ps(dx31,fscal,fjx1);
1480             fjy1             = _mm_macc_ps(dy31,fscal,fjy1);
1481             fjz1             = _mm_macc_ps(dz31,fscal,fjz1);
1482
1483             /**************************
1484              * CALCULATE INTERACTIONS *
1485              **************************/
1486
1487             /* REACTION-FIELD ELECTROSTATICS */
1488             felec            = _mm_mul_ps(qq32,_mm_msub_ps(rinv32,rinvsq32,krf2));
1489
1490             fscal            = felec;
1491
1492             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1493
1494              /* Update vectorial force */
1495             fix3             = _mm_macc_ps(dx32,fscal,fix3);
1496             fiy3             = _mm_macc_ps(dy32,fscal,fiy3);
1497             fiz3             = _mm_macc_ps(dz32,fscal,fiz3);
1498
1499             fjx2             = _mm_macc_ps(dx32,fscal,fjx2);
1500             fjy2             = _mm_macc_ps(dy32,fscal,fjy2);
1501             fjz2             = _mm_macc_ps(dz32,fscal,fjz2);
1502
1503             /**************************
1504              * CALCULATE INTERACTIONS *
1505              **************************/
1506
1507             /* REACTION-FIELD ELECTROSTATICS */
1508             felec            = _mm_mul_ps(qq33,_mm_msub_ps(rinv33,rinvsq33,krf2));
1509
1510             fscal            = felec;
1511
1512             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1513
1514              /* Update vectorial force */
1515             fix3             = _mm_macc_ps(dx33,fscal,fix3);
1516             fiy3             = _mm_macc_ps(dy33,fscal,fiy3);
1517             fiz3             = _mm_macc_ps(dz33,fscal,fiz3);
1518
1519             fjx3             = _mm_macc_ps(dx33,fscal,fjx3);
1520             fjy3             = _mm_macc_ps(dy33,fscal,fjy3);
1521             fjz3             = _mm_macc_ps(dz33,fscal,fjz3);
1522
1523             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1524             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1525             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1526             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1527
1528             gmx_mm_decrement_3rvec_4ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
1529                                                    fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1530
1531             /* Inner loop uses 270 flops */
1532         }
1533
1534         /* End of innermost loop */
1535
1536         gmx_mm_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1537                                               f+i_coord_offset+DIM,fshift+i_shift_offset);
1538
1539         /* Increment number of inner iterations */
1540         inneriter                  += j_index_end - j_index_start;
1541
1542         /* Outer loop uses 18 flops */
1543     }
1544
1545     /* Increment number of outer iterations */
1546     outeriter        += nri;
1547
1548     /* Update outer/inner flops */
1549
1550     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4W4_F,outeriter*18 + inneriter*270);
1551 }