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