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