Remove all unnecessary HAVE_CONFIG_H
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_double / nb_kernel_ElecEw_VdwNone_GeomW4P1_avx_128_fma_double.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_double 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_avx_128_fma_double.h"
48 #include "kernelutil_x86_avx_128_fma_double.h"
49
50 /*
51  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwNone_GeomW4P1_VF_avx_128_fma_double
52  * Electrostatics interaction: Ewald
53  * VdW interaction:            None
54  * Geometry:                   Water4-Particle
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecEw_VdwNone_GeomW4P1_VF_avx_128_fma_double
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 refer to j loop unrolling done with SSE double precision, e.g. for the two 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;
75     int              j_coord_offsetA,j_coord_offsetB;
76     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
77     real             rcutoff_scalar;
78     real             *shiftvec,*fshift,*x,*f;
79     __m128d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
80     int              vdwioffset1;
81     __m128d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
82     int              vdwioffset2;
83     __m128d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
84     int              vdwioffset3;
85     __m128d          ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
86     int              vdwjidx0A,vdwjidx0B;
87     __m128d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
88     __m128d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
89     __m128d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
90     __m128d          dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
91     __m128d          velec,felec,velecsum,facel,crf,krf,krf2;
92     real             *charge;
93     __m128i          ewitab;
94     __m128d          ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
95     real             *ewtab;
96     __m128d          dummy_mask,cutoff_mask;
97     __m128d          signbit   = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
98     __m128d          one     = _mm_set1_pd(1.0);
99     __m128d          two     = _mm_set1_pd(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_pd(fr->epsfac);
112     charge           = mdatoms->chargeA;
113
114     sh_ewald         = _mm_set1_pd(fr->ic->sh_ewald);
115     ewtab            = fr->ic->tabq_coul_FDV0;
116     ewtabscale       = _mm_set1_pd(fr->ic->tabq_scale);
117     ewtabhalfspace   = _mm_set1_pd(0.5/fr->ic->tabq_scale);
118
119     /* Setup water-specific parameters */
120     inr              = nlist->iinr[0];
121     iq1              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
122     iq2              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
123     iq3              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
124
125     /* Avoid stupid compiler warnings */
126     jnrA = jnrB = 0;
127     j_coord_offsetA = 0;
128     j_coord_offsetB = 0;
129
130     outeriter        = 0;
131     inneriter        = 0;
132
133     /* Start outer loop over neighborlists */
134     for(iidx=0; iidx<nri; iidx++)
135     {
136         /* Load shift vector for this list */
137         i_shift_offset   = DIM*shiftidx[iidx];
138
139         /* Load limits for loop over neighbors */
140         j_index_start    = jindex[iidx];
141         j_index_end      = jindex[iidx+1];
142
143         /* Get outer coordinate index */
144         inr              = iinr[iidx];
145         i_coord_offset   = DIM*inr;
146
147         /* Load i particle coords and add shift vector */
148         gmx_mm_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
149                                                  &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
150
151         fix1             = _mm_setzero_pd();
152         fiy1             = _mm_setzero_pd();
153         fiz1             = _mm_setzero_pd();
154         fix2             = _mm_setzero_pd();
155         fiy2             = _mm_setzero_pd();
156         fiz2             = _mm_setzero_pd();
157         fix3             = _mm_setzero_pd();
158         fiy3             = _mm_setzero_pd();
159         fiz3             = _mm_setzero_pd();
160
161         /* Reset potential sums */
162         velecsum         = _mm_setzero_pd();
163
164         /* Start inner kernel loop */
165         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
166         {
167
168             /* Get j neighbor index, and coordinate index */
169             jnrA             = jjnr[jidx];
170             jnrB             = jjnr[jidx+1];
171             j_coord_offsetA  = DIM*jnrA;
172             j_coord_offsetB  = DIM*jnrB;
173
174             /* load j atom coordinates */
175             gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
176                                               &jx0,&jy0,&jz0);
177
178             /* Calculate displacement vector */
179             dx10             = _mm_sub_pd(ix1,jx0);
180             dy10             = _mm_sub_pd(iy1,jy0);
181             dz10             = _mm_sub_pd(iz1,jz0);
182             dx20             = _mm_sub_pd(ix2,jx0);
183             dy20             = _mm_sub_pd(iy2,jy0);
184             dz20             = _mm_sub_pd(iz2,jz0);
185             dx30             = _mm_sub_pd(ix3,jx0);
186             dy30             = _mm_sub_pd(iy3,jy0);
187             dz30             = _mm_sub_pd(iz3,jz0);
188
189             /* Calculate squared distance and things based on it */
190             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
191             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
192             rsq30            = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
193
194             rinv10           = gmx_mm_invsqrt_pd(rsq10);
195             rinv20           = gmx_mm_invsqrt_pd(rsq20);
196             rinv30           = gmx_mm_invsqrt_pd(rsq30);
197
198             rinvsq10         = _mm_mul_pd(rinv10,rinv10);
199             rinvsq20         = _mm_mul_pd(rinv20,rinv20);
200             rinvsq30         = _mm_mul_pd(rinv30,rinv30);
201
202             /* Load parameters for j particles */
203             jq0              = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
204
205             fjx0             = _mm_setzero_pd();
206             fjy0             = _mm_setzero_pd();
207             fjz0             = _mm_setzero_pd();
208
209             /**************************
210              * CALCULATE INTERACTIONS *
211              **************************/
212
213             r10              = _mm_mul_pd(rsq10,rinv10);
214
215             /* Compute parameters for interactions between i and j atoms */
216             qq10             = _mm_mul_pd(iq1,jq0);
217
218             /* EWALD ELECTROSTATICS */
219
220             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
221             ewrt             = _mm_mul_pd(r10,ewtabscale);
222             ewitab           = _mm_cvttpd_epi32(ewrt);
223 #ifdef __XOP__
224             eweps            = _mm_frcz_pd(ewrt);
225 #else
226             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
227 #endif
228             twoeweps         = _mm_add_pd(eweps,eweps);
229             ewitab           = _mm_slli_epi32(ewitab,2);
230             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
231             ewtabD           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
232             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
233             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
234             ewtabFn          = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,1) +2);
235             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
236             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
237             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
238             velec            = _mm_mul_pd(qq10,_mm_sub_pd(rinv10,velec));
239             felec            = _mm_mul_pd(_mm_mul_pd(qq10,rinv10),_mm_sub_pd(rinvsq10,felec));
240
241             /* Update potential sum for this i atom from the interaction with this j atom. */
242             velecsum         = _mm_add_pd(velecsum,velec);
243
244             fscal            = felec;
245
246             /* Update vectorial force */
247             fix1             = _mm_macc_pd(dx10,fscal,fix1);
248             fiy1             = _mm_macc_pd(dy10,fscal,fiy1);
249             fiz1             = _mm_macc_pd(dz10,fscal,fiz1);
250             
251             fjx0             = _mm_macc_pd(dx10,fscal,fjx0);
252             fjy0             = _mm_macc_pd(dy10,fscal,fjy0);
253             fjz0             = _mm_macc_pd(dz10,fscal,fjz0);
254
255             /**************************
256              * CALCULATE INTERACTIONS *
257              **************************/
258
259             r20              = _mm_mul_pd(rsq20,rinv20);
260
261             /* Compute parameters for interactions between i and j atoms */
262             qq20             = _mm_mul_pd(iq2,jq0);
263
264             /* EWALD ELECTROSTATICS */
265
266             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
267             ewrt             = _mm_mul_pd(r20,ewtabscale);
268             ewitab           = _mm_cvttpd_epi32(ewrt);
269 #ifdef __XOP__
270             eweps            = _mm_frcz_pd(ewrt);
271 #else
272             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
273 #endif
274             twoeweps         = _mm_add_pd(eweps,eweps);
275             ewitab           = _mm_slli_epi32(ewitab,2);
276             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
277             ewtabD           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
278             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
279             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
280             ewtabFn          = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,1) +2);
281             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
282             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
283             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
284             velec            = _mm_mul_pd(qq20,_mm_sub_pd(rinv20,velec));
285             felec            = _mm_mul_pd(_mm_mul_pd(qq20,rinv20),_mm_sub_pd(rinvsq20,felec));
286
287             /* Update potential sum for this i atom from the interaction with this j atom. */
288             velecsum         = _mm_add_pd(velecsum,velec);
289
290             fscal            = felec;
291
292             /* Update vectorial force */
293             fix2             = _mm_macc_pd(dx20,fscal,fix2);
294             fiy2             = _mm_macc_pd(dy20,fscal,fiy2);
295             fiz2             = _mm_macc_pd(dz20,fscal,fiz2);
296             
297             fjx0             = _mm_macc_pd(dx20,fscal,fjx0);
298             fjy0             = _mm_macc_pd(dy20,fscal,fjy0);
299             fjz0             = _mm_macc_pd(dz20,fscal,fjz0);
300
301             /**************************
302              * CALCULATE INTERACTIONS *
303              **************************/
304
305             r30              = _mm_mul_pd(rsq30,rinv30);
306
307             /* Compute parameters for interactions between i and j atoms */
308             qq30             = _mm_mul_pd(iq3,jq0);
309
310             /* EWALD ELECTROSTATICS */
311
312             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
313             ewrt             = _mm_mul_pd(r30,ewtabscale);
314             ewitab           = _mm_cvttpd_epi32(ewrt);
315 #ifdef __XOP__
316             eweps            = _mm_frcz_pd(ewrt);
317 #else
318             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
319 #endif
320             twoeweps         = _mm_add_pd(eweps,eweps);
321             ewitab           = _mm_slli_epi32(ewitab,2);
322             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
323             ewtabD           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
324             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
325             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
326             ewtabFn          = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,1) +2);
327             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
328             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
329             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
330             velec            = _mm_mul_pd(qq30,_mm_sub_pd(rinv30,velec));
331             felec            = _mm_mul_pd(_mm_mul_pd(qq30,rinv30),_mm_sub_pd(rinvsq30,felec));
332
333             /* Update potential sum for this i atom from the interaction with this j atom. */
334             velecsum         = _mm_add_pd(velecsum,velec);
335
336             fscal            = felec;
337
338             /* Update vectorial force */
339             fix3             = _mm_macc_pd(dx30,fscal,fix3);
340             fiy3             = _mm_macc_pd(dy30,fscal,fiy3);
341             fiz3             = _mm_macc_pd(dz30,fscal,fiz3);
342             
343             fjx0             = _mm_macc_pd(dx30,fscal,fjx0);
344             fjy0             = _mm_macc_pd(dy30,fscal,fjy0);
345             fjz0             = _mm_macc_pd(dz30,fscal,fjz0);
346
347             gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
348
349             /* Inner loop uses 135 flops */
350         }
351
352         if(jidx<j_index_end)
353         {
354
355             jnrA             = jjnr[jidx];
356             j_coord_offsetA  = DIM*jnrA;
357
358             /* load j atom coordinates */
359             gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
360                                               &jx0,&jy0,&jz0);
361
362             /* Calculate displacement vector */
363             dx10             = _mm_sub_pd(ix1,jx0);
364             dy10             = _mm_sub_pd(iy1,jy0);
365             dz10             = _mm_sub_pd(iz1,jz0);
366             dx20             = _mm_sub_pd(ix2,jx0);
367             dy20             = _mm_sub_pd(iy2,jy0);
368             dz20             = _mm_sub_pd(iz2,jz0);
369             dx30             = _mm_sub_pd(ix3,jx0);
370             dy30             = _mm_sub_pd(iy3,jy0);
371             dz30             = _mm_sub_pd(iz3,jz0);
372
373             /* Calculate squared distance and things based on it */
374             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
375             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
376             rsq30            = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
377
378             rinv10           = gmx_mm_invsqrt_pd(rsq10);
379             rinv20           = gmx_mm_invsqrt_pd(rsq20);
380             rinv30           = gmx_mm_invsqrt_pd(rsq30);
381
382             rinvsq10         = _mm_mul_pd(rinv10,rinv10);
383             rinvsq20         = _mm_mul_pd(rinv20,rinv20);
384             rinvsq30         = _mm_mul_pd(rinv30,rinv30);
385
386             /* Load parameters for j particles */
387             jq0              = _mm_load_sd(charge+jnrA+0);
388
389             fjx0             = _mm_setzero_pd();
390             fjy0             = _mm_setzero_pd();
391             fjz0             = _mm_setzero_pd();
392
393             /**************************
394              * CALCULATE INTERACTIONS *
395              **************************/
396
397             r10              = _mm_mul_pd(rsq10,rinv10);
398
399             /* Compute parameters for interactions between i and j atoms */
400             qq10             = _mm_mul_pd(iq1,jq0);
401
402             /* EWALD ELECTROSTATICS */
403
404             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
405             ewrt             = _mm_mul_pd(r10,ewtabscale);
406             ewitab           = _mm_cvttpd_epi32(ewrt);
407 #ifdef __XOP__
408             eweps            = _mm_frcz_pd(ewrt);
409 #else
410             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
411 #endif
412             twoeweps         = _mm_add_pd(eweps,eweps);
413             ewitab           = _mm_slli_epi32(ewitab,2);
414             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
415             ewtabD           = _mm_setzero_pd();
416             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
417             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
418             ewtabFn          = _mm_setzero_pd();
419             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
420             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
421             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
422             velec            = _mm_mul_pd(qq10,_mm_sub_pd(rinv10,velec));
423             felec            = _mm_mul_pd(_mm_mul_pd(qq10,rinv10),_mm_sub_pd(rinvsq10,felec));
424
425             /* Update potential sum for this i atom from the interaction with this j atom. */
426             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
427             velecsum         = _mm_add_pd(velecsum,velec);
428
429             fscal            = felec;
430
431             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
432
433             /* Update vectorial force */
434             fix1             = _mm_macc_pd(dx10,fscal,fix1);
435             fiy1             = _mm_macc_pd(dy10,fscal,fiy1);
436             fiz1             = _mm_macc_pd(dz10,fscal,fiz1);
437             
438             fjx0             = _mm_macc_pd(dx10,fscal,fjx0);
439             fjy0             = _mm_macc_pd(dy10,fscal,fjy0);
440             fjz0             = _mm_macc_pd(dz10,fscal,fjz0);
441
442             /**************************
443              * CALCULATE INTERACTIONS *
444              **************************/
445
446             r20              = _mm_mul_pd(rsq20,rinv20);
447
448             /* Compute parameters for interactions between i and j atoms */
449             qq20             = _mm_mul_pd(iq2,jq0);
450
451             /* EWALD ELECTROSTATICS */
452
453             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
454             ewrt             = _mm_mul_pd(r20,ewtabscale);
455             ewitab           = _mm_cvttpd_epi32(ewrt);
456 #ifdef __XOP__
457             eweps            = _mm_frcz_pd(ewrt);
458 #else
459             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
460 #endif
461             twoeweps         = _mm_add_pd(eweps,eweps);
462             ewitab           = _mm_slli_epi32(ewitab,2);
463             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
464             ewtabD           = _mm_setzero_pd();
465             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
466             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
467             ewtabFn          = _mm_setzero_pd();
468             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
469             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
470             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
471             velec            = _mm_mul_pd(qq20,_mm_sub_pd(rinv20,velec));
472             felec            = _mm_mul_pd(_mm_mul_pd(qq20,rinv20),_mm_sub_pd(rinvsq20,felec));
473
474             /* Update potential sum for this i atom from the interaction with this j atom. */
475             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
476             velecsum         = _mm_add_pd(velecsum,velec);
477
478             fscal            = felec;
479
480             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
481
482             /* Update vectorial force */
483             fix2             = _mm_macc_pd(dx20,fscal,fix2);
484             fiy2             = _mm_macc_pd(dy20,fscal,fiy2);
485             fiz2             = _mm_macc_pd(dz20,fscal,fiz2);
486             
487             fjx0             = _mm_macc_pd(dx20,fscal,fjx0);
488             fjy0             = _mm_macc_pd(dy20,fscal,fjy0);
489             fjz0             = _mm_macc_pd(dz20,fscal,fjz0);
490
491             /**************************
492              * CALCULATE INTERACTIONS *
493              **************************/
494
495             r30              = _mm_mul_pd(rsq30,rinv30);
496
497             /* Compute parameters for interactions between i and j atoms */
498             qq30             = _mm_mul_pd(iq3,jq0);
499
500             /* EWALD ELECTROSTATICS */
501
502             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
503             ewrt             = _mm_mul_pd(r30,ewtabscale);
504             ewitab           = _mm_cvttpd_epi32(ewrt);
505 #ifdef __XOP__
506             eweps            = _mm_frcz_pd(ewrt);
507 #else
508             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
509 #endif
510             twoeweps         = _mm_add_pd(eweps,eweps);
511             ewitab           = _mm_slli_epi32(ewitab,2);
512             ewtabF           = _mm_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
513             ewtabD           = _mm_setzero_pd();
514             GMX_MM_TRANSPOSE2_PD(ewtabF,ewtabD);
515             ewtabV           = _mm_load_sd( ewtab + _mm_extract_epi32(ewitab,0) +2);
516             ewtabFn          = _mm_setzero_pd();
517             GMX_MM_TRANSPOSE2_PD(ewtabV,ewtabFn);
518             felec            = _mm_macc_pd(eweps,ewtabD,ewtabF);
519             velec            = _mm_nmacc_pd(_mm_mul_pd(ewtabhalfspace,eweps) ,_mm_add_pd(ewtabF,felec), ewtabV);
520             velec            = _mm_mul_pd(qq30,_mm_sub_pd(rinv30,velec));
521             felec            = _mm_mul_pd(_mm_mul_pd(qq30,rinv30),_mm_sub_pd(rinvsq30,felec));
522
523             /* Update potential sum for this i atom from the interaction with this j atom. */
524             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
525             velecsum         = _mm_add_pd(velecsum,velec);
526
527             fscal            = felec;
528
529             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
530
531             /* Update vectorial force */
532             fix3             = _mm_macc_pd(dx30,fscal,fix3);
533             fiy3             = _mm_macc_pd(dy30,fscal,fiy3);
534             fiz3             = _mm_macc_pd(dz30,fscal,fiz3);
535             
536             fjx0             = _mm_macc_pd(dx30,fscal,fjx0);
537             fjy0             = _mm_macc_pd(dy30,fscal,fjy0);
538             fjz0             = _mm_macc_pd(dz30,fscal,fjz0);
539
540             gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
541
542             /* Inner loop uses 135 flops */
543         }
544
545         /* End of innermost loop */
546
547         gmx_mm_update_iforce_3atom_swizzle_pd(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
548                                               f+i_coord_offset+DIM,fshift+i_shift_offset);
549
550         ggid                        = gid[iidx];
551         /* Update potential energies */
552         gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
553
554         /* Increment number of inner iterations */
555         inneriter                  += j_index_end - j_index_start;
556
557         /* Outer loop uses 19 flops */
558     }
559
560     /* Increment number of outer iterations */
561     outeriter        += nri;
562
563     /* Update outer/inner flops */
564
565     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_VF,outeriter*19 + inneriter*135);
566 }
567 /*
568  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwNone_GeomW4P1_F_avx_128_fma_double
569  * Electrostatics interaction: Ewald
570  * VdW interaction:            None
571  * Geometry:                   Water4-Particle
572  * Calculate force/pot:        Force
573  */
574 void
575 nb_kernel_ElecEw_VdwNone_GeomW4P1_F_avx_128_fma_double
576                     (t_nblist                    * gmx_restrict       nlist,
577                      rvec                        * gmx_restrict          xx,
578                      rvec                        * gmx_restrict          ff,
579                      t_forcerec                  * gmx_restrict          fr,
580                      t_mdatoms                   * gmx_restrict     mdatoms,
581                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
582                      t_nrnb                      * gmx_restrict        nrnb)
583 {
584     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
585      * just 0 for non-waters.
586      * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
587      * jnr indices corresponding to data put in the four positions in the SIMD register.
588      */
589     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
590     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
591     int              jnrA,jnrB;
592     int              j_coord_offsetA,j_coord_offsetB;
593     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
594     real             rcutoff_scalar;
595     real             *shiftvec,*fshift,*x,*f;
596     __m128d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
597     int              vdwioffset1;
598     __m128d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
599     int              vdwioffset2;
600     __m128d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
601     int              vdwioffset3;
602     __m128d          ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
603     int              vdwjidx0A,vdwjidx0B;
604     __m128d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
605     __m128d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
606     __m128d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
607     __m128d          dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
608     __m128d          velec,felec,velecsum,facel,crf,krf,krf2;
609     real             *charge;
610     __m128i          ewitab;
611     __m128d          ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
612     real             *ewtab;
613     __m128d          dummy_mask,cutoff_mask;
614     __m128d          signbit   = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
615     __m128d          one     = _mm_set1_pd(1.0);
616     __m128d          two     = _mm_set1_pd(2.0);
617     x                = xx[0];
618     f                = ff[0];
619
620     nri              = nlist->nri;
621     iinr             = nlist->iinr;
622     jindex           = nlist->jindex;
623     jjnr             = nlist->jjnr;
624     shiftidx         = nlist->shift;
625     gid              = nlist->gid;
626     shiftvec         = fr->shift_vec[0];
627     fshift           = fr->fshift[0];
628     facel            = _mm_set1_pd(fr->epsfac);
629     charge           = mdatoms->chargeA;
630
631     sh_ewald         = _mm_set1_pd(fr->ic->sh_ewald);
632     ewtab            = fr->ic->tabq_coul_F;
633     ewtabscale       = _mm_set1_pd(fr->ic->tabq_scale);
634     ewtabhalfspace   = _mm_set1_pd(0.5/fr->ic->tabq_scale);
635
636     /* Setup water-specific parameters */
637     inr              = nlist->iinr[0];
638     iq1              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+1]));
639     iq2              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+2]));
640     iq3              = _mm_mul_pd(facel,_mm_set1_pd(charge[inr+3]));
641
642     /* Avoid stupid compiler warnings */
643     jnrA = jnrB = 0;
644     j_coord_offsetA = 0;
645     j_coord_offsetB = 0;
646
647     outeriter        = 0;
648     inneriter        = 0;
649
650     /* Start outer loop over neighborlists */
651     for(iidx=0; iidx<nri; iidx++)
652     {
653         /* Load shift vector for this list */
654         i_shift_offset   = DIM*shiftidx[iidx];
655
656         /* Load limits for loop over neighbors */
657         j_index_start    = jindex[iidx];
658         j_index_end      = jindex[iidx+1];
659
660         /* Get outer coordinate index */
661         inr              = iinr[iidx];
662         i_coord_offset   = DIM*inr;
663
664         /* Load i particle coords and add shift vector */
665         gmx_mm_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
666                                                  &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
667
668         fix1             = _mm_setzero_pd();
669         fiy1             = _mm_setzero_pd();
670         fiz1             = _mm_setzero_pd();
671         fix2             = _mm_setzero_pd();
672         fiy2             = _mm_setzero_pd();
673         fiz2             = _mm_setzero_pd();
674         fix3             = _mm_setzero_pd();
675         fiy3             = _mm_setzero_pd();
676         fiz3             = _mm_setzero_pd();
677
678         /* Start inner kernel loop */
679         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
680         {
681
682             /* Get j neighbor index, and coordinate index */
683             jnrA             = jjnr[jidx];
684             jnrB             = jjnr[jidx+1];
685             j_coord_offsetA  = DIM*jnrA;
686             j_coord_offsetB  = DIM*jnrB;
687
688             /* load j atom coordinates */
689             gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
690                                               &jx0,&jy0,&jz0);
691
692             /* Calculate displacement vector */
693             dx10             = _mm_sub_pd(ix1,jx0);
694             dy10             = _mm_sub_pd(iy1,jy0);
695             dz10             = _mm_sub_pd(iz1,jz0);
696             dx20             = _mm_sub_pd(ix2,jx0);
697             dy20             = _mm_sub_pd(iy2,jy0);
698             dz20             = _mm_sub_pd(iz2,jz0);
699             dx30             = _mm_sub_pd(ix3,jx0);
700             dy30             = _mm_sub_pd(iy3,jy0);
701             dz30             = _mm_sub_pd(iz3,jz0);
702
703             /* Calculate squared distance and things based on it */
704             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
705             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
706             rsq30            = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
707
708             rinv10           = gmx_mm_invsqrt_pd(rsq10);
709             rinv20           = gmx_mm_invsqrt_pd(rsq20);
710             rinv30           = gmx_mm_invsqrt_pd(rsq30);
711
712             rinvsq10         = _mm_mul_pd(rinv10,rinv10);
713             rinvsq20         = _mm_mul_pd(rinv20,rinv20);
714             rinvsq30         = _mm_mul_pd(rinv30,rinv30);
715
716             /* Load parameters for j particles */
717             jq0              = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
718
719             fjx0             = _mm_setzero_pd();
720             fjy0             = _mm_setzero_pd();
721             fjz0             = _mm_setzero_pd();
722
723             /**************************
724              * CALCULATE INTERACTIONS *
725              **************************/
726
727             r10              = _mm_mul_pd(rsq10,rinv10);
728
729             /* Compute parameters for interactions between i and j atoms */
730             qq10             = _mm_mul_pd(iq1,jq0);
731
732             /* EWALD ELECTROSTATICS */
733
734             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
735             ewrt             = _mm_mul_pd(r10,ewtabscale);
736             ewitab           = _mm_cvttpd_epi32(ewrt);
737 #ifdef __XOP__
738             eweps            = _mm_frcz_pd(ewrt);
739 #else
740             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
741 #endif
742             twoeweps         = _mm_add_pd(eweps,eweps);
743             gmx_mm_load_2pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),ewtab+_mm_extract_epi32(ewitab,1),
744                                          &ewtabF,&ewtabFn);
745             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
746             felec            = _mm_mul_pd(_mm_mul_pd(qq10,rinv10),_mm_sub_pd(rinvsq10,felec));
747
748             fscal            = felec;
749
750             /* Update vectorial force */
751             fix1             = _mm_macc_pd(dx10,fscal,fix1);
752             fiy1             = _mm_macc_pd(dy10,fscal,fiy1);
753             fiz1             = _mm_macc_pd(dz10,fscal,fiz1);
754             
755             fjx0             = _mm_macc_pd(dx10,fscal,fjx0);
756             fjy0             = _mm_macc_pd(dy10,fscal,fjy0);
757             fjz0             = _mm_macc_pd(dz10,fscal,fjz0);
758
759             /**************************
760              * CALCULATE INTERACTIONS *
761              **************************/
762
763             r20              = _mm_mul_pd(rsq20,rinv20);
764
765             /* Compute parameters for interactions between i and j atoms */
766             qq20             = _mm_mul_pd(iq2,jq0);
767
768             /* EWALD ELECTROSTATICS */
769
770             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
771             ewrt             = _mm_mul_pd(r20,ewtabscale);
772             ewitab           = _mm_cvttpd_epi32(ewrt);
773 #ifdef __XOP__
774             eweps            = _mm_frcz_pd(ewrt);
775 #else
776             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
777 #endif
778             twoeweps         = _mm_add_pd(eweps,eweps);
779             gmx_mm_load_2pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),ewtab+_mm_extract_epi32(ewitab,1),
780                                          &ewtabF,&ewtabFn);
781             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
782             felec            = _mm_mul_pd(_mm_mul_pd(qq20,rinv20),_mm_sub_pd(rinvsq20,felec));
783
784             fscal            = felec;
785
786             /* Update vectorial force */
787             fix2             = _mm_macc_pd(dx20,fscal,fix2);
788             fiy2             = _mm_macc_pd(dy20,fscal,fiy2);
789             fiz2             = _mm_macc_pd(dz20,fscal,fiz2);
790             
791             fjx0             = _mm_macc_pd(dx20,fscal,fjx0);
792             fjy0             = _mm_macc_pd(dy20,fscal,fjy0);
793             fjz0             = _mm_macc_pd(dz20,fscal,fjz0);
794
795             /**************************
796              * CALCULATE INTERACTIONS *
797              **************************/
798
799             r30              = _mm_mul_pd(rsq30,rinv30);
800
801             /* Compute parameters for interactions between i and j atoms */
802             qq30             = _mm_mul_pd(iq3,jq0);
803
804             /* EWALD ELECTROSTATICS */
805
806             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
807             ewrt             = _mm_mul_pd(r30,ewtabscale);
808             ewitab           = _mm_cvttpd_epi32(ewrt);
809 #ifdef __XOP__
810             eweps            = _mm_frcz_pd(ewrt);
811 #else
812             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
813 #endif
814             twoeweps         = _mm_add_pd(eweps,eweps);
815             gmx_mm_load_2pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),ewtab+_mm_extract_epi32(ewitab,1),
816                                          &ewtabF,&ewtabFn);
817             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
818             felec            = _mm_mul_pd(_mm_mul_pd(qq30,rinv30),_mm_sub_pd(rinvsq30,felec));
819
820             fscal            = felec;
821
822             /* Update vectorial force */
823             fix3             = _mm_macc_pd(dx30,fscal,fix3);
824             fiy3             = _mm_macc_pd(dy30,fscal,fiy3);
825             fiz3             = _mm_macc_pd(dz30,fscal,fiz3);
826             
827             fjx0             = _mm_macc_pd(dx30,fscal,fjx0);
828             fjy0             = _mm_macc_pd(dy30,fscal,fjy0);
829             fjz0             = _mm_macc_pd(dz30,fscal,fjz0);
830
831             gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,fjx0,fjy0,fjz0);
832
833             /* Inner loop uses 120 flops */
834         }
835
836         if(jidx<j_index_end)
837         {
838
839             jnrA             = jjnr[jidx];
840             j_coord_offsetA  = DIM*jnrA;
841
842             /* load j atom coordinates */
843             gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
844                                               &jx0,&jy0,&jz0);
845
846             /* Calculate displacement vector */
847             dx10             = _mm_sub_pd(ix1,jx0);
848             dy10             = _mm_sub_pd(iy1,jy0);
849             dz10             = _mm_sub_pd(iz1,jz0);
850             dx20             = _mm_sub_pd(ix2,jx0);
851             dy20             = _mm_sub_pd(iy2,jy0);
852             dz20             = _mm_sub_pd(iz2,jz0);
853             dx30             = _mm_sub_pd(ix3,jx0);
854             dy30             = _mm_sub_pd(iy3,jy0);
855             dz30             = _mm_sub_pd(iz3,jz0);
856
857             /* Calculate squared distance and things based on it */
858             rsq10            = gmx_mm_calc_rsq_pd(dx10,dy10,dz10);
859             rsq20            = gmx_mm_calc_rsq_pd(dx20,dy20,dz20);
860             rsq30            = gmx_mm_calc_rsq_pd(dx30,dy30,dz30);
861
862             rinv10           = gmx_mm_invsqrt_pd(rsq10);
863             rinv20           = gmx_mm_invsqrt_pd(rsq20);
864             rinv30           = gmx_mm_invsqrt_pd(rsq30);
865
866             rinvsq10         = _mm_mul_pd(rinv10,rinv10);
867             rinvsq20         = _mm_mul_pd(rinv20,rinv20);
868             rinvsq30         = _mm_mul_pd(rinv30,rinv30);
869
870             /* Load parameters for j particles */
871             jq0              = _mm_load_sd(charge+jnrA+0);
872
873             fjx0             = _mm_setzero_pd();
874             fjy0             = _mm_setzero_pd();
875             fjz0             = _mm_setzero_pd();
876
877             /**************************
878              * CALCULATE INTERACTIONS *
879              **************************/
880
881             r10              = _mm_mul_pd(rsq10,rinv10);
882
883             /* Compute parameters for interactions between i and j atoms */
884             qq10             = _mm_mul_pd(iq1,jq0);
885
886             /* EWALD ELECTROSTATICS */
887
888             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
889             ewrt             = _mm_mul_pd(r10,ewtabscale);
890             ewitab           = _mm_cvttpd_epi32(ewrt);
891 #ifdef __XOP__
892             eweps            = _mm_frcz_pd(ewrt);
893 #else
894             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
895 #endif
896             twoeweps         = _mm_add_pd(eweps,eweps);
897             gmx_mm_load_1pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
898             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
899             felec            = _mm_mul_pd(_mm_mul_pd(qq10,rinv10),_mm_sub_pd(rinvsq10,felec));
900
901             fscal            = felec;
902
903             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
904
905             /* Update vectorial force */
906             fix1             = _mm_macc_pd(dx10,fscal,fix1);
907             fiy1             = _mm_macc_pd(dy10,fscal,fiy1);
908             fiz1             = _mm_macc_pd(dz10,fscal,fiz1);
909             
910             fjx0             = _mm_macc_pd(dx10,fscal,fjx0);
911             fjy0             = _mm_macc_pd(dy10,fscal,fjy0);
912             fjz0             = _mm_macc_pd(dz10,fscal,fjz0);
913
914             /**************************
915              * CALCULATE INTERACTIONS *
916              **************************/
917
918             r20              = _mm_mul_pd(rsq20,rinv20);
919
920             /* Compute parameters for interactions between i and j atoms */
921             qq20             = _mm_mul_pd(iq2,jq0);
922
923             /* EWALD ELECTROSTATICS */
924
925             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
926             ewrt             = _mm_mul_pd(r20,ewtabscale);
927             ewitab           = _mm_cvttpd_epi32(ewrt);
928 #ifdef __XOP__
929             eweps            = _mm_frcz_pd(ewrt);
930 #else
931             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
932 #endif
933             twoeweps         = _mm_add_pd(eweps,eweps);
934             gmx_mm_load_1pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
935             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
936             felec            = _mm_mul_pd(_mm_mul_pd(qq20,rinv20),_mm_sub_pd(rinvsq20,felec));
937
938             fscal            = felec;
939
940             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
941
942             /* Update vectorial force */
943             fix2             = _mm_macc_pd(dx20,fscal,fix2);
944             fiy2             = _mm_macc_pd(dy20,fscal,fiy2);
945             fiz2             = _mm_macc_pd(dz20,fscal,fiz2);
946             
947             fjx0             = _mm_macc_pd(dx20,fscal,fjx0);
948             fjy0             = _mm_macc_pd(dy20,fscal,fjy0);
949             fjz0             = _mm_macc_pd(dz20,fscal,fjz0);
950
951             /**************************
952              * CALCULATE INTERACTIONS *
953              **************************/
954
955             r30              = _mm_mul_pd(rsq30,rinv30);
956
957             /* Compute parameters for interactions between i and j atoms */
958             qq30             = _mm_mul_pd(iq3,jq0);
959
960             /* EWALD ELECTROSTATICS */
961
962             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
963             ewrt             = _mm_mul_pd(r30,ewtabscale);
964             ewitab           = _mm_cvttpd_epi32(ewrt);
965 #ifdef __XOP__
966             eweps            = _mm_frcz_pd(ewrt);
967 #else
968             eweps            = _mm_sub_pd(ewrt,_mm_round_pd(ewrt, _MM_FROUND_FLOOR));
969 #endif
970             twoeweps         = _mm_add_pd(eweps,eweps);
971             gmx_mm_load_1pair_swizzle_pd(ewtab+_mm_extract_epi32(ewitab,0),&ewtabF,&ewtabFn);
972             felec            = _mm_macc_pd(eweps,ewtabFn,_mm_mul_pd( _mm_sub_pd(one,eweps),ewtabF));
973             felec            = _mm_mul_pd(_mm_mul_pd(qq30,rinv30),_mm_sub_pd(rinvsq30,felec));
974
975             fscal            = felec;
976
977             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
978
979             /* Update vectorial force */
980             fix3             = _mm_macc_pd(dx30,fscal,fix3);
981             fiy3             = _mm_macc_pd(dy30,fscal,fiy3);
982             fiz3             = _mm_macc_pd(dz30,fscal,fiz3);
983             
984             fjx0             = _mm_macc_pd(dx30,fscal,fjx0);
985             fjy0             = _mm_macc_pd(dy30,fscal,fjy0);
986             fjz0             = _mm_macc_pd(dz30,fscal,fjz0);
987
988             gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,fjx0,fjy0,fjz0);
989
990             /* Inner loop uses 120 flops */
991         }
992
993         /* End of innermost loop */
994
995         gmx_mm_update_iforce_3atom_swizzle_pd(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
996                                               f+i_coord_offset+DIM,fshift+i_shift_offset);
997
998         /* Increment number of inner iterations */
999         inneriter                  += j_index_end - j_index_start;
1000
1001         /* Outer loop uses 18 flops */
1002     }
1003
1004     /* Increment number of outer iterations */
1005     outeriter        += nri;
1006
1007     /* Update outer/inner flops */
1008
1009     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W4_F,outeriter*18 + inneriter*120);
1010 }