Remove all unnecessary HAVE_CONFIG_H
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_single / nb_kernel_ElecEwSw_VdwLJSw_GeomW4P1_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 "types/simple.h"
44 #include "gromacs/math/vec.h"
45 #include "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_ElecEwSw_VdwLJSw_GeomW4P1_VF_avx_128_fma_single
52  * Electrostatics interaction: Ewald
53  * VdW interaction:            LennardJones
54  * Geometry:                   Water4-Particle
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecEwSw_VdwLJSw_GeomW4P1_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              vdwioffset1;
86     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
87     int              vdwioffset2;
88     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
89     int              vdwioffset3;
90     __m128           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
91     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
92     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
93     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
94     __m128           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
95     __m128           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
96     __m128           dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
97     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
98     real             *charge;
99     int              nvdwtype;
100     __m128           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
101     int              *vdwtype;
102     real             *vdwparam;
103     __m128           one_sixth   = _mm_set1_ps(1.0/6.0);
104     __m128           one_twelfth = _mm_set1_ps(1.0/12.0);
105     __m128i          ewitab;
106     __m128           ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
107     __m128           beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
108     real             *ewtab;
109     __m128           rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
110     real             rswitch_scalar,d_scalar;
111     __m128           dummy_mask,cutoff_mask;
112     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
113     __m128           one     = _mm_set1_ps(1.0);
114     __m128           two     = _mm_set1_ps(2.0);
115     x                = xx[0];
116     f                = ff[0];
117
118     nri              = nlist->nri;
119     iinr             = nlist->iinr;
120     jindex           = nlist->jindex;
121     jjnr             = nlist->jjnr;
122     shiftidx         = nlist->shift;
123     gid              = nlist->gid;
124     shiftvec         = fr->shift_vec[0];
125     fshift           = fr->fshift[0];
126     facel            = _mm_set1_ps(fr->epsfac);
127     charge           = mdatoms->chargeA;
128     nvdwtype         = fr->ntype;
129     vdwparam         = fr->nbfp;
130     vdwtype          = mdatoms->typeA;
131
132     sh_ewald         = _mm_set1_ps(fr->ic->sh_ewald);
133     beta             = _mm_set1_ps(fr->ic->ewaldcoeff_q);
134     beta2            = _mm_mul_ps(beta,beta);
135     beta3            = _mm_mul_ps(beta,beta2);
136     ewtab            = fr->ic->tabq_coul_FDV0;
137     ewtabscale       = _mm_set1_ps(fr->ic->tabq_scale);
138     ewtabhalfspace   = _mm_set1_ps(0.5/fr->ic->tabq_scale);
139
140     /* Setup water-specific parameters */
141     inr              = nlist->iinr[0];
142     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
143     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
144     iq3              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
145     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
146
147     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
148     rcutoff_scalar   = fr->rcoulomb;
149     rcutoff          = _mm_set1_ps(rcutoff_scalar);
150     rcutoff2         = _mm_mul_ps(rcutoff,rcutoff);
151
152     rswitch_scalar   = fr->rcoulomb_switch;
153     rswitch          = _mm_set1_ps(rswitch_scalar);
154     /* Setup switch parameters */
155     d_scalar         = rcutoff_scalar-rswitch_scalar;
156     d                = _mm_set1_ps(d_scalar);
157     swV3             = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
158     swV4             = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
159     swV5             = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
160     swF2             = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
161     swF3             = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
162     swF4             = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
163
164     /* Avoid stupid compiler warnings */
165     jnrA = jnrB = jnrC = jnrD = 0;
166     j_coord_offsetA = 0;
167     j_coord_offsetB = 0;
168     j_coord_offsetC = 0;
169     j_coord_offsetD = 0;
170
171     outeriter        = 0;
172     inneriter        = 0;
173
174     for(iidx=0;iidx<4*DIM;iidx++)
175     {
176         scratch[iidx] = 0.0;
177     }
178
179     /* Start outer loop over neighborlists */
180     for(iidx=0; iidx<nri; iidx++)
181     {
182         /* Load shift vector for this list */
183         i_shift_offset   = DIM*shiftidx[iidx];
184
185         /* Load limits for loop over neighbors */
186         j_index_start    = jindex[iidx];
187         j_index_end      = jindex[iidx+1];
188
189         /* Get outer coordinate index */
190         inr              = iinr[iidx];
191         i_coord_offset   = DIM*inr;
192
193         /* Load i particle coords and add shift vector */
194         gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
195                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
196
197         fix0             = _mm_setzero_ps();
198         fiy0             = _mm_setzero_ps();
199         fiz0             = _mm_setzero_ps();
200         fix1             = _mm_setzero_ps();
201         fiy1             = _mm_setzero_ps();
202         fiz1             = _mm_setzero_ps();
203         fix2             = _mm_setzero_ps();
204         fiy2             = _mm_setzero_ps();
205         fiz2             = _mm_setzero_ps();
206         fix3             = _mm_setzero_ps();
207         fiy3             = _mm_setzero_ps();
208         fiz3             = _mm_setzero_ps();
209
210         /* Reset potential sums */
211         velecsum         = _mm_setzero_ps();
212         vvdwsum          = _mm_setzero_ps();
213
214         /* Start inner kernel loop */
215         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
216         {
217
218             /* Get j neighbor index, and coordinate index */
219             jnrA             = jjnr[jidx];
220             jnrB             = jjnr[jidx+1];
221             jnrC             = jjnr[jidx+2];
222             jnrD             = jjnr[jidx+3];
223             j_coord_offsetA  = DIM*jnrA;
224             j_coord_offsetB  = DIM*jnrB;
225             j_coord_offsetC  = DIM*jnrC;
226             j_coord_offsetD  = DIM*jnrD;
227
228             /* load j atom coordinates */
229             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
230                                               x+j_coord_offsetC,x+j_coord_offsetD,
231                                               &jx0,&jy0,&jz0);
232
233             /* Calculate displacement vector */
234             dx00             = _mm_sub_ps(ix0,jx0);
235             dy00             = _mm_sub_ps(iy0,jy0);
236             dz00             = _mm_sub_ps(iz0,jz0);
237             dx10             = _mm_sub_ps(ix1,jx0);
238             dy10             = _mm_sub_ps(iy1,jy0);
239             dz10             = _mm_sub_ps(iz1,jz0);
240             dx20             = _mm_sub_ps(ix2,jx0);
241             dy20             = _mm_sub_ps(iy2,jy0);
242             dz20             = _mm_sub_ps(iz2,jz0);
243             dx30             = _mm_sub_ps(ix3,jx0);
244             dy30             = _mm_sub_ps(iy3,jy0);
245             dz30             = _mm_sub_ps(iz3,jz0);
246
247             /* Calculate squared distance and things based on it */
248             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
249             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
250             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
251             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
252
253             rinv00           = gmx_mm_invsqrt_ps(rsq00);
254             rinv10           = gmx_mm_invsqrt_ps(rsq10);
255             rinv20           = gmx_mm_invsqrt_ps(rsq20);
256             rinv30           = gmx_mm_invsqrt_ps(rsq30);
257
258             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
259             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
260             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
261             rinvsq30         = _mm_mul_ps(rinv30,rinv30);
262
263             /* Load parameters for j particles */
264             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
265                                                               charge+jnrC+0,charge+jnrD+0);
266             vdwjidx0A        = 2*vdwtype[jnrA+0];
267             vdwjidx0B        = 2*vdwtype[jnrB+0];
268             vdwjidx0C        = 2*vdwtype[jnrC+0];
269             vdwjidx0D        = 2*vdwtype[jnrD+0];
270
271             fjx0             = _mm_setzero_ps();
272             fjy0             = _mm_setzero_ps();
273             fjz0             = _mm_setzero_ps();
274
275             /**************************
276              * CALCULATE INTERACTIONS *
277              **************************/
278
279             if (gmx_mm_any_lt(rsq00,rcutoff2))
280             {
281
282             r00              = _mm_mul_ps(rsq00,rinv00);
283
284             /* Compute parameters for interactions between i and j atoms */
285             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
286                                          vdwparam+vdwioffset0+vdwjidx0B,
287                                          vdwparam+vdwioffset0+vdwjidx0C,
288                                          vdwparam+vdwioffset0+vdwjidx0D,
289                                          &c6_00,&c12_00);
290
291             /* LENNARD-JONES DISPERSION/REPULSION */
292
293             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
294             vvdw6            = _mm_mul_ps(c6_00,rinvsix);
295             vvdw12           = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
296             vvdw             = _mm_msub_ps(vvdw12,one_twelfth,_mm_mul_ps(vvdw6,one_sixth));
297             fvdw             = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
298
299             d                = _mm_sub_ps(r00,rswitch);
300             d                = _mm_max_ps(d,_mm_setzero_ps());
301             d2               = _mm_mul_ps(d,d);
302             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
303
304             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
305
306             /* Evaluate switch function */
307             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
308             fvdw             = _mm_msub_ps( fvdw,sw , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
309             vvdw             = _mm_mul_ps(vvdw,sw);
310             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
311
312             /* Update potential sum for this i atom from the interaction with this j atom. */
313             vvdw             = _mm_and_ps(vvdw,cutoff_mask);
314             vvdwsum          = _mm_add_ps(vvdwsum,vvdw);
315
316             fscal            = fvdw;
317
318             fscal            = _mm_and_ps(fscal,cutoff_mask);
319
320              /* Update vectorial force */
321             fix0             = _mm_macc_ps(dx00,fscal,fix0);
322             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
323             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
324
325             fjx0             = _mm_macc_ps(dx00,fscal,fjx0);
326             fjy0             = _mm_macc_ps(dy00,fscal,fjy0);
327             fjz0             = _mm_macc_ps(dz00,fscal,fjz0);
328
329             }
330
331             /**************************
332              * CALCULATE INTERACTIONS *
333              **************************/
334
335             if (gmx_mm_any_lt(rsq10,rcutoff2))
336             {
337
338             r10              = _mm_mul_ps(rsq10,rinv10);
339
340             /* Compute parameters for interactions between i and j atoms */
341             qq10             = _mm_mul_ps(iq1,jq0);
342
343             /* EWALD ELECTROSTATICS */
344
345             /* Analytical PME correction */
346             zeta2            = _mm_mul_ps(beta2,rsq10);
347             rinv3            = _mm_mul_ps(rinvsq10,rinv10);
348             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
349             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
350             felec            = _mm_mul_ps(qq10,felec);
351             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
352             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv10);
353             velec            = _mm_mul_ps(qq10,velec);
354
355             d                = _mm_sub_ps(r10,rswitch);
356             d                = _mm_max_ps(d,_mm_setzero_ps());
357             d2               = _mm_mul_ps(d,d);
358             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
359
360             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
361
362             /* Evaluate switch function */
363             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
364             felec            = _mm_msub_ps( felec,sw , _mm_mul_ps(rinv10,_mm_mul_ps(velec,dsw)) );
365             velec            = _mm_mul_ps(velec,sw);
366             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
367
368             /* Update potential sum for this i atom from the interaction with this j atom. */
369             velec            = _mm_and_ps(velec,cutoff_mask);
370             velecsum         = _mm_add_ps(velecsum,velec);
371
372             fscal            = felec;
373
374             fscal            = _mm_and_ps(fscal,cutoff_mask);
375
376              /* Update vectorial force */
377             fix1             = _mm_macc_ps(dx10,fscal,fix1);
378             fiy1             = _mm_macc_ps(dy10,fscal,fiy1);
379             fiz1             = _mm_macc_ps(dz10,fscal,fiz1);
380
381             fjx0             = _mm_macc_ps(dx10,fscal,fjx0);
382             fjy0             = _mm_macc_ps(dy10,fscal,fjy0);
383             fjz0             = _mm_macc_ps(dz10,fscal,fjz0);
384
385             }
386
387             /**************************
388              * CALCULATE INTERACTIONS *
389              **************************/
390
391             if (gmx_mm_any_lt(rsq20,rcutoff2))
392             {
393
394             r20              = _mm_mul_ps(rsq20,rinv20);
395
396             /* Compute parameters for interactions between i and j atoms */
397             qq20             = _mm_mul_ps(iq2,jq0);
398
399             /* EWALD ELECTROSTATICS */
400
401             /* Analytical PME correction */
402             zeta2            = _mm_mul_ps(beta2,rsq20);
403             rinv3            = _mm_mul_ps(rinvsq20,rinv20);
404             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
405             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
406             felec            = _mm_mul_ps(qq20,felec);
407             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
408             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv20);
409             velec            = _mm_mul_ps(qq20,velec);
410
411             d                = _mm_sub_ps(r20,rswitch);
412             d                = _mm_max_ps(d,_mm_setzero_ps());
413             d2               = _mm_mul_ps(d,d);
414             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
415
416             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
417
418             /* Evaluate switch function */
419             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
420             felec            = _mm_msub_ps( felec,sw , _mm_mul_ps(rinv20,_mm_mul_ps(velec,dsw)) );
421             velec            = _mm_mul_ps(velec,sw);
422             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
423
424             /* Update potential sum for this i atom from the interaction with this j atom. */
425             velec            = _mm_and_ps(velec,cutoff_mask);
426             velecsum         = _mm_add_ps(velecsum,velec);
427
428             fscal            = felec;
429
430             fscal            = _mm_and_ps(fscal,cutoff_mask);
431
432              /* Update vectorial force */
433             fix2             = _mm_macc_ps(dx20,fscal,fix2);
434             fiy2             = _mm_macc_ps(dy20,fscal,fiy2);
435             fiz2             = _mm_macc_ps(dz20,fscal,fiz2);
436
437             fjx0             = _mm_macc_ps(dx20,fscal,fjx0);
438             fjy0             = _mm_macc_ps(dy20,fscal,fjy0);
439             fjz0             = _mm_macc_ps(dz20,fscal,fjz0);
440
441             }
442
443             /**************************
444              * CALCULATE INTERACTIONS *
445              **************************/
446
447             if (gmx_mm_any_lt(rsq30,rcutoff2))
448             {
449
450             r30              = _mm_mul_ps(rsq30,rinv30);
451
452             /* Compute parameters for interactions between i and j atoms */
453             qq30             = _mm_mul_ps(iq3,jq0);
454
455             /* EWALD ELECTROSTATICS */
456
457             /* Analytical PME correction */
458             zeta2            = _mm_mul_ps(beta2,rsq30);
459             rinv3            = _mm_mul_ps(rinvsq30,rinv30);
460             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
461             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
462             felec            = _mm_mul_ps(qq30,felec);
463             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
464             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv30);
465             velec            = _mm_mul_ps(qq30,velec);
466
467             d                = _mm_sub_ps(r30,rswitch);
468             d                = _mm_max_ps(d,_mm_setzero_ps());
469             d2               = _mm_mul_ps(d,d);
470             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
471
472             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
473
474             /* Evaluate switch function */
475             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
476             felec            = _mm_msub_ps( felec,sw , _mm_mul_ps(rinv30,_mm_mul_ps(velec,dsw)) );
477             velec            = _mm_mul_ps(velec,sw);
478             cutoff_mask      = _mm_cmplt_ps(rsq30,rcutoff2);
479
480             /* Update potential sum for this i atom from the interaction with this j atom. */
481             velec            = _mm_and_ps(velec,cutoff_mask);
482             velecsum         = _mm_add_ps(velecsum,velec);
483
484             fscal            = felec;
485
486             fscal            = _mm_and_ps(fscal,cutoff_mask);
487
488              /* Update vectorial force */
489             fix3             = _mm_macc_ps(dx30,fscal,fix3);
490             fiy3             = _mm_macc_ps(dy30,fscal,fiy3);
491             fiz3             = _mm_macc_ps(dz30,fscal,fiz3);
492
493             fjx0             = _mm_macc_ps(dx30,fscal,fjx0);
494             fjy0             = _mm_macc_ps(dy30,fscal,fjy0);
495             fjz0             = _mm_macc_ps(dz30,fscal,fjz0);
496
497             }
498
499             fjptrA             = f+j_coord_offsetA;
500             fjptrB             = f+j_coord_offsetB;
501             fjptrC             = f+j_coord_offsetC;
502             fjptrD             = f+j_coord_offsetD;
503
504             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
505
506             /* Inner loop uses 221 flops */
507         }
508
509         if(jidx<j_index_end)
510         {
511
512             /* Get j neighbor index, and coordinate index */
513             jnrlistA         = jjnr[jidx];
514             jnrlistB         = jjnr[jidx+1];
515             jnrlistC         = jjnr[jidx+2];
516             jnrlistD         = jjnr[jidx+3];
517             /* Sign of each element will be negative for non-real atoms.
518              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
519              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
520              */
521             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
522             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
523             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
524             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
525             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
526             j_coord_offsetA  = DIM*jnrA;
527             j_coord_offsetB  = DIM*jnrB;
528             j_coord_offsetC  = DIM*jnrC;
529             j_coord_offsetD  = DIM*jnrD;
530
531             /* load j atom coordinates */
532             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
533                                               x+j_coord_offsetC,x+j_coord_offsetD,
534                                               &jx0,&jy0,&jz0);
535
536             /* Calculate displacement vector */
537             dx00             = _mm_sub_ps(ix0,jx0);
538             dy00             = _mm_sub_ps(iy0,jy0);
539             dz00             = _mm_sub_ps(iz0,jz0);
540             dx10             = _mm_sub_ps(ix1,jx0);
541             dy10             = _mm_sub_ps(iy1,jy0);
542             dz10             = _mm_sub_ps(iz1,jz0);
543             dx20             = _mm_sub_ps(ix2,jx0);
544             dy20             = _mm_sub_ps(iy2,jy0);
545             dz20             = _mm_sub_ps(iz2,jz0);
546             dx30             = _mm_sub_ps(ix3,jx0);
547             dy30             = _mm_sub_ps(iy3,jy0);
548             dz30             = _mm_sub_ps(iz3,jz0);
549
550             /* Calculate squared distance and things based on it */
551             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
552             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
553             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
554             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
555
556             rinv00           = gmx_mm_invsqrt_ps(rsq00);
557             rinv10           = gmx_mm_invsqrt_ps(rsq10);
558             rinv20           = gmx_mm_invsqrt_ps(rsq20);
559             rinv30           = gmx_mm_invsqrt_ps(rsq30);
560
561             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
562             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
563             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
564             rinvsq30         = _mm_mul_ps(rinv30,rinv30);
565
566             /* Load parameters for j particles */
567             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
568                                                               charge+jnrC+0,charge+jnrD+0);
569             vdwjidx0A        = 2*vdwtype[jnrA+0];
570             vdwjidx0B        = 2*vdwtype[jnrB+0];
571             vdwjidx0C        = 2*vdwtype[jnrC+0];
572             vdwjidx0D        = 2*vdwtype[jnrD+0];
573
574             fjx0             = _mm_setzero_ps();
575             fjy0             = _mm_setzero_ps();
576             fjz0             = _mm_setzero_ps();
577
578             /**************************
579              * CALCULATE INTERACTIONS *
580              **************************/
581
582             if (gmx_mm_any_lt(rsq00,rcutoff2))
583             {
584
585             r00              = _mm_mul_ps(rsq00,rinv00);
586             r00              = _mm_andnot_ps(dummy_mask,r00);
587
588             /* Compute parameters for interactions between i and j atoms */
589             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
590                                          vdwparam+vdwioffset0+vdwjidx0B,
591                                          vdwparam+vdwioffset0+vdwjidx0C,
592                                          vdwparam+vdwioffset0+vdwjidx0D,
593                                          &c6_00,&c12_00);
594
595             /* LENNARD-JONES DISPERSION/REPULSION */
596
597             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
598             vvdw6            = _mm_mul_ps(c6_00,rinvsix);
599             vvdw12           = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
600             vvdw             = _mm_msub_ps(vvdw12,one_twelfth,_mm_mul_ps(vvdw6,one_sixth));
601             fvdw             = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
602
603             d                = _mm_sub_ps(r00,rswitch);
604             d                = _mm_max_ps(d,_mm_setzero_ps());
605             d2               = _mm_mul_ps(d,d);
606             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
607
608             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
609
610             /* Evaluate switch function */
611             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
612             fvdw             = _mm_msub_ps( fvdw,sw , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
613             vvdw             = _mm_mul_ps(vvdw,sw);
614             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
615
616             /* Update potential sum for this i atom from the interaction with this j atom. */
617             vvdw             = _mm_and_ps(vvdw,cutoff_mask);
618             vvdw             = _mm_andnot_ps(dummy_mask,vvdw);
619             vvdwsum          = _mm_add_ps(vvdwsum,vvdw);
620
621             fscal            = fvdw;
622
623             fscal            = _mm_and_ps(fscal,cutoff_mask);
624
625             fscal            = _mm_andnot_ps(dummy_mask,fscal);
626
627              /* Update vectorial force */
628             fix0             = _mm_macc_ps(dx00,fscal,fix0);
629             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
630             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
631
632             fjx0             = _mm_macc_ps(dx00,fscal,fjx0);
633             fjy0             = _mm_macc_ps(dy00,fscal,fjy0);
634             fjz0             = _mm_macc_ps(dz00,fscal,fjz0);
635
636             }
637
638             /**************************
639              * CALCULATE INTERACTIONS *
640              **************************/
641
642             if (gmx_mm_any_lt(rsq10,rcutoff2))
643             {
644
645             r10              = _mm_mul_ps(rsq10,rinv10);
646             r10              = _mm_andnot_ps(dummy_mask,r10);
647
648             /* Compute parameters for interactions between i and j atoms */
649             qq10             = _mm_mul_ps(iq1,jq0);
650
651             /* EWALD ELECTROSTATICS */
652
653             /* Analytical PME correction */
654             zeta2            = _mm_mul_ps(beta2,rsq10);
655             rinv3            = _mm_mul_ps(rinvsq10,rinv10);
656             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
657             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
658             felec            = _mm_mul_ps(qq10,felec);
659             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
660             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv10);
661             velec            = _mm_mul_ps(qq10,velec);
662
663             d                = _mm_sub_ps(r10,rswitch);
664             d                = _mm_max_ps(d,_mm_setzero_ps());
665             d2               = _mm_mul_ps(d,d);
666             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
667
668             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
669
670             /* Evaluate switch function */
671             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
672             felec            = _mm_msub_ps( felec,sw , _mm_mul_ps(rinv10,_mm_mul_ps(velec,dsw)) );
673             velec            = _mm_mul_ps(velec,sw);
674             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
675
676             /* Update potential sum for this i atom from the interaction with this j atom. */
677             velec            = _mm_and_ps(velec,cutoff_mask);
678             velec            = _mm_andnot_ps(dummy_mask,velec);
679             velecsum         = _mm_add_ps(velecsum,velec);
680
681             fscal            = felec;
682
683             fscal            = _mm_and_ps(fscal,cutoff_mask);
684
685             fscal            = _mm_andnot_ps(dummy_mask,fscal);
686
687              /* Update vectorial force */
688             fix1             = _mm_macc_ps(dx10,fscal,fix1);
689             fiy1             = _mm_macc_ps(dy10,fscal,fiy1);
690             fiz1             = _mm_macc_ps(dz10,fscal,fiz1);
691
692             fjx0             = _mm_macc_ps(dx10,fscal,fjx0);
693             fjy0             = _mm_macc_ps(dy10,fscal,fjy0);
694             fjz0             = _mm_macc_ps(dz10,fscal,fjz0);
695
696             }
697
698             /**************************
699              * CALCULATE INTERACTIONS *
700              **************************/
701
702             if (gmx_mm_any_lt(rsq20,rcutoff2))
703             {
704
705             r20              = _mm_mul_ps(rsq20,rinv20);
706             r20              = _mm_andnot_ps(dummy_mask,r20);
707
708             /* Compute parameters for interactions between i and j atoms */
709             qq20             = _mm_mul_ps(iq2,jq0);
710
711             /* EWALD ELECTROSTATICS */
712
713             /* Analytical PME correction */
714             zeta2            = _mm_mul_ps(beta2,rsq20);
715             rinv3            = _mm_mul_ps(rinvsq20,rinv20);
716             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
717             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
718             felec            = _mm_mul_ps(qq20,felec);
719             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
720             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv20);
721             velec            = _mm_mul_ps(qq20,velec);
722
723             d                = _mm_sub_ps(r20,rswitch);
724             d                = _mm_max_ps(d,_mm_setzero_ps());
725             d2               = _mm_mul_ps(d,d);
726             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
727
728             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
729
730             /* Evaluate switch function */
731             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
732             felec            = _mm_msub_ps( felec,sw , _mm_mul_ps(rinv20,_mm_mul_ps(velec,dsw)) );
733             velec            = _mm_mul_ps(velec,sw);
734             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
735
736             /* Update potential sum for this i atom from the interaction with this j atom. */
737             velec            = _mm_and_ps(velec,cutoff_mask);
738             velec            = _mm_andnot_ps(dummy_mask,velec);
739             velecsum         = _mm_add_ps(velecsum,velec);
740
741             fscal            = felec;
742
743             fscal            = _mm_and_ps(fscal,cutoff_mask);
744
745             fscal            = _mm_andnot_ps(dummy_mask,fscal);
746
747              /* Update vectorial force */
748             fix2             = _mm_macc_ps(dx20,fscal,fix2);
749             fiy2             = _mm_macc_ps(dy20,fscal,fiy2);
750             fiz2             = _mm_macc_ps(dz20,fscal,fiz2);
751
752             fjx0             = _mm_macc_ps(dx20,fscal,fjx0);
753             fjy0             = _mm_macc_ps(dy20,fscal,fjy0);
754             fjz0             = _mm_macc_ps(dz20,fscal,fjz0);
755
756             }
757
758             /**************************
759              * CALCULATE INTERACTIONS *
760              **************************/
761
762             if (gmx_mm_any_lt(rsq30,rcutoff2))
763             {
764
765             r30              = _mm_mul_ps(rsq30,rinv30);
766             r30              = _mm_andnot_ps(dummy_mask,r30);
767
768             /* Compute parameters for interactions between i and j atoms */
769             qq30             = _mm_mul_ps(iq3,jq0);
770
771             /* EWALD ELECTROSTATICS */
772
773             /* Analytical PME correction */
774             zeta2            = _mm_mul_ps(beta2,rsq30);
775             rinv3            = _mm_mul_ps(rinvsq30,rinv30);
776             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
777             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
778             felec            = _mm_mul_ps(qq30,felec);
779             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
780             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv30);
781             velec            = _mm_mul_ps(qq30,velec);
782
783             d                = _mm_sub_ps(r30,rswitch);
784             d                = _mm_max_ps(d,_mm_setzero_ps());
785             d2               = _mm_mul_ps(d,d);
786             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
787
788             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
789
790             /* Evaluate switch function */
791             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
792             felec            = _mm_msub_ps( felec,sw , _mm_mul_ps(rinv30,_mm_mul_ps(velec,dsw)) );
793             velec            = _mm_mul_ps(velec,sw);
794             cutoff_mask      = _mm_cmplt_ps(rsq30,rcutoff2);
795
796             /* Update potential sum for this i atom from the interaction with this j atom. */
797             velec            = _mm_and_ps(velec,cutoff_mask);
798             velec            = _mm_andnot_ps(dummy_mask,velec);
799             velecsum         = _mm_add_ps(velecsum,velec);
800
801             fscal            = felec;
802
803             fscal            = _mm_and_ps(fscal,cutoff_mask);
804
805             fscal            = _mm_andnot_ps(dummy_mask,fscal);
806
807              /* Update vectorial force */
808             fix3             = _mm_macc_ps(dx30,fscal,fix3);
809             fiy3             = _mm_macc_ps(dy30,fscal,fiy3);
810             fiz3             = _mm_macc_ps(dz30,fscal,fiz3);
811
812             fjx0             = _mm_macc_ps(dx30,fscal,fjx0);
813             fjy0             = _mm_macc_ps(dy30,fscal,fjy0);
814             fjz0             = _mm_macc_ps(dz30,fscal,fjz0);
815
816             }
817
818             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
819             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
820             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
821             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
822
823             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
824
825             /* Inner loop uses 225 flops */
826         }
827
828         /* End of innermost loop */
829
830         gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
831                                               f+i_coord_offset,fshift+i_shift_offset);
832
833         ggid                        = gid[iidx];
834         /* Update potential energies */
835         gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
836         gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
837
838         /* Increment number of inner iterations */
839         inneriter                  += j_index_end - j_index_start;
840
841         /* Outer loop uses 26 flops */
842     }
843
844     /* Increment number of outer iterations */
845     outeriter        += nri;
846
847     /* Update outer/inner flops */
848
849     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*225);
850 }
851 /*
852  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSw_VdwLJSw_GeomW4P1_F_avx_128_fma_single
853  * Electrostatics interaction: Ewald
854  * VdW interaction:            LennardJones
855  * Geometry:                   Water4-Particle
856  * Calculate force/pot:        Force
857  */
858 void
859 nb_kernel_ElecEwSw_VdwLJSw_GeomW4P1_F_avx_128_fma_single
860                     (t_nblist                    * gmx_restrict       nlist,
861                      rvec                        * gmx_restrict          xx,
862                      rvec                        * gmx_restrict          ff,
863                      t_forcerec                  * gmx_restrict          fr,
864                      t_mdatoms                   * gmx_restrict     mdatoms,
865                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
866                      t_nrnb                      * gmx_restrict        nrnb)
867 {
868     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
869      * just 0 for non-waters.
870      * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
871      * jnr indices corresponding to data put in the four positions in the SIMD register.
872      */
873     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
874     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
875     int              jnrA,jnrB,jnrC,jnrD;
876     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
877     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
878     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
879     real             rcutoff_scalar;
880     real             *shiftvec,*fshift,*x,*f;
881     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
882     real             scratch[4*DIM];
883     __m128           fscal,rcutoff,rcutoff2,jidxall;
884     int              vdwioffset0;
885     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
886     int              vdwioffset1;
887     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
888     int              vdwioffset2;
889     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
890     int              vdwioffset3;
891     __m128           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
892     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
893     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
894     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
895     __m128           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
896     __m128           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
897     __m128           dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
898     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
899     real             *charge;
900     int              nvdwtype;
901     __m128           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
902     int              *vdwtype;
903     real             *vdwparam;
904     __m128           one_sixth   = _mm_set1_ps(1.0/6.0);
905     __m128           one_twelfth = _mm_set1_ps(1.0/12.0);
906     __m128i          ewitab;
907     __m128           ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
908     __m128           beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
909     real             *ewtab;
910     __m128           rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
911     real             rswitch_scalar,d_scalar;
912     __m128           dummy_mask,cutoff_mask;
913     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
914     __m128           one     = _mm_set1_ps(1.0);
915     __m128           two     = _mm_set1_ps(2.0);
916     x                = xx[0];
917     f                = ff[0];
918
919     nri              = nlist->nri;
920     iinr             = nlist->iinr;
921     jindex           = nlist->jindex;
922     jjnr             = nlist->jjnr;
923     shiftidx         = nlist->shift;
924     gid              = nlist->gid;
925     shiftvec         = fr->shift_vec[0];
926     fshift           = fr->fshift[0];
927     facel            = _mm_set1_ps(fr->epsfac);
928     charge           = mdatoms->chargeA;
929     nvdwtype         = fr->ntype;
930     vdwparam         = fr->nbfp;
931     vdwtype          = mdatoms->typeA;
932
933     sh_ewald         = _mm_set1_ps(fr->ic->sh_ewald);
934     beta             = _mm_set1_ps(fr->ic->ewaldcoeff_q);
935     beta2            = _mm_mul_ps(beta,beta);
936     beta3            = _mm_mul_ps(beta,beta2);
937     ewtab            = fr->ic->tabq_coul_FDV0;
938     ewtabscale       = _mm_set1_ps(fr->ic->tabq_scale);
939     ewtabhalfspace   = _mm_set1_ps(0.5/fr->ic->tabq_scale);
940
941     /* Setup water-specific parameters */
942     inr              = nlist->iinr[0];
943     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
944     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
945     iq3              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
946     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
947
948     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
949     rcutoff_scalar   = fr->rcoulomb;
950     rcutoff          = _mm_set1_ps(rcutoff_scalar);
951     rcutoff2         = _mm_mul_ps(rcutoff,rcutoff);
952
953     rswitch_scalar   = fr->rcoulomb_switch;
954     rswitch          = _mm_set1_ps(rswitch_scalar);
955     /* Setup switch parameters */
956     d_scalar         = rcutoff_scalar-rswitch_scalar;
957     d                = _mm_set1_ps(d_scalar);
958     swV3             = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
959     swV4             = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
960     swV5             = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
961     swF2             = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
962     swF3             = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
963     swF4             = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
964
965     /* Avoid stupid compiler warnings */
966     jnrA = jnrB = jnrC = jnrD = 0;
967     j_coord_offsetA = 0;
968     j_coord_offsetB = 0;
969     j_coord_offsetC = 0;
970     j_coord_offsetD = 0;
971
972     outeriter        = 0;
973     inneriter        = 0;
974
975     for(iidx=0;iidx<4*DIM;iidx++)
976     {
977         scratch[iidx] = 0.0;
978     }
979
980     /* Start outer loop over neighborlists */
981     for(iidx=0; iidx<nri; iidx++)
982     {
983         /* Load shift vector for this list */
984         i_shift_offset   = DIM*shiftidx[iidx];
985
986         /* Load limits for loop over neighbors */
987         j_index_start    = jindex[iidx];
988         j_index_end      = jindex[iidx+1];
989
990         /* Get outer coordinate index */
991         inr              = iinr[iidx];
992         i_coord_offset   = DIM*inr;
993
994         /* Load i particle coords and add shift vector */
995         gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
996                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
997
998         fix0             = _mm_setzero_ps();
999         fiy0             = _mm_setzero_ps();
1000         fiz0             = _mm_setzero_ps();
1001         fix1             = _mm_setzero_ps();
1002         fiy1             = _mm_setzero_ps();
1003         fiz1             = _mm_setzero_ps();
1004         fix2             = _mm_setzero_ps();
1005         fiy2             = _mm_setzero_ps();
1006         fiz2             = _mm_setzero_ps();
1007         fix3             = _mm_setzero_ps();
1008         fiy3             = _mm_setzero_ps();
1009         fiz3             = _mm_setzero_ps();
1010
1011         /* Start inner kernel loop */
1012         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
1013         {
1014
1015             /* Get j neighbor index, and coordinate index */
1016             jnrA             = jjnr[jidx];
1017             jnrB             = jjnr[jidx+1];
1018             jnrC             = jjnr[jidx+2];
1019             jnrD             = jjnr[jidx+3];
1020             j_coord_offsetA  = DIM*jnrA;
1021             j_coord_offsetB  = DIM*jnrB;
1022             j_coord_offsetC  = DIM*jnrC;
1023             j_coord_offsetD  = DIM*jnrD;
1024
1025             /* load j atom coordinates */
1026             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1027                                               x+j_coord_offsetC,x+j_coord_offsetD,
1028                                               &jx0,&jy0,&jz0);
1029
1030             /* Calculate displacement vector */
1031             dx00             = _mm_sub_ps(ix0,jx0);
1032             dy00             = _mm_sub_ps(iy0,jy0);
1033             dz00             = _mm_sub_ps(iz0,jz0);
1034             dx10             = _mm_sub_ps(ix1,jx0);
1035             dy10             = _mm_sub_ps(iy1,jy0);
1036             dz10             = _mm_sub_ps(iz1,jz0);
1037             dx20             = _mm_sub_ps(ix2,jx0);
1038             dy20             = _mm_sub_ps(iy2,jy0);
1039             dz20             = _mm_sub_ps(iz2,jz0);
1040             dx30             = _mm_sub_ps(ix3,jx0);
1041             dy30             = _mm_sub_ps(iy3,jy0);
1042             dz30             = _mm_sub_ps(iz3,jz0);
1043
1044             /* Calculate squared distance and things based on it */
1045             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
1046             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
1047             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
1048             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
1049
1050             rinv00           = gmx_mm_invsqrt_ps(rsq00);
1051             rinv10           = gmx_mm_invsqrt_ps(rsq10);
1052             rinv20           = gmx_mm_invsqrt_ps(rsq20);
1053             rinv30           = gmx_mm_invsqrt_ps(rsq30);
1054
1055             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
1056             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
1057             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
1058             rinvsq30         = _mm_mul_ps(rinv30,rinv30);
1059
1060             /* Load parameters for j particles */
1061             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1062                                                               charge+jnrC+0,charge+jnrD+0);
1063             vdwjidx0A        = 2*vdwtype[jnrA+0];
1064             vdwjidx0B        = 2*vdwtype[jnrB+0];
1065             vdwjidx0C        = 2*vdwtype[jnrC+0];
1066             vdwjidx0D        = 2*vdwtype[jnrD+0];
1067
1068             fjx0             = _mm_setzero_ps();
1069             fjy0             = _mm_setzero_ps();
1070             fjz0             = _mm_setzero_ps();
1071
1072             /**************************
1073              * CALCULATE INTERACTIONS *
1074              **************************/
1075
1076             if (gmx_mm_any_lt(rsq00,rcutoff2))
1077             {
1078
1079             r00              = _mm_mul_ps(rsq00,rinv00);
1080
1081             /* Compute parameters for interactions between i and j atoms */
1082             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
1083                                          vdwparam+vdwioffset0+vdwjidx0B,
1084                                          vdwparam+vdwioffset0+vdwjidx0C,
1085                                          vdwparam+vdwioffset0+vdwjidx0D,
1086                                          &c6_00,&c12_00);
1087
1088             /* LENNARD-JONES DISPERSION/REPULSION */
1089
1090             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
1091             vvdw6            = _mm_mul_ps(c6_00,rinvsix);
1092             vvdw12           = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
1093             vvdw             = _mm_msub_ps(vvdw12,one_twelfth,_mm_mul_ps(vvdw6,one_sixth));
1094             fvdw             = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
1095
1096             d                = _mm_sub_ps(r00,rswitch);
1097             d                = _mm_max_ps(d,_mm_setzero_ps());
1098             d2               = _mm_mul_ps(d,d);
1099             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
1100
1101             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
1102
1103             /* Evaluate switch function */
1104             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1105             fvdw             = _mm_msub_ps( fvdw,sw , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
1106             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
1107
1108             fscal            = fvdw;
1109
1110             fscal            = _mm_and_ps(fscal,cutoff_mask);
1111
1112              /* Update vectorial force */
1113             fix0             = _mm_macc_ps(dx00,fscal,fix0);
1114             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
1115             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
1116
1117             fjx0             = _mm_macc_ps(dx00,fscal,fjx0);
1118             fjy0             = _mm_macc_ps(dy00,fscal,fjy0);
1119             fjz0             = _mm_macc_ps(dz00,fscal,fjz0);
1120
1121             }
1122
1123             /**************************
1124              * CALCULATE INTERACTIONS *
1125              **************************/
1126
1127             if (gmx_mm_any_lt(rsq10,rcutoff2))
1128             {
1129
1130             r10              = _mm_mul_ps(rsq10,rinv10);
1131
1132             /* Compute parameters for interactions between i and j atoms */
1133             qq10             = _mm_mul_ps(iq1,jq0);
1134
1135             /* EWALD ELECTROSTATICS */
1136
1137             /* Analytical PME correction */
1138             zeta2            = _mm_mul_ps(beta2,rsq10);
1139             rinv3            = _mm_mul_ps(rinvsq10,rinv10);
1140             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
1141             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
1142             felec            = _mm_mul_ps(qq10,felec);
1143             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
1144             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv10);
1145             velec            = _mm_mul_ps(qq10,velec);
1146
1147             d                = _mm_sub_ps(r10,rswitch);
1148             d                = _mm_max_ps(d,_mm_setzero_ps());
1149             d2               = _mm_mul_ps(d,d);
1150             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
1151
1152             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
1153
1154             /* Evaluate switch function */
1155             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1156             felec            = _mm_msub_ps( felec,sw , _mm_mul_ps(rinv10,_mm_mul_ps(velec,dsw)) );
1157             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
1158
1159             fscal            = felec;
1160
1161             fscal            = _mm_and_ps(fscal,cutoff_mask);
1162
1163              /* Update vectorial force */
1164             fix1             = _mm_macc_ps(dx10,fscal,fix1);
1165             fiy1             = _mm_macc_ps(dy10,fscal,fiy1);
1166             fiz1             = _mm_macc_ps(dz10,fscal,fiz1);
1167
1168             fjx0             = _mm_macc_ps(dx10,fscal,fjx0);
1169             fjy0             = _mm_macc_ps(dy10,fscal,fjy0);
1170             fjz0             = _mm_macc_ps(dz10,fscal,fjz0);
1171
1172             }
1173
1174             /**************************
1175              * CALCULATE INTERACTIONS *
1176              **************************/
1177
1178             if (gmx_mm_any_lt(rsq20,rcutoff2))
1179             {
1180
1181             r20              = _mm_mul_ps(rsq20,rinv20);
1182
1183             /* Compute parameters for interactions between i and j atoms */
1184             qq20             = _mm_mul_ps(iq2,jq0);
1185
1186             /* EWALD ELECTROSTATICS */
1187
1188             /* Analytical PME correction */
1189             zeta2            = _mm_mul_ps(beta2,rsq20);
1190             rinv3            = _mm_mul_ps(rinvsq20,rinv20);
1191             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
1192             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
1193             felec            = _mm_mul_ps(qq20,felec);
1194             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
1195             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv20);
1196             velec            = _mm_mul_ps(qq20,velec);
1197
1198             d                = _mm_sub_ps(r20,rswitch);
1199             d                = _mm_max_ps(d,_mm_setzero_ps());
1200             d2               = _mm_mul_ps(d,d);
1201             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
1202
1203             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
1204
1205             /* Evaluate switch function */
1206             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1207             felec            = _mm_msub_ps( felec,sw , _mm_mul_ps(rinv20,_mm_mul_ps(velec,dsw)) );
1208             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
1209
1210             fscal            = felec;
1211
1212             fscal            = _mm_and_ps(fscal,cutoff_mask);
1213
1214              /* Update vectorial force */
1215             fix2             = _mm_macc_ps(dx20,fscal,fix2);
1216             fiy2             = _mm_macc_ps(dy20,fscal,fiy2);
1217             fiz2             = _mm_macc_ps(dz20,fscal,fiz2);
1218
1219             fjx0             = _mm_macc_ps(dx20,fscal,fjx0);
1220             fjy0             = _mm_macc_ps(dy20,fscal,fjy0);
1221             fjz0             = _mm_macc_ps(dz20,fscal,fjz0);
1222
1223             }
1224
1225             /**************************
1226              * CALCULATE INTERACTIONS *
1227              **************************/
1228
1229             if (gmx_mm_any_lt(rsq30,rcutoff2))
1230             {
1231
1232             r30              = _mm_mul_ps(rsq30,rinv30);
1233
1234             /* Compute parameters for interactions between i and j atoms */
1235             qq30             = _mm_mul_ps(iq3,jq0);
1236
1237             /* EWALD ELECTROSTATICS */
1238
1239             /* Analytical PME correction */
1240             zeta2            = _mm_mul_ps(beta2,rsq30);
1241             rinv3            = _mm_mul_ps(rinvsq30,rinv30);
1242             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
1243             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
1244             felec            = _mm_mul_ps(qq30,felec);
1245             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
1246             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv30);
1247             velec            = _mm_mul_ps(qq30,velec);
1248
1249             d                = _mm_sub_ps(r30,rswitch);
1250             d                = _mm_max_ps(d,_mm_setzero_ps());
1251             d2               = _mm_mul_ps(d,d);
1252             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
1253
1254             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
1255
1256             /* Evaluate switch function */
1257             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1258             felec            = _mm_msub_ps( felec,sw , _mm_mul_ps(rinv30,_mm_mul_ps(velec,dsw)) );
1259             cutoff_mask      = _mm_cmplt_ps(rsq30,rcutoff2);
1260
1261             fscal            = felec;
1262
1263             fscal            = _mm_and_ps(fscal,cutoff_mask);
1264
1265              /* Update vectorial force */
1266             fix3             = _mm_macc_ps(dx30,fscal,fix3);
1267             fiy3             = _mm_macc_ps(dy30,fscal,fiy3);
1268             fiz3             = _mm_macc_ps(dz30,fscal,fiz3);
1269
1270             fjx0             = _mm_macc_ps(dx30,fscal,fjx0);
1271             fjy0             = _mm_macc_ps(dy30,fscal,fjy0);
1272             fjz0             = _mm_macc_ps(dz30,fscal,fjz0);
1273
1274             }
1275
1276             fjptrA             = f+j_coord_offsetA;
1277             fjptrB             = f+j_coord_offsetB;
1278             fjptrC             = f+j_coord_offsetC;
1279             fjptrD             = f+j_coord_offsetD;
1280
1281             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1282
1283             /* Inner loop uses 209 flops */
1284         }
1285
1286         if(jidx<j_index_end)
1287         {
1288
1289             /* Get j neighbor index, and coordinate index */
1290             jnrlistA         = jjnr[jidx];
1291             jnrlistB         = jjnr[jidx+1];
1292             jnrlistC         = jjnr[jidx+2];
1293             jnrlistD         = jjnr[jidx+3];
1294             /* Sign of each element will be negative for non-real atoms.
1295              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1296              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1297              */
1298             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
1299             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
1300             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
1301             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
1302             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
1303             j_coord_offsetA  = DIM*jnrA;
1304             j_coord_offsetB  = DIM*jnrB;
1305             j_coord_offsetC  = DIM*jnrC;
1306             j_coord_offsetD  = DIM*jnrD;
1307
1308             /* load j atom coordinates */
1309             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1310                                               x+j_coord_offsetC,x+j_coord_offsetD,
1311                                               &jx0,&jy0,&jz0);
1312
1313             /* Calculate displacement vector */
1314             dx00             = _mm_sub_ps(ix0,jx0);
1315             dy00             = _mm_sub_ps(iy0,jy0);
1316             dz00             = _mm_sub_ps(iz0,jz0);
1317             dx10             = _mm_sub_ps(ix1,jx0);
1318             dy10             = _mm_sub_ps(iy1,jy0);
1319             dz10             = _mm_sub_ps(iz1,jz0);
1320             dx20             = _mm_sub_ps(ix2,jx0);
1321             dy20             = _mm_sub_ps(iy2,jy0);
1322             dz20             = _mm_sub_ps(iz2,jz0);
1323             dx30             = _mm_sub_ps(ix3,jx0);
1324             dy30             = _mm_sub_ps(iy3,jy0);
1325             dz30             = _mm_sub_ps(iz3,jz0);
1326
1327             /* Calculate squared distance and things based on it */
1328             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
1329             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
1330             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
1331             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
1332
1333             rinv00           = gmx_mm_invsqrt_ps(rsq00);
1334             rinv10           = gmx_mm_invsqrt_ps(rsq10);
1335             rinv20           = gmx_mm_invsqrt_ps(rsq20);
1336             rinv30           = gmx_mm_invsqrt_ps(rsq30);
1337
1338             rinvsq00         = _mm_mul_ps(rinv00,rinv00);
1339             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
1340             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
1341             rinvsq30         = _mm_mul_ps(rinv30,rinv30);
1342
1343             /* Load parameters for j particles */
1344             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
1345                                                               charge+jnrC+0,charge+jnrD+0);
1346             vdwjidx0A        = 2*vdwtype[jnrA+0];
1347             vdwjidx0B        = 2*vdwtype[jnrB+0];
1348             vdwjidx0C        = 2*vdwtype[jnrC+0];
1349             vdwjidx0D        = 2*vdwtype[jnrD+0];
1350
1351             fjx0             = _mm_setzero_ps();
1352             fjy0             = _mm_setzero_ps();
1353             fjz0             = _mm_setzero_ps();
1354
1355             /**************************
1356              * CALCULATE INTERACTIONS *
1357              **************************/
1358
1359             if (gmx_mm_any_lt(rsq00,rcutoff2))
1360             {
1361
1362             r00              = _mm_mul_ps(rsq00,rinv00);
1363             r00              = _mm_andnot_ps(dummy_mask,r00);
1364
1365             /* Compute parameters for interactions between i and j atoms */
1366             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
1367                                          vdwparam+vdwioffset0+vdwjidx0B,
1368                                          vdwparam+vdwioffset0+vdwjidx0C,
1369                                          vdwparam+vdwioffset0+vdwjidx0D,
1370                                          &c6_00,&c12_00);
1371
1372             /* LENNARD-JONES DISPERSION/REPULSION */
1373
1374             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
1375             vvdw6            = _mm_mul_ps(c6_00,rinvsix);
1376             vvdw12           = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
1377             vvdw             = _mm_msub_ps(vvdw12,one_twelfth,_mm_mul_ps(vvdw6,one_sixth));
1378             fvdw             = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
1379
1380             d                = _mm_sub_ps(r00,rswitch);
1381             d                = _mm_max_ps(d,_mm_setzero_ps());
1382             d2               = _mm_mul_ps(d,d);
1383             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
1384
1385             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
1386
1387             /* Evaluate switch function */
1388             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1389             fvdw             = _mm_msub_ps( fvdw,sw , _mm_mul_ps(rinv00,_mm_mul_ps(vvdw,dsw)) );
1390             cutoff_mask      = _mm_cmplt_ps(rsq00,rcutoff2);
1391
1392             fscal            = fvdw;
1393
1394             fscal            = _mm_and_ps(fscal,cutoff_mask);
1395
1396             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1397
1398              /* Update vectorial force */
1399             fix0             = _mm_macc_ps(dx00,fscal,fix0);
1400             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
1401             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
1402
1403             fjx0             = _mm_macc_ps(dx00,fscal,fjx0);
1404             fjy0             = _mm_macc_ps(dy00,fscal,fjy0);
1405             fjz0             = _mm_macc_ps(dz00,fscal,fjz0);
1406
1407             }
1408
1409             /**************************
1410              * CALCULATE INTERACTIONS *
1411              **************************/
1412
1413             if (gmx_mm_any_lt(rsq10,rcutoff2))
1414             {
1415
1416             r10              = _mm_mul_ps(rsq10,rinv10);
1417             r10              = _mm_andnot_ps(dummy_mask,r10);
1418
1419             /* Compute parameters for interactions between i and j atoms */
1420             qq10             = _mm_mul_ps(iq1,jq0);
1421
1422             /* EWALD ELECTROSTATICS */
1423
1424             /* Analytical PME correction */
1425             zeta2            = _mm_mul_ps(beta2,rsq10);
1426             rinv3            = _mm_mul_ps(rinvsq10,rinv10);
1427             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
1428             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
1429             felec            = _mm_mul_ps(qq10,felec);
1430             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
1431             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv10);
1432             velec            = _mm_mul_ps(qq10,velec);
1433
1434             d                = _mm_sub_ps(r10,rswitch);
1435             d                = _mm_max_ps(d,_mm_setzero_ps());
1436             d2               = _mm_mul_ps(d,d);
1437             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
1438
1439             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
1440
1441             /* Evaluate switch function */
1442             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1443             felec            = _mm_msub_ps( felec,sw , _mm_mul_ps(rinv10,_mm_mul_ps(velec,dsw)) );
1444             cutoff_mask      = _mm_cmplt_ps(rsq10,rcutoff2);
1445
1446             fscal            = felec;
1447
1448             fscal            = _mm_and_ps(fscal,cutoff_mask);
1449
1450             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1451
1452              /* Update vectorial force */
1453             fix1             = _mm_macc_ps(dx10,fscal,fix1);
1454             fiy1             = _mm_macc_ps(dy10,fscal,fiy1);
1455             fiz1             = _mm_macc_ps(dz10,fscal,fiz1);
1456
1457             fjx0             = _mm_macc_ps(dx10,fscal,fjx0);
1458             fjy0             = _mm_macc_ps(dy10,fscal,fjy0);
1459             fjz0             = _mm_macc_ps(dz10,fscal,fjz0);
1460
1461             }
1462
1463             /**************************
1464              * CALCULATE INTERACTIONS *
1465              **************************/
1466
1467             if (gmx_mm_any_lt(rsq20,rcutoff2))
1468             {
1469
1470             r20              = _mm_mul_ps(rsq20,rinv20);
1471             r20              = _mm_andnot_ps(dummy_mask,r20);
1472
1473             /* Compute parameters for interactions between i and j atoms */
1474             qq20             = _mm_mul_ps(iq2,jq0);
1475
1476             /* EWALD ELECTROSTATICS */
1477
1478             /* Analytical PME correction */
1479             zeta2            = _mm_mul_ps(beta2,rsq20);
1480             rinv3            = _mm_mul_ps(rinvsq20,rinv20);
1481             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
1482             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
1483             felec            = _mm_mul_ps(qq20,felec);
1484             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
1485             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv20);
1486             velec            = _mm_mul_ps(qq20,velec);
1487
1488             d                = _mm_sub_ps(r20,rswitch);
1489             d                = _mm_max_ps(d,_mm_setzero_ps());
1490             d2               = _mm_mul_ps(d,d);
1491             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
1492
1493             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
1494
1495             /* Evaluate switch function */
1496             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1497             felec            = _mm_msub_ps( felec,sw , _mm_mul_ps(rinv20,_mm_mul_ps(velec,dsw)) );
1498             cutoff_mask      = _mm_cmplt_ps(rsq20,rcutoff2);
1499
1500             fscal            = felec;
1501
1502             fscal            = _mm_and_ps(fscal,cutoff_mask);
1503
1504             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1505
1506              /* Update vectorial force */
1507             fix2             = _mm_macc_ps(dx20,fscal,fix2);
1508             fiy2             = _mm_macc_ps(dy20,fscal,fiy2);
1509             fiz2             = _mm_macc_ps(dz20,fscal,fiz2);
1510
1511             fjx0             = _mm_macc_ps(dx20,fscal,fjx0);
1512             fjy0             = _mm_macc_ps(dy20,fscal,fjy0);
1513             fjz0             = _mm_macc_ps(dz20,fscal,fjz0);
1514
1515             }
1516
1517             /**************************
1518              * CALCULATE INTERACTIONS *
1519              **************************/
1520
1521             if (gmx_mm_any_lt(rsq30,rcutoff2))
1522             {
1523
1524             r30              = _mm_mul_ps(rsq30,rinv30);
1525             r30              = _mm_andnot_ps(dummy_mask,r30);
1526
1527             /* Compute parameters for interactions between i and j atoms */
1528             qq30             = _mm_mul_ps(iq3,jq0);
1529
1530             /* EWALD ELECTROSTATICS */
1531
1532             /* Analytical PME correction */
1533             zeta2            = _mm_mul_ps(beta2,rsq30);
1534             rinv3            = _mm_mul_ps(rinvsq30,rinv30);
1535             pmecorrF         = gmx_mm_pmecorrF_ps(zeta2);
1536             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
1537             felec            = _mm_mul_ps(qq30,felec);
1538             pmecorrV         = gmx_mm_pmecorrV_ps(zeta2);
1539             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv30);
1540             velec            = _mm_mul_ps(qq30,velec);
1541
1542             d                = _mm_sub_ps(r30,rswitch);
1543             d                = _mm_max_ps(d,_mm_setzero_ps());
1544             d2               = _mm_mul_ps(d,d);
1545             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
1546
1547             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
1548
1549             /* Evaluate switch function */
1550             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
1551             felec            = _mm_msub_ps( felec,sw , _mm_mul_ps(rinv30,_mm_mul_ps(velec,dsw)) );
1552             cutoff_mask      = _mm_cmplt_ps(rsq30,rcutoff2);
1553
1554             fscal            = felec;
1555
1556             fscal            = _mm_and_ps(fscal,cutoff_mask);
1557
1558             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1559
1560              /* Update vectorial force */
1561             fix3             = _mm_macc_ps(dx30,fscal,fix3);
1562             fiy3             = _mm_macc_ps(dy30,fscal,fiy3);
1563             fiz3             = _mm_macc_ps(dz30,fscal,fiz3);
1564
1565             fjx0             = _mm_macc_ps(dx30,fscal,fjx0);
1566             fjy0             = _mm_macc_ps(dy30,fscal,fjy0);
1567             fjz0             = _mm_macc_ps(dz30,fscal,fjz0);
1568
1569             }
1570
1571             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1572             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1573             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1574             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1575
1576             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1577
1578             /* Inner loop uses 213 flops */
1579         }
1580
1581         /* End of innermost loop */
1582
1583         gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1584                                               f+i_coord_offset,fshift+i_shift_offset);
1585
1586         /* Increment number of inner iterations */
1587         inneriter                  += j_index_end - j_index_start;
1588
1589         /* Outer loop uses 24 flops */
1590     }
1591
1592     /* Increment number of outer iterations */
1593     outeriter        += nri;
1594
1595     /* Update outer/inner flops */
1596
1597     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*213);
1598 }