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