Remove all unnecessary HAVE_CONFIG_H
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_double / nb_kernel_ElecCoul_VdwLJ_GeomP1P1_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_ElecCoul_VdwLJ_GeomP1P1_VF_avx_128_fma_double
52  * Electrostatics interaction: Coulomb
53  * VdW interaction:            LennardJones
54  * Geometry:                   Particle-Particle
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecCoul_VdwLJ_GeomP1P1_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              vdwioffset0;
81     __m128d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
82     int              vdwjidx0A,vdwjidx0B;
83     __m128d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
84     __m128d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
85     __m128d          velec,felec,velecsum,facel,crf,krf,krf2;
86     real             *charge;
87     int              nvdwtype;
88     __m128d          rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
89     int              *vdwtype;
90     real             *vdwparam;
91     __m128d          one_sixth   = _mm_set1_pd(1.0/6.0);
92     __m128d          one_twelfth = _mm_set1_pd(1.0/12.0);
93     __m128d          dummy_mask,cutoff_mask;
94     __m128d          signbit   = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
95     __m128d          one     = _mm_set1_pd(1.0);
96     __m128d          two     = _mm_set1_pd(2.0);
97     x                = xx[0];
98     f                = ff[0];
99
100     nri              = nlist->nri;
101     iinr             = nlist->iinr;
102     jindex           = nlist->jindex;
103     jjnr             = nlist->jjnr;
104     shiftidx         = nlist->shift;
105     gid              = nlist->gid;
106     shiftvec         = fr->shift_vec[0];
107     fshift           = fr->fshift[0];
108     facel            = _mm_set1_pd(fr->epsfac);
109     charge           = mdatoms->chargeA;
110     nvdwtype         = fr->ntype;
111     vdwparam         = fr->nbfp;
112     vdwtype          = mdatoms->typeA;
113
114     /* Avoid stupid compiler warnings */
115     jnrA = jnrB = 0;
116     j_coord_offsetA = 0;
117     j_coord_offsetB = 0;
118
119     outeriter        = 0;
120     inneriter        = 0;
121
122     /* Start outer loop over neighborlists */
123     for(iidx=0; iidx<nri; iidx++)
124     {
125         /* Load shift vector for this list */
126         i_shift_offset   = DIM*shiftidx[iidx];
127
128         /* Load limits for loop over neighbors */
129         j_index_start    = jindex[iidx];
130         j_index_end      = jindex[iidx+1];
131
132         /* Get outer coordinate index */
133         inr              = iinr[iidx];
134         i_coord_offset   = DIM*inr;
135
136         /* Load i particle coords and add shift vector */
137         gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
138
139         fix0             = _mm_setzero_pd();
140         fiy0             = _mm_setzero_pd();
141         fiz0             = _mm_setzero_pd();
142
143         /* Load parameters for i particles */
144         iq0              = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
145         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
146
147         /* Reset potential sums */
148         velecsum         = _mm_setzero_pd();
149         vvdwsum          = _mm_setzero_pd();
150
151         /* Start inner kernel loop */
152         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
153         {
154
155             /* Get j neighbor index, and coordinate index */
156             jnrA             = jjnr[jidx];
157             jnrB             = jjnr[jidx+1];
158             j_coord_offsetA  = DIM*jnrA;
159             j_coord_offsetB  = DIM*jnrB;
160
161             /* load j atom coordinates */
162             gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
163                                               &jx0,&jy0,&jz0);
164
165             /* Calculate displacement vector */
166             dx00             = _mm_sub_pd(ix0,jx0);
167             dy00             = _mm_sub_pd(iy0,jy0);
168             dz00             = _mm_sub_pd(iz0,jz0);
169
170             /* Calculate squared distance and things based on it */
171             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
172
173             rinv00           = gmx_mm_invsqrt_pd(rsq00);
174
175             rinvsq00         = _mm_mul_pd(rinv00,rinv00);
176
177             /* Load parameters for j particles */
178             jq0              = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
179             vdwjidx0A        = 2*vdwtype[jnrA+0];
180             vdwjidx0B        = 2*vdwtype[jnrB+0];
181
182             /**************************
183              * CALCULATE INTERACTIONS *
184              **************************/
185
186             /* Compute parameters for interactions between i and j atoms */
187             qq00             = _mm_mul_pd(iq0,jq0);
188             gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
189                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
190
191             /* COULOMB ELECTROSTATICS */
192             velec            = _mm_mul_pd(qq00,rinv00);
193             felec            = _mm_mul_pd(velec,rinvsq00);
194
195             /* LENNARD-JONES DISPERSION/REPULSION */
196
197             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
198             vvdw6            = _mm_mul_pd(c6_00,rinvsix);
199             vvdw12           = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
200             vvdw             = _mm_msub_pd( vvdw12,one_twelfth, _mm_mul_pd(vvdw6,one_sixth) );
201             fvdw             = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
202
203             /* Update potential sum for this i atom from the interaction with this j atom. */
204             velecsum         = _mm_add_pd(velecsum,velec);
205             vvdwsum          = _mm_add_pd(vvdwsum,vvdw);
206
207             fscal            = _mm_add_pd(felec,fvdw);
208
209             /* Update vectorial force */
210             fix0             = _mm_macc_pd(dx00,fscal,fix0);
211             fiy0             = _mm_macc_pd(dy00,fscal,fiy0);
212             fiz0             = _mm_macc_pd(dz00,fscal,fiz0);
213             
214             gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
215                                                    _mm_mul_pd(dx00,fscal),
216                                                    _mm_mul_pd(dy00,fscal),
217                                                    _mm_mul_pd(dz00,fscal));
218
219             /* Inner loop uses 43 flops */
220         }
221
222         if(jidx<j_index_end)
223         {
224
225             jnrA             = jjnr[jidx];
226             j_coord_offsetA  = DIM*jnrA;
227
228             /* load j atom coordinates */
229             gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
230                                               &jx0,&jy0,&jz0);
231
232             /* Calculate displacement vector */
233             dx00             = _mm_sub_pd(ix0,jx0);
234             dy00             = _mm_sub_pd(iy0,jy0);
235             dz00             = _mm_sub_pd(iz0,jz0);
236
237             /* Calculate squared distance and things based on it */
238             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
239
240             rinv00           = gmx_mm_invsqrt_pd(rsq00);
241
242             rinvsq00         = _mm_mul_pd(rinv00,rinv00);
243
244             /* Load parameters for j particles */
245             jq0              = _mm_load_sd(charge+jnrA+0);
246             vdwjidx0A        = 2*vdwtype[jnrA+0];
247
248             /**************************
249              * CALCULATE INTERACTIONS *
250              **************************/
251
252             /* Compute parameters for interactions between i and j atoms */
253             qq00             = _mm_mul_pd(iq0,jq0);
254             gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
255
256             /* COULOMB ELECTROSTATICS */
257             velec            = _mm_mul_pd(qq00,rinv00);
258             felec            = _mm_mul_pd(velec,rinvsq00);
259
260             /* LENNARD-JONES DISPERSION/REPULSION */
261
262             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
263             vvdw6            = _mm_mul_pd(c6_00,rinvsix);
264             vvdw12           = _mm_mul_pd(c12_00,_mm_mul_pd(rinvsix,rinvsix));
265             vvdw             = _mm_msub_pd( vvdw12,one_twelfth, _mm_mul_pd(vvdw6,one_sixth) );
266             fvdw             = _mm_mul_pd(_mm_sub_pd(vvdw12,vvdw6),rinvsq00);
267
268             /* Update potential sum for this i atom from the interaction with this j atom. */
269             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
270             velecsum         = _mm_add_pd(velecsum,velec);
271             vvdw             = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
272             vvdwsum          = _mm_add_pd(vvdwsum,vvdw);
273
274             fscal            = _mm_add_pd(felec,fvdw);
275
276             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
277
278             /* Update vectorial force */
279             fix0             = _mm_macc_pd(dx00,fscal,fix0);
280             fiy0             = _mm_macc_pd(dy00,fscal,fiy0);
281             fiz0             = _mm_macc_pd(dz00,fscal,fiz0);
282             
283             gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
284                                                    _mm_mul_pd(dx00,fscal),
285                                                    _mm_mul_pd(dy00,fscal),
286                                                    _mm_mul_pd(dz00,fscal));
287
288             /* Inner loop uses 43 flops */
289         }
290
291         /* End of innermost loop */
292
293         gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
294                                               f+i_coord_offset,fshift+i_shift_offset);
295
296         ggid                        = gid[iidx];
297         /* Update potential energies */
298         gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
299         gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
300
301         /* Increment number of inner iterations */
302         inneriter                  += j_index_end - j_index_start;
303
304         /* Outer loop uses 9 flops */
305     }
306
307     /* Increment number of outer iterations */
308     outeriter        += nri;
309
310     /* Update outer/inner flops */
311
312     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*43);
313 }
314 /*
315  * Gromacs nonbonded kernel:   nb_kernel_ElecCoul_VdwLJ_GeomP1P1_F_avx_128_fma_double
316  * Electrostatics interaction: Coulomb
317  * VdW interaction:            LennardJones
318  * Geometry:                   Particle-Particle
319  * Calculate force/pot:        Force
320  */
321 void
322 nb_kernel_ElecCoul_VdwLJ_GeomP1P1_F_avx_128_fma_double
323                     (t_nblist                    * gmx_restrict       nlist,
324                      rvec                        * gmx_restrict          xx,
325                      rvec                        * gmx_restrict          ff,
326                      t_forcerec                  * gmx_restrict          fr,
327                      t_mdatoms                   * gmx_restrict     mdatoms,
328                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
329                      t_nrnb                      * gmx_restrict        nrnb)
330 {
331     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
332      * just 0 for non-waters.
333      * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
334      * jnr indices corresponding to data put in the four positions in the SIMD register.
335      */
336     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
337     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
338     int              jnrA,jnrB;
339     int              j_coord_offsetA,j_coord_offsetB;
340     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
341     real             rcutoff_scalar;
342     real             *shiftvec,*fshift,*x,*f;
343     __m128d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
344     int              vdwioffset0;
345     __m128d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
346     int              vdwjidx0A,vdwjidx0B;
347     __m128d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
348     __m128d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
349     __m128d          velec,felec,velecsum,facel,crf,krf,krf2;
350     real             *charge;
351     int              nvdwtype;
352     __m128d          rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
353     int              *vdwtype;
354     real             *vdwparam;
355     __m128d          one_sixth   = _mm_set1_pd(1.0/6.0);
356     __m128d          one_twelfth = _mm_set1_pd(1.0/12.0);
357     __m128d          dummy_mask,cutoff_mask;
358     __m128d          signbit   = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
359     __m128d          one     = _mm_set1_pd(1.0);
360     __m128d          two     = _mm_set1_pd(2.0);
361     x                = xx[0];
362     f                = ff[0];
363
364     nri              = nlist->nri;
365     iinr             = nlist->iinr;
366     jindex           = nlist->jindex;
367     jjnr             = nlist->jjnr;
368     shiftidx         = nlist->shift;
369     gid              = nlist->gid;
370     shiftvec         = fr->shift_vec[0];
371     fshift           = fr->fshift[0];
372     facel            = _mm_set1_pd(fr->epsfac);
373     charge           = mdatoms->chargeA;
374     nvdwtype         = fr->ntype;
375     vdwparam         = fr->nbfp;
376     vdwtype          = mdatoms->typeA;
377
378     /* Avoid stupid compiler warnings */
379     jnrA = jnrB = 0;
380     j_coord_offsetA = 0;
381     j_coord_offsetB = 0;
382
383     outeriter        = 0;
384     inneriter        = 0;
385
386     /* Start outer loop over neighborlists */
387     for(iidx=0; iidx<nri; iidx++)
388     {
389         /* Load shift vector for this list */
390         i_shift_offset   = DIM*shiftidx[iidx];
391
392         /* Load limits for loop over neighbors */
393         j_index_start    = jindex[iidx];
394         j_index_end      = jindex[iidx+1];
395
396         /* Get outer coordinate index */
397         inr              = iinr[iidx];
398         i_coord_offset   = DIM*inr;
399
400         /* Load i particle coords and add shift vector */
401         gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
402
403         fix0             = _mm_setzero_pd();
404         fiy0             = _mm_setzero_pd();
405         fiz0             = _mm_setzero_pd();
406
407         /* Load parameters for i particles */
408         iq0              = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
409         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
410
411         /* Start inner kernel loop */
412         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
413         {
414
415             /* Get j neighbor index, and coordinate index */
416             jnrA             = jjnr[jidx];
417             jnrB             = jjnr[jidx+1];
418             j_coord_offsetA  = DIM*jnrA;
419             j_coord_offsetB  = DIM*jnrB;
420
421             /* load j atom coordinates */
422             gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
423                                               &jx0,&jy0,&jz0);
424
425             /* Calculate displacement vector */
426             dx00             = _mm_sub_pd(ix0,jx0);
427             dy00             = _mm_sub_pd(iy0,jy0);
428             dz00             = _mm_sub_pd(iz0,jz0);
429
430             /* Calculate squared distance and things based on it */
431             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
432
433             rinv00           = gmx_mm_invsqrt_pd(rsq00);
434
435             rinvsq00         = _mm_mul_pd(rinv00,rinv00);
436
437             /* Load parameters for j particles */
438             jq0              = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
439             vdwjidx0A        = 2*vdwtype[jnrA+0];
440             vdwjidx0B        = 2*vdwtype[jnrB+0];
441
442             /**************************
443              * CALCULATE INTERACTIONS *
444              **************************/
445
446             /* Compute parameters for interactions between i and j atoms */
447             qq00             = _mm_mul_pd(iq0,jq0);
448             gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
449                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
450
451             /* COULOMB ELECTROSTATICS */
452             velec            = _mm_mul_pd(qq00,rinv00);
453             felec            = _mm_mul_pd(velec,rinvsq00);
454
455             /* LENNARD-JONES DISPERSION/REPULSION */
456
457             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
458             fvdw             = _mm_mul_pd(_mm_msub_pd(c12_00,rinvsix,c6_00),_mm_mul_pd(rinvsix,rinvsq00));
459
460             fscal            = _mm_add_pd(felec,fvdw);
461
462             /* Update vectorial force */
463             fix0             = _mm_macc_pd(dx00,fscal,fix0);
464             fiy0             = _mm_macc_pd(dy00,fscal,fiy0);
465             fiz0             = _mm_macc_pd(dz00,fscal,fiz0);
466             
467             gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,
468                                                    _mm_mul_pd(dx00,fscal),
469                                                    _mm_mul_pd(dy00,fscal),
470                                                    _mm_mul_pd(dz00,fscal));
471
472             /* Inner loop uses 37 flops */
473         }
474
475         if(jidx<j_index_end)
476         {
477
478             jnrA             = jjnr[jidx];
479             j_coord_offsetA  = DIM*jnrA;
480
481             /* load j atom coordinates */
482             gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
483                                               &jx0,&jy0,&jz0);
484
485             /* Calculate displacement vector */
486             dx00             = _mm_sub_pd(ix0,jx0);
487             dy00             = _mm_sub_pd(iy0,jy0);
488             dz00             = _mm_sub_pd(iz0,jz0);
489
490             /* Calculate squared distance and things based on it */
491             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
492
493             rinv00           = gmx_mm_invsqrt_pd(rsq00);
494
495             rinvsq00         = _mm_mul_pd(rinv00,rinv00);
496
497             /* Load parameters for j particles */
498             jq0              = _mm_load_sd(charge+jnrA+0);
499             vdwjidx0A        = 2*vdwtype[jnrA+0];
500
501             /**************************
502              * CALCULATE INTERACTIONS *
503              **************************/
504
505             /* Compute parameters for interactions between i and j atoms */
506             qq00             = _mm_mul_pd(iq0,jq0);
507             gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
508
509             /* COULOMB ELECTROSTATICS */
510             velec            = _mm_mul_pd(qq00,rinv00);
511             felec            = _mm_mul_pd(velec,rinvsq00);
512
513             /* LENNARD-JONES DISPERSION/REPULSION */
514
515             rinvsix          = _mm_mul_pd(_mm_mul_pd(rinvsq00,rinvsq00),rinvsq00);
516             fvdw             = _mm_mul_pd(_mm_msub_pd(c12_00,rinvsix,c6_00),_mm_mul_pd(rinvsix,rinvsq00));
517
518             fscal            = _mm_add_pd(felec,fvdw);
519
520             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
521
522             /* Update vectorial force */
523             fix0             = _mm_macc_pd(dx00,fscal,fix0);
524             fiy0             = _mm_macc_pd(dy00,fscal,fiy0);
525             fiz0             = _mm_macc_pd(dz00,fscal,fiz0);
526             
527             gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,
528                                                    _mm_mul_pd(dx00,fscal),
529                                                    _mm_mul_pd(dy00,fscal),
530                                                    _mm_mul_pd(dz00,fscal));
531
532             /* Inner loop uses 37 flops */
533         }
534
535         /* End of innermost loop */
536
537         gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
538                                               f+i_coord_offset,fshift+i_shift_offset);
539
540         /* Increment number of inner iterations */
541         inneriter                  += j_index_end - j_index_start;
542
543         /* Outer loop uses 7 flops */
544     }
545
546     /* Increment number of outer iterations */
547     outeriter        += nri;
548
549     /* Update outer/inner flops */
550
551     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*37);
552 }