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