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