Compile nonbonded kernels as C++
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_double / nb_kernel_ElecCoul_VdwLJ_GeomW3W3_avx_256_double.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2017,2018, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS avx_256_double kernel generator.
37  */
38 #include "gmxpre.h"
39
40 #include "config.h"
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "gromacs/gmxlib/nrnb.h"
46
47 #include "kernelutil_x86_avx_256_double.h"
48
49 /*
50  * Gromacs nonbonded kernel:   nb_kernel_ElecCoul_VdwLJ_GeomW3W3_VF_avx_256_double
51  * Electrostatics interaction: Coulomb
52  * VdW interaction:            LennardJones
53  * Geometry:                   Water3-Water3
54  * Calculate force/pot:        PotentialAndForce
55  */
56 void
57 nb_kernel_ElecCoul_VdwLJ_GeomW3W3_VF_avx_256_double
58                     (t_nblist                    * gmx_restrict       nlist,
59                      rvec                        * gmx_restrict          xx,
60                      rvec                        * gmx_restrict          ff,
61                      struct t_forcerec           * gmx_restrict          fr,
62                      t_mdatoms                   * gmx_restrict     mdatoms,
63                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
64                      t_nrnb                      * gmx_restrict        nrnb)
65 {
66     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
67      * just 0 for non-waters.
68      * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
69      * jnr indices corresponding to data put in the four positions in the SIMD register.
70      */
71     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
72     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
73     int              jnrA,jnrB,jnrC,jnrD;
74     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
75     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
76     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
77     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
78     real             rcutoff_scalar;
79     real             *shiftvec,*fshift,*x,*f;
80     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
81     real             scratch[4*DIM];
82     __m256d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
83     real *           vdwioffsetptr0;
84     __m256d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
85     real *           vdwioffsetptr1;
86     __m256d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
87     real *           vdwioffsetptr2;
88     __m256d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
89     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
90     __m256d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
91     int              vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D;
92     __m256d          jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
93     int              vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D;
94     __m256d          jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
95     __m256d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
96     __m256d          dx01,dy01,dz01,rsq01,rinv01,rinvsq01,r01,qq01,c6_01,c12_01;
97     __m256d          dx02,dy02,dz02,rsq02,rinv02,rinvsq02,r02,qq02,c6_02,c12_02;
98     __m256d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
99     __m256d          dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
100     __m256d          dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
101     __m256d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
102     __m256d          dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
103     __m256d          dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
104     __m256d          velec,felec,velecsum,facel,crf,krf,krf2;
105     real             *charge;
106     int              nvdwtype;
107     __m256d          rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
108     int              *vdwtype;
109     real             *vdwparam;
110     __m256d          one_sixth   = _mm256_set1_pd(1.0/6.0);
111     __m256d          one_twelfth = _mm256_set1_pd(1.0/12.0);
112     __m256d          dummy_mask,cutoff_mask;
113     __m128           tmpmask0,tmpmask1;
114     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
115     __m256d          one     = _mm256_set1_pd(1.0);
116     __m256d          two     = _mm256_set1_pd(2.0);
117     x                = xx[0];
118     f                = ff[0];
119
120     nri              = nlist->nri;
121     iinr             = nlist->iinr;
122     jindex           = nlist->jindex;
123     jjnr             = nlist->jjnr;
124     shiftidx         = nlist->shift;
125     gid              = nlist->gid;
126     shiftvec         = fr->shift_vec[0];
127     fshift           = fr->fshift[0];
128     facel            = _mm256_set1_pd(fr->ic->epsfac);
129     charge           = mdatoms->chargeA;
130     nvdwtype         = fr->ntype;
131     vdwparam         = fr->nbfp;
132     vdwtype          = mdatoms->typeA;
133
134     /* Setup water-specific parameters */
135     inr              = nlist->iinr[0];
136     iq0              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
137     iq1              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+1]));
138     iq2              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+2]));
139     vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
140
141     jq0              = _mm256_set1_pd(charge[inr+0]);
142     jq1              = _mm256_set1_pd(charge[inr+1]);
143     jq2              = _mm256_set1_pd(charge[inr+2]);
144     vdwjidx0A        = 2*vdwtype[inr+0];
145     qq00             = _mm256_mul_pd(iq0,jq0);
146     c6_00            = _mm256_set1_pd(vdwioffsetptr0[vdwjidx0A]);
147     c12_00           = _mm256_set1_pd(vdwioffsetptr0[vdwjidx0A+1]);
148     qq01             = _mm256_mul_pd(iq0,jq1);
149     qq02             = _mm256_mul_pd(iq0,jq2);
150     qq10             = _mm256_mul_pd(iq1,jq0);
151     qq11             = _mm256_mul_pd(iq1,jq1);
152     qq12             = _mm256_mul_pd(iq1,jq2);
153     qq20             = _mm256_mul_pd(iq2,jq0);
154     qq21             = _mm256_mul_pd(iq2,jq1);
155     qq22             = _mm256_mul_pd(iq2,jq2);
156
157     /* Avoid stupid compiler warnings */
158     jnrA = jnrB = jnrC = jnrD = 0;
159     j_coord_offsetA = 0;
160     j_coord_offsetB = 0;
161     j_coord_offsetC = 0;
162     j_coord_offsetD = 0;
163
164     outeriter        = 0;
165     inneriter        = 0;
166
167     for(iidx=0;iidx<4*DIM;iidx++)
168     {
169         scratch[iidx] = 0.0;
170     }
171
172     /* Start outer loop over neighborlists */
173     for(iidx=0; iidx<nri; iidx++)
174     {
175         /* Load shift vector for this list */
176         i_shift_offset   = DIM*shiftidx[iidx];
177
178         /* Load limits for loop over neighbors */
179         j_index_start    = jindex[iidx];
180         j_index_end      = jindex[iidx+1];
181
182         /* Get outer coordinate index */
183         inr              = iinr[iidx];
184         i_coord_offset   = DIM*inr;
185
186         /* Load i particle coords and add shift vector */
187         gmx_mm256_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
188                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
189
190         fix0             = _mm256_setzero_pd();
191         fiy0             = _mm256_setzero_pd();
192         fiz0             = _mm256_setzero_pd();
193         fix1             = _mm256_setzero_pd();
194         fiy1             = _mm256_setzero_pd();
195         fiz1             = _mm256_setzero_pd();
196         fix2             = _mm256_setzero_pd();
197         fiy2             = _mm256_setzero_pd();
198         fiz2             = _mm256_setzero_pd();
199
200         /* Reset potential sums */
201         velecsum         = _mm256_setzero_pd();
202         vvdwsum          = _mm256_setzero_pd();
203
204         /* Start inner kernel loop */
205         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
206         {
207
208             /* Get j neighbor index, and coordinate index */
209             jnrA             = jjnr[jidx];
210             jnrB             = jjnr[jidx+1];
211             jnrC             = jjnr[jidx+2];
212             jnrD             = jjnr[jidx+3];
213             j_coord_offsetA  = DIM*jnrA;
214             j_coord_offsetB  = DIM*jnrB;
215             j_coord_offsetC  = DIM*jnrC;
216             j_coord_offsetD  = DIM*jnrD;
217
218             /* load j atom coordinates */
219             gmx_mm256_load_3rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
220                                                  x+j_coord_offsetC,x+j_coord_offsetD,
221                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
222
223             /* Calculate displacement vector */
224             dx00             = _mm256_sub_pd(ix0,jx0);
225             dy00             = _mm256_sub_pd(iy0,jy0);
226             dz00             = _mm256_sub_pd(iz0,jz0);
227             dx01             = _mm256_sub_pd(ix0,jx1);
228             dy01             = _mm256_sub_pd(iy0,jy1);
229             dz01             = _mm256_sub_pd(iz0,jz1);
230             dx02             = _mm256_sub_pd(ix0,jx2);
231             dy02             = _mm256_sub_pd(iy0,jy2);
232             dz02             = _mm256_sub_pd(iz0,jz2);
233             dx10             = _mm256_sub_pd(ix1,jx0);
234             dy10             = _mm256_sub_pd(iy1,jy0);
235             dz10             = _mm256_sub_pd(iz1,jz0);
236             dx11             = _mm256_sub_pd(ix1,jx1);
237             dy11             = _mm256_sub_pd(iy1,jy1);
238             dz11             = _mm256_sub_pd(iz1,jz1);
239             dx12             = _mm256_sub_pd(ix1,jx2);
240             dy12             = _mm256_sub_pd(iy1,jy2);
241             dz12             = _mm256_sub_pd(iz1,jz2);
242             dx20             = _mm256_sub_pd(ix2,jx0);
243             dy20             = _mm256_sub_pd(iy2,jy0);
244             dz20             = _mm256_sub_pd(iz2,jz0);
245             dx21             = _mm256_sub_pd(ix2,jx1);
246             dy21             = _mm256_sub_pd(iy2,jy1);
247             dz21             = _mm256_sub_pd(iz2,jz1);
248             dx22             = _mm256_sub_pd(ix2,jx2);
249             dy22             = _mm256_sub_pd(iy2,jy2);
250             dz22             = _mm256_sub_pd(iz2,jz2);
251
252             /* Calculate squared distance and things based on it */
253             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
254             rsq01            = gmx_mm256_calc_rsq_pd(dx01,dy01,dz01);
255             rsq02            = gmx_mm256_calc_rsq_pd(dx02,dy02,dz02);
256             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
257             rsq11            = gmx_mm256_calc_rsq_pd(dx11,dy11,dz11);
258             rsq12            = gmx_mm256_calc_rsq_pd(dx12,dy12,dz12);
259             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
260             rsq21            = gmx_mm256_calc_rsq_pd(dx21,dy21,dz21);
261             rsq22            = gmx_mm256_calc_rsq_pd(dx22,dy22,dz22);
262
263             rinv00           = avx256_invsqrt_d(rsq00);
264             rinv01           = avx256_invsqrt_d(rsq01);
265             rinv02           = avx256_invsqrt_d(rsq02);
266             rinv10           = avx256_invsqrt_d(rsq10);
267             rinv11           = avx256_invsqrt_d(rsq11);
268             rinv12           = avx256_invsqrt_d(rsq12);
269             rinv20           = avx256_invsqrt_d(rsq20);
270             rinv21           = avx256_invsqrt_d(rsq21);
271             rinv22           = avx256_invsqrt_d(rsq22);
272
273             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
274             rinvsq01         = _mm256_mul_pd(rinv01,rinv01);
275             rinvsq02         = _mm256_mul_pd(rinv02,rinv02);
276             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
277             rinvsq11         = _mm256_mul_pd(rinv11,rinv11);
278             rinvsq12         = _mm256_mul_pd(rinv12,rinv12);
279             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
280             rinvsq21         = _mm256_mul_pd(rinv21,rinv21);
281             rinvsq22         = _mm256_mul_pd(rinv22,rinv22);
282
283             fjx0             = _mm256_setzero_pd();
284             fjy0             = _mm256_setzero_pd();
285             fjz0             = _mm256_setzero_pd();
286             fjx1             = _mm256_setzero_pd();
287             fjy1             = _mm256_setzero_pd();
288             fjz1             = _mm256_setzero_pd();
289             fjx2             = _mm256_setzero_pd();
290             fjy2             = _mm256_setzero_pd();
291             fjz2             = _mm256_setzero_pd();
292
293             /**************************
294              * CALCULATE INTERACTIONS *
295              **************************/
296
297             /* COULOMB ELECTROSTATICS */
298             velec            = _mm256_mul_pd(qq00,rinv00);
299             felec            = _mm256_mul_pd(velec,rinvsq00);
300
301             /* LENNARD-JONES DISPERSION/REPULSION */
302
303             rinvsix          = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
304             vvdw6            = _mm256_mul_pd(c6_00,rinvsix);
305             vvdw12           = _mm256_mul_pd(c12_00,_mm256_mul_pd(rinvsix,rinvsix));
306             vvdw             = _mm256_sub_pd( _mm256_mul_pd(vvdw12,one_twelfth) , _mm256_mul_pd(vvdw6,one_sixth) );
307             fvdw             = _mm256_mul_pd(_mm256_sub_pd(vvdw12,vvdw6),rinvsq00);
308
309             /* Update potential sum for this i atom from the interaction with this j atom. */
310             velecsum         = _mm256_add_pd(velecsum,velec);
311             vvdwsum          = _mm256_add_pd(vvdwsum,vvdw);
312
313             fscal            = _mm256_add_pd(felec,fvdw);
314
315             /* Calculate temporary vectorial force */
316             tx               = _mm256_mul_pd(fscal,dx00);
317             ty               = _mm256_mul_pd(fscal,dy00);
318             tz               = _mm256_mul_pd(fscal,dz00);
319
320             /* Update vectorial force */
321             fix0             = _mm256_add_pd(fix0,tx);
322             fiy0             = _mm256_add_pd(fiy0,ty);
323             fiz0             = _mm256_add_pd(fiz0,tz);
324
325             fjx0             = _mm256_add_pd(fjx0,tx);
326             fjy0             = _mm256_add_pd(fjy0,ty);
327             fjz0             = _mm256_add_pd(fjz0,tz);
328
329             /**************************
330              * CALCULATE INTERACTIONS *
331              **************************/
332
333             /* COULOMB ELECTROSTATICS */
334             velec            = _mm256_mul_pd(qq01,rinv01);
335             felec            = _mm256_mul_pd(velec,rinvsq01);
336
337             /* Update potential sum for this i atom from the interaction with this j atom. */
338             velecsum         = _mm256_add_pd(velecsum,velec);
339
340             fscal            = felec;
341
342             /* Calculate temporary vectorial force */
343             tx               = _mm256_mul_pd(fscal,dx01);
344             ty               = _mm256_mul_pd(fscal,dy01);
345             tz               = _mm256_mul_pd(fscal,dz01);
346
347             /* Update vectorial force */
348             fix0             = _mm256_add_pd(fix0,tx);
349             fiy0             = _mm256_add_pd(fiy0,ty);
350             fiz0             = _mm256_add_pd(fiz0,tz);
351
352             fjx1             = _mm256_add_pd(fjx1,tx);
353             fjy1             = _mm256_add_pd(fjy1,ty);
354             fjz1             = _mm256_add_pd(fjz1,tz);
355
356             /**************************
357              * CALCULATE INTERACTIONS *
358              **************************/
359
360             /* COULOMB ELECTROSTATICS */
361             velec            = _mm256_mul_pd(qq02,rinv02);
362             felec            = _mm256_mul_pd(velec,rinvsq02);
363
364             /* Update potential sum for this i atom from the interaction with this j atom. */
365             velecsum         = _mm256_add_pd(velecsum,velec);
366
367             fscal            = felec;
368
369             /* Calculate temporary vectorial force */
370             tx               = _mm256_mul_pd(fscal,dx02);
371             ty               = _mm256_mul_pd(fscal,dy02);
372             tz               = _mm256_mul_pd(fscal,dz02);
373
374             /* Update vectorial force */
375             fix0             = _mm256_add_pd(fix0,tx);
376             fiy0             = _mm256_add_pd(fiy0,ty);
377             fiz0             = _mm256_add_pd(fiz0,tz);
378
379             fjx2             = _mm256_add_pd(fjx2,tx);
380             fjy2             = _mm256_add_pd(fjy2,ty);
381             fjz2             = _mm256_add_pd(fjz2,tz);
382
383             /**************************
384              * CALCULATE INTERACTIONS *
385              **************************/
386
387             /* COULOMB ELECTROSTATICS */
388             velec            = _mm256_mul_pd(qq10,rinv10);
389             felec            = _mm256_mul_pd(velec,rinvsq10);
390
391             /* Update potential sum for this i atom from the interaction with this j atom. */
392             velecsum         = _mm256_add_pd(velecsum,velec);
393
394             fscal            = felec;
395
396             /* Calculate temporary vectorial force */
397             tx               = _mm256_mul_pd(fscal,dx10);
398             ty               = _mm256_mul_pd(fscal,dy10);
399             tz               = _mm256_mul_pd(fscal,dz10);
400
401             /* Update vectorial force */
402             fix1             = _mm256_add_pd(fix1,tx);
403             fiy1             = _mm256_add_pd(fiy1,ty);
404             fiz1             = _mm256_add_pd(fiz1,tz);
405
406             fjx0             = _mm256_add_pd(fjx0,tx);
407             fjy0             = _mm256_add_pd(fjy0,ty);
408             fjz0             = _mm256_add_pd(fjz0,tz);
409
410             /**************************
411              * CALCULATE INTERACTIONS *
412              **************************/
413
414             /* COULOMB ELECTROSTATICS */
415             velec            = _mm256_mul_pd(qq11,rinv11);
416             felec            = _mm256_mul_pd(velec,rinvsq11);
417
418             /* Update potential sum for this i atom from the interaction with this j atom. */
419             velecsum         = _mm256_add_pd(velecsum,velec);
420
421             fscal            = felec;
422
423             /* Calculate temporary vectorial force */
424             tx               = _mm256_mul_pd(fscal,dx11);
425             ty               = _mm256_mul_pd(fscal,dy11);
426             tz               = _mm256_mul_pd(fscal,dz11);
427
428             /* Update vectorial force */
429             fix1             = _mm256_add_pd(fix1,tx);
430             fiy1             = _mm256_add_pd(fiy1,ty);
431             fiz1             = _mm256_add_pd(fiz1,tz);
432
433             fjx1             = _mm256_add_pd(fjx1,tx);
434             fjy1             = _mm256_add_pd(fjy1,ty);
435             fjz1             = _mm256_add_pd(fjz1,tz);
436
437             /**************************
438              * CALCULATE INTERACTIONS *
439              **************************/
440
441             /* COULOMB ELECTROSTATICS */
442             velec            = _mm256_mul_pd(qq12,rinv12);
443             felec            = _mm256_mul_pd(velec,rinvsq12);
444
445             /* Update potential sum for this i atom from the interaction with this j atom. */
446             velecsum         = _mm256_add_pd(velecsum,velec);
447
448             fscal            = felec;
449
450             /* Calculate temporary vectorial force */
451             tx               = _mm256_mul_pd(fscal,dx12);
452             ty               = _mm256_mul_pd(fscal,dy12);
453             tz               = _mm256_mul_pd(fscal,dz12);
454
455             /* Update vectorial force */
456             fix1             = _mm256_add_pd(fix1,tx);
457             fiy1             = _mm256_add_pd(fiy1,ty);
458             fiz1             = _mm256_add_pd(fiz1,tz);
459
460             fjx2             = _mm256_add_pd(fjx2,tx);
461             fjy2             = _mm256_add_pd(fjy2,ty);
462             fjz2             = _mm256_add_pd(fjz2,tz);
463
464             /**************************
465              * CALCULATE INTERACTIONS *
466              **************************/
467
468             /* COULOMB ELECTROSTATICS */
469             velec            = _mm256_mul_pd(qq20,rinv20);
470             felec            = _mm256_mul_pd(velec,rinvsq20);
471
472             /* Update potential sum for this i atom from the interaction with this j atom. */
473             velecsum         = _mm256_add_pd(velecsum,velec);
474
475             fscal            = felec;
476
477             /* Calculate temporary vectorial force */
478             tx               = _mm256_mul_pd(fscal,dx20);
479             ty               = _mm256_mul_pd(fscal,dy20);
480             tz               = _mm256_mul_pd(fscal,dz20);
481
482             /* Update vectorial force */
483             fix2             = _mm256_add_pd(fix2,tx);
484             fiy2             = _mm256_add_pd(fiy2,ty);
485             fiz2             = _mm256_add_pd(fiz2,tz);
486
487             fjx0             = _mm256_add_pd(fjx0,tx);
488             fjy0             = _mm256_add_pd(fjy0,ty);
489             fjz0             = _mm256_add_pd(fjz0,tz);
490
491             /**************************
492              * CALCULATE INTERACTIONS *
493              **************************/
494
495             /* COULOMB ELECTROSTATICS */
496             velec            = _mm256_mul_pd(qq21,rinv21);
497             felec            = _mm256_mul_pd(velec,rinvsq21);
498
499             /* Update potential sum for this i atom from the interaction with this j atom. */
500             velecsum         = _mm256_add_pd(velecsum,velec);
501
502             fscal            = felec;
503
504             /* Calculate temporary vectorial force */
505             tx               = _mm256_mul_pd(fscal,dx21);
506             ty               = _mm256_mul_pd(fscal,dy21);
507             tz               = _mm256_mul_pd(fscal,dz21);
508
509             /* Update vectorial force */
510             fix2             = _mm256_add_pd(fix2,tx);
511             fiy2             = _mm256_add_pd(fiy2,ty);
512             fiz2             = _mm256_add_pd(fiz2,tz);
513
514             fjx1             = _mm256_add_pd(fjx1,tx);
515             fjy1             = _mm256_add_pd(fjy1,ty);
516             fjz1             = _mm256_add_pd(fjz1,tz);
517
518             /**************************
519              * CALCULATE INTERACTIONS *
520              **************************/
521
522             /* COULOMB ELECTROSTATICS */
523             velec            = _mm256_mul_pd(qq22,rinv22);
524             felec            = _mm256_mul_pd(velec,rinvsq22);
525
526             /* Update potential sum for this i atom from the interaction with this j atom. */
527             velecsum         = _mm256_add_pd(velecsum,velec);
528
529             fscal            = felec;
530
531             /* Calculate temporary vectorial force */
532             tx               = _mm256_mul_pd(fscal,dx22);
533             ty               = _mm256_mul_pd(fscal,dy22);
534             tz               = _mm256_mul_pd(fscal,dz22);
535
536             /* Update vectorial force */
537             fix2             = _mm256_add_pd(fix2,tx);
538             fiy2             = _mm256_add_pd(fiy2,ty);
539             fiz2             = _mm256_add_pd(fiz2,tz);
540
541             fjx2             = _mm256_add_pd(fjx2,tx);
542             fjy2             = _mm256_add_pd(fjy2,ty);
543             fjz2             = _mm256_add_pd(fjz2,tz);
544
545             fjptrA             = f+j_coord_offsetA;
546             fjptrB             = f+j_coord_offsetB;
547             fjptrC             = f+j_coord_offsetC;
548             fjptrD             = f+j_coord_offsetD;
549
550             gmx_mm256_decrement_3rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
551                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
552
553             /* Inner loop uses 255 flops */
554         }
555
556         if(jidx<j_index_end)
557         {
558
559             /* Get j neighbor index, and coordinate index */
560             jnrlistA         = jjnr[jidx];
561             jnrlistB         = jjnr[jidx+1];
562             jnrlistC         = jjnr[jidx+2];
563             jnrlistD         = jjnr[jidx+3];
564             /* Sign of each element will be negative for non-real atoms.
565              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
566              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
567              */
568             tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
569
570             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
571             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
572             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
573
574             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
575             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
576             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
577             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
578             j_coord_offsetA  = DIM*jnrA;
579             j_coord_offsetB  = DIM*jnrB;
580             j_coord_offsetC  = DIM*jnrC;
581             j_coord_offsetD  = DIM*jnrD;
582
583             /* load j atom coordinates */
584             gmx_mm256_load_3rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
585                                                  x+j_coord_offsetC,x+j_coord_offsetD,
586                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
587
588             /* Calculate displacement vector */
589             dx00             = _mm256_sub_pd(ix0,jx0);
590             dy00             = _mm256_sub_pd(iy0,jy0);
591             dz00             = _mm256_sub_pd(iz0,jz0);
592             dx01             = _mm256_sub_pd(ix0,jx1);
593             dy01             = _mm256_sub_pd(iy0,jy1);
594             dz01             = _mm256_sub_pd(iz0,jz1);
595             dx02             = _mm256_sub_pd(ix0,jx2);
596             dy02             = _mm256_sub_pd(iy0,jy2);
597             dz02             = _mm256_sub_pd(iz0,jz2);
598             dx10             = _mm256_sub_pd(ix1,jx0);
599             dy10             = _mm256_sub_pd(iy1,jy0);
600             dz10             = _mm256_sub_pd(iz1,jz0);
601             dx11             = _mm256_sub_pd(ix1,jx1);
602             dy11             = _mm256_sub_pd(iy1,jy1);
603             dz11             = _mm256_sub_pd(iz1,jz1);
604             dx12             = _mm256_sub_pd(ix1,jx2);
605             dy12             = _mm256_sub_pd(iy1,jy2);
606             dz12             = _mm256_sub_pd(iz1,jz2);
607             dx20             = _mm256_sub_pd(ix2,jx0);
608             dy20             = _mm256_sub_pd(iy2,jy0);
609             dz20             = _mm256_sub_pd(iz2,jz0);
610             dx21             = _mm256_sub_pd(ix2,jx1);
611             dy21             = _mm256_sub_pd(iy2,jy1);
612             dz21             = _mm256_sub_pd(iz2,jz1);
613             dx22             = _mm256_sub_pd(ix2,jx2);
614             dy22             = _mm256_sub_pd(iy2,jy2);
615             dz22             = _mm256_sub_pd(iz2,jz2);
616
617             /* Calculate squared distance and things based on it */
618             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
619             rsq01            = gmx_mm256_calc_rsq_pd(dx01,dy01,dz01);
620             rsq02            = gmx_mm256_calc_rsq_pd(dx02,dy02,dz02);
621             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
622             rsq11            = gmx_mm256_calc_rsq_pd(dx11,dy11,dz11);
623             rsq12            = gmx_mm256_calc_rsq_pd(dx12,dy12,dz12);
624             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
625             rsq21            = gmx_mm256_calc_rsq_pd(dx21,dy21,dz21);
626             rsq22            = gmx_mm256_calc_rsq_pd(dx22,dy22,dz22);
627
628             rinv00           = avx256_invsqrt_d(rsq00);
629             rinv01           = avx256_invsqrt_d(rsq01);
630             rinv02           = avx256_invsqrt_d(rsq02);
631             rinv10           = avx256_invsqrt_d(rsq10);
632             rinv11           = avx256_invsqrt_d(rsq11);
633             rinv12           = avx256_invsqrt_d(rsq12);
634             rinv20           = avx256_invsqrt_d(rsq20);
635             rinv21           = avx256_invsqrt_d(rsq21);
636             rinv22           = avx256_invsqrt_d(rsq22);
637
638             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
639             rinvsq01         = _mm256_mul_pd(rinv01,rinv01);
640             rinvsq02         = _mm256_mul_pd(rinv02,rinv02);
641             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
642             rinvsq11         = _mm256_mul_pd(rinv11,rinv11);
643             rinvsq12         = _mm256_mul_pd(rinv12,rinv12);
644             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
645             rinvsq21         = _mm256_mul_pd(rinv21,rinv21);
646             rinvsq22         = _mm256_mul_pd(rinv22,rinv22);
647
648             fjx0             = _mm256_setzero_pd();
649             fjy0             = _mm256_setzero_pd();
650             fjz0             = _mm256_setzero_pd();
651             fjx1             = _mm256_setzero_pd();
652             fjy1             = _mm256_setzero_pd();
653             fjz1             = _mm256_setzero_pd();
654             fjx2             = _mm256_setzero_pd();
655             fjy2             = _mm256_setzero_pd();
656             fjz2             = _mm256_setzero_pd();
657
658             /**************************
659              * CALCULATE INTERACTIONS *
660              **************************/
661
662             /* COULOMB ELECTROSTATICS */
663             velec            = _mm256_mul_pd(qq00,rinv00);
664             felec            = _mm256_mul_pd(velec,rinvsq00);
665
666             /* LENNARD-JONES DISPERSION/REPULSION */
667
668             rinvsix          = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
669             vvdw6            = _mm256_mul_pd(c6_00,rinvsix);
670             vvdw12           = _mm256_mul_pd(c12_00,_mm256_mul_pd(rinvsix,rinvsix));
671             vvdw             = _mm256_sub_pd( _mm256_mul_pd(vvdw12,one_twelfth) , _mm256_mul_pd(vvdw6,one_sixth) );
672             fvdw             = _mm256_mul_pd(_mm256_sub_pd(vvdw12,vvdw6),rinvsq00);
673
674             /* Update potential sum for this i atom from the interaction with this j atom. */
675             velec            = _mm256_andnot_pd(dummy_mask,velec);
676             velecsum         = _mm256_add_pd(velecsum,velec);
677             vvdw             = _mm256_andnot_pd(dummy_mask,vvdw);
678             vvdwsum          = _mm256_add_pd(vvdwsum,vvdw);
679
680             fscal            = _mm256_add_pd(felec,fvdw);
681
682             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
683
684             /* Calculate temporary vectorial force */
685             tx               = _mm256_mul_pd(fscal,dx00);
686             ty               = _mm256_mul_pd(fscal,dy00);
687             tz               = _mm256_mul_pd(fscal,dz00);
688
689             /* Update vectorial force */
690             fix0             = _mm256_add_pd(fix0,tx);
691             fiy0             = _mm256_add_pd(fiy0,ty);
692             fiz0             = _mm256_add_pd(fiz0,tz);
693
694             fjx0             = _mm256_add_pd(fjx0,tx);
695             fjy0             = _mm256_add_pd(fjy0,ty);
696             fjz0             = _mm256_add_pd(fjz0,tz);
697
698             /**************************
699              * CALCULATE INTERACTIONS *
700              **************************/
701
702             /* COULOMB ELECTROSTATICS */
703             velec            = _mm256_mul_pd(qq01,rinv01);
704             felec            = _mm256_mul_pd(velec,rinvsq01);
705
706             /* Update potential sum for this i atom from the interaction with this j atom. */
707             velec            = _mm256_andnot_pd(dummy_mask,velec);
708             velecsum         = _mm256_add_pd(velecsum,velec);
709
710             fscal            = felec;
711
712             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
713
714             /* Calculate temporary vectorial force */
715             tx               = _mm256_mul_pd(fscal,dx01);
716             ty               = _mm256_mul_pd(fscal,dy01);
717             tz               = _mm256_mul_pd(fscal,dz01);
718
719             /* Update vectorial force */
720             fix0             = _mm256_add_pd(fix0,tx);
721             fiy0             = _mm256_add_pd(fiy0,ty);
722             fiz0             = _mm256_add_pd(fiz0,tz);
723
724             fjx1             = _mm256_add_pd(fjx1,tx);
725             fjy1             = _mm256_add_pd(fjy1,ty);
726             fjz1             = _mm256_add_pd(fjz1,tz);
727
728             /**************************
729              * CALCULATE INTERACTIONS *
730              **************************/
731
732             /* COULOMB ELECTROSTATICS */
733             velec            = _mm256_mul_pd(qq02,rinv02);
734             felec            = _mm256_mul_pd(velec,rinvsq02);
735
736             /* Update potential sum for this i atom from the interaction with this j atom. */
737             velec            = _mm256_andnot_pd(dummy_mask,velec);
738             velecsum         = _mm256_add_pd(velecsum,velec);
739
740             fscal            = felec;
741
742             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
743
744             /* Calculate temporary vectorial force */
745             tx               = _mm256_mul_pd(fscal,dx02);
746             ty               = _mm256_mul_pd(fscal,dy02);
747             tz               = _mm256_mul_pd(fscal,dz02);
748
749             /* Update vectorial force */
750             fix0             = _mm256_add_pd(fix0,tx);
751             fiy0             = _mm256_add_pd(fiy0,ty);
752             fiz0             = _mm256_add_pd(fiz0,tz);
753
754             fjx2             = _mm256_add_pd(fjx2,tx);
755             fjy2             = _mm256_add_pd(fjy2,ty);
756             fjz2             = _mm256_add_pd(fjz2,tz);
757
758             /**************************
759              * CALCULATE INTERACTIONS *
760              **************************/
761
762             /* COULOMB ELECTROSTATICS */
763             velec            = _mm256_mul_pd(qq10,rinv10);
764             felec            = _mm256_mul_pd(velec,rinvsq10);
765
766             /* Update potential sum for this i atom from the interaction with this j atom. */
767             velec            = _mm256_andnot_pd(dummy_mask,velec);
768             velecsum         = _mm256_add_pd(velecsum,velec);
769
770             fscal            = felec;
771
772             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
773
774             /* Calculate temporary vectorial force */
775             tx               = _mm256_mul_pd(fscal,dx10);
776             ty               = _mm256_mul_pd(fscal,dy10);
777             tz               = _mm256_mul_pd(fscal,dz10);
778
779             /* Update vectorial force */
780             fix1             = _mm256_add_pd(fix1,tx);
781             fiy1             = _mm256_add_pd(fiy1,ty);
782             fiz1             = _mm256_add_pd(fiz1,tz);
783
784             fjx0             = _mm256_add_pd(fjx0,tx);
785             fjy0             = _mm256_add_pd(fjy0,ty);
786             fjz0             = _mm256_add_pd(fjz0,tz);
787
788             /**************************
789              * CALCULATE INTERACTIONS *
790              **************************/
791
792             /* COULOMB ELECTROSTATICS */
793             velec            = _mm256_mul_pd(qq11,rinv11);
794             felec            = _mm256_mul_pd(velec,rinvsq11);
795
796             /* Update potential sum for this i atom from the interaction with this j atom. */
797             velec            = _mm256_andnot_pd(dummy_mask,velec);
798             velecsum         = _mm256_add_pd(velecsum,velec);
799
800             fscal            = felec;
801
802             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
803
804             /* Calculate temporary vectorial force */
805             tx               = _mm256_mul_pd(fscal,dx11);
806             ty               = _mm256_mul_pd(fscal,dy11);
807             tz               = _mm256_mul_pd(fscal,dz11);
808
809             /* Update vectorial force */
810             fix1             = _mm256_add_pd(fix1,tx);
811             fiy1             = _mm256_add_pd(fiy1,ty);
812             fiz1             = _mm256_add_pd(fiz1,tz);
813
814             fjx1             = _mm256_add_pd(fjx1,tx);
815             fjy1             = _mm256_add_pd(fjy1,ty);
816             fjz1             = _mm256_add_pd(fjz1,tz);
817
818             /**************************
819              * CALCULATE INTERACTIONS *
820              **************************/
821
822             /* COULOMB ELECTROSTATICS */
823             velec            = _mm256_mul_pd(qq12,rinv12);
824             felec            = _mm256_mul_pd(velec,rinvsq12);
825
826             /* Update potential sum for this i atom from the interaction with this j atom. */
827             velec            = _mm256_andnot_pd(dummy_mask,velec);
828             velecsum         = _mm256_add_pd(velecsum,velec);
829
830             fscal            = felec;
831
832             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
833
834             /* Calculate temporary vectorial force */
835             tx               = _mm256_mul_pd(fscal,dx12);
836             ty               = _mm256_mul_pd(fscal,dy12);
837             tz               = _mm256_mul_pd(fscal,dz12);
838
839             /* Update vectorial force */
840             fix1             = _mm256_add_pd(fix1,tx);
841             fiy1             = _mm256_add_pd(fiy1,ty);
842             fiz1             = _mm256_add_pd(fiz1,tz);
843
844             fjx2             = _mm256_add_pd(fjx2,tx);
845             fjy2             = _mm256_add_pd(fjy2,ty);
846             fjz2             = _mm256_add_pd(fjz2,tz);
847
848             /**************************
849              * CALCULATE INTERACTIONS *
850              **************************/
851
852             /* COULOMB ELECTROSTATICS */
853             velec            = _mm256_mul_pd(qq20,rinv20);
854             felec            = _mm256_mul_pd(velec,rinvsq20);
855
856             /* Update potential sum for this i atom from the interaction with this j atom. */
857             velec            = _mm256_andnot_pd(dummy_mask,velec);
858             velecsum         = _mm256_add_pd(velecsum,velec);
859
860             fscal            = felec;
861
862             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
863
864             /* Calculate temporary vectorial force */
865             tx               = _mm256_mul_pd(fscal,dx20);
866             ty               = _mm256_mul_pd(fscal,dy20);
867             tz               = _mm256_mul_pd(fscal,dz20);
868
869             /* Update vectorial force */
870             fix2             = _mm256_add_pd(fix2,tx);
871             fiy2             = _mm256_add_pd(fiy2,ty);
872             fiz2             = _mm256_add_pd(fiz2,tz);
873
874             fjx0             = _mm256_add_pd(fjx0,tx);
875             fjy0             = _mm256_add_pd(fjy0,ty);
876             fjz0             = _mm256_add_pd(fjz0,tz);
877
878             /**************************
879              * CALCULATE INTERACTIONS *
880              **************************/
881
882             /* COULOMB ELECTROSTATICS */
883             velec            = _mm256_mul_pd(qq21,rinv21);
884             felec            = _mm256_mul_pd(velec,rinvsq21);
885
886             /* Update potential sum for this i atom from the interaction with this j atom. */
887             velec            = _mm256_andnot_pd(dummy_mask,velec);
888             velecsum         = _mm256_add_pd(velecsum,velec);
889
890             fscal            = felec;
891
892             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
893
894             /* Calculate temporary vectorial force */
895             tx               = _mm256_mul_pd(fscal,dx21);
896             ty               = _mm256_mul_pd(fscal,dy21);
897             tz               = _mm256_mul_pd(fscal,dz21);
898
899             /* Update vectorial force */
900             fix2             = _mm256_add_pd(fix2,tx);
901             fiy2             = _mm256_add_pd(fiy2,ty);
902             fiz2             = _mm256_add_pd(fiz2,tz);
903
904             fjx1             = _mm256_add_pd(fjx1,tx);
905             fjy1             = _mm256_add_pd(fjy1,ty);
906             fjz1             = _mm256_add_pd(fjz1,tz);
907
908             /**************************
909              * CALCULATE INTERACTIONS *
910              **************************/
911
912             /* COULOMB ELECTROSTATICS */
913             velec            = _mm256_mul_pd(qq22,rinv22);
914             felec            = _mm256_mul_pd(velec,rinvsq22);
915
916             /* Update potential sum for this i atom from the interaction with this j atom. */
917             velec            = _mm256_andnot_pd(dummy_mask,velec);
918             velecsum         = _mm256_add_pd(velecsum,velec);
919
920             fscal            = felec;
921
922             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
923
924             /* Calculate temporary vectorial force */
925             tx               = _mm256_mul_pd(fscal,dx22);
926             ty               = _mm256_mul_pd(fscal,dy22);
927             tz               = _mm256_mul_pd(fscal,dz22);
928
929             /* Update vectorial force */
930             fix2             = _mm256_add_pd(fix2,tx);
931             fiy2             = _mm256_add_pd(fiy2,ty);
932             fiz2             = _mm256_add_pd(fiz2,tz);
933
934             fjx2             = _mm256_add_pd(fjx2,tx);
935             fjy2             = _mm256_add_pd(fjy2,ty);
936             fjz2             = _mm256_add_pd(fjz2,tz);
937
938             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
939             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
940             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
941             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
942
943             gmx_mm256_decrement_3rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
944                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
945
946             /* Inner loop uses 255 flops */
947         }
948
949         /* End of innermost loop */
950
951         gmx_mm256_update_iforce_3atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
952                                                  f+i_coord_offset,fshift+i_shift_offset);
953
954         ggid                        = gid[iidx];
955         /* Update potential energies */
956         gmx_mm256_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
957         gmx_mm256_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
958
959         /* Increment number of inner iterations */
960         inneriter                  += j_index_end - j_index_start;
961
962         /* Outer loop uses 20 flops */
963     }
964
965     /* Increment number of outer iterations */
966     outeriter        += nri;
967
968     /* Update outer/inner flops */
969
970     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3W3_VF,outeriter*20 + inneriter*255);
971 }
972 /*
973  * Gromacs nonbonded kernel:   nb_kernel_ElecCoul_VdwLJ_GeomW3W3_F_avx_256_double
974  * Electrostatics interaction: Coulomb
975  * VdW interaction:            LennardJones
976  * Geometry:                   Water3-Water3
977  * Calculate force/pot:        Force
978  */
979 void
980 nb_kernel_ElecCoul_VdwLJ_GeomW3W3_F_avx_256_double
981                     (t_nblist                    * gmx_restrict       nlist,
982                      rvec                        * gmx_restrict          xx,
983                      rvec                        * gmx_restrict          ff,
984                      struct t_forcerec           * gmx_restrict          fr,
985                      t_mdatoms                   * gmx_restrict     mdatoms,
986                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
987                      t_nrnb                      * gmx_restrict        nrnb)
988 {
989     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
990      * just 0 for non-waters.
991      * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
992      * jnr indices corresponding to data put in the four positions in the SIMD register.
993      */
994     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
995     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
996     int              jnrA,jnrB,jnrC,jnrD;
997     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
998     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
999     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
1000     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
1001     real             rcutoff_scalar;
1002     real             *shiftvec,*fshift,*x,*f;
1003     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
1004     real             scratch[4*DIM];
1005     __m256d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
1006     real *           vdwioffsetptr0;
1007     __m256d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
1008     real *           vdwioffsetptr1;
1009     __m256d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
1010     real *           vdwioffsetptr2;
1011     __m256d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
1012     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
1013     __m256d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
1014     int              vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D;
1015     __m256d          jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
1016     int              vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D;
1017     __m256d          jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
1018     __m256d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
1019     __m256d          dx01,dy01,dz01,rsq01,rinv01,rinvsq01,r01,qq01,c6_01,c12_01;
1020     __m256d          dx02,dy02,dz02,rsq02,rinv02,rinvsq02,r02,qq02,c6_02,c12_02;
1021     __m256d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
1022     __m256d          dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
1023     __m256d          dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
1024     __m256d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
1025     __m256d          dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
1026     __m256d          dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
1027     __m256d          velec,felec,velecsum,facel,crf,krf,krf2;
1028     real             *charge;
1029     int              nvdwtype;
1030     __m256d          rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
1031     int              *vdwtype;
1032     real             *vdwparam;
1033     __m256d          one_sixth   = _mm256_set1_pd(1.0/6.0);
1034     __m256d          one_twelfth = _mm256_set1_pd(1.0/12.0);
1035     __m256d          dummy_mask,cutoff_mask;
1036     __m128           tmpmask0,tmpmask1;
1037     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
1038     __m256d          one     = _mm256_set1_pd(1.0);
1039     __m256d          two     = _mm256_set1_pd(2.0);
1040     x                = xx[0];
1041     f                = ff[0];
1042
1043     nri              = nlist->nri;
1044     iinr             = nlist->iinr;
1045     jindex           = nlist->jindex;
1046     jjnr             = nlist->jjnr;
1047     shiftidx         = nlist->shift;
1048     gid              = nlist->gid;
1049     shiftvec         = fr->shift_vec[0];
1050     fshift           = fr->fshift[0];
1051     facel            = _mm256_set1_pd(fr->ic->epsfac);
1052     charge           = mdatoms->chargeA;
1053     nvdwtype         = fr->ntype;
1054     vdwparam         = fr->nbfp;
1055     vdwtype          = mdatoms->typeA;
1056
1057     /* Setup water-specific parameters */
1058     inr              = nlist->iinr[0];
1059     iq0              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
1060     iq1              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+1]));
1061     iq2              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+2]));
1062     vdwioffsetptr0   = vdwparam+2*nvdwtype*vdwtype[inr+0];
1063
1064     jq0              = _mm256_set1_pd(charge[inr+0]);
1065     jq1              = _mm256_set1_pd(charge[inr+1]);
1066     jq2              = _mm256_set1_pd(charge[inr+2]);
1067     vdwjidx0A        = 2*vdwtype[inr+0];
1068     qq00             = _mm256_mul_pd(iq0,jq0);
1069     c6_00            = _mm256_set1_pd(vdwioffsetptr0[vdwjidx0A]);
1070     c12_00           = _mm256_set1_pd(vdwioffsetptr0[vdwjidx0A+1]);
1071     qq01             = _mm256_mul_pd(iq0,jq1);
1072     qq02             = _mm256_mul_pd(iq0,jq2);
1073     qq10             = _mm256_mul_pd(iq1,jq0);
1074     qq11             = _mm256_mul_pd(iq1,jq1);
1075     qq12             = _mm256_mul_pd(iq1,jq2);
1076     qq20             = _mm256_mul_pd(iq2,jq0);
1077     qq21             = _mm256_mul_pd(iq2,jq1);
1078     qq22             = _mm256_mul_pd(iq2,jq2);
1079
1080     /* Avoid stupid compiler warnings */
1081     jnrA = jnrB = jnrC = jnrD = 0;
1082     j_coord_offsetA = 0;
1083     j_coord_offsetB = 0;
1084     j_coord_offsetC = 0;
1085     j_coord_offsetD = 0;
1086
1087     outeriter        = 0;
1088     inneriter        = 0;
1089
1090     for(iidx=0;iidx<4*DIM;iidx++)
1091     {
1092         scratch[iidx] = 0.0;
1093     }
1094
1095     /* Start outer loop over neighborlists */
1096     for(iidx=0; iidx<nri; iidx++)
1097     {
1098         /* Load shift vector for this list */
1099         i_shift_offset   = DIM*shiftidx[iidx];
1100
1101         /* Load limits for loop over neighbors */
1102         j_index_start    = jindex[iidx];
1103         j_index_end      = jindex[iidx+1];
1104
1105         /* Get outer coordinate index */
1106         inr              = iinr[iidx];
1107         i_coord_offset   = DIM*inr;
1108
1109         /* Load i particle coords and add shift vector */
1110         gmx_mm256_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
1111                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
1112
1113         fix0             = _mm256_setzero_pd();
1114         fiy0             = _mm256_setzero_pd();
1115         fiz0             = _mm256_setzero_pd();
1116         fix1             = _mm256_setzero_pd();
1117         fiy1             = _mm256_setzero_pd();
1118         fiz1             = _mm256_setzero_pd();
1119         fix2             = _mm256_setzero_pd();
1120         fiy2             = _mm256_setzero_pd();
1121         fiz2             = _mm256_setzero_pd();
1122
1123         /* Start inner kernel loop */
1124         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
1125         {
1126
1127             /* Get j neighbor index, and coordinate index */
1128             jnrA             = jjnr[jidx];
1129             jnrB             = jjnr[jidx+1];
1130             jnrC             = jjnr[jidx+2];
1131             jnrD             = jjnr[jidx+3];
1132             j_coord_offsetA  = DIM*jnrA;
1133             j_coord_offsetB  = DIM*jnrB;
1134             j_coord_offsetC  = DIM*jnrC;
1135             j_coord_offsetD  = DIM*jnrD;
1136
1137             /* load j atom coordinates */
1138             gmx_mm256_load_3rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
1139                                                  x+j_coord_offsetC,x+j_coord_offsetD,
1140                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
1141
1142             /* Calculate displacement vector */
1143             dx00             = _mm256_sub_pd(ix0,jx0);
1144             dy00             = _mm256_sub_pd(iy0,jy0);
1145             dz00             = _mm256_sub_pd(iz0,jz0);
1146             dx01             = _mm256_sub_pd(ix0,jx1);
1147             dy01             = _mm256_sub_pd(iy0,jy1);
1148             dz01             = _mm256_sub_pd(iz0,jz1);
1149             dx02             = _mm256_sub_pd(ix0,jx2);
1150             dy02             = _mm256_sub_pd(iy0,jy2);
1151             dz02             = _mm256_sub_pd(iz0,jz2);
1152             dx10             = _mm256_sub_pd(ix1,jx0);
1153             dy10             = _mm256_sub_pd(iy1,jy0);
1154             dz10             = _mm256_sub_pd(iz1,jz0);
1155             dx11             = _mm256_sub_pd(ix1,jx1);
1156             dy11             = _mm256_sub_pd(iy1,jy1);
1157             dz11             = _mm256_sub_pd(iz1,jz1);
1158             dx12             = _mm256_sub_pd(ix1,jx2);
1159             dy12             = _mm256_sub_pd(iy1,jy2);
1160             dz12             = _mm256_sub_pd(iz1,jz2);
1161             dx20             = _mm256_sub_pd(ix2,jx0);
1162             dy20             = _mm256_sub_pd(iy2,jy0);
1163             dz20             = _mm256_sub_pd(iz2,jz0);
1164             dx21             = _mm256_sub_pd(ix2,jx1);
1165             dy21             = _mm256_sub_pd(iy2,jy1);
1166             dz21             = _mm256_sub_pd(iz2,jz1);
1167             dx22             = _mm256_sub_pd(ix2,jx2);
1168             dy22             = _mm256_sub_pd(iy2,jy2);
1169             dz22             = _mm256_sub_pd(iz2,jz2);
1170
1171             /* Calculate squared distance and things based on it */
1172             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
1173             rsq01            = gmx_mm256_calc_rsq_pd(dx01,dy01,dz01);
1174             rsq02            = gmx_mm256_calc_rsq_pd(dx02,dy02,dz02);
1175             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
1176             rsq11            = gmx_mm256_calc_rsq_pd(dx11,dy11,dz11);
1177             rsq12            = gmx_mm256_calc_rsq_pd(dx12,dy12,dz12);
1178             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
1179             rsq21            = gmx_mm256_calc_rsq_pd(dx21,dy21,dz21);
1180             rsq22            = gmx_mm256_calc_rsq_pd(dx22,dy22,dz22);
1181
1182             rinv00           = avx256_invsqrt_d(rsq00);
1183             rinv01           = avx256_invsqrt_d(rsq01);
1184             rinv02           = avx256_invsqrt_d(rsq02);
1185             rinv10           = avx256_invsqrt_d(rsq10);
1186             rinv11           = avx256_invsqrt_d(rsq11);
1187             rinv12           = avx256_invsqrt_d(rsq12);
1188             rinv20           = avx256_invsqrt_d(rsq20);
1189             rinv21           = avx256_invsqrt_d(rsq21);
1190             rinv22           = avx256_invsqrt_d(rsq22);
1191
1192             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
1193             rinvsq01         = _mm256_mul_pd(rinv01,rinv01);
1194             rinvsq02         = _mm256_mul_pd(rinv02,rinv02);
1195             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
1196             rinvsq11         = _mm256_mul_pd(rinv11,rinv11);
1197             rinvsq12         = _mm256_mul_pd(rinv12,rinv12);
1198             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
1199             rinvsq21         = _mm256_mul_pd(rinv21,rinv21);
1200             rinvsq22         = _mm256_mul_pd(rinv22,rinv22);
1201
1202             fjx0             = _mm256_setzero_pd();
1203             fjy0             = _mm256_setzero_pd();
1204             fjz0             = _mm256_setzero_pd();
1205             fjx1             = _mm256_setzero_pd();
1206             fjy1             = _mm256_setzero_pd();
1207             fjz1             = _mm256_setzero_pd();
1208             fjx2             = _mm256_setzero_pd();
1209             fjy2             = _mm256_setzero_pd();
1210             fjz2             = _mm256_setzero_pd();
1211
1212             /**************************
1213              * CALCULATE INTERACTIONS *
1214              **************************/
1215
1216             /* COULOMB ELECTROSTATICS */
1217             velec            = _mm256_mul_pd(qq00,rinv00);
1218             felec            = _mm256_mul_pd(velec,rinvsq00);
1219
1220             /* LENNARD-JONES DISPERSION/REPULSION */
1221
1222             rinvsix          = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
1223             fvdw             = _mm256_mul_pd(_mm256_sub_pd(_mm256_mul_pd(c12_00,rinvsix),c6_00),_mm256_mul_pd(rinvsix,rinvsq00));
1224
1225             fscal            = _mm256_add_pd(felec,fvdw);
1226
1227             /* Calculate temporary vectorial force */
1228             tx               = _mm256_mul_pd(fscal,dx00);
1229             ty               = _mm256_mul_pd(fscal,dy00);
1230             tz               = _mm256_mul_pd(fscal,dz00);
1231
1232             /* Update vectorial force */
1233             fix0             = _mm256_add_pd(fix0,tx);
1234             fiy0             = _mm256_add_pd(fiy0,ty);
1235             fiz0             = _mm256_add_pd(fiz0,tz);
1236
1237             fjx0             = _mm256_add_pd(fjx0,tx);
1238             fjy0             = _mm256_add_pd(fjy0,ty);
1239             fjz0             = _mm256_add_pd(fjz0,tz);
1240
1241             /**************************
1242              * CALCULATE INTERACTIONS *
1243              **************************/
1244
1245             /* COULOMB ELECTROSTATICS */
1246             velec            = _mm256_mul_pd(qq01,rinv01);
1247             felec            = _mm256_mul_pd(velec,rinvsq01);
1248
1249             fscal            = felec;
1250
1251             /* Calculate temporary vectorial force */
1252             tx               = _mm256_mul_pd(fscal,dx01);
1253             ty               = _mm256_mul_pd(fscal,dy01);
1254             tz               = _mm256_mul_pd(fscal,dz01);
1255
1256             /* Update vectorial force */
1257             fix0             = _mm256_add_pd(fix0,tx);
1258             fiy0             = _mm256_add_pd(fiy0,ty);
1259             fiz0             = _mm256_add_pd(fiz0,tz);
1260
1261             fjx1             = _mm256_add_pd(fjx1,tx);
1262             fjy1             = _mm256_add_pd(fjy1,ty);
1263             fjz1             = _mm256_add_pd(fjz1,tz);
1264
1265             /**************************
1266              * CALCULATE INTERACTIONS *
1267              **************************/
1268
1269             /* COULOMB ELECTROSTATICS */
1270             velec            = _mm256_mul_pd(qq02,rinv02);
1271             felec            = _mm256_mul_pd(velec,rinvsq02);
1272
1273             fscal            = felec;
1274
1275             /* Calculate temporary vectorial force */
1276             tx               = _mm256_mul_pd(fscal,dx02);
1277             ty               = _mm256_mul_pd(fscal,dy02);
1278             tz               = _mm256_mul_pd(fscal,dz02);
1279
1280             /* Update vectorial force */
1281             fix0             = _mm256_add_pd(fix0,tx);
1282             fiy0             = _mm256_add_pd(fiy0,ty);
1283             fiz0             = _mm256_add_pd(fiz0,tz);
1284
1285             fjx2             = _mm256_add_pd(fjx2,tx);
1286             fjy2             = _mm256_add_pd(fjy2,ty);
1287             fjz2             = _mm256_add_pd(fjz2,tz);
1288
1289             /**************************
1290              * CALCULATE INTERACTIONS *
1291              **************************/
1292
1293             /* COULOMB ELECTROSTATICS */
1294             velec            = _mm256_mul_pd(qq10,rinv10);
1295             felec            = _mm256_mul_pd(velec,rinvsq10);
1296
1297             fscal            = felec;
1298
1299             /* Calculate temporary vectorial force */
1300             tx               = _mm256_mul_pd(fscal,dx10);
1301             ty               = _mm256_mul_pd(fscal,dy10);
1302             tz               = _mm256_mul_pd(fscal,dz10);
1303
1304             /* Update vectorial force */
1305             fix1             = _mm256_add_pd(fix1,tx);
1306             fiy1             = _mm256_add_pd(fiy1,ty);
1307             fiz1             = _mm256_add_pd(fiz1,tz);
1308
1309             fjx0             = _mm256_add_pd(fjx0,tx);
1310             fjy0             = _mm256_add_pd(fjy0,ty);
1311             fjz0             = _mm256_add_pd(fjz0,tz);
1312
1313             /**************************
1314              * CALCULATE INTERACTIONS *
1315              **************************/
1316
1317             /* COULOMB ELECTROSTATICS */
1318             velec            = _mm256_mul_pd(qq11,rinv11);
1319             felec            = _mm256_mul_pd(velec,rinvsq11);
1320
1321             fscal            = felec;
1322
1323             /* Calculate temporary vectorial force */
1324             tx               = _mm256_mul_pd(fscal,dx11);
1325             ty               = _mm256_mul_pd(fscal,dy11);
1326             tz               = _mm256_mul_pd(fscal,dz11);
1327
1328             /* Update vectorial force */
1329             fix1             = _mm256_add_pd(fix1,tx);
1330             fiy1             = _mm256_add_pd(fiy1,ty);
1331             fiz1             = _mm256_add_pd(fiz1,tz);
1332
1333             fjx1             = _mm256_add_pd(fjx1,tx);
1334             fjy1             = _mm256_add_pd(fjy1,ty);
1335             fjz1             = _mm256_add_pd(fjz1,tz);
1336
1337             /**************************
1338              * CALCULATE INTERACTIONS *
1339              **************************/
1340
1341             /* COULOMB ELECTROSTATICS */
1342             velec            = _mm256_mul_pd(qq12,rinv12);
1343             felec            = _mm256_mul_pd(velec,rinvsq12);
1344
1345             fscal            = felec;
1346
1347             /* Calculate temporary vectorial force */
1348             tx               = _mm256_mul_pd(fscal,dx12);
1349             ty               = _mm256_mul_pd(fscal,dy12);
1350             tz               = _mm256_mul_pd(fscal,dz12);
1351
1352             /* Update vectorial force */
1353             fix1             = _mm256_add_pd(fix1,tx);
1354             fiy1             = _mm256_add_pd(fiy1,ty);
1355             fiz1             = _mm256_add_pd(fiz1,tz);
1356
1357             fjx2             = _mm256_add_pd(fjx2,tx);
1358             fjy2             = _mm256_add_pd(fjy2,ty);
1359             fjz2             = _mm256_add_pd(fjz2,tz);
1360
1361             /**************************
1362              * CALCULATE INTERACTIONS *
1363              **************************/
1364
1365             /* COULOMB ELECTROSTATICS */
1366             velec            = _mm256_mul_pd(qq20,rinv20);
1367             felec            = _mm256_mul_pd(velec,rinvsq20);
1368
1369             fscal            = felec;
1370
1371             /* Calculate temporary vectorial force */
1372             tx               = _mm256_mul_pd(fscal,dx20);
1373             ty               = _mm256_mul_pd(fscal,dy20);
1374             tz               = _mm256_mul_pd(fscal,dz20);
1375
1376             /* Update vectorial force */
1377             fix2             = _mm256_add_pd(fix2,tx);
1378             fiy2             = _mm256_add_pd(fiy2,ty);
1379             fiz2             = _mm256_add_pd(fiz2,tz);
1380
1381             fjx0             = _mm256_add_pd(fjx0,tx);
1382             fjy0             = _mm256_add_pd(fjy0,ty);
1383             fjz0             = _mm256_add_pd(fjz0,tz);
1384
1385             /**************************
1386              * CALCULATE INTERACTIONS *
1387              **************************/
1388
1389             /* COULOMB ELECTROSTATICS */
1390             velec            = _mm256_mul_pd(qq21,rinv21);
1391             felec            = _mm256_mul_pd(velec,rinvsq21);
1392
1393             fscal            = felec;
1394
1395             /* Calculate temporary vectorial force */
1396             tx               = _mm256_mul_pd(fscal,dx21);
1397             ty               = _mm256_mul_pd(fscal,dy21);
1398             tz               = _mm256_mul_pd(fscal,dz21);
1399
1400             /* Update vectorial force */
1401             fix2             = _mm256_add_pd(fix2,tx);
1402             fiy2             = _mm256_add_pd(fiy2,ty);
1403             fiz2             = _mm256_add_pd(fiz2,tz);
1404
1405             fjx1             = _mm256_add_pd(fjx1,tx);
1406             fjy1             = _mm256_add_pd(fjy1,ty);
1407             fjz1             = _mm256_add_pd(fjz1,tz);
1408
1409             /**************************
1410              * CALCULATE INTERACTIONS *
1411              **************************/
1412
1413             /* COULOMB ELECTROSTATICS */
1414             velec            = _mm256_mul_pd(qq22,rinv22);
1415             felec            = _mm256_mul_pd(velec,rinvsq22);
1416
1417             fscal            = felec;
1418
1419             /* Calculate temporary vectorial force */
1420             tx               = _mm256_mul_pd(fscal,dx22);
1421             ty               = _mm256_mul_pd(fscal,dy22);
1422             tz               = _mm256_mul_pd(fscal,dz22);
1423
1424             /* Update vectorial force */
1425             fix2             = _mm256_add_pd(fix2,tx);
1426             fiy2             = _mm256_add_pd(fiy2,ty);
1427             fiz2             = _mm256_add_pd(fiz2,tz);
1428
1429             fjx2             = _mm256_add_pd(fjx2,tx);
1430             fjy2             = _mm256_add_pd(fjy2,ty);
1431             fjz2             = _mm256_add_pd(fjz2,tz);
1432
1433             fjptrA             = f+j_coord_offsetA;
1434             fjptrB             = f+j_coord_offsetB;
1435             fjptrC             = f+j_coord_offsetC;
1436             fjptrD             = f+j_coord_offsetD;
1437
1438             gmx_mm256_decrement_3rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
1439                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1440
1441             /* Inner loop uses 241 flops */
1442         }
1443
1444         if(jidx<j_index_end)
1445         {
1446
1447             /* Get j neighbor index, and coordinate index */
1448             jnrlistA         = jjnr[jidx];
1449             jnrlistB         = jjnr[jidx+1];
1450             jnrlistC         = jjnr[jidx+2];
1451             jnrlistD         = jjnr[jidx+3];
1452             /* Sign of each element will be negative for non-real atoms.
1453              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1454              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
1455              */
1456             tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
1457
1458             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
1459             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
1460             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
1461
1462             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
1463             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
1464             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
1465             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
1466             j_coord_offsetA  = DIM*jnrA;
1467             j_coord_offsetB  = DIM*jnrB;
1468             j_coord_offsetC  = DIM*jnrC;
1469             j_coord_offsetD  = DIM*jnrD;
1470
1471             /* load j atom coordinates */
1472             gmx_mm256_load_3rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
1473                                                  x+j_coord_offsetC,x+j_coord_offsetD,
1474                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
1475
1476             /* Calculate displacement vector */
1477             dx00             = _mm256_sub_pd(ix0,jx0);
1478             dy00             = _mm256_sub_pd(iy0,jy0);
1479             dz00             = _mm256_sub_pd(iz0,jz0);
1480             dx01             = _mm256_sub_pd(ix0,jx1);
1481             dy01             = _mm256_sub_pd(iy0,jy1);
1482             dz01             = _mm256_sub_pd(iz0,jz1);
1483             dx02             = _mm256_sub_pd(ix0,jx2);
1484             dy02             = _mm256_sub_pd(iy0,jy2);
1485             dz02             = _mm256_sub_pd(iz0,jz2);
1486             dx10             = _mm256_sub_pd(ix1,jx0);
1487             dy10             = _mm256_sub_pd(iy1,jy0);
1488             dz10             = _mm256_sub_pd(iz1,jz0);
1489             dx11             = _mm256_sub_pd(ix1,jx1);
1490             dy11             = _mm256_sub_pd(iy1,jy1);
1491             dz11             = _mm256_sub_pd(iz1,jz1);
1492             dx12             = _mm256_sub_pd(ix1,jx2);
1493             dy12             = _mm256_sub_pd(iy1,jy2);
1494             dz12             = _mm256_sub_pd(iz1,jz2);
1495             dx20             = _mm256_sub_pd(ix2,jx0);
1496             dy20             = _mm256_sub_pd(iy2,jy0);
1497             dz20             = _mm256_sub_pd(iz2,jz0);
1498             dx21             = _mm256_sub_pd(ix2,jx1);
1499             dy21             = _mm256_sub_pd(iy2,jy1);
1500             dz21             = _mm256_sub_pd(iz2,jz1);
1501             dx22             = _mm256_sub_pd(ix2,jx2);
1502             dy22             = _mm256_sub_pd(iy2,jy2);
1503             dz22             = _mm256_sub_pd(iz2,jz2);
1504
1505             /* Calculate squared distance and things based on it */
1506             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
1507             rsq01            = gmx_mm256_calc_rsq_pd(dx01,dy01,dz01);
1508             rsq02            = gmx_mm256_calc_rsq_pd(dx02,dy02,dz02);
1509             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
1510             rsq11            = gmx_mm256_calc_rsq_pd(dx11,dy11,dz11);
1511             rsq12            = gmx_mm256_calc_rsq_pd(dx12,dy12,dz12);
1512             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
1513             rsq21            = gmx_mm256_calc_rsq_pd(dx21,dy21,dz21);
1514             rsq22            = gmx_mm256_calc_rsq_pd(dx22,dy22,dz22);
1515
1516             rinv00           = avx256_invsqrt_d(rsq00);
1517             rinv01           = avx256_invsqrt_d(rsq01);
1518             rinv02           = avx256_invsqrt_d(rsq02);
1519             rinv10           = avx256_invsqrt_d(rsq10);
1520             rinv11           = avx256_invsqrt_d(rsq11);
1521             rinv12           = avx256_invsqrt_d(rsq12);
1522             rinv20           = avx256_invsqrt_d(rsq20);
1523             rinv21           = avx256_invsqrt_d(rsq21);
1524             rinv22           = avx256_invsqrt_d(rsq22);
1525
1526             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
1527             rinvsq01         = _mm256_mul_pd(rinv01,rinv01);
1528             rinvsq02         = _mm256_mul_pd(rinv02,rinv02);
1529             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
1530             rinvsq11         = _mm256_mul_pd(rinv11,rinv11);
1531             rinvsq12         = _mm256_mul_pd(rinv12,rinv12);
1532             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
1533             rinvsq21         = _mm256_mul_pd(rinv21,rinv21);
1534             rinvsq22         = _mm256_mul_pd(rinv22,rinv22);
1535
1536             fjx0             = _mm256_setzero_pd();
1537             fjy0             = _mm256_setzero_pd();
1538             fjz0             = _mm256_setzero_pd();
1539             fjx1             = _mm256_setzero_pd();
1540             fjy1             = _mm256_setzero_pd();
1541             fjz1             = _mm256_setzero_pd();
1542             fjx2             = _mm256_setzero_pd();
1543             fjy2             = _mm256_setzero_pd();
1544             fjz2             = _mm256_setzero_pd();
1545
1546             /**************************
1547              * CALCULATE INTERACTIONS *
1548              **************************/
1549
1550             /* COULOMB ELECTROSTATICS */
1551             velec            = _mm256_mul_pd(qq00,rinv00);
1552             felec            = _mm256_mul_pd(velec,rinvsq00);
1553
1554             /* LENNARD-JONES DISPERSION/REPULSION */
1555
1556             rinvsix          = _mm256_mul_pd(_mm256_mul_pd(rinvsq00,rinvsq00),rinvsq00);
1557             fvdw             = _mm256_mul_pd(_mm256_sub_pd(_mm256_mul_pd(c12_00,rinvsix),c6_00),_mm256_mul_pd(rinvsix,rinvsq00));
1558
1559             fscal            = _mm256_add_pd(felec,fvdw);
1560
1561             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1562
1563             /* Calculate temporary vectorial force */
1564             tx               = _mm256_mul_pd(fscal,dx00);
1565             ty               = _mm256_mul_pd(fscal,dy00);
1566             tz               = _mm256_mul_pd(fscal,dz00);
1567
1568             /* Update vectorial force */
1569             fix0             = _mm256_add_pd(fix0,tx);
1570             fiy0             = _mm256_add_pd(fiy0,ty);
1571             fiz0             = _mm256_add_pd(fiz0,tz);
1572
1573             fjx0             = _mm256_add_pd(fjx0,tx);
1574             fjy0             = _mm256_add_pd(fjy0,ty);
1575             fjz0             = _mm256_add_pd(fjz0,tz);
1576
1577             /**************************
1578              * CALCULATE INTERACTIONS *
1579              **************************/
1580
1581             /* COULOMB ELECTROSTATICS */
1582             velec            = _mm256_mul_pd(qq01,rinv01);
1583             felec            = _mm256_mul_pd(velec,rinvsq01);
1584
1585             fscal            = felec;
1586
1587             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1588
1589             /* Calculate temporary vectorial force */
1590             tx               = _mm256_mul_pd(fscal,dx01);
1591             ty               = _mm256_mul_pd(fscal,dy01);
1592             tz               = _mm256_mul_pd(fscal,dz01);
1593
1594             /* Update vectorial force */
1595             fix0             = _mm256_add_pd(fix0,tx);
1596             fiy0             = _mm256_add_pd(fiy0,ty);
1597             fiz0             = _mm256_add_pd(fiz0,tz);
1598
1599             fjx1             = _mm256_add_pd(fjx1,tx);
1600             fjy1             = _mm256_add_pd(fjy1,ty);
1601             fjz1             = _mm256_add_pd(fjz1,tz);
1602
1603             /**************************
1604              * CALCULATE INTERACTIONS *
1605              **************************/
1606
1607             /* COULOMB ELECTROSTATICS */
1608             velec            = _mm256_mul_pd(qq02,rinv02);
1609             felec            = _mm256_mul_pd(velec,rinvsq02);
1610
1611             fscal            = felec;
1612
1613             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1614
1615             /* Calculate temporary vectorial force */
1616             tx               = _mm256_mul_pd(fscal,dx02);
1617             ty               = _mm256_mul_pd(fscal,dy02);
1618             tz               = _mm256_mul_pd(fscal,dz02);
1619
1620             /* Update vectorial force */
1621             fix0             = _mm256_add_pd(fix0,tx);
1622             fiy0             = _mm256_add_pd(fiy0,ty);
1623             fiz0             = _mm256_add_pd(fiz0,tz);
1624
1625             fjx2             = _mm256_add_pd(fjx2,tx);
1626             fjy2             = _mm256_add_pd(fjy2,ty);
1627             fjz2             = _mm256_add_pd(fjz2,tz);
1628
1629             /**************************
1630              * CALCULATE INTERACTIONS *
1631              **************************/
1632
1633             /* COULOMB ELECTROSTATICS */
1634             velec            = _mm256_mul_pd(qq10,rinv10);
1635             felec            = _mm256_mul_pd(velec,rinvsq10);
1636
1637             fscal            = felec;
1638
1639             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1640
1641             /* Calculate temporary vectorial force */
1642             tx               = _mm256_mul_pd(fscal,dx10);
1643             ty               = _mm256_mul_pd(fscal,dy10);
1644             tz               = _mm256_mul_pd(fscal,dz10);
1645
1646             /* Update vectorial force */
1647             fix1             = _mm256_add_pd(fix1,tx);
1648             fiy1             = _mm256_add_pd(fiy1,ty);
1649             fiz1             = _mm256_add_pd(fiz1,tz);
1650
1651             fjx0             = _mm256_add_pd(fjx0,tx);
1652             fjy0             = _mm256_add_pd(fjy0,ty);
1653             fjz0             = _mm256_add_pd(fjz0,tz);
1654
1655             /**************************
1656              * CALCULATE INTERACTIONS *
1657              **************************/
1658
1659             /* COULOMB ELECTROSTATICS */
1660             velec            = _mm256_mul_pd(qq11,rinv11);
1661             felec            = _mm256_mul_pd(velec,rinvsq11);
1662
1663             fscal            = felec;
1664
1665             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1666
1667             /* Calculate temporary vectorial force */
1668             tx               = _mm256_mul_pd(fscal,dx11);
1669             ty               = _mm256_mul_pd(fscal,dy11);
1670             tz               = _mm256_mul_pd(fscal,dz11);
1671
1672             /* Update vectorial force */
1673             fix1             = _mm256_add_pd(fix1,tx);
1674             fiy1             = _mm256_add_pd(fiy1,ty);
1675             fiz1             = _mm256_add_pd(fiz1,tz);
1676
1677             fjx1             = _mm256_add_pd(fjx1,tx);
1678             fjy1             = _mm256_add_pd(fjy1,ty);
1679             fjz1             = _mm256_add_pd(fjz1,tz);
1680
1681             /**************************
1682              * CALCULATE INTERACTIONS *
1683              **************************/
1684
1685             /* COULOMB ELECTROSTATICS */
1686             velec            = _mm256_mul_pd(qq12,rinv12);
1687             felec            = _mm256_mul_pd(velec,rinvsq12);
1688
1689             fscal            = felec;
1690
1691             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1692
1693             /* Calculate temporary vectorial force */
1694             tx               = _mm256_mul_pd(fscal,dx12);
1695             ty               = _mm256_mul_pd(fscal,dy12);
1696             tz               = _mm256_mul_pd(fscal,dz12);
1697
1698             /* Update vectorial force */
1699             fix1             = _mm256_add_pd(fix1,tx);
1700             fiy1             = _mm256_add_pd(fiy1,ty);
1701             fiz1             = _mm256_add_pd(fiz1,tz);
1702
1703             fjx2             = _mm256_add_pd(fjx2,tx);
1704             fjy2             = _mm256_add_pd(fjy2,ty);
1705             fjz2             = _mm256_add_pd(fjz2,tz);
1706
1707             /**************************
1708              * CALCULATE INTERACTIONS *
1709              **************************/
1710
1711             /* COULOMB ELECTROSTATICS */
1712             velec            = _mm256_mul_pd(qq20,rinv20);
1713             felec            = _mm256_mul_pd(velec,rinvsq20);
1714
1715             fscal            = felec;
1716
1717             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1718
1719             /* Calculate temporary vectorial force */
1720             tx               = _mm256_mul_pd(fscal,dx20);
1721             ty               = _mm256_mul_pd(fscal,dy20);
1722             tz               = _mm256_mul_pd(fscal,dz20);
1723
1724             /* Update vectorial force */
1725             fix2             = _mm256_add_pd(fix2,tx);
1726             fiy2             = _mm256_add_pd(fiy2,ty);
1727             fiz2             = _mm256_add_pd(fiz2,tz);
1728
1729             fjx0             = _mm256_add_pd(fjx0,tx);
1730             fjy0             = _mm256_add_pd(fjy0,ty);
1731             fjz0             = _mm256_add_pd(fjz0,tz);
1732
1733             /**************************
1734              * CALCULATE INTERACTIONS *
1735              **************************/
1736
1737             /* COULOMB ELECTROSTATICS */
1738             velec            = _mm256_mul_pd(qq21,rinv21);
1739             felec            = _mm256_mul_pd(velec,rinvsq21);
1740
1741             fscal            = felec;
1742
1743             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1744
1745             /* Calculate temporary vectorial force */
1746             tx               = _mm256_mul_pd(fscal,dx21);
1747             ty               = _mm256_mul_pd(fscal,dy21);
1748             tz               = _mm256_mul_pd(fscal,dz21);
1749
1750             /* Update vectorial force */
1751             fix2             = _mm256_add_pd(fix2,tx);
1752             fiy2             = _mm256_add_pd(fiy2,ty);
1753             fiz2             = _mm256_add_pd(fiz2,tz);
1754
1755             fjx1             = _mm256_add_pd(fjx1,tx);
1756             fjy1             = _mm256_add_pd(fjy1,ty);
1757             fjz1             = _mm256_add_pd(fjz1,tz);
1758
1759             /**************************
1760              * CALCULATE INTERACTIONS *
1761              **************************/
1762
1763             /* COULOMB ELECTROSTATICS */
1764             velec            = _mm256_mul_pd(qq22,rinv22);
1765             felec            = _mm256_mul_pd(velec,rinvsq22);
1766
1767             fscal            = felec;
1768
1769             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1770
1771             /* Calculate temporary vectorial force */
1772             tx               = _mm256_mul_pd(fscal,dx22);
1773             ty               = _mm256_mul_pd(fscal,dy22);
1774             tz               = _mm256_mul_pd(fscal,dz22);
1775
1776             /* Update vectorial force */
1777             fix2             = _mm256_add_pd(fix2,tx);
1778             fiy2             = _mm256_add_pd(fiy2,ty);
1779             fiz2             = _mm256_add_pd(fiz2,tz);
1780
1781             fjx2             = _mm256_add_pd(fjx2,tx);
1782             fjy2             = _mm256_add_pd(fjy2,ty);
1783             fjz2             = _mm256_add_pd(fjz2,tz);
1784
1785             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1786             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1787             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1788             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1789
1790             gmx_mm256_decrement_3rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
1791                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1792
1793             /* Inner loop uses 241 flops */
1794         }
1795
1796         /* End of innermost loop */
1797
1798         gmx_mm256_update_iforce_3atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1799                                                  f+i_coord_offset,fshift+i_shift_offset);
1800
1801         /* Increment number of inner iterations */
1802         inneriter                  += j_index_end - j_index_start;
1803
1804         /* Outer loop uses 18 flops */
1805     }
1806
1807     /* Increment number of outer iterations */
1808     outeriter        += nri;
1809
1810     /* Update outer/inner flops */
1811
1812     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW_W3W3_F,outeriter*18 + inneriter*241);
1813 }