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