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