de209ba8535c8555b100636fcaeb9a102ff20ee4
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_single / nb_kernel_ElecRFCut_VdwLJSh_GeomW4W4_avx_256_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_256_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_256_single.h"
48 #include "kernelutil_x86_avx_256_single.h"
49
50 /*
51  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwLJSh_GeomW4W4_VF_avx_256_single
52  * Electrostatics interaction: ReactionField
53  * VdW interaction:            LennardJones
54  * Geometry:                   Water4-Water4
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecRFCut_VdwLJSh_GeomW4W4_VF_avx_256_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,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight 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              jnrE,jnrF,jnrG,jnrH;
76     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
77     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
78     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
79     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
80     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
81     real             rcutoff_scalar;
82     real             *shiftvec,*fshift,*x,*f;
83     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
84     real             scratch[4*DIM];
85     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
86     real *           vdwioffsetptr0;
87     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
88     real *           vdwioffsetptr1;
89     __m256           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
90     real *           vdwioffsetptr2;
91     __m256           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
92     real *           vdwioffsetptr3;
93     __m256           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
94     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
95     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
96     int              vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D,vdwjidx1E,vdwjidx1F,vdwjidx1G,vdwjidx1H;
97     __m256           jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
98     int              vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D,vdwjidx2E,vdwjidx2F,vdwjidx2G,vdwjidx2H;
99     __m256           jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
100     int              vdwjidx3A,vdwjidx3B,vdwjidx3C,vdwjidx3D,vdwjidx3E,vdwjidx3F,vdwjidx3G,vdwjidx3H;
101     __m256           jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
102     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
103     __m256           dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
104     __m256           dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
105     __m256           dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
106     __m256           dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
107     __m256           dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
108     __m256           dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
109     __m256           dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
110     __m256           dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
111     __m256           dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
112     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
113     real             *charge;
114     int              nvdwtype;
115     __m256           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
116     int              *vdwtype;
117     real             *vdwparam;
118     __m256           one_sixth   = _mm256_set1_ps(1.0/6.0);
119     __m256           one_twelfth = _mm256_set1_ps(1.0/12.0);
120     __m256           dummy_mask,cutoff_mask;
121     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
122     __m256           one     = _mm256_set1_ps(1.0);
123     __m256           two     = _mm256_set1_ps(2.0);
124     x                = xx[0];
125     f                = ff[0];
126
127     nri              = nlist->nri;
128     iinr             = nlist->iinr;
129     jindex           = nlist->jindex;
130     jjnr             = nlist->jjnr;
131     shiftidx         = nlist->shift;
132     gid              = nlist->gid;
133     shiftvec         = fr->shift_vec[0];
134     fshift           = fr->fshift[0];
135     facel            = _mm256_set1_ps(fr->epsfac);
136     charge           = mdatoms->chargeA;
137     krf              = _mm256_set1_ps(fr->ic->k_rf);
138     krf2             = _mm256_set1_ps(fr->ic->k_rf*2.0);
139     crf              = _mm256_set1_ps(fr->ic->c_rf);
140     nvdwtype         = fr->ntype;
141     vdwparam         = fr->nbfp;
142     vdwtype          = mdatoms->typeA;
143
144     /* Setup water-specific parameters */
145     inr              = nlist->iinr[0];
146     iq1              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
147     iq2              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
148     iq3              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+3]));
149     vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
150
151     jq1              = _mm256_set1_ps(charge[inr+1]);
152     jq2              = _mm256_set1_ps(charge[inr+2]);
153     jq3              = _mm256_set1_ps(charge[inr+3]);
154     vdwjidx0A        = 2*vdwtype[inr+0];
155     c6_00            = _mm256_set1_ps(vdwioffsetptr0[vdwjidx0A]);
156     c12_00           = _mm256_set1_ps(vdwioffsetptr0[vdwjidx0A+1]);
157     qq11             = _mm256_mul_ps(iq1,jq1);
158     qq12             = _mm256_mul_ps(iq1,jq2);
159     qq13             = _mm256_mul_ps(iq1,jq3);
160     qq21             = _mm256_mul_ps(iq2,jq1);
161     qq22             = _mm256_mul_ps(iq2,jq2);
162     qq23             = _mm256_mul_ps(iq2,jq3);
163     qq31             = _mm256_mul_ps(iq3,jq1);
164     qq32             = _mm256_mul_ps(iq3,jq2);
165     qq33             = _mm256_mul_ps(iq3,jq3);
166
167     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
168     rcutoff_scalar   = fr->rcoulomb;
169     rcutoff          = _mm256_set1_ps(rcutoff_scalar);
170     rcutoff2         = _mm256_mul_ps(rcutoff,rcutoff);
171
172     sh_vdw_invrcut6  = _mm256_set1_ps(fr->ic->sh_invrc6);
173     rvdw             = _mm256_set1_ps(fr->rvdw);
174
175     /* Avoid stupid compiler warnings */
176     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
177     j_coord_offsetA = 0;
178     j_coord_offsetB = 0;
179     j_coord_offsetC = 0;
180     j_coord_offsetD = 0;
181     j_coord_offsetE = 0;
182     j_coord_offsetF = 0;
183     j_coord_offsetG = 0;
184     j_coord_offsetH = 0;
185
186     outeriter        = 0;
187     inneriter        = 0;
188
189     for(iidx=0;iidx<4*DIM;iidx++)
190     {
191         scratch[iidx] = 0.0;
192     }
193
194     /* Start outer loop over neighborlists */
195     for(iidx=0; iidx<nri; iidx++)
196     {
197         /* Load shift vector for this list */
198         i_shift_offset   = DIM*shiftidx[iidx];
199
200         /* Load limits for loop over neighbors */
201         j_index_start    = jindex[iidx];
202         j_index_end      = jindex[iidx+1];
203
204         /* Get outer coordinate index */
205         inr              = iinr[iidx];
206         i_coord_offset   = DIM*inr;
207
208         /* Load i particle coords and add shift vector */
209         gmx_mm256_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
210                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
211
212         fix0             = _mm256_setzero_ps();
213         fiy0             = _mm256_setzero_ps();
214         fiz0             = _mm256_setzero_ps();
215         fix1             = _mm256_setzero_ps();
216         fiy1             = _mm256_setzero_ps();
217         fiz1             = _mm256_setzero_ps();
218         fix2             = _mm256_setzero_ps();
219         fiy2             = _mm256_setzero_ps();
220         fiz2             = _mm256_setzero_ps();
221         fix3             = _mm256_setzero_ps();
222         fiy3             = _mm256_setzero_ps();
223         fiz3             = _mm256_setzero_ps();
224
225         /* Reset potential sums */
226         velecsum         = _mm256_setzero_ps();
227         vvdwsum          = _mm256_setzero_ps();
228
229         /* Start inner kernel loop */
230         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
231         {
232
233             /* Get j neighbor index, and coordinate index */
234             jnrA             = jjnr[jidx];
235             jnrB             = jjnr[jidx+1];
236             jnrC             = jjnr[jidx+2];
237             jnrD             = jjnr[jidx+3];
238             jnrE             = jjnr[jidx+4];
239             jnrF             = jjnr[jidx+5];
240             jnrG             = jjnr[jidx+6];
241             jnrH             = jjnr[jidx+7];
242             j_coord_offsetA  = DIM*jnrA;
243             j_coord_offsetB  = DIM*jnrB;
244             j_coord_offsetC  = DIM*jnrC;
245             j_coord_offsetD  = DIM*jnrD;
246             j_coord_offsetE  = DIM*jnrE;
247             j_coord_offsetF  = DIM*jnrF;
248             j_coord_offsetG  = DIM*jnrG;
249             j_coord_offsetH  = DIM*jnrH;
250
251             /* load j atom coordinates */
252             gmx_mm256_load_4rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
253                                                  x+j_coord_offsetC,x+j_coord_offsetD,
254                                                  x+j_coord_offsetE,x+j_coord_offsetF,
255                                                  x+j_coord_offsetG,x+j_coord_offsetH,
256                                                  &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
257                                                  &jy2,&jz2,&jx3,&jy3,&jz3);
258
259             /* Calculate displacement vector */
260             dx00             = _mm256_sub_ps(ix0,jx0);
261             dy00             = _mm256_sub_ps(iy0,jy0);
262             dz00             = _mm256_sub_ps(iz0,jz0);
263             dx11             = _mm256_sub_ps(ix1,jx1);
264             dy11             = _mm256_sub_ps(iy1,jy1);
265             dz11             = _mm256_sub_ps(iz1,jz1);
266             dx12             = _mm256_sub_ps(ix1,jx2);
267             dy12             = _mm256_sub_ps(iy1,jy2);
268             dz12             = _mm256_sub_ps(iz1,jz2);
269             dx13             = _mm256_sub_ps(ix1,jx3);
270             dy13             = _mm256_sub_ps(iy1,jy3);
271             dz13             = _mm256_sub_ps(iz1,jz3);
272             dx21             = _mm256_sub_ps(ix2,jx1);
273             dy21             = _mm256_sub_ps(iy2,jy1);
274             dz21             = _mm256_sub_ps(iz2,jz1);
275             dx22             = _mm256_sub_ps(ix2,jx2);
276             dy22             = _mm256_sub_ps(iy2,jy2);
277             dz22             = _mm256_sub_ps(iz2,jz2);
278             dx23             = _mm256_sub_ps(ix2,jx3);
279             dy23             = _mm256_sub_ps(iy2,jy3);
280             dz23             = _mm256_sub_ps(iz2,jz3);
281             dx31             = _mm256_sub_ps(ix3,jx1);
282             dy31             = _mm256_sub_ps(iy3,jy1);
283             dz31             = _mm256_sub_ps(iz3,jz1);
284             dx32             = _mm256_sub_ps(ix3,jx2);
285             dy32             = _mm256_sub_ps(iy3,jy2);
286             dz32             = _mm256_sub_ps(iz3,jz2);
287             dx33             = _mm256_sub_ps(ix3,jx3);
288             dy33             = _mm256_sub_ps(iy3,jy3);
289             dz33             = _mm256_sub_ps(iz3,jz3);
290
291             /* Calculate squared distance and things based on it */
292             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
293             rsq11            = gmx_mm256_calc_rsq_ps(dx11,dy11,dz11);
294             rsq12            = gmx_mm256_calc_rsq_ps(dx12,dy12,dz12);
295             rsq13            = gmx_mm256_calc_rsq_ps(dx13,dy13,dz13);
296             rsq21            = gmx_mm256_calc_rsq_ps(dx21,dy21,dz21);
297             rsq22            = gmx_mm256_calc_rsq_ps(dx22,dy22,dz22);
298             rsq23            = gmx_mm256_calc_rsq_ps(dx23,dy23,dz23);
299             rsq31            = gmx_mm256_calc_rsq_ps(dx31,dy31,dz31);
300             rsq32            = gmx_mm256_calc_rsq_ps(dx32,dy32,dz32);
301             rsq33            = gmx_mm256_calc_rsq_ps(dx33,dy33,dz33);
302
303             rinv11           = gmx_mm256_invsqrt_ps(rsq11);
304             rinv12           = gmx_mm256_invsqrt_ps(rsq12);
305             rinv13           = gmx_mm256_invsqrt_ps(rsq13);
306             rinv21           = gmx_mm256_invsqrt_ps(rsq21);
307             rinv22           = gmx_mm256_invsqrt_ps(rsq22);
308             rinv23           = gmx_mm256_invsqrt_ps(rsq23);
309             rinv31           = gmx_mm256_invsqrt_ps(rsq31);
310             rinv32           = gmx_mm256_invsqrt_ps(rsq32);
311             rinv33           = gmx_mm256_invsqrt_ps(rsq33);
312
313             rinvsq00         = gmx_mm256_inv_ps(rsq00);
314             rinvsq11         = _mm256_mul_ps(rinv11,rinv11);
315             rinvsq12         = _mm256_mul_ps(rinv12,rinv12);
316             rinvsq13         = _mm256_mul_ps(rinv13,rinv13);
317             rinvsq21         = _mm256_mul_ps(rinv21,rinv21);
318             rinvsq22         = _mm256_mul_ps(rinv22,rinv22);
319             rinvsq23         = _mm256_mul_ps(rinv23,rinv23);
320             rinvsq31         = _mm256_mul_ps(rinv31,rinv31);
321             rinvsq32         = _mm256_mul_ps(rinv32,rinv32);
322             rinvsq33         = _mm256_mul_ps(rinv33,rinv33);
323
324             fjx0             = _mm256_setzero_ps();
325             fjy0             = _mm256_setzero_ps();
326             fjz0             = _mm256_setzero_ps();
327             fjx1             = _mm256_setzero_ps();
328             fjy1             = _mm256_setzero_ps();
329             fjz1             = _mm256_setzero_ps();
330             fjx2             = _mm256_setzero_ps();
331             fjy2             = _mm256_setzero_ps();
332             fjz2             = _mm256_setzero_ps();
333             fjx3             = _mm256_setzero_ps();
334             fjy3             = _mm256_setzero_ps();
335             fjz3             = _mm256_setzero_ps();
336
337             /**************************
338              * CALCULATE INTERACTIONS *
339              **************************/
340
341             if (gmx_mm256_any_lt(rsq00,rcutoff2))
342             {
343
344             /* LENNARD-JONES DISPERSION/REPULSION */
345
346             rinvsix          = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
347             vvdw6            = _mm256_mul_ps(c6_00,rinvsix);
348             vvdw12           = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
349             vvdw             = _mm256_sub_ps(_mm256_mul_ps( _mm256_sub_ps(vvdw12 , _mm256_mul_ps(c12_00,_mm256_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
350                                           _mm256_mul_ps( _mm256_sub_ps(vvdw6,_mm256_mul_ps(c6_00,sh_vdw_invrcut6)),one_sixth));
351             fvdw             = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
352
353             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
354
355             /* Update potential sum for this i atom from the interaction with this j atom. */
356             vvdw             = _mm256_and_ps(vvdw,cutoff_mask);
357             vvdwsum          = _mm256_add_ps(vvdwsum,vvdw);
358
359             fscal            = fvdw;
360
361             fscal            = _mm256_and_ps(fscal,cutoff_mask);
362
363             /* Calculate temporary vectorial force */
364             tx               = _mm256_mul_ps(fscal,dx00);
365             ty               = _mm256_mul_ps(fscal,dy00);
366             tz               = _mm256_mul_ps(fscal,dz00);
367
368             /* Update vectorial force */
369             fix0             = _mm256_add_ps(fix0,tx);
370             fiy0             = _mm256_add_ps(fiy0,ty);
371             fiz0             = _mm256_add_ps(fiz0,tz);
372
373             fjx0             = _mm256_add_ps(fjx0,tx);
374             fjy0             = _mm256_add_ps(fjy0,ty);
375             fjz0             = _mm256_add_ps(fjz0,tz);
376
377             }
378
379             /**************************
380              * CALCULATE INTERACTIONS *
381              **************************/
382
383             if (gmx_mm256_any_lt(rsq11,rcutoff2))
384             {
385
386             /* REACTION-FIELD ELECTROSTATICS */
387             velec            = _mm256_mul_ps(qq11,_mm256_sub_ps(_mm256_add_ps(rinv11,_mm256_mul_ps(krf,rsq11)),crf));
388             felec            = _mm256_mul_ps(qq11,_mm256_sub_ps(_mm256_mul_ps(rinv11,rinvsq11),krf2));
389
390             cutoff_mask      = _mm256_cmp_ps(rsq11,rcutoff2,_CMP_LT_OQ);
391
392             /* Update potential sum for this i atom from the interaction with this j atom. */
393             velec            = _mm256_and_ps(velec,cutoff_mask);
394             velecsum         = _mm256_add_ps(velecsum,velec);
395
396             fscal            = felec;
397
398             fscal            = _mm256_and_ps(fscal,cutoff_mask);
399
400             /* Calculate temporary vectorial force */
401             tx               = _mm256_mul_ps(fscal,dx11);
402             ty               = _mm256_mul_ps(fscal,dy11);
403             tz               = _mm256_mul_ps(fscal,dz11);
404
405             /* Update vectorial force */
406             fix1             = _mm256_add_ps(fix1,tx);
407             fiy1             = _mm256_add_ps(fiy1,ty);
408             fiz1             = _mm256_add_ps(fiz1,tz);
409
410             fjx1             = _mm256_add_ps(fjx1,tx);
411             fjy1             = _mm256_add_ps(fjy1,ty);
412             fjz1             = _mm256_add_ps(fjz1,tz);
413
414             }
415
416             /**************************
417              * CALCULATE INTERACTIONS *
418              **************************/
419
420             if (gmx_mm256_any_lt(rsq12,rcutoff2))
421             {
422
423             /* REACTION-FIELD ELECTROSTATICS */
424             velec            = _mm256_mul_ps(qq12,_mm256_sub_ps(_mm256_add_ps(rinv12,_mm256_mul_ps(krf,rsq12)),crf));
425             felec            = _mm256_mul_ps(qq12,_mm256_sub_ps(_mm256_mul_ps(rinv12,rinvsq12),krf2));
426
427             cutoff_mask      = _mm256_cmp_ps(rsq12,rcutoff2,_CMP_LT_OQ);
428
429             /* Update potential sum for this i atom from the interaction with this j atom. */
430             velec            = _mm256_and_ps(velec,cutoff_mask);
431             velecsum         = _mm256_add_ps(velecsum,velec);
432
433             fscal            = felec;
434
435             fscal            = _mm256_and_ps(fscal,cutoff_mask);
436
437             /* Calculate temporary vectorial force */
438             tx               = _mm256_mul_ps(fscal,dx12);
439             ty               = _mm256_mul_ps(fscal,dy12);
440             tz               = _mm256_mul_ps(fscal,dz12);
441
442             /* Update vectorial force */
443             fix1             = _mm256_add_ps(fix1,tx);
444             fiy1             = _mm256_add_ps(fiy1,ty);
445             fiz1             = _mm256_add_ps(fiz1,tz);
446
447             fjx2             = _mm256_add_ps(fjx2,tx);
448             fjy2             = _mm256_add_ps(fjy2,ty);
449             fjz2             = _mm256_add_ps(fjz2,tz);
450
451             }
452
453             /**************************
454              * CALCULATE INTERACTIONS *
455              **************************/
456
457             if (gmx_mm256_any_lt(rsq13,rcutoff2))
458             {
459
460             /* REACTION-FIELD ELECTROSTATICS */
461             velec            = _mm256_mul_ps(qq13,_mm256_sub_ps(_mm256_add_ps(rinv13,_mm256_mul_ps(krf,rsq13)),crf));
462             felec            = _mm256_mul_ps(qq13,_mm256_sub_ps(_mm256_mul_ps(rinv13,rinvsq13),krf2));
463
464             cutoff_mask      = _mm256_cmp_ps(rsq13,rcutoff2,_CMP_LT_OQ);
465
466             /* Update potential sum for this i atom from the interaction with this j atom. */
467             velec            = _mm256_and_ps(velec,cutoff_mask);
468             velecsum         = _mm256_add_ps(velecsum,velec);
469
470             fscal            = felec;
471
472             fscal            = _mm256_and_ps(fscal,cutoff_mask);
473
474             /* Calculate temporary vectorial force */
475             tx               = _mm256_mul_ps(fscal,dx13);
476             ty               = _mm256_mul_ps(fscal,dy13);
477             tz               = _mm256_mul_ps(fscal,dz13);
478
479             /* Update vectorial force */
480             fix1             = _mm256_add_ps(fix1,tx);
481             fiy1             = _mm256_add_ps(fiy1,ty);
482             fiz1             = _mm256_add_ps(fiz1,tz);
483
484             fjx3             = _mm256_add_ps(fjx3,tx);
485             fjy3             = _mm256_add_ps(fjy3,ty);
486             fjz3             = _mm256_add_ps(fjz3,tz);
487
488             }
489
490             /**************************
491              * CALCULATE INTERACTIONS *
492              **************************/
493
494             if (gmx_mm256_any_lt(rsq21,rcutoff2))
495             {
496
497             /* REACTION-FIELD ELECTROSTATICS */
498             velec            = _mm256_mul_ps(qq21,_mm256_sub_ps(_mm256_add_ps(rinv21,_mm256_mul_ps(krf,rsq21)),crf));
499             felec            = _mm256_mul_ps(qq21,_mm256_sub_ps(_mm256_mul_ps(rinv21,rinvsq21),krf2));
500
501             cutoff_mask      = _mm256_cmp_ps(rsq21,rcutoff2,_CMP_LT_OQ);
502
503             /* Update potential sum for this i atom from the interaction with this j atom. */
504             velec            = _mm256_and_ps(velec,cutoff_mask);
505             velecsum         = _mm256_add_ps(velecsum,velec);
506
507             fscal            = felec;
508
509             fscal            = _mm256_and_ps(fscal,cutoff_mask);
510
511             /* Calculate temporary vectorial force */
512             tx               = _mm256_mul_ps(fscal,dx21);
513             ty               = _mm256_mul_ps(fscal,dy21);
514             tz               = _mm256_mul_ps(fscal,dz21);
515
516             /* Update vectorial force */
517             fix2             = _mm256_add_ps(fix2,tx);
518             fiy2             = _mm256_add_ps(fiy2,ty);
519             fiz2             = _mm256_add_ps(fiz2,tz);
520
521             fjx1             = _mm256_add_ps(fjx1,tx);
522             fjy1             = _mm256_add_ps(fjy1,ty);
523             fjz1             = _mm256_add_ps(fjz1,tz);
524
525             }
526
527             /**************************
528              * CALCULATE INTERACTIONS *
529              **************************/
530
531             if (gmx_mm256_any_lt(rsq22,rcutoff2))
532             {
533
534             /* REACTION-FIELD ELECTROSTATICS */
535             velec            = _mm256_mul_ps(qq22,_mm256_sub_ps(_mm256_add_ps(rinv22,_mm256_mul_ps(krf,rsq22)),crf));
536             felec            = _mm256_mul_ps(qq22,_mm256_sub_ps(_mm256_mul_ps(rinv22,rinvsq22),krf2));
537
538             cutoff_mask      = _mm256_cmp_ps(rsq22,rcutoff2,_CMP_LT_OQ);
539
540             /* Update potential sum for this i atom from the interaction with this j atom. */
541             velec            = _mm256_and_ps(velec,cutoff_mask);
542             velecsum         = _mm256_add_ps(velecsum,velec);
543
544             fscal            = felec;
545
546             fscal            = _mm256_and_ps(fscal,cutoff_mask);
547
548             /* Calculate temporary vectorial force */
549             tx               = _mm256_mul_ps(fscal,dx22);
550             ty               = _mm256_mul_ps(fscal,dy22);
551             tz               = _mm256_mul_ps(fscal,dz22);
552
553             /* Update vectorial force */
554             fix2             = _mm256_add_ps(fix2,tx);
555             fiy2             = _mm256_add_ps(fiy2,ty);
556             fiz2             = _mm256_add_ps(fiz2,tz);
557
558             fjx2             = _mm256_add_ps(fjx2,tx);
559             fjy2             = _mm256_add_ps(fjy2,ty);
560             fjz2             = _mm256_add_ps(fjz2,tz);
561
562             }
563
564             /**************************
565              * CALCULATE INTERACTIONS *
566              **************************/
567
568             if (gmx_mm256_any_lt(rsq23,rcutoff2))
569             {
570
571             /* REACTION-FIELD ELECTROSTATICS */
572             velec            = _mm256_mul_ps(qq23,_mm256_sub_ps(_mm256_add_ps(rinv23,_mm256_mul_ps(krf,rsq23)),crf));
573             felec            = _mm256_mul_ps(qq23,_mm256_sub_ps(_mm256_mul_ps(rinv23,rinvsq23),krf2));
574
575             cutoff_mask      = _mm256_cmp_ps(rsq23,rcutoff2,_CMP_LT_OQ);
576
577             /* Update potential sum for this i atom from the interaction with this j atom. */
578             velec            = _mm256_and_ps(velec,cutoff_mask);
579             velecsum         = _mm256_add_ps(velecsum,velec);
580
581             fscal            = felec;
582
583             fscal            = _mm256_and_ps(fscal,cutoff_mask);
584
585             /* Calculate temporary vectorial force */
586             tx               = _mm256_mul_ps(fscal,dx23);
587             ty               = _mm256_mul_ps(fscal,dy23);
588             tz               = _mm256_mul_ps(fscal,dz23);
589
590             /* Update vectorial force */
591             fix2             = _mm256_add_ps(fix2,tx);
592             fiy2             = _mm256_add_ps(fiy2,ty);
593             fiz2             = _mm256_add_ps(fiz2,tz);
594
595             fjx3             = _mm256_add_ps(fjx3,tx);
596             fjy3             = _mm256_add_ps(fjy3,ty);
597             fjz3             = _mm256_add_ps(fjz3,tz);
598
599             }
600
601             /**************************
602              * CALCULATE INTERACTIONS *
603              **************************/
604
605             if (gmx_mm256_any_lt(rsq31,rcutoff2))
606             {
607
608             /* REACTION-FIELD ELECTROSTATICS */
609             velec            = _mm256_mul_ps(qq31,_mm256_sub_ps(_mm256_add_ps(rinv31,_mm256_mul_ps(krf,rsq31)),crf));
610             felec            = _mm256_mul_ps(qq31,_mm256_sub_ps(_mm256_mul_ps(rinv31,rinvsq31),krf2));
611
612             cutoff_mask      = _mm256_cmp_ps(rsq31,rcutoff2,_CMP_LT_OQ);
613
614             /* Update potential sum for this i atom from the interaction with this j atom. */
615             velec            = _mm256_and_ps(velec,cutoff_mask);
616             velecsum         = _mm256_add_ps(velecsum,velec);
617
618             fscal            = felec;
619
620             fscal            = _mm256_and_ps(fscal,cutoff_mask);
621
622             /* Calculate temporary vectorial force */
623             tx               = _mm256_mul_ps(fscal,dx31);
624             ty               = _mm256_mul_ps(fscal,dy31);
625             tz               = _mm256_mul_ps(fscal,dz31);
626
627             /* Update vectorial force */
628             fix3             = _mm256_add_ps(fix3,tx);
629             fiy3             = _mm256_add_ps(fiy3,ty);
630             fiz3             = _mm256_add_ps(fiz3,tz);
631
632             fjx1             = _mm256_add_ps(fjx1,tx);
633             fjy1             = _mm256_add_ps(fjy1,ty);
634             fjz1             = _mm256_add_ps(fjz1,tz);
635
636             }
637
638             /**************************
639              * CALCULATE INTERACTIONS *
640              **************************/
641
642             if (gmx_mm256_any_lt(rsq32,rcutoff2))
643             {
644
645             /* REACTION-FIELD ELECTROSTATICS */
646             velec            = _mm256_mul_ps(qq32,_mm256_sub_ps(_mm256_add_ps(rinv32,_mm256_mul_ps(krf,rsq32)),crf));
647             felec            = _mm256_mul_ps(qq32,_mm256_sub_ps(_mm256_mul_ps(rinv32,rinvsq32),krf2));
648
649             cutoff_mask      = _mm256_cmp_ps(rsq32,rcutoff2,_CMP_LT_OQ);
650
651             /* Update potential sum for this i atom from the interaction with this j atom. */
652             velec            = _mm256_and_ps(velec,cutoff_mask);
653             velecsum         = _mm256_add_ps(velecsum,velec);
654
655             fscal            = felec;
656
657             fscal            = _mm256_and_ps(fscal,cutoff_mask);
658
659             /* Calculate temporary vectorial force */
660             tx               = _mm256_mul_ps(fscal,dx32);
661             ty               = _mm256_mul_ps(fscal,dy32);
662             tz               = _mm256_mul_ps(fscal,dz32);
663
664             /* Update vectorial force */
665             fix3             = _mm256_add_ps(fix3,tx);
666             fiy3             = _mm256_add_ps(fiy3,ty);
667             fiz3             = _mm256_add_ps(fiz3,tz);
668
669             fjx2             = _mm256_add_ps(fjx2,tx);
670             fjy2             = _mm256_add_ps(fjy2,ty);
671             fjz2             = _mm256_add_ps(fjz2,tz);
672
673             }
674
675             /**************************
676              * CALCULATE INTERACTIONS *
677              **************************/
678
679             if (gmx_mm256_any_lt(rsq33,rcutoff2))
680             {
681
682             /* REACTION-FIELD ELECTROSTATICS */
683             velec            = _mm256_mul_ps(qq33,_mm256_sub_ps(_mm256_add_ps(rinv33,_mm256_mul_ps(krf,rsq33)),crf));
684             felec            = _mm256_mul_ps(qq33,_mm256_sub_ps(_mm256_mul_ps(rinv33,rinvsq33),krf2));
685
686             cutoff_mask      = _mm256_cmp_ps(rsq33,rcutoff2,_CMP_LT_OQ);
687
688             /* Update potential sum for this i atom from the interaction with this j atom. */
689             velec            = _mm256_and_ps(velec,cutoff_mask);
690             velecsum         = _mm256_add_ps(velecsum,velec);
691
692             fscal            = felec;
693
694             fscal            = _mm256_and_ps(fscal,cutoff_mask);
695
696             /* Calculate temporary vectorial force */
697             tx               = _mm256_mul_ps(fscal,dx33);
698             ty               = _mm256_mul_ps(fscal,dy33);
699             tz               = _mm256_mul_ps(fscal,dz33);
700
701             /* Update vectorial force */
702             fix3             = _mm256_add_ps(fix3,tx);
703             fiy3             = _mm256_add_ps(fiy3,ty);
704             fiz3             = _mm256_add_ps(fiz3,tz);
705
706             fjx3             = _mm256_add_ps(fjx3,tx);
707             fjy3             = _mm256_add_ps(fjy3,ty);
708             fjz3             = _mm256_add_ps(fjz3,tz);
709
710             }
711
712             fjptrA             = f+j_coord_offsetA;
713             fjptrB             = f+j_coord_offsetB;
714             fjptrC             = f+j_coord_offsetC;
715             fjptrD             = f+j_coord_offsetD;
716             fjptrE             = f+j_coord_offsetE;
717             fjptrF             = f+j_coord_offsetF;
718             fjptrG             = f+j_coord_offsetG;
719             fjptrH             = f+j_coord_offsetH;
720
721             gmx_mm256_decrement_4rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
722                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
723                                                       fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
724
725             /* Inner loop uses 368 flops */
726         }
727
728         if(jidx<j_index_end)
729         {
730
731             /* Get j neighbor index, and coordinate index */
732             jnrlistA         = jjnr[jidx];
733             jnrlistB         = jjnr[jidx+1];
734             jnrlistC         = jjnr[jidx+2];
735             jnrlistD         = jjnr[jidx+3];
736             jnrlistE         = jjnr[jidx+4];
737             jnrlistF         = jjnr[jidx+5];
738             jnrlistG         = jjnr[jidx+6];
739             jnrlistH         = jjnr[jidx+7];
740             /* Sign of each element will be negative for non-real atoms.
741              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
742              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
743              */
744             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
745                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
746                                             
747             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
748             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
749             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
750             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
751             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
752             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
753             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
754             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
755             j_coord_offsetA  = DIM*jnrA;
756             j_coord_offsetB  = DIM*jnrB;
757             j_coord_offsetC  = DIM*jnrC;
758             j_coord_offsetD  = DIM*jnrD;
759             j_coord_offsetE  = DIM*jnrE;
760             j_coord_offsetF  = DIM*jnrF;
761             j_coord_offsetG  = DIM*jnrG;
762             j_coord_offsetH  = DIM*jnrH;
763
764             /* load j atom coordinates */
765             gmx_mm256_load_4rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
766                                                  x+j_coord_offsetC,x+j_coord_offsetD,
767                                                  x+j_coord_offsetE,x+j_coord_offsetF,
768                                                  x+j_coord_offsetG,x+j_coord_offsetH,
769                                                  &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
770                                                  &jy2,&jz2,&jx3,&jy3,&jz3);
771
772             /* Calculate displacement vector */
773             dx00             = _mm256_sub_ps(ix0,jx0);
774             dy00             = _mm256_sub_ps(iy0,jy0);
775             dz00             = _mm256_sub_ps(iz0,jz0);
776             dx11             = _mm256_sub_ps(ix1,jx1);
777             dy11             = _mm256_sub_ps(iy1,jy1);
778             dz11             = _mm256_sub_ps(iz1,jz1);
779             dx12             = _mm256_sub_ps(ix1,jx2);
780             dy12             = _mm256_sub_ps(iy1,jy2);
781             dz12             = _mm256_sub_ps(iz1,jz2);
782             dx13             = _mm256_sub_ps(ix1,jx3);
783             dy13             = _mm256_sub_ps(iy1,jy3);
784             dz13             = _mm256_sub_ps(iz1,jz3);
785             dx21             = _mm256_sub_ps(ix2,jx1);
786             dy21             = _mm256_sub_ps(iy2,jy1);
787             dz21             = _mm256_sub_ps(iz2,jz1);
788             dx22             = _mm256_sub_ps(ix2,jx2);
789             dy22             = _mm256_sub_ps(iy2,jy2);
790             dz22             = _mm256_sub_ps(iz2,jz2);
791             dx23             = _mm256_sub_ps(ix2,jx3);
792             dy23             = _mm256_sub_ps(iy2,jy3);
793             dz23             = _mm256_sub_ps(iz2,jz3);
794             dx31             = _mm256_sub_ps(ix3,jx1);
795             dy31             = _mm256_sub_ps(iy3,jy1);
796             dz31             = _mm256_sub_ps(iz3,jz1);
797             dx32             = _mm256_sub_ps(ix3,jx2);
798             dy32             = _mm256_sub_ps(iy3,jy2);
799             dz32             = _mm256_sub_ps(iz3,jz2);
800             dx33             = _mm256_sub_ps(ix3,jx3);
801             dy33             = _mm256_sub_ps(iy3,jy3);
802             dz33             = _mm256_sub_ps(iz3,jz3);
803
804             /* Calculate squared distance and things based on it */
805             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
806             rsq11            = gmx_mm256_calc_rsq_ps(dx11,dy11,dz11);
807             rsq12            = gmx_mm256_calc_rsq_ps(dx12,dy12,dz12);
808             rsq13            = gmx_mm256_calc_rsq_ps(dx13,dy13,dz13);
809             rsq21            = gmx_mm256_calc_rsq_ps(dx21,dy21,dz21);
810             rsq22            = gmx_mm256_calc_rsq_ps(dx22,dy22,dz22);
811             rsq23            = gmx_mm256_calc_rsq_ps(dx23,dy23,dz23);
812             rsq31            = gmx_mm256_calc_rsq_ps(dx31,dy31,dz31);
813             rsq32            = gmx_mm256_calc_rsq_ps(dx32,dy32,dz32);
814             rsq33            = gmx_mm256_calc_rsq_ps(dx33,dy33,dz33);
815
816             rinv11           = gmx_mm256_invsqrt_ps(rsq11);
817             rinv12           = gmx_mm256_invsqrt_ps(rsq12);
818             rinv13           = gmx_mm256_invsqrt_ps(rsq13);
819             rinv21           = gmx_mm256_invsqrt_ps(rsq21);
820             rinv22           = gmx_mm256_invsqrt_ps(rsq22);
821             rinv23           = gmx_mm256_invsqrt_ps(rsq23);
822             rinv31           = gmx_mm256_invsqrt_ps(rsq31);
823             rinv32           = gmx_mm256_invsqrt_ps(rsq32);
824             rinv33           = gmx_mm256_invsqrt_ps(rsq33);
825
826             rinvsq00         = gmx_mm256_inv_ps(rsq00);
827             rinvsq11         = _mm256_mul_ps(rinv11,rinv11);
828             rinvsq12         = _mm256_mul_ps(rinv12,rinv12);
829             rinvsq13         = _mm256_mul_ps(rinv13,rinv13);
830             rinvsq21         = _mm256_mul_ps(rinv21,rinv21);
831             rinvsq22         = _mm256_mul_ps(rinv22,rinv22);
832             rinvsq23         = _mm256_mul_ps(rinv23,rinv23);
833             rinvsq31         = _mm256_mul_ps(rinv31,rinv31);
834             rinvsq32         = _mm256_mul_ps(rinv32,rinv32);
835             rinvsq33         = _mm256_mul_ps(rinv33,rinv33);
836
837             fjx0             = _mm256_setzero_ps();
838             fjy0             = _mm256_setzero_ps();
839             fjz0             = _mm256_setzero_ps();
840             fjx1             = _mm256_setzero_ps();
841             fjy1             = _mm256_setzero_ps();
842             fjz1             = _mm256_setzero_ps();
843             fjx2             = _mm256_setzero_ps();
844             fjy2             = _mm256_setzero_ps();
845             fjz2             = _mm256_setzero_ps();
846             fjx3             = _mm256_setzero_ps();
847             fjy3             = _mm256_setzero_ps();
848             fjz3             = _mm256_setzero_ps();
849
850             /**************************
851              * CALCULATE INTERACTIONS *
852              **************************/
853
854             if (gmx_mm256_any_lt(rsq00,rcutoff2))
855             {
856
857             /* LENNARD-JONES DISPERSION/REPULSION */
858
859             rinvsix          = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
860             vvdw6            = _mm256_mul_ps(c6_00,rinvsix);
861             vvdw12           = _mm256_mul_ps(c12_00,_mm256_mul_ps(rinvsix,rinvsix));
862             vvdw             = _mm256_sub_ps(_mm256_mul_ps( _mm256_sub_ps(vvdw12 , _mm256_mul_ps(c12_00,_mm256_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
863                                           _mm256_mul_ps( _mm256_sub_ps(vvdw6,_mm256_mul_ps(c6_00,sh_vdw_invrcut6)),one_sixth));
864             fvdw             = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq00);
865
866             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
867
868             /* Update potential sum for this i atom from the interaction with this j atom. */
869             vvdw             = _mm256_and_ps(vvdw,cutoff_mask);
870             vvdw             = _mm256_andnot_ps(dummy_mask,vvdw);
871             vvdwsum          = _mm256_add_ps(vvdwsum,vvdw);
872
873             fscal            = fvdw;
874
875             fscal            = _mm256_and_ps(fscal,cutoff_mask);
876
877             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
878
879             /* Calculate temporary vectorial force */
880             tx               = _mm256_mul_ps(fscal,dx00);
881             ty               = _mm256_mul_ps(fscal,dy00);
882             tz               = _mm256_mul_ps(fscal,dz00);
883
884             /* Update vectorial force */
885             fix0             = _mm256_add_ps(fix0,tx);
886             fiy0             = _mm256_add_ps(fiy0,ty);
887             fiz0             = _mm256_add_ps(fiz0,tz);
888
889             fjx0             = _mm256_add_ps(fjx0,tx);
890             fjy0             = _mm256_add_ps(fjy0,ty);
891             fjz0             = _mm256_add_ps(fjz0,tz);
892
893             }
894
895             /**************************
896              * CALCULATE INTERACTIONS *
897              **************************/
898
899             if (gmx_mm256_any_lt(rsq11,rcutoff2))
900             {
901
902             /* REACTION-FIELD ELECTROSTATICS */
903             velec            = _mm256_mul_ps(qq11,_mm256_sub_ps(_mm256_add_ps(rinv11,_mm256_mul_ps(krf,rsq11)),crf));
904             felec            = _mm256_mul_ps(qq11,_mm256_sub_ps(_mm256_mul_ps(rinv11,rinvsq11),krf2));
905
906             cutoff_mask      = _mm256_cmp_ps(rsq11,rcutoff2,_CMP_LT_OQ);
907
908             /* Update potential sum for this i atom from the interaction with this j atom. */
909             velec            = _mm256_and_ps(velec,cutoff_mask);
910             velec            = _mm256_andnot_ps(dummy_mask,velec);
911             velecsum         = _mm256_add_ps(velecsum,velec);
912
913             fscal            = felec;
914
915             fscal            = _mm256_and_ps(fscal,cutoff_mask);
916
917             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
918
919             /* Calculate temporary vectorial force */
920             tx               = _mm256_mul_ps(fscal,dx11);
921             ty               = _mm256_mul_ps(fscal,dy11);
922             tz               = _mm256_mul_ps(fscal,dz11);
923
924             /* Update vectorial force */
925             fix1             = _mm256_add_ps(fix1,tx);
926             fiy1             = _mm256_add_ps(fiy1,ty);
927             fiz1             = _mm256_add_ps(fiz1,tz);
928
929             fjx1             = _mm256_add_ps(fjx1,tx);
930             fjy1             = _mm256_add_ps(fjy1,ty);
931             fjz1             = _mm256_add_ps(fjz1,tz);
932
933             }
934
935             /**************************
936              * CALCULATE INTERACTIONS *
937              **************************/
938
939             if (gmx_mm256_any_lt(rsq12,rcutoff2))
940             {
941
942             /* REACTION-FIELD ELECTROSTATICS */
943             velec            = _mm256_mul_ps(qq12,_mm256_sub_ps(_mm256_add_ps(rinv12,_mm256_mul_ps(krf,rsq12)),crf));
944             felec            = _mm256_mul_ps(qq12,_mm256_sub_ps(_mm256_mul_ps(rinv12,rinvsq12),krf2));
945
946             cutoff_mask      = _mm256_cmp_ps(rsq12,rcutoff2,_CMP_LT_OQ);
947
948             /* Update potential sum for this i atom from the interaction with this j atom. */
949             velec            = _mm256_and_ps(velec,cutoff_mask);
950             velec            = _mm256_andnot_ps(dummy_mask,velec);
951             velecsum         = _mm256_add_ps(velecsum,velec);
952
953             fscal            = felec;
954
955             fscal            = _mm256_and_ps(fscal,cutoff_mask);
956
957             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
958
959             /* Calculate temporary vectorial force */
960             tx               = _mm256_mul_ps(fscal,dx12);
961             ty               = _mm256_mul_ps(fscal,dy12);
962             tz               = _mm256_mul_ps(fscal,dz12);
963
964             /* Update vectorial force */
965             fix1             = _mm256_add_ps(fix1,tx);
966             fiy1             = _mm256_add_ps(fiy1,ty);
967             fiz1             = _mm256_add_ps(fiz1,tz);
968
969             fjx2             = _mm256_add_ps(fjx2,tx);
970             fjy2             = _mm256_add_ps(fjy2,ty);
971             fjz2             = _mm256_add_ps(fjz2,tz);
972
973             }
974
975             /**************************
976              * CALCULATE INTERACTIONS *
977              **************************/
978
979             if (gmx_mm256_any_lt(rsq13,rcutoff2))
980             {
981
982             /* REACTION-FIELD ELECTROSTATICS */
983             velec            = _mm256_mul_ps(qq13,_mm256_sub_ps(_mm256_add_ps(rinv13,_mm256_mul_ps(krf,rsq13)),crf));
984             felec            = _mm256_mul_ps(qq13,_mm256_sub_ps(_mm256_mul_ps(rinv13,rinvsq13),krf2));
985
986             cutoff_mask      = _mm256_cmp_ps(rsq13,rcutoff2,_CMP_LT_OQ);
987
988             /* Update potential sum for this i atom from the interaction with this j atom. */
989             velec            = _mm256_and_ps(velec,cutoff_mask);
990             velec            = _mm256_andnot_ps(dummy_mask,velec);
991             velecsum         = _mm256_add_ps(velecsum,velec);
992
993             fscal            = felec;
994
995             fscal            = _mm256_and_ps(fscal,cutoff_mask);
996
997             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
998
999             /* Calculate temporary vectorial force */
1000             tx               = _mm256_mul_ps(fscal,dx13);
1001             ty               = _mm256_mul_ps(fscal,dy13);
1002             tz               = _mm256_mul_ps(fscal,dz13);
1003
1004             /* Update vectorial force */
1005             fix1             = _mm256_add_ps(fix1,tx);
1006             fiy1             = _mm256_add_ps(fiy1,ty);
1007             fiz1             = _mm256_add_ps(fiz1,tz);
1008
1009             fjx3             = _mm256_add_ps(fjx3,tx);
1010             fjy3             = _mm256_add_ps(fjy3,ty);
1011             fjz3             = _mm256_add_ps(fjz3,tz);
1012
1013             }
1014
1015             /**************************
1016              * CALCULATE INTERACTIONS *
1017              **************************/
1018
1019             if (gmx_mm256_any_lt(rsq21,rcutoff2))
1020             {
1021
1022             /* REACTION-FIELD ELECTROSTATICS */
1023             velec            = _mm256_mul_ps(qq21,_mm256_sub_ps(_mm256_add_ps(rinv21,_mm256_mul_ps(krf,rsq21)),crf));
1024             felec            = _mm256_mul_ps(qq21,_mm256_sub_ps(_mm256_mul_ps(rinv21,rinvsq21),krf2));
1025
1026             cutoff_mask      = _mm256_cmp_ps(rsq21,rcutoff2,_CMP_LT_OQ);
1027
1028             /* Update potential sum for this i atom from the interaction with this j atom. */
1029             velec            = _mm256_and_ps(velec,cutoff_mask);
1030             velec            = _mm256_andnot_ps(dummy_mask,velec);
1031             velecsum         = _mm256_add_ps(velecsum,velec);
1032
1033             fscal            = felec;
1034
1035             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1036
1037             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1038
1039             /* Calculate temporary vectorial force */
1040             tx               = _mm256_mul_ps(fscal,dx21);
1041             ty               = _mm256_mul_ps(fscal,dy21);
1042             tz               = _mm256_mul_ps(fscal,dz21);
1043
1044             /* Update vectorial force */
1045             fix2             = _mm256_add_ps(fix2,tx);
1046             fiy2             = _mm256_add_ps(fiy2,ty);
1047             fiz2             = _mm256_add_ps(fiz2,tz);
1048
1049             fjx1             = _mm256_add_ps(fjx1,tx);
1050             fjy1             = _mm256_add_ps(fjy1,ty);
1051             fjz1             = _mm256_add_ps(fjz1,tz);
1052
1053             }
1054
1055             /**************************
1056              * CALCULATE INTERACTIONS *
1057              **************************/
1058
1059             if (gmx_mm256_any_lt(rsq22,rcutoff2))
1060             {
1061
1062             /* REACTION-FIELD ELECTROSTATICS */
1063             velec            = _mm256_mul_ps(qq22,_mm256_sub_ps(_mm256_add_ps(rinv22,_mm256_mul_ps(krf,rsq22)),crf));
1064             felec            = _mm256_mul_ps(qq22,_mm256_sub_ps(_mm256_mul_ps(rinv22,rinvsq22),krf2));
1065
1066             cutoff_mask      = _mm256_cmp_ps(rsq22,rcutoff2,_CMP_LT_OQ);
1067
1068             /* Update potential sum for this i atom from the interaction with this j atom. */
1069             velec            = _mm256_and_ps(velec,cutoff_mask);
1070             velec            = _mm256_andnot_ps(dummy_mask,velec);
1071             velecsum         = _mm256_add_ps(velecsum,velec);
1072
1073             fscal            = felec;
1074
1075             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1076
1077             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1078
1079             /* Calculate temporary vectorial force */
1080             tx               = _mm256_mul_ps(fscal,dx22);
1081             ty               = _mm256_mul_ps(fscal,dy22);
1082             tz               = _mm256_mul_ps(fscal,dz22);
1083
1084             /* Update vectorial force */
1085             fix2             = _mm256_add_ps(fix2,tx);
1086             fiy2             = _mm256_add_ps(fiy2,ty);
1087             fiz2             = _mm256_add_ps(fiz2,tz);
1088
1089             fjx2             = _mm256_add_ps(fjx2,tx);
1090             fjy2             = _mm256_add_ps(fjy2,ty);
1091             fjz2             = _mm256_add_ps(fjz2,tz);
1092
1093             }
1094
1095             /**************************
1096              * CALCULATE INTERACTIONS *
1097              **************************/
1098
1099             if (gmx_mm256_any_lt(rsq23,rcutoff2))
1100             {
1101
1102             /* REACTION-FIELD ELECTROSTATICS */
1103             velec            = _mm256_mul_ps(qq23,_mm256_sub_ps(_mm256_add_ps(rinv23,_mm256_mul_ps(krf,rsq23)),crf));
1104             felec            = _mm256_mul_ps(qq23,_mm256_sub_ps(_mm256_mul_ps(rinv23,rinvsq23),krf2));
1105
1106             cutoff_mask      = _mm256_cmp_ps(rsq23,rcutoff2,_CMP_LT_OQ);
1107
1108             /* Update potential sum for this i atom from the interaction with this j atom. */
1109             velec            = _mm256_and_ps(velec,cutoff_mask);
1110             velec            = _mm256_andnot_ps(dummy_mask,velec);
1111             velecsum         = _mm256_add_ps(velecsum,velec);
1112
1113             fscal            = felec;
1114
1115             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1116
1117             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1118
1119             /* Calculate temporary vectorial force */
1120             tx               = _mm256_mul_ps(fscal,dx23);
1121             ty               = _mm256_mul_ps(fscal,dy23);
1122             tz               = _mm256_mul_ps(fscal,dz23);
1123
1124             /* Update vectorial force */
1125             fix2             = _mm256_add_ps(fix2,tx);
1126             fiy2             = _mm256_add_ps(fiy2,ty);
1127             fiz2             = _mm256_add_ps(fiz2,tz);
1128
1129             fjx3             = _mm256_add_ps(fjx3,tx);
1130             fjy3             = _mm256_add_ps(fjy3,ty);
1131             fjz3             = _mm256_add_ps(fjz3,tz);
1132
1133             }
1134
1135             /**************************
1136              * CALCULATE INTERACTIONS *
1137              **************************/
1138
1139             if (gmx_mm256_any_lt(rsq31,rcutoff2))
1140             {
1141
1142             /* REACTION-FIELD ELECTROSTATICS */
1143             velec            = _mm256_mul_ps(qq31,_mm256_sub_ps(_mm256_add_ps(rinv31,_mm256_mul_ps(krf,rsq31)),crf));
1144             felec            = _mm256_mul_ps(qq31,_mm256_sub_ps(_mm256_mul_ps(rinv31,rinvsq31),krf2));
1145
1146             cutoff_mask      = _mm256_cmp_ps(rsq31,rcutoff2,_CMP_LT_OQ);
1147
1148             /* Update potential sum for this i atom from the interaction with this j atom. */
1149             velec            = _mm256_and_ps(velec,cutoff_mask);
1150             velec            = _mm256_andnot_ps(dummy_mask,velec);
1151             velecsum         = _mm256_add_ps(velecsum,velec);
1152
1153             fscal            = felec;
1154
1155             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1156
1157             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1158
1159             /* Calculate temporary vectorial force */
1160             tx               = _mm256_mul_ps(fscal,dx31);
1161             ty               = _mm256_mul_ps(fscal,dy31);
1162             tz               = _mm256_mul_ps(fscal,dz31);
1163
1164             /* Update vectorial force */
1165             fix3             = _mm256_add_ps(fix3,tx);
1166             fiy3             = _mm256_add_ps(fiy3,ty);
1167             fiz3             = _mm256_add_ps(fiz3,tz);
1168
1169             fjx1             = _mm256_add_ps(fjx1,tx);
1170             fjy1             = _mm256_add_ps(fjy1,ty);
1171             fjz1             = _mm256_add_ps(fjz1,tz);
1172
1173             }
1174
1175             /**************************
1176              * CALCULATE INTERACTIONS *
1177              **************************/
1178
1179             if (gmx_mm256_any_lt(rsq32,rcutoff2))
1180             {
1181
1182             /* REACTION-FIELD ELECTROSTATICS */
1183             velec            = _mm256_mul_ps(qq32,_mm256_sub_ps(_mm256_add_ps(rinv32,_mm256_mul_ps(krf,rsq32)),crf));
1184             felec            = _mm256_mul_ps(qq32,_mm256_sub_ps(_mm256_mul_ps(rinv32,rinvsq32),krf2));
1185
1186             cutoff_mask      = _mm256_cmp_ps(rsq32,rcutoff2,_CMP_LT_OQ);
1187
1188             /* Update potential sum for this i atom from the interaction with this j atom. */
1189             velec            = _mm256_and_ps(velec,cutoff_mask);
1190             velec            = _mm256_andnot_ps(dummy_mask,velec);
1191             velecsum         = _mm256_add_ps(velecsum,velec);
1192
1193             fscal            = felec;
1194
1195             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1196
1197             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1198
1199             /* Calculate temporary vectorial force */
1200             tx               = _mm256_mul_ps(fscal,dx32);
1201             ty               = _mm256_mul_ps(fscal,dy32);
1202             tz               = _mm256_mul_ps(fscal,dz32);
1203
1204             /* Update vectorial force */
1205             fix3             = _mm256_add_ps(fix3,tx);
1206             fiy3             = _mm256_add_ps(fiy3,ty);
1207             fiz3             = _mm256_add_ps(fiz3,tz);
1208
1209             fjx2             = _mm256_add_ps(fjx2,tx);
1210             fjy2             = _mm256_add_ps(fjy2,ty);
1211             fjz2             = _mm256_add_ps(fjz2,tz);
1212
1213             }
1214
1215             /**************************
1216              * CALCULATE INTERACTIONS *
1217              **************************/
1218
1219             if (gmx_mm256_any_lt(rsq33,rcutoff2))
1220             {
1221
1222             /* REACTION-FIELD ELECTROSTATICS */
1223             velec            = _mm256_mul_ps(qq33,_mm256_sub_ps(_mm256_add_ps(rinv33,_mm256_mul_ps(krf,rsq33)),crf));
1224             felec            = _mm256_mul_ps(qq33,_mm256_sub_ps(_mm256_mul_ps(rinv33,rinvsq33),krf2));
1225
1226             cutoff_mask      = _mm256_cmp_ps(rsq33,rcutoff2,_CMP_LT_OQ);
1227
1228             /* Update potential sum for this i atom from the interaction with this j atom. */
1229             velec            = _mm256_and_ps(velec,cutoff_mask);
1230             velec            = _mm256_andnot_ps(dummy_mask,velec);
1231             velecsum         = _mm256_add_ps(velecsum,velec);
1232
1233             fscal            = felec;
1234
1235             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1236
1237             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1238
1239             /* Calculate temporary vectorial force */
1240             tx               = _mm256_mul_ps(fscal,dx33);
1241             ty               = _mm256_mul_ps(fscal,dy33);
1242             tz               = _mm256_mul_ps(fscal,dz33);
1243
1244             /* Update vectorial force */
1245             fix3             = _mm256_add_ps(fix3,tx);
1246             fiy3             = _mm256_add_ps(fiy3,ty);
1247             fiz3             = _mm256_add_ps(fiz3,tz);
1248
1249             fjx3             = _mm256_add_ps(fjx3,tx);
1250             fjy3             = _mm256_add_ps(fjy3,ty);
1251             fjz3             = _mm256_add_ps(fjz3,tz);
1252
1253             }
1254
1255             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1256             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1257             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1258             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1259             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1260             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1261             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1262             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1263
1264             gmx_mm256_decrement_4rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
1265                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
1266                                                       fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1267
1268             /* Inner loop uses 368 flops */
1269         }
1270
1271         /* End of innermost loop */
1272
1273         gmx_mm256_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1274                                                  f+i_coord_offset,fshift+i_shift_offset);
1275
1276         ggid                        = gid[iidx];
1277         /* Update potential energies */
1278         gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
1279         gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
1280
1281         /* Increment number of inner iterations */
1282         inneriter                  += j_index_end - j_index_start;
1283
1284         /* Outer loop uses 26 flops */
1285     }
1286
1287     /* Increment number of outer iterations */
1288     outeriter        += nri;
1289
1290     /* Update outer/inner flops */
1291
1292     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4W4_VF,outeriter*26 + inneriter*368);
1293 }
1294 /*
1295  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwLJSh_GeomW4W4_F_avx_256_single
1296  * Electrostatics interaction: ReactionField
1297  * VdW interaction:            LennardJones
1298  * Geometry:                   Water4-Water4
1299  * Calculate force/pot:        Force
1300  */
1301 void
1302 nb_kernel_ElecRFCut_VdwLJSh_GeomW4W4_F_avx_256_single
1303                     (t_nblist                    * gmx_restrict       nlist,
1304                      rvec                        * gmx_restrict          xx,
1305                      rvec                        * gmx_restrict          ff,
1306                      t_forcerec                  * gmx_restrict          fr,
1307                      t_mdatoms                   * gmx_restrict     mdatoms,
1308                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
1309                      t_nrnb                      * gmx_restrict        nrnb)
1310 {
1311     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
1312      * just 0 for non-waters.
1313      * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
1314      * jnr indices corresponding to data put in the four positions in the SIMD register.
1315      */
1316     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
1317     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
1318     int              jnrA,jnrB,jnrC,jnrD;
1319     int              jnrE,jnrF,jnrG,jnrH;
1320     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
1321     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
1322     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
1323     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
1324     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
1325     real             rcutoff_scalar;
1326     real             *shiftvec,*fshift,*x,*f;
1327     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
1328     real             scratch[4*DIM];
1329     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
1330     real *           vdwioffsetptr0;
1331     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
1332     real *           vdwioffsetptr1;
1333     __m256           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
1334     real *           vdwioffsetptr2;
1335     __m256           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
1336     real *           vdwioffsetptr3;
1337     __m256           ix3,iy3,iz3,fix3,fiy3,fiz3,iq3,isai3;
1338     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
1339     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
1340     int              vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D,vdwjidx1E,vdwjidx1F,vdwjidx1G,vdwjidx1H;
1341     __m256           jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
1342     int              vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D,vdwjidx2E,vdwjidx2F,vdwjidx2G,vdwjidx2H;
1343     __m256           jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
1344     int              vdwjidx3A,vdwjidx3B,vdwjidx3C,vdwjidx3D,vdwjidx3E,vdwjidx3F,vdwjidx3G,vdwjidx3H;
1345     __m256           jx3,jy3,jz3,fjx3,fjy3,fjz3,jq3,isaj3;
1346     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
1347     __m256           dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
1348     __m256           dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
1349     __m256           dx13,dy13,dz13,rsq13,rinv13,rinvsq13,r13,qq13,c6_13,c12_13;
1350     __m256           dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
1351     __m256           dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
1352     __m256           dx23,dy23,dz23,rsq23,rinv23,rinvsq23,r23,qq23,c6_23,c12_23;
1353     __m256           dx31,dy31,dz31,rsq31,rinv31,rinvsq31,r31,qq31,c6_31,c12_31;
1354     __m256           dx32,dy32,dz32,rsq32,rinv32,rinvsq32,r32,qq32,c6_32,c12_32;
1355     __m256           dx33,dy33,dz33,rsq33,rinv33,rinvsq33,r33,qq33,c6_33,c12_33;
1356     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
1357     real             *charge;
1358     int              nvdwtype;
1359     __m256           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
1360     int              *vdwtype;
1361     real             *vdwparam;
1362     __m256           one_sixth   = _mm256_set1_ps(1.0/6.0);
1363     __m256           one_twelfth = _mm256_set1_ps(1.0/12.0);
1364     __m256           dummy_mask,cutoff_mask;
1365     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
1366     __m256           one     = _mm256_set1_ps(1.0);
1367     __m256           two     = _mm256_set1_ps(2.0);
1368     x                = xx[0];
1369     f                = ff[0];
1370
1371     nri              = nlist->nri;
1372     iinr             = nlist->iinr;
1373     jindex           = nlist->jindex;
1374     jjnr             = nlist->jjnr;
1375     shiftidx         = nlist->shift;
1376     gid              = nlist->gid;
1377     shiftvec         = fr->shift_vec[0];
1378     fshift           = fr->fshift[0];
1379     facel            = _mm256_set1_ps(fr->epsfac);
1380     charge           = mdatoms->chargeA;
1381     krf              = _mm256_set1_ps(fr->ic->k_rf);
1382     krf2             = _mm256_set1_ps(fr->ic->k_rf*2.0);
1383     crf              = _mm256_set1_ps(fr->ic->c_rf);
1384     nvdwtype         = fr->ntype;
1385     vdwparam         = fr->nbfp;
1386     vdwtype          = mdatoms->typeA;
1387
1388     /* Setup water-specific parameters */
1389     inr              = nlist->iinr[0];
1390     iq1              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
1391     iq2              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
1392     iq3              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+3]));
1393     vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
1394
1395     jq1              = _mm256_set1_ps(charge[inr+1]);
1396     jq2              = _mm256_set1_ps(charge[inr+2]);
1397     jq3              = _mm256_set1_ps(charge[inr+3]);
1398     vdwjidx0A        = 2*vdwtype[inr+0];
1399     c6_00            = _mm256_set1_ps(vdwioffsetptr0[vdwjidx0A]);
1400     c12_00           = _mm256_set1_ps(vdwioffsetptr0[vdwjidx0A+1]);
1401     qq11             = _mm256_mul_ps(iq1,jq1);
1402     qq12             = _mm256_mul_ps(iq1,jq2);
1403     qq13             = _mm256_mul_ps(iq1,jq3);
1404     qq21             = _mm256_mul_ps(iq2,jq1);
1405     qq22             = _mm256_mul_ps(iq2,jq2);
1406     qq23             = _mm256_mul_ps(iq2,jq3);
1407     qq31             = _mm256_mul_ps(iq3,jq1);
1408     qq32             = _mm256_mul_ps(iq3,jq2);
1409     qq33             = _mm256_mul_ps(iq3,jq3);
1410
1411     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
1412     rcutoff_scalar   = fr->rcoulomb;
1413     rcutoff          = _mm256_set1_ps(rcutoff_scalar);
1414     rcutoff2         = _mm256_mul_ps(rcutoff,rcutoff);
1415
1416     sh_vdw_invrcut6  = _mm256_set1_ps(fr->ic->sh_invrc6);
1417     rvdw             = _mm256_set1_ps(fr->rvdw);
1418
1419     /* Avoid stupid compiler warnings */
1420     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
1421     j_coord_offsetA = 0;
1422     j_coord_offsetB = 0;
1423     j_coord_offsetC = 0;
1424     j_coord_offsetD = 0;
1425     j_coord_offsetE = 0;
1426     j_coord_offsetF = 0;
1427     j_coord_offsetG = 0;
1428     j_coord_offsetH = 0;
1429
1430     outeriter        = 0;
1431     inneriter        = 0;
1432
1433     for(iidx=0;iidx<4*DIM;iidx++)
1434     {
1435         scratch[iidx] = 0.0;
1436     }
1437
1438     /* Start outer loop over neighborlists */
1439     for(iidx=0; iidx<nri; iidx++)
1440     {
1441         /* Load shift vector for this list */
1442         i_shift_offset   = DIM*shiftidx[iidx];
1443
1444         /* Load limits for loop over neighbors */
1445         j_index_start    = jindex[iidx];
1446         j_index_end      = jindex[iidx+1];
1447
1448         /* Get outer coordinate index */
1449         inr              = iinr[iidx];
1450         i_coord_offset   = DIM*inr;
1451
1452         /* Load i particle coords and add shift vector */
1453         gmx_mm256_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
1454                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
1455
1456         fix0             = _mm256_setzero_ps();
1457         fiy0             = _mm256_setzero_ps();
1458         fiz0             = _mm256_setzero_ps();
1459         fix1             = _mm256_setzero_ps();
1460         fiy1             = _mm256_setzero_ps();
1461         fiz1             = _mm256_setzero_ps();
1462         fix2             = _mm256_setzero_ps();
1463         fiy2             = _mm256_setzero_ps();
1464         fiz2             = _mm256_setzero_ps();
1465         fix3             = _mm256_setzero_ps();
1466         fiy3             = _mm256_setzero_ps();
1467         fiz3             = _mm256_setzero_ps();
1468
1469         /* Start inner kernel loop */
1470         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
1471         {
1472
1473             /* Get j neighbor index, and coordinate index */
1474             jnrA             = jjnr[jidx];
1475             jnrB             = jjnr[jidx+1];
1476             jnrC             = jjnr[jidx+2];
1477             jnrD             = jjnr[jidx+3];
1478             jnrE             = jjnr[jidx+4];
1479             jnrF             = jjnr[jidx+5];
1480             jnrG             = jjnr[jidx+6];
1481             jnrH             = jjnr[jidx+7];
1482             j_coord_offsetA  = DIM*jnrA;
1483             j_coord_offsetB  = DIM*jnrB;
1484             j_coord_offsetC  = DIM*jnrC;
1485             j_coord_offsetD  = DIM*jnrD;
1486             j_coord_offsetE  = DIM*jnrE;
1487             j_coord_offsetF  = DIM*jnrF;
1488             j_coord_offsetG  = DIM*jnrG;
1489             j_coord_offsetH  = DIM*jnrH;
1490
1491             /* load j atom coordinates */
1492             gmx_mm256_load_4rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1493                                                  x+j_coord_offsetC,x+j_coord_offsetD,
1494                                                  x+j_coord_offsetE,x+j_coord_offsetF,
1495                                                  x+j_coord_offsetG,x+j_coord_offsetH,
1496                                                  &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
1497                                                  &jy2,&jz2,&jx3,&jy3,&jz3);
1498
1499             /* Calculate displacement vector */
1500             dx00             = _mm256_sub_ps(ix0,jx0);
1501             dy00             = _mm256_sub_ps(iy0,jy0);
1502             dz00             = _mm256_sub_ps(iz0,jz0);
1503             dx11             = _mm256_sub_ps(ix1,jx1);
1504             dy11             = _mm256_sub_ps(iy1,jy1);
1505             dz11             = _mm256_sub_ps(iz1,jz1);
1506             dx12             = _mm256_sub_ps(ix1,jx2);
1507             dy12             = _mm256_sub_ps(iy1,jy2);
1508             dz12             = _mm256_sub_ps(iz1,jz2);
1509             dx13             = _mm256_sub_ps(ix1,jx3);
1510             dy13             = _mm256_sub_ps(iy1,jy3);
1511             dz13             = _mm256_sub_ps(iz1,jz3);
1512             dx21             = _mm256_sub_ps(ix2,jx1);
1513             dy21             = _mm256_sub_ps(iy2,jy1);
1514             dz21             = _mm256_sub_ps(iz2,jz1);
1515             dx22             = _mm256_sub_ps(ix2,jx2);
1516             dy22             = _mm256_sub_ps(iy2,jy2);
1517             dz22             = _mm256_sub_ps(iz2,jz2);
1518             dx23             = _mm256_sub_ps(ix2,jx3);
1519             dy23             = _mm256_sub_ps(iy2,jy3);
1520             dz23             = _mm256_sub_ps(iz2,jz3);
1521             dx31             = _mm256_sub_ps(ix3,jx1);
1522             dy31             = _mm256_sub_ps(iy3,jy1);
1523             dz31             = _mm256_sub_ps(iz3,jz1);
1524             dx32             = _mm256_sub_ps(ix3,jx2);
1525             dy32             = _mm256_sub_ps(iy3,jy2);
1526             dz32             = _mm256_sub_ps(iz3,jz2);
1527             dx33             = _mm256_sub_ps(ix3,jx3);
1528             dy33             = _mm256_sub_ps(iy3,jy3);
1529             dz33             = _mm256_sub_ps(iz3,jz3);
1530
1531             /* Calculate squared distance and things based on it */
1532             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
1533             rsq11            = gmx_mm256_calc_rsq_ps(dx11,dy11,dz11);
1534             rsq12            = gmx_mm256_calc_rsq_ps(dx12,dy12,dz12);
1535             rsq13            = gmx_mm256_calc_rsq_ps(dx13,dy13,dz13);
1536             rsq21            = gmx_mm256_calc_rsq_ps(dx21,dy21,dz21);
1537             rsq22            = gmx_mm256_calc_rsq_ps(dx22,dy22,dz22);
1538             rsq23            = gmx_mm256_calc_rsq_ps(dx23,dy23,dz23);
1539             rsq31            = gmx_mm256_calc_rsq_ps(dx31,dy31,dz31);
1540             rsq32            = gmx_mm256_calc_rsq_ps(dx32,dy32,dz32);
1541             rsq33            = gmx_mm256_calc_rsq_ps(dx33,dy33,dz33);
1542
1543             rinv11           = gmx_mm256_invsqrt_ps(rsq11);
1544             rinv12           = gmx_mm256_invsqrt_ps(rsq12);
1545             rinv13           = gmx_mm256_invsqrt_ps(rsq13);
1546             rinv21           = gmx_mm256_invsqrt_ps(rsq21);
1547             rinv22           = gmx_mm256_invsqrt_ps(rsq22);
1548             rinv23           = gmx_mm256_invsqrt_ps(rsq23);
1549             rinv31           = gmx_mm256_invsqrt_ps(rsq31);
1550             rinv32           = gmx_mm256_invsqrt_ps(rsq32);
1551             rinv33           = gmx_mm256_invsqrt_ps(rsq33);
1552
1553             rinvsq00         = gmx_mm256_inv_ps(rsq00);
1554             rinvsq11         = _mm256_mul_ps(rinv11,rinv11);
1555             rinvsq12         = _mm256_mul_ps(rinv12,rinv12);
1556             rinvsq13         = _mm256_mul_ps(rinv13,rinv13);
1557             rinvsq21         = _mm256_mul_ps(rinv21,rinv21);
1558             rinvsq22         = _mm256_mul_ps(rinv22,rinv22);
1559             rinvsq23         = _mm256_mul_ps(rinv23,rinv23);
1560             rinvsq31         = _mm256_mul_ps(rinv31,rinv31);
1561             rinvsq32         = _mm256_mul_ps(rinv32,rinv32);
1562             rinvsq33         = _mm256_mul_ps(rinv33,rinv33);
1563
1564             fjx0             = _mm256_setzero_ps();
1565             fjy0             = _mm256_setzero_ps();
1566             fjz0             = _mm256_setzero_ps();
1567             fjx1             = _mm256_setzero_ps();
1568             fjy1             = _mm256_setzero_ps();
1569             fjz1             = _mm256_setzero_ps();
1570             fjx2             = _mm256_setzero_ps();
1571             fjy2             = _mm256_setzero_ps();
1572             fjz2             = _mm256_setzero_ps();
1573             fjx3             = _mm256_setzero_ps();
1574             fjy3             = _mm256_setzero_ps();
1575             fjz3             = _mm256_setzero_ps();
1576
1577             /**************************
1578              * CALCULATE INTERACTIONS *
1579              **************************/
1580
1581             if (gmx_mm256_any_lt(rsq00,rcutoff2))
1582             {
1583
1584             /* LENNARD-JONES DISPERSION/REPULSION */
1585
1586             rinvsix          = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
1587             fvdw             = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_00,rinvsix),c6_00),_mm256_mul_ps(rinvsix,rinvsq00));
1588
1589             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
1590
1591             fscal            = fvdw;
1592
1593             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1594
1595             /* Calculate temporary vectorial force */
1596             tx               = _mm256_mul_ps(fscal,dx00);
1597             ty               = _mm256_mul_ps(fscal,dy00);
1598             tz               = _mm256_mul_ps(fscal,dz00);
1599
1600             /* Update vectorial force */
1601             fix0             = _mm256_add_ps(fix0,tx);
1602             fiy0             = _mm256_add_ps(fiy0,ty);
1603             fiz0             = _mm256_add_ps(fiz0,tz);
1604
1605             fjx0             = _mm256_add_ps(fjx0,tx);
1606             fjy0             = _mm256_add_ps(fjy0,ty);
1607             fjz0             = _mm256_add_ps(fjz0,tz);
1608
1609             }
1610
1611             /**************************
1612              * CALCULATE INTERACTIONS *
1613              **************************/
1614
1615             if (gmx_mm256_any_lt(rsq11,rcutoff2))
1616             {
1617
1618             /* REACTION-FIELD ELECTROSTATICS */
1619             felec            = _mm256_mul_ps(qq11,_mm256_sub_ps(_mm256_mul_ps(rinv11,rinvsq11),krf2));
1620
1621             cutoff_mask      = _mm256_cmp_ps(rsq11,rcutoff2,_CMP_LT_OQ);
1622
1623             fscal            = felec;
1624
1625             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1626
1627             /* Calculate temporary vectorial force */
1628             tx               = _mm256_mul_ps(fscal,dx11);
1629             ty               = _mm256_mul_ps(fscal,dy11);
1630             tz               = _mm256_mul_ps(fscal,dz11);
1631
1632             /* Update vectorial force */
1633             fix1             = _mm256_add_ps(fix1,tx);
1634             fiy1             = _mm256_add_ps(fiy1,ty);
1635             fiz1             = _mm256_add_ps(fiz1,tz);
1636
1637             fjx1             = _mm256_add_ps(fjx1,tx);
1638             fjy1             = _mm256_add_ps(fjy1,ty);
1639             fjz1             = _mm256_add_ps(fjz1,tz);
1640
1641             }
1642
1643             /**************************
1644              * CALCULATE INTERACTIONS *
1645              **************************/
1646
1647             if (gmx_mm256_any_lt(rsq12,rcutoff2))
1648             {
1649
1650             /* REACTION-FIELD ELECTROSTATICS */
1651             felec            = _mm256_mul_ps(qq12,_mm256_sub_ps(_mm256_mul_ps(rinv12,rinvsq12),krf2));
1652
1653             cutoff_mask      = _mm256_cmp_ps(rsq12,rcutoff2,_CMP_LT_OQ);
1654
1655             fscal            = felec;
1656
1657             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1658
1659             /* Calculate temporary vectorial force */
1660             tx               = _mm256_mul_ps(fscal,dx12);
1661             ty               = _mm256_mul_ps(fscal,dy12);
1662             tz               = _mm256_mul_ps(fscal,dz12);
1663
1664             /* Update vectorial force */
1665             fix1             = _mm256_add_ps(fix1,tx);
1666             fiy1             = _mm256_add_ps(fiy1,ty);
1667             fiz1             = _mm256_add_ps(fiz1,tz);
1668
1669             fjx2             = _mm256_add_ps(fjx2,tx);
1670             fjy2             = _mm256_add_ps(fjy2,ty);
1671             fjz2             = _mm256_add_ps(fjz2,tz);
1672
1673             }
1674
1675             /**************************
1676              * CALCULATE INTERACTIONS *
1677              **************************/
1678
1679             if (gmx_mm256_any_lt(rsq13,rcutoff2))
1680             {
1681
1682             /* REACTION-FIELD ELECTROSTATICS */
1683             felec            = _mm256_mul_ps(qq13,_mm256_sub_ps(_mm256_mul_ps(rinv13,rinvsq13),krf2));
1684
1685             cutoff_mask      = _mm256_cmp_ps(rsq13,rcutoff2,_CMP_LT_OQ);
1686
1687             fscal            = felec;
1688
1689             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1690
1691             /* Calculate temporary vectorial force */
1692             tx               = _mm256_mul_ps(fscal,dx13);
1693             ty               = _mm256_mul_ps(fscal,dy13);
1694             tz               = _mm256_mul_ps(fscal,dz13);
1695
1696             /* Update vectorial force */
1697             fix1             = _mm256_add_ps(fix1,tx);
1698             fiy1             = _mm256_add_ps(fiy1,ty);
1699             fiz1             = _mm256_add_ps(fiz1,tz);
1700
1701             fjx3             = _mm256_add_ps(fjx3,tx);
1702             fjy3             = _mm256_add_ps(fjy3,ty);
1703             fjz3             = _mm256_add_ps(fjz3,tz);
1704
1705             }
1706
1707             /**************************
1708              * CALCULATE INTERACTIONS *
1709              **************************/
1710
1711             if (gmx_mm256_any_lt(rsq21,rcutoff2))
1712             {
1713
1714             /* REACTION-FIELD ELECTROSTATICS */
1715             felec            = _mm256_mul_ps(qq21,_mm256_sub_ps(_mm256_mul_ps(rinv21,rinvsq21),krf2));
1716
1717             cutoff_mask      = _mm256_cmp_ps(rsq21,rcutoff2,_CMP_LT_OQ);
1718
1719             fscal            = felec;
1720
1721             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1722
1723             /* Calculate temporary vectorial force */
1724             tx               = _mm256_mul_ps(fscal,dx21);
1725             ty               = _mm256_mul_ps(fscal,dy21);
1726             tz               = _mm256_mul_ps(fscal,dz21);
1727
1728             /* Update vectorial force */
1729             fix2             = _mm256_add_ps(fix2,tx);
1730             fiy2             = _mm256_add_ps(fiy2,ty);
1731             fiz2             = _mm256_add_ps(fiz2,tz);
1732
1733             fjx1             = _mm256_add_ps(fjx1,tx);
1734             fjy1             = _mm256_add_ps(fjy1,ty);
1735             fjz1             = _mm256_add_ps(fjz1,tz);
1736
1737             }
1738
1739             /**************************
1740              * CALCULATE INTERACTIONS *
1741              **************************/
1742
1743             if (gmx_mm256_any_lt(rsq22,rcutoff2))
1744             {
1745
1746             /* REACTION-FIELD ELECTROSTATICS */
1747             felec            = _mm256_mul_ps(qq22,_mm256_sub_ps(_mm256_mul_ps(rinv22,rinvsq22),krf2));
1748
1749             cutoff_mask      = _mm256_cmp_ps(rsq22,rcutoff2,_CMP_LT_OQ);
1750
1751             fscal            = felec;
1752
1753             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1754
1755             /* Calculate temporary vectorial force */
1756             tx               = _mm256_mul_ps(fscal,dx22);
1757             ty               = _mm256_mul_ps(fscal,dy22);
1758             tz               = _mm256_mul_ps(fscal,dz22);
1759
1760             /* Update vectorial force */
1761             fix2             = _mm256_add_ps(fix2,tx);
1762             fiy2             = _mm256_add_ps(fiy2,ty);
1763             fiz2             = _mm256_add_ps(fiz2,tz);
1764
1765             fjx2             = _mm256_add_ps(fjx2,tx);
1766             fjy2             = _mm256_add_ps(fjy2,ty);
1767             fjz2             = _mm256_add_ps(fjz2,tz);
1768
1769             }
1770
1771             /**************************
1772              * CALCULATE INTERACTIONS *
1773              **************************/
1774
1775             if (gmx_mm256_any_lt(rsq23,rcutoff2))
1776             {
1777
1778             /* REACTION-FIELD ELECTROSTATICS */
1779             felec            = _mm256_mul_ps(qq23,_mm256_sub_ps(_mm256_mul_ps(rinv23,rinvsq23),krf2));
1780
1781             cutoff_mask      = _mm256_cmp_ps(rsq23,rcutoff2,_CMP_LT_OQ);
1782
1783             fscal            = felec;
1784
1785             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1786
1787             /* Calculate temporary vectorial force */
1788             tx               = _mm256_mul_ps(fscal,dx23);
1789             ty               = _mm256_mul_ps(fscal,dy23);
1790             tz               = _mm256_mul_ps(fscal,dz23);
1791
1792             /* Update vectorial force */
1793             fix2             = _mm256_add_ps(fix2,tx);
1794             fiy2             = _mm256_add_ps(fiy2,ty);
1795             fiz2             = _mm256_add_ps(fiz2,tz);
1796
1797             fjx3             = _mm256_add_ps(fjx3,tx);
1798             fjy3             = _mm256_add_ps(fjy3,ty);
1799             fjz3             = _mm256_add_ps(fjz3,tz);
1800
1801             }
1802
1803             /**************************
1804              * CALCULATE INTERACTIONS *
1805              **************************/
1806
1807             if (gmx_mm256_any_lt(rsq31,rcutoff2))
1808             {
1809
1810             /* REACTION-FIELD ELECTROSTATICS */
1811             felec            = _mm256_mul_ps(qq31,_mm256_sub_ps(_mm256_mul_ps(rinv31,rinvsq31),krf2));
1812
1813             cutoff_mask      = _mm256_cmp_ps(rsq31,rcutoff2,_CMP_LT_OQ);
1814
1815             fscal            = felec;
1816
1817             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1818
1819             /* Calculate temporary vectorial force */
1820             tx               = _mm256_mul_ps(fscal,dx31);
1821             ty               = _mm256_mul_ps(fscal,dy31);
1822             tz               = _mm256_mul_ps(fscal,dz31);
1823
1824             /* Update vectorial force */
1825             fix3             = _mm256_add_ps(fix3,tx);
1826             fiy3             = _mm256_add_ps(fiy3,ty);
1827             fiz3             = _mm256_add_ps(fiz3,tz);
1828
1829             fjx1             = _mm256_add_ps(fjx1,tx);
1830             fjy1             = _mm256_add_ps(fjy1,ty);
1831             fjz1             = _mm256_add_ps(fjz1,tz);
1832
1833             }
1834
1835             /**************************
1836              * CALCULATE INTERACTIONS *
1837              **************************/
1838
1839             if (gmx_mm256_any_lt(rsq32,rcutoff2))
1840             {
1841
1842             /* REACTION-FIELD ELECTROSTATICS */
1843             felec            = _mm256_mul_ps(qq32,_mm256_sub_ps(_mm256_mul_ps(rinv32,rinvsq32),krf2));
1844
1845             cutoff_mask      = _mm256_cmp_ps(rsq32,rcutoff2,_CMP_LT_OQ);
1846
1847             fscal            = felec;
1848
1849             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1850
1851             /* Calculate temporary vectorial force */
1852             tx               = _mm256_mul_ps(fscal,dx32);
1853             ty               = _mm256_mul_ps(fscal,dy32);
1854             tz               = _mm256_mul_ps(fscal,dz32);
1855
1856             /* Update vectorial force */
1857             fix3             = _mm256_add_ps(fix3,tx);
1858             fiy3             = _mm256_add_ps(fiy3,ty);
1859             fiz3             = _mm256_add_ps(fiz3,tz);
1860
1861             fjx2             = _mm256_add_ps(fjx2,tx);
1862             fjy2             = _mm256_add_ps(fjy2,ty);
1863             fjz2             = _mm256_add_ps(fjz2,tz);
1864
1865             }
1866
1867             /**************************
1868              * CALCULATE INTERACTIONS *
1869              **************************/
1870
1871             if (gmx_mm256_any_lt(rsq33,rcutoff2))
1872             {
1873
1874             /* REACTION-FIELD ELECTROSTATICS */
1875             felec            = _mm256_mul_ps(qq33,_mm256_sub_ps(_mm256_mul_ps(rinv33,rinvsq33),krf2));
1876
1877             cutoff_mask      = _mm256_cmp_ps(rsq33,rcutoff2,_CMP_LT_OQ);
1878
1879             fscal            = felec;
1880
1881             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1882
1883             /* Calculate temporary vectorial force */
1884             tx               = _mm256_mul_ps(fscal,dx33);
1885             ty               = _mm256_mul_ps(fscal,dy33);
1886             tz               = _mm256_mul_ps(fscal,dz33);
1887
1888             /* Update vectorial force */
1889             fix3             = _mm256_add_ps(fix3,tx);
1890             fiy3             = _mm256_add_ps(fiy3,ty);
1891             fiz3             = _mm256_add_ps(fiz3,tz);
1892
1893             fjx3             = _mm256_add_ps(fjx3,tx);
1894             fjy3             = _mm256_add_ps(fjy3,ty);
1895             fjz3             = _mm256_add_ps(fjz3,tz);
1896
1897             }
1898
1899             fjptrA             = f+j_coord_offsetA;
1900             fjptrB             = f+j_coord_offsetB;
1901             fjptrC             = f+j_coord_offsetC;
1902             fjptrD             = f+j_coord_offsetD;
1903             fjptrE             = f+j_coord_offsetE;
1904             fjptrF             = f+j_coord_offsetF;
1905             fjptrG             = f+j_coord_offsetG;
1906             fjptrH             = f+j_coord_offsetH;
1907
1908             gmx_mm256_decrement_4rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
1909                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
1910                                                       fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1911
1912             /* Inner loop uses 303 flops */
1913         }
1914
1915         if(jidx<j_index_end)
1916         {
1917
1918             /* Get j neighbor index, and coordinate index */
1919             jnrlistA         = jjnr[jidx];
1920             jnrlistB         = jjnr[jidx+1];
1921             jnrlistC         = jjnr[jidx+2];
1922             jnrlistD         = jjnr[jidx+3];
1923             jnrlistE         = jjnr[jidx+4];
1924             jnrlistF         = jjnr[jidx+5];
1925             jnrlistG         = jjnr[jidx+6];
1926             jnrlistH         = jjnr[jidx+7];
1927             /* Sign of each element will be negative for non-real atoms.
1928              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1929              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
1930              */
1931             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
1932                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
1933                                             
1934             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
1935             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
1936             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
1937             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
1938             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
1939             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
1940             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
1941             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
1942             j_coord_offsetA  = DIM*jnrA;
1943             j_coord_offsetB  = DIM*jnrB;
1944             j_coord_offsetC  = DIM*jnrC;
1945             j_coord_offsetD  = DIM*jnrD;
1946             j_coord_offsetE  = DIM*jnrE;
1947             j_coord_offsetF  = DIM*jnrF;
1948             j_coord_offsetG  = DIM*jnrG;
1949             j_coord_offsetH  = DIM*jnrH;
1950
1951             /* load j atom coordinates */
1952             gmx_mm256_load_4rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
1953                                                  x+j_coord_offsetC,x+j_coord_offsetD,
1954                                                  x+j_coord_offsetE,x+j_coord_offsetF,
1955                                                  x+j_coord_offsetG,x+j_coord_offsetH,
1956                                                  &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
1957                                                  &jy2,&jz2,&jx3,&jy3,&jz3);
1958
1959             /* Calculate displacement vector */
1960             dx00             = _mm256_sub_ps(ix0,jx0);
1961             dy00             = _mm256_sub_ps(iy0,jy0);
1962             dz00             = _mm256_sub_ps(iz0,jz0);
1963             dx11             = _mm256_sub_ps(ix1,jx1);
1964             dy11             = _mm256_sub_ps(iy1,jy1);
1965             dz11             = _mm256_sub_ps(iz1,jz1);
1966             dx12             = _mm256_sub_ps(ix1,jx2);
1967             dy12             = _mm256_sub_ps(iy1,jy2);
1968             dz12             = _mm256_sub_ps(iz1,jz2);
1969             dx13             = _mm256_sub_ps(ix1,jx3);
1970             dy13             = _mm256_sub_ps(iy1,jy3);
1971             dz13             = _mm256_sub_ps(iz1,jz3);
1972             dx21             = _mm256_sub_ps(ix2,jx1);
1973             dy21             = _mm256_sub_ps(iy2,jy1);
1974             dz21             = _mm256_sub_ps(iz2,jz1);
1975             dx22             = _mm256_sub_ps(ix2,jx2);
1976             dy22             = _mm256_sub_ps(iy2,jy2);
1977             dz22             = _mm256_sub_ps(iz2,jz2);
1978             dx23             = _mm256_sub_ps(ix2,jx3);
1979             dy23             = _mm256_sub_ps(iy2,jy3);
1980             dz23             = _mm256_sub_ps(iz2,jz3);
1981             dx31             = _mm256_sub_ps(ix3,jx1);
1982             dy31             = _mm256_sub_ps(iy3,jy1);
1983             dz31             = _mm256_sub_ps(iz3,jz1);
1984             dx32             = _mm256_sub_ps(ix3,jx2);
1985             dy32             = _mm256_sub_ps(iy3,jy2);
1986             dz32             = _mm256_sub_ps(iz3,jz2);
1987             dx33             = _mm256_sub_ps(ix3,jx3);
1988             dy33             = _mm256_sub_ps(iy3,jy3);
1989             dz33             = _mm256_sub_ps(iz3,jz3);
1990
1991             /* Calculate squared distance and things based on it */
1992             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
1993             rsq11            = gmx_mm256_calc_rsq_ps(dx11,dy11,dz11);
1994             rsq12            = gmx_mm256_calc_rsq_ps(dx12,dy12,dz12);
1995             rsq13            = gmx_mm256_calc_rsq_ps(dx13,dy13,dz13);
1996             rsq21            = gmx_mm256_calc_rsq_ps(dx21,dy21,dz21);
1997             rsq22            = gmx_mm256_calc_rsq_ps(dx22,dy22,dz22);
1998             rsq23            = gmx_mm256_calc_rsq_ps(dx23,dy23,dz23);
1999             rsq31            = gmx_mm256_calc_rsq_ps(dx31,dy31,dz31);
2000             rsq32            = gmx_mm256_calc_rsq_ps(dx32,dy32,dz32);
2001             rsq33            = gmx_mm256_calc_rsq_ps(dx33,dy33,dz33);
2002
2003             rinv11           = gmx_mm256_invsqrt_ps(rsq11);
2004             rinv12           = gmx_mm256_invsqrt_ps(rsq12);
2005             rinv13           = gmx_mm256_invsqrt_ps(rsq13);
2006             rinv21           = gmx_mm256_invsqrt_ps(rsq21);
2007             rinv22           = gmx_mm256_invsqrt_ps(rsq22);
2008             rinv23           = gmx_mm256_invsqrt_ps(rsq23);
2009             rinv31           = gmx_mm256_invsqrt_ps(rsq31);
2010             rinv32           = gmx_mm256_invsqrt_ps(rsq32);
2011             rinv33           = gmx_mm256_invsqrt_ps(rsq33);
2012
2013             rinvsq00         = gmx_mm256_inv_ps(rsq00);
2014             rinvsq11         = _mm256_mul_ps(rinv11,rinv11);
2015             rinvsq12         = _mm256_mul_ps(rinv12,rinv12);
2016             rinvsq13         = _mm256_mul_ps(rinv13,rinv13);
2017             rinvsq21         = _mm256_mul_ps(rinv21,rinv21);
2018             rinvsq22         = _mm256_mul_ps(rinv22,rinv22);
2019             rinvsq23         = _mm256_mul_ps(rinv23,rinv23);
2020             rinvsq31         = _mm256_mul_ps(rinv31,rinv31);
2021             rinvsq32         = _mm256_mul_ps(rinv32,rinv32);
2022             rinvsq33         = _mm256_mul_ps(rinv33,rinv33);
2023
2024             fjx0             = _mm256_setzero_ps();
2025             fjy0             = _mm256_setzero_ps();
2026             fjz0             = _mm256_setzero_ps();
2027             fjx1             = _mm256_setzero_ps();
2028             fjy1             = _mm256_setzero_ps();
2029             fjz1             = _mm256_setzero_ps();
2030             fjx2             = _mm256_setzero_ps();
2031             fjy2             = _mm256_setzero_ps();
2032             fjz2             = _mm256_setzero_ps();
2033             fjx3             = _mm256_setzero_ps();
2034             fjy3             = _mm256_setzero_ps();
2035             fjz3             = _mm256_setzero_ps();
2036
2037             /**************************
2038              * CALCULATE INTERACTIONS *
2039              **************************/
2040
2041             if (gmx_mm256_any_lt(rsq00,rcutoff2))
2042             {
2043
2044             /* LENNARD-JONES DISPERSION/REPULSION */
2045
2046             rinvsix          = _mm256_mul_ps(_mm256_mul_ps(rinvsq00,rinvsq00),rinvsq00);
2047             fvdw             = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_00,rinvsix),c6_00),_mm256_mul_ps(rinvsix,rinvsq00));
2048
2049             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
2050
2051             fscal            = fvdw;
2052
2053             fscal            = _mm256_and_ps(fscal,cutoff_mask);
2054
2055             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
2056
2057             /* Calculate temporary vectorial force */
2058             tx               = _mm256_mul_ps(fscal,dx00);
2059             ty               = _mm256_mul_ps(fscal,dy00);
2060             tz               = _mm256_mul_ps(fscal,dz00);
2061
2062             /* Update vectorial force */
2063             fix0             = _mm256_add_ps(fix0,tx);
2064             fiy0             = _mm256_add_ps(fiy0,ty);
2065             fiz0             = _mm256_add_ps(fiz0,tz);
2066
2067             fjx0             = _mm256_add_ps(fjx0,tx);
2068             fjy0             = _mm256_add_ps(fjy0,ty);
2069             fjz0             = _mm256_add_ps(fjz0,tz);
2070
2071             }
2072
2073             /**************************
2074              * CALCULATE INTERACTIONS *
2075              **************************/
2076
2077             if (gmx_mm256_any_lt(rsq11,rcutoff2))
2078             {
2079
2080             /* REACTION-FIELD ELECTROSTATICS */
2081             felec            = _mm256_mul_ps(qq11,_mm256_sub_ps(_mm256_mul_ps(rinv11,rinvsq11),krf2));
2082
2083             cutoff_mask      = _mm256_cmp_ps(rsq11,rcutoff2,_CMP_LT_OQ);
2084
2085             fscal            = felec;
2086
2087             fscal            = _mm256_and_ps(fscal,cutoff_mask);
2088
2089             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
2090
2091             /* Calculate temporary vectorial force */
2092             tx               = _mm256_mul_ps(fscal,dx11);
2093             ty               = _mm256_mul_ps(fscal,dy11);
2094             tz               = _mm256_mul_ps(fscal,dz11);
2095
2096             /* Update vectorial force */
2097             fix1             = _mm256_add_ps(fix1,tx);
2098             fiy1             = _mm256_add_ps(fiy1,ty);
2099             fiz1             = _mm256_add_ps(fiz1,tz);
2100
2101             fjx1             = _mm256_add_ps(fjx1,tx);
2102             fjy1             = _mm256_add_ps(fjy1,ty);
2103             fjz1             = _mm256_add_ps(fjz1,tz);
2104
2105             }
2106
2107             /**************************
2108              * CALCULATE INTERACTIONS *
2109              **************************/
2110
2111             if (gmx_mm256_any_lt(rsq12,rcutoff2))
2112             {
2113
2114             /* REACTION-FIELD ELECTROSTATICS */
2115             felec            = _mm256_mul_ps(qq12,_mm256_sub_ps(_mm256_mul_ps(rinv12,rinvsq12),krf2));
2116
2117             cutoff_mask      = _mm256_cmp_ps(rsq12,rcutoff2,_CMP_LT_OQ);
2118
2119             fscal            = felec;
2120
2121             fscal            = _mm256_and_ps(fscal,cutoff_mask);
2122
2123             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
2124
2125             /* Calculate temporary vectorial force */
2126             tx               = _mm256_mul_ps(fscal,dx12);
2127             ty               = _mm256_mul_ps(fscal,dy12);
2128             tz               = _mm256_mul_ps(fscal,dz12);
2129
2130             /* Update vectorial force */
2131             fix1             = _mm256_add_ps(fix1,tx);
2132             fiy1             = _mm256_add_ps(fiy1,ty);
2133             fiz1             = _mm256_add_ps(fiz1,tz);
2134
2135             fjx2             = _mm256_add_ps(fjx2,tx);
2136             fjy2             = _mm256_add_ps(fjy2,ty);
2137             fjz2             = _mm256_add_ps(fjz2,tz);
2138
2139             }
2140
2141             /**************************
2142              * CALCULATE INTERACTIONS *
2143              **************************/
2144
2145             if (gmx_mm256_any_lt(rsq13,rcutoff2))
2146             {
2147
2148             /* REACTION-FIELD ELECTROSTATICS */
2149             felec            = _mm256_mul_ps(qq13,_mm256_sub_ps(_mm256_mul_ps(rinv13,rinvsq13),krf2));
2150
2151             cutoff_mask      = _mm256_cmp_ps(rsq13,rcutoff2,_CMP_LT_OQ);
2152
2153             fscal            = felec;
2154
2155             fscal            = _mm256_and_ps(fscal,cutoff_mask);
2156
2157             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
2158
2159             /* Calculate temporary vectorial force */
2160             tx               = _mm256_mul_ps(fscal,dx13);
2161             ty               = _mm256_mul_ps(fscal,dy13);
2162             tz               = _mm256_mul_ps(fscal,dz13);
2163
2164             /* Update vectorial force */
2165             fix1             = _mm256_add_ps(fix1,tx);
2166             fiy1             = _mm256_add_ps(fiy1,ty);
2167             fiz1             = _mm256_add_ps(fiz1,tz);
2168
2169             fjx3             = _mm256_add_ps(fjx3,tx);
2170             fjy3             = _mm256_add_ps(fjy3,ty);
2171             fjz3             = _mm256_add_ps(fjz3,tz);
2172
2173             }
2174
2175             /**************************
2176              * CALCULATE INTERACTIONS *
2177              **************************/
2178
2179             if (gmx_mm256_any_lt(rsq21,rcutoff2))
2180             {
2181
2182             /* REACTION-FIELD ELECTROSTATICS */
2183             felec            = _mm256_mul_ps(qq21,_mm256_sub_ps(_mm256_mul_ps(rinv21,rinvsq21),krf2));
2184
2185             cutoff_mask      = _mm256_cmp_ps(rsq21,rcutoff2,_CMP_LT_OQ);
2186
2187             fscal            = felec;
2188
2189             fscal            = _mm256_and_ps(fscal,cutoff_mask);
2190
2191             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
2192
2193             /* Calculate temporary vectorial force */
2194             tx               = _mm256_mul_ps(fscal,dx21);
2195             ty               = _mm256_mul_ps(fscal,dy21);
2196             tz               = _mm256_mul_ps(fscal,dz21);
2197
2198             /* Update vectorial force */
2199             fix2             = _mm256_add_ps(fix2,tx);
2200             fiy2             = _mm256_add_ps(fiy2,ty);
2201             fiz2             = _mm256_add_ps(fiz2,tz);
2202
2203             fjx1             = _mm256_add_ps(fjx1,tx);
2204             fjy1             = _mm256_add_ps(fjy1,ty);
2205             fjz1             = _mm256_add_ps(fjz1,tz);
2206
2207             }
2208
2209             /**************************
2210              * CALCULATE INTERACTIONS *
2211              **************************/
2212
2213             if (gmx_mm256_any_lt(rsq22,rcutoff2))
2214             {
2215
2216             /* REACTION-FIELD ELECTROSTATICS */
2217             felec            = _mm256_mul_ps(qq22,_mm256_sub_ps(_mm256_mul_ps(rinv22,rinvsq22),krf2));
2218
2219             cutoff_mask      = _mm256_cmp_ps(rsq22,rcutoff2,_CMP_LT_OQ);
2220
2221             fscal            = felec;
2222
2223             fscal            = _mm256_and_ps(fscal,cutoff_mask);
2224
2225             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
2226
2227             /* Calculate temporary vectorial force */
2228             tx               = _mm256_mul_ps(fscal,dx22);
2229             ty               = _mm256_mul_ps(fscal,dy22);
2230             tz               = _mm256_mul_ps(fscal,dz22);
2231
2232             /* Update vectorial force */
2233             fix2             = _mm256_add_ps(fix2,tx);
2234             fiy2             = _mm256_add_ps(fiy2,ty);
2235             fiz2             = _mm256_add_ps(fiz2,tz);
2236
2237             fjx2             = _mm256_add_ps(fjx2,tx);
2238             fjy2             = _mm256_add_ps(fjy2,ty);
2239             fjz2             = _mm256_add_ps(fjz2,tz);
2240
2241             }
2242
2243             /**************************
2244              * CALCULATE INTERACTIONS *
2245              **************************/
2246
2247             if (gmx_mm256_any_lt(rsq23,rcutoff2))
2248             {
2249
2250             /* REACTION-FIELD ELECTROSTATICS */
2251             felec            = _mm256_mul_ps(qq23,_mm256_sub_ps(_mm256_mul_ps(rinv23,rinvsq23),krf2));
2252
2253             cutoff_mask      = _mm256_cmp_ps(rsq23,rcutoff2,_CMP_LT_OQ);
2254
2255             fscal            = felec;
2256
2257             fscal            = _mm256_and_ps(fscal,cutoff_mask);
2258
2259             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
2260
2261             /* Calculate temporary vectorial force */
2262             tx               = _mm256_mul_ps(fscal,dx23);
2263             ty               = _mm256_mul_ps(fscal,dy23);
2264             tz               = _mm256_mul_ps(fscal,dz23);
2265
2266             /* Update vectorial force */
2267             fix2             = _mm256_add_ps(fix2,tx);
2268             fiy2             = _mm256_add_ps(fiy2,ty);
2269             fiz2             = _mm256_add_ps(fiz2,tz);
2270
2271             fjx3             = _mm256_add_ps(fjx3,tx);
2272             fjy3             = _mm256_add_ps(fjy3,ty);
2273             fjz3             = _mm256_add_ps(fjz3,tz);
2274
2275             }
2276
2277             /**************************
2278              * CALCULATE INTERACTIONS *
2279              **************************/
2280
2281             if (gmx_mm256_any_lt(rsq31,rcutoff2))
2282             {
2283
2284             /* REACTION-FIELD ELECTROSTATICS */
2285             felec            = _mm256_mul_ps(qq31,_mm256_sub_ps(_mm256_mul_ps(rinv31,rinvsq31),krf2));
2286
2287             cutoff_mask      = _mm256_cmp_ps(rsq31,rcutoff2,_CMP_LT_OQ);
2288
2289             fscal            = felec;
2290
2291             fscal            = _mm256_and_ps(fscal,cutoff_mask);
2292
2293             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
2294
2295             /* Calculate temporary vectorial force */
2296             tx               = _mm256_mul_ps(fscal,dx31);
2297             ty               = _mm256_mul_ps(fscal,dy31);
2298             tz               = _mm256_mul_ps(fscal,dz31);
2299
2300             /* Update vectorial force */
2301             fix3             = _mm256_add_ps(fix3,tx);
2302             fiy3             = _mm256_add_ps(fiy3,ty);
2303             fiz3             = _mm256_add_ps(fiz3,tz);
2304
2305             fjx1             = _mm256_add_ps(fjx1,tx);
2306             fjy1             = _mm256_add_ps(fjy1,ty);
2307             fjz1             = _mm256_add_ps(fjz1,tz);
2308
2309             }
2310
2311             /**************************
2312              * CALCULATE INTERACTIONS *
2313              **************************/
2314
2315             if (gmx_mm256_any_lt(rsq32,rcutoff2))
2316             {
2317
2318             /* REACTION-FIELD ELECTROSTATICS */
2319             felec            = _mm256_mul_ps(qq32,_mm256_sub_ps(_mm256_mul_ps(rinv32,rinvsq32),krf2));
2320
2321             cutoff_mask      = _mm256_cmp_ps(rsq32,rcutoff2,_CMP_LT_OQ);
2322
2323             fscal            = felec;
2324
2325             fscal            = _mm256_and_ps(fscal,cutoff_mask);
2326
2327             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
2328
2329             /* Calculate temporary vectorial force */
2330             tx               = _mm256_mul_ps(fscal,dx32);
2331             ty               = _mm256_mul_ps(fscal,dy32);
2332             tz               = _mm256_mul_ps(fscal,dz32);
2333
2334             /* Update vectorial force */
2335             fix3             = _mm256_add_ps(fix3,tx);
2336             fiy3             = _mm256_add_ps(fiy3,ty);
2337             fiz3             = _mm256_add_ps(fiz3,tz);
2338
2339             fjx2             = _mm256_add_ps(fjx2,tx);
2340             fjy2             = _mm256_add_ps(fjy2,ty);
2341             fjz2             = _mm256_add_ps(fjz2,tz);
2342
2343             }
2344
2345             /**************************
2346              * CALCULATE INTERACTIONS *
2347              **************************/
2348
2349             if (gmx_mm256_any_lt(rsq33,rcutoff2))
2350             {
2351
2352             /* REACTION-FIELD ELECTROSTATICS */
2353             felec            = _mm256_mul_ps(qq33,_mm256_sub_ps(_mm256_mul_ps(rinv33,rinvsq33),krf2));
2354
2355             cutoff_mask      = _mm256_cmp_ps(rsq33,rcutoff2,_CMP_LT_OQ);
2356
2357             fscal            = felec;
2358
2359             fscal            = _mm256_and_ps(fscal,cutoff_mask);
2360
2361             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
2362
2363             /* Calculate temporary vectorial force */
2364             tx               = _mm256_mul_ps(fscal,dx33);
2365             ty               = _mm256_mul_ps(fscal,dy33);
2366             tz               = _mm256_mul_ps(fscal,dz33);
2367
2368             /* Update vectorial force */
2369             fix3             = _mm256_add_ps(fix3,tx);
2370             fiy3             = _mm256_add_ps(fiy3,ty);
2371             fiz3             = _mm256_add_ps(fiz3,tz);
2372
2373             fjx3             = _mm256_add_ps(fjx3,tx);
2374             fjy3             = _mm256_add_ps(fjy3,ty);
2375             fjz3             = _mm256_add_ps(fjz3,tz);
2376
2377             }
2378
2379             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
2380             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
2381             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
2382             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
2383             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
2384             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
2385             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
2386             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
2387
2388             gmx_mm256_decrement_4rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
2389                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
2390                                                       fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
2391
2392             /* Inner loop uses 303 flops */
2393         }
2394
2395         /* End of innermost loop */
2396
2397         gmx_mm256_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
2398                                                  f+i_coord_offset,fshift+i_shift_offset);
2399
2400         /* Increment number of inner iterations */
2401         inneriter                  += j_index_end - j_index_start;
2402
2403         /* Outer loop uses 24 flops */
2404     }
2405
2406     /* Increment number of outer iterations */
2407     outeriter        += nri;
2408
2409     /* Update outer/inner flops */
2410
2411     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W4W4_F,outeriter*24 + inneriter*303);
2412 }