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