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