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