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