995dbac81f08132826249ab807173d8918746571
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_sse4_1_double / nb_kernel_ElecCSTab_VdwCSTab_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_ElecCSTab_VdwCSTab_GeomP1P1_VF_sse4_1_double
52  * Electrostatics interaction: CubicSplineTable
53  * VdW interaction:            CubicSplineTable
54  * Geometry:                   Particle-Particle
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecCSTab_VdwCSTab_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     __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     __m128i          vfitab;
94     __m128i          ifour       = _mm_set1_epi32(4);
95     __m128d          rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
96     real             *vftab;
97     __m128d          dummy_mask,cutoff_mask;
98     __m128d          signbit   = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
99     __m128d          one     = _mm_set1_pd(1.0);
100     __m128d          two     = _mm_set1_pd(2.0);
101     x                = xx[0];
102     f                = ff[0];
103
104     nri              = nlist->nri;
105     iinr             = nlist->iinr;
106     jindex           = nlist->jindex;
107     jjnr             = nlist->jjnr;
108     shiftidx         = nlist->shift;
109     gid              = nlist->gid;
110     shiftvec         = fr->shift_vec[0];
111     fshift           = fr->fshift[0];
112     facel            = _mm_set1_pd(fr->epsfac);
113     charge           = mdatoms->chargeA;
114     nvdwtype         = fr->ntype;
115     vdwparam         = fr->nbfp;
116     vdwtype          = mdatoms->typeA;
117
118     vftab            = kernel_data->table_elec_vdw->data;
119     vftabscale       = _mm_set1_pd(kernel_data->table_elec_vdw->scale);
120
121     /* Avoid stupid compiler warnings */
122     jnrA = jnrB = 0;
123     j_coord_offsetA = 0;
124     j_coord_offsetB = 0;
125
126     outeriter        = 0;
127     inneriter        = 0;
128
129     /* Start outer loop over neighborlists */
130     for(iidx=0; iidx<nri; iidx++)
131     {
132         /* Load shift vector for this list */
133         i_shift_offset   = DIM*shiftidx[iidx];
134
135         /* Load limits for loop over neighbors */
136         j_index_start    = jindex[iidx];
137         j_index_end      = jindex[iidx+1];
138
139         /* Get outer coordinate index */
140         inr              = iinr[iidx];
141         i_coord_offset   = DIM*inr;
142
143         /* Load i particle coords and add shift vector */
144         gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
145
146         fix0             = _mm_setzero_pd();
147         fiy0             = _mm_setzero_pd();
148         fiz0             = _mm_setzero_pd();
149
150         /* Load parameters for i particles */
151         iq0              = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
152         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
153
154         /* Reset potential sums */
155         velecsum         = _mm_setzero_pd();
156         vvdwsum          = _mm_setzero_pd();
157
158         /* Start inner kernel loop */
159         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
160         {
161
162             /* Get j neighbor index, and coordinate index */
163             jnrA             = jjnr[jidx];
164             jnrB             = jjnr[jidx+1];
165             j_coord_offsetA  = DIM*jnrA;
166             j_coord_offsetB  = DIM*jnrB;
167
168             /* load j atom coordinates */
169             gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
170                                               &jx0,&jy0,&jz0);
171
172             /* Calculate displacement vector */
173             dx00             = _mm_sub_pd(ix0,jx0);
174             dy00             = _mm_sub_pd(iy0,jy0);
175             dz00             = _mm_sub_pd(iz0,jz0);
176
177             /* Calculate squared distance and things based on it */
178             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
179
180             rinv00           = gmx_mm_invsqrt_pd(rsq00);
181
182             /* Load parameters for j particles */
183             jq0              = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
184             vdwjidx0A        = 2*vdwtype[jnrA+0];
185             vdwjidx0B        = 2*vdwtype[jnrB+0];
186
187             /**************************
188              * CALCULATE INTERACTIONS *
189              **************************/
190
191             r00              = _mm_mul_pd(rsq00,rinv00);
192
193             /* Compute parameters for interactions between i and j atoms */
194             qq00             = _mm_mul_pd(iq0,jq0);
195             gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
196                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
197
198             /* Calculate table index by multiplying r with table scale and truncate to integer */
199             rt               = _mm_mul_pd(r00,vftabscale);
200             vfitab           = _mm_cvttpd_epi32(rt);
201             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
202             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
203
204             /* CUBIC SPLINE TABLE ELECTROSTATICS */
205             Y                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
206             F                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
207             GMX_MM_TRANSPOSE2_PD(Y,F);
208             G                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
209             H                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
210             GMX_MM_TRANSPOSE2_PD(G,H);
211             Heps             = _mm_mul_pd(vfeps,H);
212             Fp               = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
213             VV               = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
214             velec            = _mm_mul_pd(qq00,VV);
215             FF               = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
216             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
217
218             /* CUBIC SPLINE TABLE DISPERSION */
219             vfitab           = _mm_add_epi32(vfitab,ifour);
220             Y                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
221             F                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
222             GMX_MM_TRANSPOSE2_PD(Y,F);
223             G                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
224             H                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
225             GMX_MM_TRANSPOSE2_PD(G,H);
226             Heps             = _mm_mul_pd(vfeps,H);
227             Fp               = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
228             VV               = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
229             vvdw6            = _mm_mul_pd(c6_00,VV);
230             FF               = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
231             fvdw6            = _mm_mul_pd(c6_00,FF);
232
233             /* CUBIC SPLINE TABLE REPULSION */
234             vfitab           = _mm_add_epi32(vfitab,ifour);
235             Y                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
236             F                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
237             GMX_MM_TRANSPOSE2_PD(Y,F);
238             G                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
239             H                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
240             GMX_MM_TRANSPOSE2_PD(G,H);
241             Heps             = _mm_mul_pd(vfeps,H);
242             Fp               = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
243             VV               = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
244             vvdw12           = _mm_mul_pd(c12_00,VV);
245             FF               = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
246             fvdw12           = _mm_mul_pd(c12_00,FF);
247             vvdw             = _mm_add_pd(vvdw12,vvdw6);
248             fvdw             = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
249
250             /* Update potential sum for this i atom from the interaction with this j atom. */
251             velecsum         = _mm_add_pd(velecsum,velec);
252             vvdwsum          = _mm_add_pd(vvdwsum,vvdw);
253
254             fscal            = _mm_add_pd(felec,fvdw);
255
256             /* Calculate temporary vectorial force */
257             tx               = _mm_mul_pd(fscal,dx00);
258             ty               = _mm_mul_pd(fscal,dy00);
259             tz               = _mm_mul_pd(fscal,dz00);
260
261             /* Update vectorial force */
262             fix0             = _mm_add_pd(fix0,tx);
263             fiy0             = _mm_add_pd(fiy0,ty);
264             fiz0             = _mm_add_pd(fiz0,tz);
265
266             gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
267
268             /* Inner loop uses 73 flops */
269         }
270
271         if(jidx<j_index_end)
272         {
273
274             jnrA             = jjnr[jidx];
275             j_coord_offsetA  = DIM*jnrA;
276
277             /* load j atom coordinates */
278             gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
279                                               &jx0,&jy0,&jz0);
280
281             /* Calculate displacement vector */
282             dx00             = _mm_sub_pd(ix0,jx0);
283             dy00             = _mm_sub_pd(iy0,jy0);
284             dz00             = _mm_sub_pd(iz0,jz0);
285
286             /* Calculate squared distance and things based on it */
287             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
288
289             rinv00           = gmx_mm_invsqrt_pd(rsq00);
290
291             /* Load parameters for j particles */
292             jq0              = _mm_load_sd(charge+jnrA+0);
293             vdwjidx0A        = 2*vdwtype[jnrA+0];
294
295             /**************************
296              * CALCULATE INTERACTIONS *
297              **************************/
298
299             r00              = _mm_mul_pd(rsq00,rinv00);
300
301             /* Compute parameters for interactions between i and j atoms */
302             qq00             = _mm_mul_pd(iq0,jq0);
303             gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
304
305             /* Calculate table index by multiplying r with table scale and truncate to integer */
306             rt               = _mm_mul_pd(r00,vftabscale);
307             vfitab           = _mm_cvttpd_epi32(rt);
308             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
309             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
310
311             /* CUBIC SPLINE TABLE ELECTROSTATICS */
312             Y                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
313             F                = _mm_setzero_pd();
314             GMX_MM_TRANSPOSE2_PD(Y,F);
315             G                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
316             H                = _mm_setzero_pd();
317             GMX_MM_TRANSPOSE2_PD(G,H);
318             Heps             = _mm_mul_pd(vfeps,H);
319             Fp               = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
320             VV               = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
321             velec            = _mm_mul_pd(qq00,VV);
322             FF               = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
323             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
324
325             /* CUBIC SPLINE TABLE DISPERSION */
326             vfitab           = _mm_add_epi32(vfitab,ifour);
327             Y                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
328             F                = _mm_setzero_pd();
329             GMX_MM_TRANSPOSE2_PD(Y,F);
330             G                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
331             H                = _mm_setzero_pd();
332             GMX_MM_TRANSPOSE2_PD(G,H);
333             Heps             = _mm_mul_pd(vfeps,H);
334             Fp               = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
335             VV               = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
336             vvdw6            = _mm_mul_pd(c6_00,VV);
337             FF               = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
338             fvdw6            = _mm_mul_pd(c6_00,FF);
339
340             /* CUBIC SPLINE TABLE REPULSION */
341             vfitab           = _mm_add_epi32(vfitab,ifour);
342             Y                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
343             F                = _mm_setzero_pd();
344             GMX_MM_TRANSPOSE2_PD(Y,F);
345             G                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
346             H                = _mm_setzero_pd();
347             GMX_MM_TRANSPOSE2_PD(G,H);
348             Heps             = _mm_mul_pd(vfeps,H);
349             Fp               = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
350             VV               = _mm_add_pd(Y,_mm_mul_pd(vfeps,Fp));
351             vvdw12           = _mm_mul_pd(c12_00,VV);
352             FF               = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
353             fvdw12           = _mm_mul_pd(c12_00,FF);
354             vvdw             = _mm_add_pd(vvdw12,vvdw6);
355             fvdw             = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
356
357             /* Update potential sum for this i atom from the interaction with this j atom. */
358             velec            = _mm_unpacklo_pd(velec,_mm_setzero_pd());
359             velecsum         = _mm_add_pd(velecsum,velec);
360             vvdw             = _mm_unpacklo_pd(vvdw,_mm_setzero_pd());
361             vvdwsum          = _mm_add_pd(vvdwsum,vvdw);
362
363             fscal            = _mm_add_pd(felec,fvdw);
364
365             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
366
367             /* Calculate temporary vectorial force */
368             tx               = _mm_mul_pd(fscal,dx00);
369             ty               = _mm_mul_pd(fscal,dy00);
370             tz               = _mm_mul_pd(fscal,dz00);
371
372             /* Update vectorial force */
373             fix0             = _mm_add_pd(fix0,tx);
374             fiy0             = _mm_add_pd(fiy0,ty);
375             fiz0             = _mm_add_pd(fiz0,tz);
376
377             gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
378
379             /* Inner loop uses 73 flops */
380         }
381
382         /* End of innermost loop */
383
384         gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
385                                               f+i_coord_offset,fshift+i_shift_offset);
386
387         ggid                        = gid[iidx];
388         /* Update potential energies */
389         gmx_mm_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
390         gmx_mm_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
391
392         /* Increment number of inner iterations */
393         inneriter                  += j_index_end - j_index_start;
394
395         /* Outer loop uses 9 flops */
396     }
397
398     /* Increment number of outer iterations */
399     outeriter        += nri;
400
401     /* Update outer/inner flops */
402
403     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_VF,outeriter*9 + inneriter*73);
404 }
405 /*
406  * Gromacs nonbonded kernel:   nb_kernel_ElecCSTab_VdwCSTab_GeomP1P1_F_sse4_1_double
407  * Electrostatics interaction: CubicSplineTable
408  * VdW interaction:            CubicSplineTable
409  * Geometry:                   Particle-Particle
410  * Calculate force/pot:        Force
411  */
412 void
413 nb_kernel_ElecCSTab_VdwCSTab_GeomP1P1_F_sse4_1_double
414                     (t_nblist                    * gmx_restrict       nlist,
415                      rvec                        * gmx_restrict          xx,
416                      rvec                        * gmx_restrict          ff,
417                      t_forcerec                  * gmx_restrict          fr,
418                      t_mdatoms                   * gmx_restrict     mdatoms,
419                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
420                      t_nrnb                      * gmx_restrict        nrnb)
421 {
422     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
423      * just 0 for non-waters.
424      * Suffixes A,B refer to j loop unrolling done with SSE double precision, e.g. for the two different
425      * jnr indices corresponding to data put in the four positions in the SIMD register.
426      */
427     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
428     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
429     int              jnrA,jnrB;
430     int              j_coord_offsetA,j_coord_offsetB;
431     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
432     real             rcutoff_scalar;
433     real             *shiftvec,*fshift,*x,*f;
434     __m128d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
435     int              vdwioffset0;
436     __m128d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
437     int              vdwjidx0A,vdwjidx0B;
438     __m128d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
439     __m128d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
440     __m128d          velec,felec,velecsum,facel,crf,krf,krf2;
441     real             *charge;
442     int              nvdwtype;
443     __m128d          rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
444     int              *vdwtype;
445     real             *vdwparam;
446     __m128d          one_sixth   = _mm_set1_pd(1.0/6.0);
447     __m128d          one_twelfth = _mm_set1_pd(1.0/12.0);
448     __m128i          vfitab;
449     __m128i          ifour       = _mm_set1_epi32(4);
450     __m128d          rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
451     real             *vftab;
452     __m128d          dummy_mask,cutoff_mask;
453     __m128d          signbit   = gmx_mm_castsi128_pd( _mm_set_epi32(0x80000000,0x00000000,0x80000000,0x00000000) );
454     __m128d          one     = _mm_set1_pd(1.0);
455     __m128d          two     = _mm_set1_pd(2.0);
456     x                = xx[0];
457     f                = ff[0];
458
459     nri              = nlist->nri;
460     iinr             = nlist->iinr;
461     jindex           = nlist->jindex;
462     jjnr             = nlist->jjnr;
463     shiftidx         = nlist->shift;
464     gid              = nlist->gid;
465     shiftvec         = fr->shift_vec[0];
466     fshift           = fr->fshift[0];
467     facel            = _mm_set1_pd(fr->epsfac);
468     charge           = mdatoms->chargeA;
469     nvdwtype         = fr->ntype;
470     vdwparam         = fr->nbfp;
471     vdwtype          = mdatoms->typeA;
472
473     vftab            = kernel_data->table_elec_vdw->data;
474     vftabscale       = _mm_set1_pd(kernel_data->table_elec_vdw->scale);
475
476     /* Avoid stupid compiler warnings */
477     jnrA = jnrB = 0;
478     j_coord_offsetA = 0;
479     j_coord_offsetB = 0;
480
481     outeriter        = 0;
482     inneriter        = 0;
483
484     /* Start outer loop over neighborlists */
485     for(iidx=0; iidx<nri; iidx++)
486     {
487         /* Load shift vector for this list */
488         i_shift_offset   = DIM*shiftidx[iidx];
489
490         /* Load limits for loop over neighbors */
491         j_index_start    = jindex[iidx];
492         j_index_end      = jindex[iidx+1];
493
494         /* Get outer coordinate index */
495         inr              = iinr[iidx];
496         i_coord_offset   = DIM*inr;
497
498         /* Load i particle coords and add shift vector */
499         gmx_mm_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
500
501         fix0             = _mm_setzero_pd();
502         fiy0             = _mm_setzero_pd();
503         fiz0             = _mm_setzero_pd();
504
505         /* Load parameters for i particles */
506         iq0              = _mm_mul_pd(facel,_mm_load1_pd(charge+inr+0));
507         vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
508
509         /* Start inner kernel loop */
510         for(jidx=j_index_start; jidx<j_index_end-1; jidx+=2)
511         {
512
513             /* Get j neighbor index, and coordinate index */
514             jnrA             = jjnr[jidx];
515             jnrB             = jjnr[jidx+1];
516             j_coord_offsetA  = DIM*jnrA;
517             j_coord_offsetB  = DIM*jnrB;
518
519             /* load j atom coordinates */
520             gmx_mm_load_1rvec_2ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
521                                               &jx0,&jy0,&jz0);
522
523             /* Calculate displacement vector */
524             dx00             = _mm_sub_pd(ix0,jx0);
525             dy00             = _mm_sub_pd(iy0,jy0);
526             dz00             = _mm_sub_pd(iz0,jz0);
527
528             /* Calculate squared distance and things based on it */
529             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
530
531             rinv00           = gmx_mm_invsqrt_pd(rsq00);
532
533             /* Load parameters for j particles */
534             jq0              = gmx_mm_load_2real_swizzle_pd(charge+jnrA+0,charge+jnrB+0);
535             vdwjidx0A        = 2*vdwtype[jnrA+0];
536             vdwjidx0B        = 2*vdwtype[jnrB+0];
537
538             /**************************
539              * CALCULATE INTERACTIONS *
540              **************************/
541
542             r00              = _mm_mul_pd(rsq00,rinv00);
543
544             /* Compute parameters for interactions between i and j atoms */
545             qq00             = _mm_mul_pd(iq0,jq0);
546             gmx_mm_load_2pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,
547                                          vdwparam+vdwioffset0+vdwjidx0B,&c6_00,&c12_00);
548
549             /* Calculate table index by multiplying r with table scale and truncate to integer */
550             rt               = _mm_mul_pd(r00,vftabscale);
551             vfitab           = _mm_cvttpd_epi32(rt);
552             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
553             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
554
555             /* CUBIC SPLINE TABLE ELECTROSTATICS */
556             Y                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
557             F                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
558             GMX_MM_TRANSPOSE2_PD(Y,F);
559             G                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
560             H                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
561             GMX_MM_TRANSPOSE2_PD(G,H);
562             Heps             = _mm_mul_pd(vfeps,H);
563             Fp               = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
564             FF               = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
565             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
566
567             /* CUBIC SPLINE TABLE DISPERSION */
568             vfitab           = _mm_add_epi32(vfitab,ifour);
569             Y                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
570             F                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
571             GMX_MM_TRANSPOSE2_PD(Y,F);
572             G                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
573             H                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
574             GMX_MM_TRANSPOSE2_PD(G,H);
575             Heps             = _mm_mul_pd(vfeps,H);
576             Fp               = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
577             FF               = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
578             fvdw6            = _mm_mul_pd(c6_00,FF);
579
580             /* CUBIC SPLINE TABLE REPULSION */
581             vfitab           = _mm_add_epi32(vfitab,ifour);
582             Y                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
583             F                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) );
584             GMX_MM_TRANSPOSE2_PD(Y,F);
585             G                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
586             H                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,1) +2);
587             GMX_MM_TRANSPOSE2_PD(G,H);
588             Heps             = _mm_mul_pd(vfeps,H);
589             Fp               = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
590             FF               = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
591             fvdw12           = _mm_mul_pd(c12_00,FF);
592             fvdw             = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
593
594             fscal            = _mm_add_pd(felec,fvdw);
595
596             /* Calculate temporary vectorial force */
597             tx               = _mm_mul_pd(fscal,dx00);
598             ty               = _mm_mul_pd(fscal,dy00);
599             tz               = _mm_mul_pd(fscal,dz00);
600
601             /* Update vectorial force */
602             fix0             = _mm_add_pd(fix0,tx);
603             fiy0             = _mm_add_pd(fiy0,ty);
604             fiz0             = _mm_add_pd(fiz0,tz);
605
606             gmx_mm_decrement_1rvec_2ptr_swizzle_pd(f+j_coord_offsetA,f+j_coord_offsetB,tx,ty,tz);
607
608             /* Inner loop uses 61 flops */
609         }
610
611         if(jidx<j_index_end)
612         {
613
614             jnrA             = jjnr[jidx];
615             j_coord_offsetA  = DIM*jnrA;
616
617             /* load j atom coordinates */
618             gmx_mm_load_1rvec_1ptr_swizzle_pd(x+j_coord_offsetA,
619                                               &jx0,&jy0,&jz0);
620
621             /* Calculate displacement vector */
622             dx00             = _mm_sub_pd(ix0,jx0);
623             dy00             = _mm_sub_pd(iy0,jy0);
624             dz00             = _mm_sub_pd(iz0,jz0);
625
626             /* Calculate squared distance and things based on it */
627             rsq00            = gmx_mm_calc_rsq_pd(dx00,dy00,dz00);
628
629             rinv00           = gmx_mm_invsqrt_pd(rsq00);
630
631             /* Load parameters for j particles */
632             jq0              = _mm_load_sd(charge+jnrA+0);
633             vdwjidx0A        = 2*vdwtype[jnrA+0];
634
635             /**************************
636              * CALCULATE INTERACTIONS *
637              **************************/
638
639             r00              = _mm_mul_pd(rsq00,rinv00);
640
641             /* Compute parameters for interactions between i and j atoms */
642             qq00             = _mm_mul_pd(iq0,jq0);
643             gmx_mm_load_1pair_swizzle_pd(vdwparam+vdwioffset0+vdwjidx0A,&c6_00,&c12_00);
644
645             /* Calculate table index by multiplying r with table scale and truncate to integer */
646             rt               = _mm_mul_pd(r00,vftabscale);
647             vfitab           = _mm_cvttpd_epi32(rt);
648             vfeps            = _mm_sub_pd(rt,_mm_round_pd(rt, _MM_FROUND_FLOOR));
649             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
650
651             /* CUBIC SPLINE TABLE ELECTROSTATICS */
652             Y                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
653             F                = _mm_setzero_pd();
654             GMX_MM_TRANSPOSE2_PD(Y,F);
655             G                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
656             H                = _mm_setzero_pd();
657             GMX_MM_TRANSPOSE2_PD(G,H);
658             Heps             = _mm_mul_pd(vfeps,H);
659             Fp               = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
660             FF               = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
661             felec            = _mm_xor_pd(signbit,_mm_mul_pd(_mm_mul_pd(qq00,FF),_mm_mul_pd(vftabscale,rinv00)));
662
663             /* CUBIC SPLINE TABLE DISPERSION */
664             vfitab           = _mm_add_epi32(vfitab,ifour);
665             Y                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
666             F                = _mm_setzero_pd();
667             GMX_MM_TRANSPOSE2_PD(Y,F);
668             G                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
669             H                = _mm_setzero_pd();
670             GMX_MM_TRANSPOSE2_PD(G,H);
671             Heps             = _mm_mul_pd(vfeps,H);
672             Fp               = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
673             FF               = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
674             fvdw6            = _mm_mul_pd(c6_00,FF);
675
676             /* CUBIC SPLINE TABLE REPULSION */
677             vfitab           = _mm_add_epi32(vfitab,ifour);
678             Y                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) );
679             F                = _mm_setzero_pd();
680             GMX_MM_TRANSPOSE2_PD(Y,F);
681             G                = _mm_load_pd( vftab + gmx_mm_extract_epi32(vfitab,0) +2);
682             H                = _mm_setzero_pd();
683             GMX_MM_TRANSPOSE2_PD(G,H);
684             Heps             = _mm_mul_pd(vfeps,H);
685             Fp               = _mm_add_pd(F,_mm_mul_pd(vfeps,_mm_add_pd(G,Heps)));
686             FF               = _mm_add_pd(Fp,_mm_mul_pd(vfeps,_mm_add_pd(G,_mm_add_pd(Heps,Heps))));
687             fvdw12           = _mm_mul_pd(c12_00,FF);
688             fvdw             = _mm_xor_pd(signbit,_mm_mul_pd(_mm_add_pd(fvdw6,fvdw12),_mm_mul_pd(vftabscale,rinv00)));
689
690             fscal            = _mm_add_pd(felec,fvdw);
691
692             fscal            = _mm_unpacklo_pd(fscal,_mm_setzero_pd());
693
694             /* Calculate temporary vectorial force */
695             tx               = _mm_mul_pd(fscal,dx00);
696             ty               = _mm_mul_pd(fscal,dy00);
697             tz               = _mm_mul_pd(fscal,dz00);
698
699             /* Update vectorial force */
700             fix0             = _mm_add_pd(fix0,tx);
701             fiy0             = _mm_add_pd(fiy0,ty);
702             fiz0             = _mm_add_pd(fiz0,tz);
703
704             gmx_mm_decrement_1rvec_1ptr_swizzle_pd(f+j_coord_offsetA,tx,ty,tz);
705
706             /* Inner loop uses 61 flops */
707         }
708
709         /* End of innermost loop */
710
711         gmx_mm_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
712                                               f+i_coord_offset,fshift+i_shift_offset);
713
714         /* Increment number of inner iterations */
715         inneriter                  += j_index_end - j_index_start;
716
717         /* Outer loop uses 7 flops */
718     }
719
720     /* Increment number of outer iterations */
721     outeriter        += nri;
722
723     /* Update outer/inner flops */
724
725     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_F,outeriter*7 + inneriter*61);
726 }