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