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