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