Added option to gmx nmeig to print ZPE.
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_single / nb_kernel_ElecCoul_VdwLJ_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,2015,2017, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS avx_128_fma_single kernel generator.
37  */
38 #include "gmxpre.h"
39
40 #include "config.h"
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "gromacs/gmxlib/nrnb.h"
46
47 #include "kernelutil_x86_avx_128_fma_single.h"
48
49 /*
50  * Gromacs nonbonded kernel:   nb_kernel_ElecCoul_VdwLJ_GeomW4P1_VF_avx_128_fma_single
51  * Electrostatics interaction: Coulomb
52  * VdW interaction:            LennardJones
53  * Geometry:                   Water4-Particle
54  * Calculate force/pot:        PotentialAndForce
55  */
56 void
57 nb_kernel_ElecCoul_VdwLJ_GeomW4P1_VF_avx_128_fma_single
58                     (t_nblist                    * gmx_restrict       nlist,
59                      rvec                        * gmx_restrict          xx,
60                      rvec                        * gmx_restrict          ff,
61                      struct t_forcerec           * gmx_restrict          fr,
62                      t_mdatoms                   * gmx_restrict     mdatoms,
63                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
64                      t_nrnb                      * gmx_restrict        nrnb)
65 {
66     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
67      * just 0 for non-waters.
68      * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
69      * jnr indices corresponding to data put in the four positions in the SIMD register.
70      */
71     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
72     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
73     int              jnrA,jnrB,jnrC,jnrD;
74     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
75     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
76     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
77     real             rcutoff_scalar;
78     real             *shiftvec,*fshift,*x,*f;
79     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
80     real             scratch[4*DIM];
81     __m128           fscal,rcutoff,rcutoff2,jidxall;
82     int              vdwioffset0;
83     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
84     int              vdwioffset1;
85     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
86     int              vdwioffset2;
87     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
88     int              vdwioffset3;
89     __m128           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
90     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
91     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
92     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
93     __m128           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
94     __m128           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
95     __m128           dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
96     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
97     real             *charge;
98     int              nvdwtype;
99     __m128           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
100     int              *vdwtype;
101     real             *vdwparam;
102     __m128           one_sixth   = _mm_set1_ps(1.0/6.0);
103     __m128           one_twelfth = _mm_set1_ps(1.0/12.0);
104     __m128           dummy_mask,cutoff_mask;
105     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
106     __m128           one     = _mm_set1_ps(1.0);
107     __m128           two     = _mm_set1_ps(2.0);
108     x                = xx[0];
109     f                = ff[0];
110
111     nri              = nlist->nri;
112     iinr             = nlist->iinr;
113     jindex           = nlist->jindex;
114     jjnr             = nlist->jjnr;
115     shiftidx         = nlist->shift;
116     gid              = nlist->gid;
117     shiftvec         = fr->shift_vec[0];
118     fshift           = fr->fshift[0];
119     facel            = _mm_set1_ps(fr->ic->epsfac);
120     charge           = mdatoms->chargeA;
121     nvdwtype         = fr->ntype;
122     vdwparam         = fr->nbfp;
123     vdwtype          = mdatoms->typeA;
124
125     /* Setup water-specific parameters */
126     inr              = nlist->iinr[0];
127     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
128     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
129     iq3              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
130     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
131
132     /* Avoid stupid compiler warnings */
133     jnrA = jnrB = jnrC = jnrD = 0;
134     j_coord_offsetA = 0;
135     j_coord_offsetB = 0;
136     j_coord_offsetC = 0;
137     j_coord_offsetD = 0;
138
139     outeriter        = 0;
140     inneriter        = 0;
141
142     for(iidx=0;iidx<4*DIM;iidx++)
143     {
144         scratch[iidx] = 0.0;
145     }
146
147     /* Start outer loop over neighborlists */
148     for(iidx=0; iidx<nri; iidx++)
149     {
150         /* Load shift vector for this list */
151         i_shift_offset   = DIM*shiftidx[iidx];
152
153         /* Load limits for loop over neighbors */
154         j_index_start    = jindex[iidx];
155         j_index_end      = jindex[iidx+1];
156
157         /* Get outer coordinate index */
158         inr              = iinr[iidx];
159         i_coord_offset   = DIM*inr;
160
161         /* Load i particle coords and add shift vector */
162         gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
163                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
164
165         fix0             = _mm_setzero_ps();
166         fiy0             = _mm_setzero_ps();
167         fiz0             = _mm_setzero_ps();
168         fix1             = _mm_setzero_ps();
169         fiy1             = _mm_setzero_ps();
170         fiz1             = _mm_setzero_ps();
171         fix2             = _mm_setzero_ps();
172         fiy2             = _mm_setzero_ps();
173         fiz2             = _mm_setzero_ps();
174         fix3             = _mm_setzero_ps();
175         fiy3             = _mm_setzero_ps();
176         fiz3             = _mm_setzero_ps();
177
178         /* Reset potential sums */
179         velecsum         = _mm_setzero_ps();
180         vvdwsum          = _mm_setzero_ps();
181
182         /* Start inner kernel loop */
183         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
184         {
185
186             /* Get j neighbor index, and coordinate index */
187             jnrA             = jjnr[jidx];
188             jnrB             = jjnr[jidx+1];
189             jnrC             = jjnr[jidx+2];
190             jnrD             = jjnr[jidx+3];
191             j_coord_offsetA  = DIM*jnrA;
192             j_coord_offsetB  = DIM*jnrB;
193             j_coord_offsetC  = DIM*jnrC;
194             j_coord_offsetD  = DIM*jnrD;
195
196             /* load j atom coordinates */
197             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
198                                               x+j_coord_offsetC,x+j_coord_offsetD,
199                                               &jx0,&jy0,&jz0);
200
201             /* Calculate displacement vector */
202             dx00             = _mm_sub_ps(ix0,jx0);
203             dy00             = _mm_sub_ps(iy0,jy0);
204             dz00             = _mm_sub_ps(iz0,jz0);
205             dx10             = _mm_sub_ps(ix1,jx0);
206             dy10             = _mm_sub_ps(iy1,jy0);
207             dz10             = _mm_sub_ps(iz1,jz0);
208             dx20             = _mm_sub_ps(ix2,jx0);
209             dy20             = _mm_sub_ps(iy2,jy0);
210             dz20             = _mm_sub_ps(iz2,jz0);
211             dx30             = _mm_sub_ps(ix3,jx0);
212             dy30             = _mm_sub_ps(iy3,jy0);
213             dz30             = _mm_sub_ps(iz3,jz0);
214
215             /* Calculate squared distance and things based on it */
216             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
217             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
218             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
219             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
220
221             rinv10           = avx128fma_invsqrt_f(rsq10);
222             rinv20           = avx128fma_invsqrt_f(rsq20);
223             rinv30           = avx128fma_invsqrt_f(rsq30);
224
225             rinvsq00         = avx128fma_inv_f(rsq00);
226             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
227             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
228             rinvsq30         = _mm_mul_ps(rinv30,rinv30);
229
230             /* Load parameters for j particles */
231             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
232                                                               charge+jnrC+0,charge+jnrD+0);
233             vdwjidx0A        = 2*vdwtype[jnrA+0];
234             vdwjidx0B        = 2*vdwtype[jnrB+0];
235             vdwjidx0C        = 2*vdwtype[jnrC+0];
236             vdwjidx0D        = 2*vdwtype[jnrD+0];
237
238             fjx0             = _mm_setzero_ps();
239             fjy0             = _mm_setzero_ps();
240             fjz0             = _mm_setzero_ps();
241
242             /**************************
243              * CALCULATE INTERACTIONS *
244              **************************/
245
246             /* Compute parameters for interactions between i and j atoms */
247             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
248                                          vdwparam+vdwioffset0+vdwjidx0B,
249                                          vdwparam+vdwioffset0+vdwjidx0C,
250                                          vdwparam+vdwioffset0+vdwjidx0D,
251                                          &c6_00,&c12_00);
252
253             /* LENNARD-JONES DISPERSION/REPULSION */
254
255             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
256             vvdw6            = _mm_mul_ps(c6_00,rinvsix);
257             vvdw12           = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
258             vvdw             = _mm_msub_ps(vvdw12,one_twelfth,_mm_mul_ps(vvdw6,one_sixth));
259             fvdw             = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
260
261             /* Update potential sum for this i atom from the interaction with this j atom. */
262             vvdwsum          = _mm_add_ps(vvdwsum,vvdw);
263
264             fscal            = fvdw;
265
266              /* Update vectorial force */
267             fix0             = _mm_macc_ps(dx00,fscal,fix0);
268             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
269             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
270
271             fjx0             = _mm_macc_ps(dx00,fscal,fjx0);
272             fjy0             = _mm_macc_ps(dy00,fscal,fjy0);
273             fjz0             = _mm_macc_ps(dz00,fscal,fjz0);
274
275             /**************************
276              * CALCULATE INTERACTIONS *
277              **************************/
278
279             /* Compute parameters for interactions between i and j atoms */
280             qq10             = _mm_mul_ps(iq1,jq0);
281
282             /* COULOMB ELECTROSTATICS */
283             velec            = _mm_mul_ps(qq10,rinv10);
284             felec            = _mm_mul_ps(velec,rinvsq10);
285
286             /* Update potential sum for this i atom from the interaction with this j atom. */
287             velecsum         = _mm_add_ps(velecsum,velec);
288
289             fscal            = felec;
290
291              /* Update vectorial force */
292             fix1             = _mm_macc_ps(dx10,fscal,fix1);
293             fiy1             = _mm_macc_ps(dy10,fscal,fiy1);
294             fiz1             = _mm_macc_ps(dz10,fscal,fiz1);
295
296             fjx0             = _mm_macc_ps(dx10,fscal,fjx0);
297             fjy0             = _mm_macc_ps(dy10,fscal,fjy0);
298             fjz0             = _mm_macc_ps(dz10,fscal,fjz0);
299
300             /**************************
301              * CALCULATE INTERACTIONS *
302              **************************/
303
304             /* Compute parameters for interactions between i and j atoms */
305             qq20             = _mm_mul_ps(iq2,jq0);
306
307             /* COULOMB ELECTROSTATICS */
308             velec            = _mm_mul_ps(qq20,rinv20);
309             felec            = _mm_mul_ps(velec,rinvsq20);
310
311             /* Update potential sum for this i atom from the interaction with this j atom. */
312             velecsum         = _mm_add_ps(velecsum,velec);
313
314             fscal            = felec;
315
316              /* Update vectorial force */
317             fix2             = _mm_macc_ps(dx20,fscal,fix2);
318             fiy2             = _mm_macc_ps(dy20,fscal,fiy2);
319             fiz2             = _mm_macc_ps(dz20,fscal,fiz2);
320
321             fjx0             = _mm_macc_ps(dx20,fscal,fjx0);
322             fjy0             = _mm_macc_ps(dy20,fscal,fjy0);
323             fjz0             = _mm_macc_ps(dz20,fscal,fjz0);
324
325             /**************************
326              * CALCULATE INTERACTIONS *
327              **************************/
328
329             /* Compute parameters for interactions between i and j atoms */
330             qq30             = _mm_mul_ps(iq3,jq0);
331
332             /* COULOMB ELECTROSTATICS */
333             velec            = _mm_mul_ps(qq30,rinv30);
334             felec            = _mm_mul_ps(velec,rinvsq30);
335
336             /* Update potential sum for this i atom from the interaction with this j atom. */
337             velecsum         = _mm_add_ps(velecsum,velec);
338
339             fscal            = felec;
340
341              /* Update vectorial force */
342             fix3             = _mm_macc_ps(dx30,fscal,fix3);
343             fiy3             = _mm_macc_ps(dy30,fscal,fiy3);
344             fiz3             = _mm_macc_ps(dz30,fscal,fiz3);
345
346             fjx0             = _mm_macc_ps(dx30,fscal,fjx0);
347             fjy0             = _mm_macc_ps(dy30,fscal,fjy0);
348             fjz0             = _mm_macc_ps(dz30,fscal,fjz0);
349
350             fjptrA             = f+j_coord_offsetA;
351             fjptrB             = f+j_coord_offsetB;
352             fjptrC             = f+j_coord_offsetC;
353             fjptrD             = f+j_coord_offsetD;
354
355             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
356
357             /* Inner loop uses 128 flops */
358         }
359
360         if(jidx<j_index_end)
361         {
362
363             /* Get j neighbor index, and coordinate index */
364             jnrlistA         = jjnr[jidx];
365             jnrlistB         = jjnr[jidx+1];
366             jnrlistC         = jjnr[jidx+2];
367             jnrlistD         = jjnr[jidx+3];
368             /* Sign of each element will be negative for non-real atoms.
369              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
370              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
371              */
372             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
373             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
374             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
375             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
376             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
377             j_coord_offsetA  = DIM*jnrA;
378             j_coord_offsetB  = DIM*jnrB;
379             j_coord_offsetC  = DIM*jnrC;
380             j_coord_offsetD  = DIM*jnrD;
381
382             /* load j atom coordinates */
383             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
384                                               x+j_coord_offsetC,x+j_coord_offsetD,
385                                               &jx0,&jy0,&jz0);
386
387             /* Calculate displacement vector */
388             dx00             = _mm_sub_ps(ix0,jx0);
389             dy00             = _mm_sub_ps(iy0,jy0);
390             dz00             = _mm_sub_ps(iz0,jz0);
391             dx10             = _mm_sub_ps(ix1,jx0);
392             dy10             = _mm_sub_ps(iy1,jy0);
393             dz10             = _mm_sub_ps(iz1,jz0);
394             dx20             = _mm_sub_ps(ix2,jx0);
395             dy20             = _mm_sub_ps(iy2,jy0);
396             dz20             = _mm_sub_ps(iz2,jz0);
397             dx30             = _mm_sub_ps(ix3,jx0);
398             dy30             = _mm_sub_ps(iy3,jy0);
399             dz30             = _mm_sub_ps(iz3,jz0);
400
401             /* Calculate squared distance and things based on it */
402             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
403             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
404             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
405             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
406
407             rinv10           = avx128fma_invsqrt_f(rsq10);
408             rinv20           = avx128fma_invsqrt_f(rsq20);
409             rinv30           = avx128fma_invsqrt_f(rsq30);
410
411             rinvsq00         = avx128fma_inv_f(rsq00);
412             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
413             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
414             rinvsq30         = _mm_mul_ps(rinv30,rinv30);
415
416             /* Load parameters for j particles */
417             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
418                                                               charge+jnrC+0,charge+jnrD+0);
419             vdwjidx0A        = 2*vdwtype[jnrA+0];
420             vdwjidx0B        = 2*vdwtype[jnrB+0];
421             vdwjidx0C        = 2*vdwtype[jnrC+0];
422             vdwjidx0D        = 2*vdwtype[jnrD+0];
423
424             fjx0             = _mm_setzero_ps();
425             fjy0             = _mm_setzero_ps();
426             fjz0             = _mm_setzero_ps();
427
428             /**************************
429              * CALCULATE INTERACTIONS *
430              **************************/
431
432             /* Compute parameters for interactions between i and j atoms */
433             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
434                                          vdwparam+vdwioffset0+vdwjidx0B,
435                                          vdwparam+vdwioffset0+vdwjidx0C,
436                                          vdwparam+vdwioffset0+vdwjidx0D,
437                                          &c6_00,&c12_00);
438
439             /* LENNARD-JONES DISPERSION/REPULSION */
440
441             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
442             vvdw6            = _mm_mul_ps(c6_00,rinvsix);
443             vvdw12           = _mm_mul_ps(c12_00,_mm_mul_ps(rinvsix,rinvsix));
444             vvdw             = _mm_msub_ps(vvdw12,one_twelfth,_mm_mul_ps(vvdw6,one_sixth));
445             fvdw             = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq00);
446
447             /* Update potential sum for this i atom from the interaction with this j atom. */
448             vvdw             = _mm_andnot_ps(dummy_mask,vvdw);
449             vvdwsum          = _mm_add_ps(vvdwsum,vvdw);
450
451             fscal            = fvdw;
452
453             fscal            = _mm_andnot_ps(dummy_mask,fscal);
454
455              /* Update vectorial force */
456             fix0             = _mm_macc_ps(dx00,fscal,fix0);
457             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
458             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
459
460             fjx0             = _mm_macc_ps(dx00,fscal,fjx0);
461             fjy0             = _mm_macc_ps(dy00,fscal,fjy0);
462             fjz0             = _mm_macc_ps(dz00,fscal,fjz0);
463
464             /**************************
465              * CALCULATE INTERACTIONS *
466              **************************/
467
468             /* Compute parameters for interactions between i and j atoms */
469             qq10             = _mm_mul_ps(iq1,jq0);
470
471             /* COULOMB ELECTROSTATICS */
472             velec            = _mm_mul_ps(qq10,rinv10);
473             felec            = _mm_mul_ps(velec,rinvsq10);
474
475             /* Update potential sum for this i atom from the interaction with this j atom. */
476             velec            = _mm_andnot_ps(dummy_mask,velec);
477             velecsum         = _mm_add_ps(velecsum,velec);
478
479             fscal            = felec;
480
481             fscal            = _mm_andnot_ps(dummy_mask,fscal);
482
483              /* Update vectorial force */
484             fix1             = _mm_macc_ps(dx10,fscal,fix1);
485             fiy1             = _mm_macc_ps(dy10,fscal,fiy1);
486             fiz1             = _mm_macc_ps(dz10,fscal,fiz1);
487
488             fjx0             = _mm_macc_ps(dx10,fscal,fjx0);
489             fjy0             = _mm_macc_ps(dy10,fscal,fjy0);
490             fjz0             = _mm_macc_ps(dz10,fscal,fjz0);
491
492             /**************************
493              * CALCULATE INTERACTIONS *
494              **************************/
495
496             /* Compute parameters for interactions between i and j atoms */
497             qq20             = _mm_mul_ps(iq2,jq0);
498
499             /* COULOMB ELECTROSTATICS */
500             velec            = _mm_mul_ps(qq20,rinv20);
501             felec            = _mm_mul_ps(velec,rinvsq20);
502
503             /* Update potential sum for this i atom from the interaction with this j atom. */
504             velec            = _mm_andnot_ps(dummy_mask,velec);
505             velecsum         = _mm_add_ps(velecsum,velec);
506
507             fscal            = felec;
508
509             fscal            = _mm_andnot_ps(dummy_mask,fscal);
510
511              /* Update vectorial force */
512             fix2             = _mm_macc_ps(dx20,fscal,fix2);
513             fiy2             = _mm_macc_ps(dy20,fscal,fiy2);
514             fiz2             = _mm_macc_ps(dz20,fscal,fiz2);
515
516             fjx0             = _mm_macc_ps(dx20,fscal,fjx0);
517             fjy0             = _mm_macc_ps(dy20,fscal,fjy0);
518             fjz0             = _mm_macc_ps(dz20,fscal,fjz0);
519
520             /**************************
521              * CALCULATE INTERACTIONS *
522              **************************/
523
524             /* Compute parameters for interactions between i and j atoms */
525             qq30             = _mm_mul_ps(iq3,jq0);
526
527             /* COULOMB ELECTROSTATICS */
528             velec            = _mm_mul_ps(qq30,rinv30);
529             felec            = _mm_mul_ps(velec,rinvsq30);
530
531             /* Update potential sum for this i atom from the interaction with this j atom. */
532             velec            = _mm_andnot_ps(dummy_mask,velec);
533             velecsum         = _mm_add_ps(velecsum,velec);
534
535             fscal            = felec;
536
537             fscal            = _mm_andnot_ps(dummy_mask,fscal);
538
539              /* Update vectorial force */
540             fix3             = _mm_macc_ps(dx30,fscal,fix3);
541             fiy3             = _mm_macc_ps(dy30,fscal,fiy3);
542             fiz3             = _mm_macc_ps(dz30,fscal,fiz3);
543
544             fjx0             = _mm_macc_ps(dx30,fscal,fjx0);
545             fjy0             = _mm_macc_ps(dy30,fscal,fjy0);
546             fjz0             = _mm_macc_ps(dz30,fscal,fjz0);
547
548             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
549             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
550             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
551             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
552
553             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
554
555             /* Inner loop uses 128 flops */
556         }
557
558         /* End of innermost loop */
559
560         gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
561                                               f+i_coord_offset,fshift+i_shift_offset);
562
563         ggid                        = gid[iidx];
564         /* Update potential energies */
565         gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
566         gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
567
568         /* Increment number of inner iterations */
569         inneriter                  += j_index_end - j_index_start;
570
571         /* Outer loop uses 26 flops */
572     }
573
574     /* Increment number of outer iterations */
575     outeriter        += nri;
576
577     /* Update outer/inner flops */
578
579     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_VF,outeriter*26 + inneriter*128);
580 }
581 /*
582  * Gromacs nonbonded kernel:   nb_kernel_ElecCoul_VdwLJ_GeomW4P1_F_avx_128_fma_single
583  * Electrostatics interaction: Coulomb
584  * VdW interaction:            LennardJones
585  * Geometry:                   Water4-Particle
586  * Calculate force/pot:        Force
587  */
588 void
589 nb_kernel_ElecCoul_VdwLJ_GeomW4P1_F_avx_128_fma_single
590                     (t_nblist                    * gmx_restrict       nlist,
591                      rvec                        * gmx_restrict          xx,
592                      rvec                        * gmx_restrict          ff,
593                      struct t_forcerec           * gmx_restrict          fr,
594                      t_mdatoms                   * gmx_restrict     mdatoms,
595                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
596                      t_nrnb                      * gmx_restrict        nrnb)
597 {
598     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
599      * just 0 for non-waters.
600      * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
601      * jnr indices corresponding to data put in the four positions in the SIMD register.
602      */
603     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
604     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
605     int              jnrA,jnrB,jnrC,jnrD;
606     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
607     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
608     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
609     real             rcutoff_scalar;
610     real             *shiftvec,*fshift,*x,*f;
611     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
612     real             scratch[4*DIM];
613     __m128           fscal,rcutoff,rcutoff2,jidxall;
614     int              vdwioffset0;
615     __m128           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
616     int              vdwioffset1;
617     __m128           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
618     int              vdwioffset2;
619     __m128           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
620     int              vdwioffset3;
621     __m128           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
622     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
623     __m128           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
624     __m128           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
625     __m128           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
626     __m128           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
627     __m128           dx30,dy30,dz30,rsq30,rinv30,rinvsq30,r30,qq30,c6_30,c12_30;
628     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
629     real             *charge;
630     int              nvdwtype;
631     __m128           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
632     int              *vdwtype;
633     real             *vdwparam;
634     __m128           one_sixth   = _mm_set1_ps(1.0/6.0);
635     __m128           one_twelfth = _mm_set1_ps(1.0/12.0);
636     __m128           dummy_mask,cutoff_mask;
637     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
638     __m128           one     = _mm_set1_ps(1.0);
639     __m128           two     = _mm_set1_ps(2.0);
640     x                = xx[0];
641     f                = ff[0];
642
643     nri              = nlist->nri;
644     iinr             = nlist->iinr;
645     jindex           = nlist->jindex;
646     jjnr             = nlist->jjnr;
647     shiftidx         = nlist->shift;
648     gid              = nlist->gid;
649     shiftvec         = fr->shift_vec[0];
650     fshift           = fr->fshift[0];
651     facel            = _mm_set1_ps(fr->ic->epsfac);
652     charge           = mdatoms->chargeA;
653     nvdwtype         = fr->ntype;
654     vdwparam         = fr->nbfp;
655     vdwtype          = mdatoms->typeA;
656
657     /* Setup water-specific parameters */
658     inr              = nlist->iinr[0];
659     iq1              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+1]));
660     iq2              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+2]));
661     iq3              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+3]));
662     vdwioffset0      = 2*nvdwtype*vdwtype[inr+0];
663
664     /* Avoid stupid compiler warnings */
665     jnrA = jnrB = jnrC = jnrD = 0;
666     j_coord_offsetA = 0;
667     j_coord_offsetB = 0;
668     j_coord_offsetC = 0;
669     j_coord_offsetD = 0;
670
671     outeriter        = 0;
672     inneriter        = 0;
673
674     for(iidx=0;iidx<4*DIM;iidx++)
675     {
676         scratch[iidx] = 0.0;
677     }
678
679     /* Start outer loop over neighborlists */
680     for(iidx=0; iidx<nri; iidx++)
681     {
682         /* Load shift vector for this list */
683         i_shift_offset   = DIM*shiftidx[iidx];
684
685         /* Load limits for loop over neighbors */
686         j_index_start    = jindex[iidx];
687         j_index_end      = jindex[iidx+1];
688
689         /* Get outer coordinate index */
690         inr              = iinr[iidx];
691         i_coord_offset   = DIM*inr;
692
693         /* Load i particle coords and add shift vector */
694         gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
695                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
696
697         fix0             = _mm_setzero_ps();
698         fiy0             = _mm_setzero_ps();
699         fiz0             = _mm_setzero_ps();
700         fix1             = _mm_setzero_ps();
701         fiy1             = _mm_setzero_ps();
702         fiz1             = _mm_setzero_ps();
703         fix2             = _mm_setzero_ps();
704         fiy2             = _mm_setzero_ps();
705         fiz2             = _mm_setzero_ps();
706         fix3             = _mm_setzero_ps();
707         fiy3             = _mm_setzero_ps();
708         fiz3             = _mm_setzero_ps();
709
710         /* Start inner kernel loop */
711         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
712         {
713
714             /* Get j neighbor index, and coordinate index */
715             jnrA             = jjnr[jidx];
716             jnrB             = jjnr[jidx+1];
717             jnrC             = jjnr[jidx+2];
718             jnrD             = jjnr[jidx+3];
719             j_coord_offsetA  = DIM*jnrA;
720             j_coord_offsetB  = DIM*jnrB;
721             j_coord_offsetC  = DIM*jnrC;
722             j_coord_offsetD  = DIM*jnrD;
723
724             /* load j atom coordinates */
725             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
726                                               x+j_coord_offsetC,x+j_coord_offsetD,
727                                               &jx0,&jy0,&jz0);
728
729             /* Calculate displacement vector */
730             dx00             = _mm_sub_ps(ix0,jx0);
731             dy00             = _mm_sub_ps(iy0,jy0);
732             dz00             = _mm_sub_ps(iz0,jz0);
733             dx10             = _mm_sub_ps(ix1,jx0);
734             dy10             = _mm_sub_ps(iy1,jy0);
735             dz10             = _mm_sub_ps(iz1,jz0);
736             dx20             = _mm_sub_ps(ix2,jx0);
737             dy20             = _mm_sub_ps(iy2,jy0);
738             dz20             = _mm_sub_ps(iz2,jz0);
739             dx30             = _mm_sub_ps(ix3,jx0);
740             dy30             = _mm_sub_ps(iy3,jy0);
741             dz30             = _mm_sub_ps(iz3,jz0);
742
743             /* Calculate squared distance and things based on it */
744             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
745             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
746             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
747             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
748
749             rinv10           = avx128fma_invsqrt_f(rsq10);
750             rinv20           = avx128fma_invsqrt_f(rsq20);
751             rinv30           = avx128fma_invsqrt_f(rsq30);
752
753             rinvsq00         = avx128fma_inv_f(rsq00);
754             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
755             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
756             rinvsq30         = _mm_mul_ps(rinv30,rinv30);
757
758             /* Load parameters for j particles */
759             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
760                                                               charge+jnrC+0,charge+jnrD+0);
761             vdwjidx0A        = 2*vdwtype[jnrA+0];
762             vdwjidx0B        = 2*vdwtype[jnrB+0];
763             vdwjidx0C        = 2*vdwtype[jnrC+0];
764             vdwjidx0D        = 2*vdwtype[jnrD+0];
765
766             fjx0             = _mm_setzero_ps();
767             fjy0             = _mm_setzero_ps();
768             fjz0             = _mm_setzero_ps();
769
770             /**************************
771              * CALCULATE INTERACTIONS *
772              **************************/
773
774             /* Compute parameters for interactions between i and j atoms */
775             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
776                                          vdwparam+vdwioffset0+vdwjidx0B,
777                                          vdwparam+vdwioffset0+vdwjidx0C,
778                                          vdwparam+vdwioffset0+vdwjidx0D,
779                                          &c6_00,&c12_00);
780
781             /* LENNARD-JONES DISPERSION/REPULSION */
782
783             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
784             fvdw             = _mm_mul_ps(_mm_msub_ps(c12_00,rinvsix,c6_00),_mm_mul_ps(rinvsix,rinvsq00));
785
786             fscal            = fvdw;
787
788              /* Update vectorial force */
789             fix0             = _mm_macc_ps(dx00,fscal,fix0);
790             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
791             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
792
793             fjx0             = _mm_macc_ps(dx00,fscal,fjx0);
794             fjy0             = _mm_macc_ps(dy00,fscal,fjy0);
795             fjz0             = _mm_macc_ps(dz00,fscal,fjz0);
796
797             /**************************
798              * CALCULATE INTERACTIONS *
799              **************************/
800
801             /* Compute parameters for interactions between i and j atoms */
802             qq10             = _mm_mul_ps(iq1,jq0);
803
804             /* COULOMB ELECTROSTATICS */
805             velec            = _mm_mul_ps(qq10,rinv10);
806             felec            = _mm_mul_ps(velec,rinvsq10);
807
808             fscal            = felec;
809
810              /* Update vectorial force */
811             fix1             = _mm_macc_ps(dx10,fscal,fix1);
812             fiy1             = _mm_macc_ps(dy10,fscal,fiy1);
813             fiz1             = _mm_macc_ps(dz10,fscal,fiz1);
814
815             fjx0             = _mm_macc_ps(dx10,fscal,fjx0);
816             fjy0             = _mm_macc_ps(dy10,fscal,fjy0);
817             fjz0             = _mm_macc_ps(dz10,fscal,fjz0);
818
819             /**************************
820              * CALCULATE INTERACTIONS *
821              **************************/
822
823             /* Compute parameters for interactions between i and j atoms */
824             qq20             = _mm_mul_ps(iq2,jq0);
825
826             /* COULOMB ELECTROSTATICS */
827             velec            = _mm_mul_ps(qq20,rinv20);
828             felec            = _mm_mul_ps(velec,rinvsq20);
829
830             fscal            = felec;
831
832              /* Update vectorial force */
833             fix2             = _mm_macc_ps(dx20,fscal,fix2);
834             fiy2             = _mm_macc_ps(dy20,fscal,fiy2);
835             fiz2             = _mm_macc_ps(dz20,fscal,fiz2);
836
837             fjx0             = _mm_macc_ps(dx20,fscal,fjx0);
838             fjy0             = _mm_macc_ps(dy20,fscal,fjy0);
839             fjz0             = _mm_macc_ps(dz20,fscal,fjz0);
840
841             /**************************
842              * CALCULATE INTERACTIONS *
843              **************************/
844
845             /* Compute parameters for interactions between i and j atoms */
846             qq30             = _mm_mul_ps(iq3,jq0);
847
848             /* COULOMB ELECTROSTATICS */
849             velec            = _mm_mul_ps(qq30,rinv30);
850             felec            = _mm_mul_ps(velec,rinvsq30);
851
852             fscal            = felec;
853
854              /* Update vectorial force */
855             fix3             = _mm_macc_ps(dx30,fscal,fix3);
856             fiy3             = _mm_macc_ps(dy30,fscal,fiy3);
857             fiz3             = _mm_macc_ps(dz30,fscal,fiz3);
858
859             fjx0             = _mm_macc_ps(dx30,fscal,fjx0);
860             fjy0             = _mm_macc_ps(dy30,fscal,fjy0);
861             fjz0             = _mm_macc_ps(dz30,fscal,fjz0);
862
863             fjptrA             = f+j_coord_offsetA;
864             fjptrB             = f+j_coord_offsetB;
865             fjptrC             = f+j_coord_offsetC;
866             fjptrD             = f+j_coord_offsetD;
867
868             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
869
870             /* Inner loop uses 120 flops */
871         }
872
873         if(jidx<j_index_end)
874         {
875
876             /* Get j neighbor index, and coordinate index */
877             jnrlistA         = jjnr[jidx];
878             jnrlistB         = jjnr[jidx+1];
879             jnrlistC         = jjnr[jidx+2];
880             jnrlistD         = jjnr[jidx+3];
881             /* Sign of each element will be negative for non-real atoms.
882              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
883              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
884              */
885             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
886             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
887             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
888             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
889             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
890             j_coord_offsetA  = DIM*jnrA;
891             j_coord_offsetB  = DIM*jnrB;
892             j_coord_offsetC  = DIM*jnrC;
893             j_coord_offsetD  = DIM*jnrD;
894
895             /* load j atom coordinates */
896             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
897                                               x+j_coord_offsetC,x+j_coord_offsetD,
898                                               &jx0,&jy0,&jz0);
899
900             /* Calculate displacement vector */
901             dx00             = _mm_sub_ps(ix0,jx0);
902             dy00             = _mm_sub_ps(iy0,jy0);
903             dz00             = _mm_sub_ps(iz0,jz0);
904             dx10             = _mm_sub_ps(ix1,jx0);
905             dy10             = _mm_sub_ps(iy1,jy0);
906             dz10             = _mm_sub_ps(iz1,jz0);
907             dx20             = _mm_sub_ps(ix2,jx0);
908             dy20             = _mm_sub_ps(iy2,jy0);
909             dz20             = _mm_sub_ps(iz2,jz0);
910             dx30             = _mm_sub_ps(ix3,jx0);
911             dy30             = _mm_sub_ps(iy3,jy0);
912             dz30             = _mm_sub_ps(iz3,jz0);
913
914             /* Calculate squared distance and things based on it */
915             rsq00            = gmx_mm_calc_rsq_ps(dx00,dy00,dz00);
916             rsq10            = gmx_mm_calc_rsq_ps(dx10,dy10,dz10);
917             rsq20            = gmx_mm_calc_rsq_ps(dx20,dy20,dz20);
918             rsq30            = gmx_mm_calc_rsq_ps(dx30,dy30,dz30);
919
920             rinv10           = avx128fma_invsqrt_f(rsq10);
921             rinv20           = avx128fma_invsqrt_f(rsq20);
922             rinv30           = avx128fma_invsqrt_f(rsq30);
923
924             rinvsq00         = avx128fma_inv_f(rsq00);
925             rinvsq10         = _mm_mul_ps(rinv10,rinv10);
926             rinvsq20         = _mm_mul_ps(rinv20,rinv20);
927             rinvsq30         = _mm_mul_ps(rinv30,rinv30);
928
929             /* Load parameters for j particles */
930             jq0              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
931                                                               charge+jnrC+0,charge+jnrD+0);
932             vdwjidx0A        = 2*vdwtype[jnrA+0];
933             vdwjidx0B        = 2*vdwtype[jnrB+0];
934             vdwjidx0C        = 2*vdwtype[jnrC+0];
935             vdwjidx0D        = 2*vdwtype[jnrD+0];
936
937             fjx0             = _mm_setzero_ps();
938             fjy0             = _mm_setzero_ps();
939             fjz0             = _mm_setzero_ps();
940
941             /**************************
942              * CALCULATE INTERACTIONS *
943              **************************/
944
945             /* Compute parameters for interactions between i and j atoms */
946             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset0+vdwjidx0A,
947                                          vdwparam+vdwioffset0+vdwjidx0B,
948                                          vdwparam+vdwioffset0+vdwjidx0C,
949                                          vdwparam+vdwioffset0+vdwjidx0D,
950                                          &c6_00,&c12_00);
951
952             /* LENNARD-JONES DISPERSION/REPULSION */
953
954             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq00,rinvsq00),rinvsq00);
955             fvdw             = _mm_mul_ps(_mm_msub_ps(c12_00,rinvsix,c6_00),_mm_mul_ps(rinvsix,rinvsq00));
956
957             fscal            = fvdw;
958
959             fscal            = _mm_andnot_ps(dummy_mask,fscal);
960
961              /* Update vectorial force */
962             fix0             = _mm_macc_ps(dx00,fscal,fix0);
963             fiy0             = _mm_macc_ps(dy00,fscal,fiy0);
964             fiz0             = _mm_macc_ps(dz00,fscal,fiz0);
965
966             fjx0             = _mm_macc_ps(dx00,fscal,fjx0);
967             fjy0             = _mm_macc_ps(dy00,fscal,fjy0);
968             fjz0             = _mm_macc_ps(dz00,fscal,fjz0);
969
970             /**************************
971              * CALCULATE INTERACTIONS *
972              **************************/
973
974             /* Compute parameters for interactions between i and j atoms */
975             qq10             = _mm_mul_ps(iq1,jq0);
976
977             /* COULOMB ELECTROSTATICS */
978             velec            = _mm_mul_ps(qq10,rinv10);
979             felec            = _mm_mul_ps(velec,rinvsq10);
980
981             fscal            = felec;
982
983             fscal            = _mm_andnot_ps(dummy_mask,fscal);
984
985              /* Update vectorial force */
986             fix1             = _mm_macc_ps(dx10,fscal,fix1);
987             fiy1             = _mm_macc_ps(dy10,fscal,fiy1);
988             fiz1             = _mm_macc_ps(dz10,fscal,fiz1);
989
990             fjx0             = _mm_macc_ps(dx10,fscal,fjx0);
991             fjy0             = _mm_macc_ps(dy10,fscal,fjy0);
992             fjz0             = _mm_macc_ps(dz10,fscal,fjz0);
993
994             /**************************
995              * CALCULATE INTERACTIONS *
996              **************************/
997
998             /* Compute parameters for interactions between i and j atoms */
999             qq20             = _mm_mul_ps(iq2,jq0);
1000
1001             /* COULOMB ELECTROSTATICS */
1002             velec            = _mm_mul_ps(qq20,rinv20);
1003             felec            = _mm_mul_ps(velec,rinvsq20);
1004
1005             fscal            = felec;
1006
1007             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1008
1009              /* Update vectorial force */
1010             fix2             = _mm_macc_ps(dx20,fscal,fix2);
1011             fiy2             = _mm_macc_ps(dy20,fscal,fiy2);
1012             fiz2             = _mm_macc_ps(dz20,fscal,fiz2);
1013
1014             fjx0             = _mm_macc_ps(dx20,fscal,fjx0);
1015             fjy0             = _mm_macc_ps(dy20,fscal,fjy0);
1016             fjz0             = _mm_macc_ps(dz20,fscal,fjz0);
1017
1018             /**************************
1019              * CALCULATE INTERACTIONS *
1020              **************************/
1021
1022             /* Compute parameters for interactions between i and j atoms */
1023             qq30             = _mm_mul_ps(iq3,jq0);
1024
1025             /* COULOMB ELECTROSTATICS */
1026             velec            = _mm_mul_ps(qq30,rinv30);
1027             felec            = _mm_mul_ps(velec,rinvsq30);
1028
1029             fscal            = felec;
1030
1031             fscal            = _mm_andnot_ps(dummy_mask,fscal);
1032
1033              /* Update vectorial force */
1034             fix3             = _mm_macc_ps(dx30,fscal,fix3);
1035             fiy3             = _mm_macc_ps(dy30,fscal,fiy3);
1036             fiz3             = _mm_macc_ps(dz30,fscal,fiz3);
1037
1038             fjx0             = _mm_macc_ps(dx30,fscal,fjx0);
1039             fjy0             = _mm_macc_ps(dy30,fscal,fjy0);
1040             fjz0             = _mm_macc_ps(dz30,fscal,fjz0);
1041
1042             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1043             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1044             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1045             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1046
1047             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1048
1049             /* Inner loop uses 120 flops */
1050         }
1051
1052         /* End of innermost loop */
1053
1054         gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1055                                               f+i_coord_offset,fshift+i_shift_offset);
1056
1057         /* Increment number of inner iterations */
1058         inneriter                  += j_index_end - j_index_start;
1059
1060         /* Outer loop uses 24 flops */
1061     }
1062
1063     /* Increment number of outer iterations */
1064     outeriter        += nri;
1065
1066     /* Update outer/inner flops */
1067
1068     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4_F,outeriter*24 + inneriter*120);
1069 }