Use full path for legacyheaders
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_double / nb_kernel_ElecEw_VdwNone_GeomW3P1_avx_256_double.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_double 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_double.h"
48 #include "kernelutil_x86_avx_256_double.h"
49
50 /*
51  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwNone_GeomW3P1_VF_avx_256_double
52  * Electrostatics interaction: Ewald
53  * VdW interaction:            None
54  * Geometry:                   Water3-Particle
55  * Calculate force/pot:        PotentialAndForce
56  */
57 void
58 nb_kernel_ElecEw_VdwNone_GeomW3P1_VF_avx_256_double
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 refer to j loop unrolling done with AVX, e.g. for the four 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              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
76     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
77     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
78     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
79     real             rcutoff_scalar;
80     real             *shiftvec,*fshift,*x,*f;
81     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
82     real             scratch[4*DIM];
83     __m256d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
84     real *           vdwioffsetptr0;
85     __m256d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
86     real *           vdwioffsetptr1;
87     __m256d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
88     real *           vdwioffsetptr2;
89     __m256d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
90     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
91     __m256d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
92     __m256d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
93     __m256d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
94     __m256d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
95     __m256d          velec,felec,velecsum,facel,crf,krf,krf2;
96     real             *charge;
97     __m128i          ewitab;
98     __m256d          ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
99     __m256d          beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
100     real             *ewtab;
101     __m256d          dummy_mask,cutoff_mask;
102     __m128           tmpmask0,tmpmask1;
103     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
104     __m256d          one     = _mm256_set1_pd(1.0);
105     __m256d          two     = _mm256_set1_pd(2.0);
106     x                = xx[0];
107     f                = ff[0];
108
109     nri              = nlist->nri;
110     iinr             = nlist->iinr;
111     jindex           = nlist->jindex;
112     jjnr             = nlist->jjnr;
113     shiftidx         = nlist->shift;
114     gid              = nlist->gid;
115     shiftvec         = fr->shift_vec[0];
116     fshift           = fr->fshift[0];
117     facel            = _mm256_set1_pd(fr->epsfac);
118     charge           = mdatoms->chargeA;
119
120     sh_ewald         = _mm256_set1_pd(fr->ic->sh_ewald);
121     beta             = _mm256_set1_pd(fr->ic->ewaldcoeff_q);
122     beta2            = _mm256_mul_pd(beta,beta);
123     beta3            = _mm256_mul_pd(beta,beta2);
124
125     ewtab            = fr->ic->tabq_coul_FDV0;
126     ewtabscale       = _mm256_set1_pd(fr->ic->tabq_scale);
127     ewtabhalfspace   = _mm256_set1_pd(0.5/fr->ic->tabq_scale);
128
129     /* Setup water-specific parameters */
130     inr              = nlist->iinr[0];
131     iq0              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
132     iq1              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+1]));
133     iq2              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+2]));
134
135     /* Avoid stupid compiler warnings */
136     jnrA = jnrB = jnrC = jnrD = 0;
137     j_coord_offsetA = 0;
138     j_coord_offsetB = 0;
139     j_coord_offsetC = 0;
140     j_coord_offsetD = 0;
141
142     outeriter        = 0;
143     inneriter        = 0;
144
145     for(iidx=0;iidx<4*DIM;iidx++)
146     {
147         scratch[iidx] = 0.0;
148     }
149
150     /* Start outer loop over neighborlists */
151     for(iidx=0; iidx<nri; iidx++)
152     {
153         /* Load shift vector for this list */
154         i_shift_offset   = DIM*shiftidx[iidx];
155
156         /* Load limits for loop over neighbors */
157         j_index_start    = jindex[iidx];
158         j_index_end      = jindex[iidx+1];
159
160         /* Get outer coordinate index */
161         inr              = iinr[iidx];
162         i_coord_offset   = DIM*inr;
163
164         /* Load i particle coords and add shift vector */
165         gmx_mm256_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
166                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
167
168         fix0             = _mm256_setzero_pd();
169         fiy0             = _mm256_setzero_pd();
170         fiz0             = _mm256_setzero_pd();
171         fix1             = _mm256_setzero_pd();
172         fiy1             = _mm256_setzero_pd();
173         fiz1             = _mm256_setzero_pd();
174         fix2             = _mm256_setzero_pd();
175         fiy2             = _mm256_setzero_pd();
176         fiz2             = _mm256_setzero_pd();
177
178         /* Reset potential sums */
179         velecsum         = _mm256_setzero_pd();
180
181         /* Start inner kernel loop */
182         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
183         {
184
185             /* Get j neighbor index, and coordinate index */
186             jnrA             = jjnr[jidx];
187             jnrB             = jjnr[jidx+1];
188             jnrC             = jjnr[jidx+2];
189             jnrD             = jjnr[jidx+3];
190             j_coord_offsetA  = DIM*jnrA;
191             j_coord_offsetB  = DIM*jnrB;
192             j_coord_offsetC  = DIM*jnrC;
193             j_coord_offsetD  = DIM*jnrD;
194
195             /* load j atom coordinates */
196             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
197                                                  x+j_coord_offsetC,x+j_coord_offsetD,
198                                                  &jx0,&jy0,&jz0);
199
200             /* Calculate displacement vector */
201             dx00             = _mm256_sub_pd(ix0,jx0);
202             dy00             = _mm256_sub_pd(iy0,jy0);
203             dz00             = _mm256_sub_pd(iz0,jz0);
204             dx10             = _mm256_sub_pd(ix1,jx0);
205             dy10             = _mm256_sub_pd(iy1,jy0);
206             dz10             = _mm256_sub_pd(iz1,jz0);
207             dx20             = _mm256_sub_pd(ix2,jx0);
208             dy20             = _mm256_sub_pd(iy2,jy0);
209             dz20             = _mm256_sub_pd(iz2,jz0);
210
211             /* Calculate squared distance and things based on it */
212             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
213             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
214             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
215
216             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
217             rinv10           = gmx_mm256_invsqrt_pd(rsq10);
218             rinv20           = gmx_mm256_invsqrt_pd(rsq20);
219
220             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
221             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
222             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
223
224             /* Load parameters for j particles */
225             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
226                                                                  charge+jnrC+0,charge+jnrD+0);
227
228             fjx0             = _mm256_setzero_pd();
229             fjy0             = _mm256_setzero_pd();
230             fjz0             = _mm256_setzero_pd();
231
232             /**************************
233              * CALCULATE INTERACTIONS *
234              **************************/
235
236             r00              = _mm256_mul_pd(rsq00,rinv00);
237
238             /* Compute parameters for interactions between i and j atoms */
239             qq00             = _mm256_mul_pd(iq0,jq0);
240
241             /* EWALD ELECTROSTATICS */
242
243             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
244             ewrt             = _mm256_mul_pd(r00,ewtabscale);
245             ewitab           = _mm256_cvttpd_epi32(ewrt);
246             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
247             ewitab           = _mm_slli_epi32(ewitab,2);
248             ewtabF           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
249             ewtabD           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
250             ewtabV           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,2) );
251             ewtabFn          = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,3) );
252             GMX_MM256_FULLTRANSPOSE4_PD(ewtabF,ewtabD,ewtabV,ewtabFn);
253             felec            = _mm256_add_pd(ewtabF,_mm256_mul_pd(eweps,ewtabD));
254             velec            = _mm256_sub_pd(ewtabV,_mm256_mul_pd(_mm256_mul_pd(ewtabhalfspace,eweps),_mm256_add_pd(ewtabF,felec)));
255             velec            = _mm256_mul_pd(qq00,_mm256_sub_pd(rinv00,velec));
256             felec            = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
257
258             /* Update potential sum for this i atom from the interaction with this j atom. */
259             velecsum         = _mm256_add_pd(velecsum,velec);
260
261             fscal            = felec;
262
263             /* Calculate temporary vectorial force */
264             tx               = _mm256_mul_pd(fscal,dx00);
265             ty               = _mm256_mul_pd(fscal,dy00);
266             tz               = _mm256_mul_pd(fscal,dz00);
267
268             /* Update vectorial force */
269             fix0             = _mm256_add_pd(fix0,tx);
270             fiy0             = _mm256_add_pd(fiy0,ty);
271             fiz0             = _mm256_add_pd(fiz0,tz);
272
273             fjx0             = _mm256_add_pd(fjx0,tx);
274             fjy0             = _mm256_add_pd(fjy0,ty);
275             fjz0             = _mm256_add_pd(fjz0,tz);
276
277             /**************************
278              * CALCULATE INTERACTIONS *
279              **************************/
280
281             r10              = _mm256_mul_pd(rsq10,rinv10);
282
283             /* Compute parameters for interactions between i and j atoms */
284             qq10             = _mm256_mul_pd(iq1,jq0);
285
286             /* EWALD ELECTROSTATICS */
287
288             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
289             ewrt             = _mm256_mul_pd(r10,ewtabscale);
290             ewitab           = _mm256_cvttpd_epi32(ewrt);
291             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
292             ewitab           = _mm_slli_epi32(ewitab,2);
293             ewtabF           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
294             ewtabD           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
295             ewtabV           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,2) );
296             ewtabFn          = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,3) );
297             GMX_MM256_FULLTRANSPOSE4_PD(ewtabF,ewtabD,ewtabV,ewtabFn);
298             felec            = _mm256_add_pd(ewtabF,_mm256_mul_pd(eweps,ewtabD));
299             velec            = _mm256_sub_pd(ewtabV,_mm256_mul_pd(_mm256_mul_pd(ewtabhalfspace,eweps),_mm256_add_pd(ewtabF,felec)));
300             velec            = _mm256_mul_pd(qq10,_mm256_sub_pd(rinv10,velec));
301             felec            = _mm256_mul_pd(_mm256_mul_pd(qq10,rinv10),_mm256_sub_pd(rinvsq10,felec));
302
303             /* Update potential sum for this i atom from the interaction with this j atom. */
304             velecsum         = _mm256_add_pd(velecsum,velec);
305
306             fscal            = felec;
307
308             /* Calculate temporary vectorial force */
309             tx               = _mm256_mul_pd(fscal,dx10);
310             ty               = _mm256_mul_pd(fscal,dy10);
311             tz               = _mm256_mul_pd(fscal,dz10);
312
313             /* Update vectorial force */
314             fix1             = _mm256_add_pd(fix1,tx);
315             fiy1             = _mm256_add_pd(fiy1,ty);
316             fiz1             = _mm256_add_pd(fiz1,tz);
317
318             fjx0             = _mm256_add_pd(fjx0,tx);
319             fjy0             = _mm256_add_pd(fjy0,ty);
320             fjz0             = _mm256_add_pd(fjz0,tz);
321
322             /**************************
323              * CALCULATE INTERACTIONS *
324              **************************/
325
326             r20              = _mm256_mul_pd(rsq20,rinv20);
327
328             /* Compute parameters for interactions between i and j atoms */
329             qq20             = _mm256_mul_pd(iq2,jq0);
330
331             /* EWALD ELECTROSTATICS */
332
333             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
334             ewrt             = _mm256_mul_pd(r20,ewtabscale);
335             ewitab           = _mm256_cvttpd_epi32(ewrt);
336             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
337             ewitab           = _mm_slli_epi32(ewitab,2);
338             ewtabF           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
339             ewtabD           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
340             ewtabV           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,2) );
341             ewtabFn          = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,3) );
342             GMX_MM256_FULLTRANSPOSE4_PD(ewtabF,ewtabD,ewtabV,ewtabFn);
343             felec            = _mm256_add_pd(ewtabF,_mm256_mul_pd(eweps,ewtabD));
344             velec            = _mm256_sub_pd(ewtabV,_mm256_mul_pd(_mm256_mul_pd(ewtabhalfspace,eweps),_mm256_add_pd(ewtabF,felec)));
345             velec            = _mm256_mul_pd(qq20,_mm256_sub_pd(rinv20,velec));
346             felec            = _mm256_mul_pd(_mm256_mul_pd(qq20,rinv20),_mm256_sub_pd(rinvsq20,felec));
347
348             /* Update potential sum for this i atom from the interaction with this j atom. */
349             velecsum         = _mm256_add_pd(velecsum,velec);
350
351             fscal            = felec;
352
353             /* Calculate temporary vectorial force */
354             tx               = _mm256_mul_pd(fscal,dx20);
355             ty               = _mm256_mul_pd(fscal,dy20);
356             tz               = _mm256_mul_pd(fscal,dz20);
357
358             /* Update vectorial force */
359             fix2             = _mm256_add_pd(fix2,tx);
360             fiy2             = _mm256_add_pd(fiy2,ty);
361             fiz2             = _mm256_add_pd(fiz2,tz);
362
363             fjx0             = _mm256_add_pd(fjx0,tx);
364             fjy0             = _mm256_add_pd(fjy0,ty);
365             fjz0             = _mm256_add_pd(fjz0,tz);
366
367             fjptrA             = f+j_coord_offsetA;
368             fjptrB             = f+j_coord_offsetB;
369             fjptrC             = f+j_coord_offsetC;
370             fjptrD             = f+j_coord_offsetD;
371
372             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
373
374             /* Inner loop uses 126 flops */
375         }
376
377         if(jidx<j_index_end)
378         {
379
380             /* Get j neighbor index, and coordinate index */
381             jnrlistA         = jjnr[jidx];
382             jnrlistB         = jjnr[jidx+1];
383             jnrlistC         = jjnr[jidx+2];
384             jnrlistD         = jjnr[jidx+3];
385             /* Sign of each element will be negative for non-real atoms.
386              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
387              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
388              */
389             tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
390
391             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
392             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
393             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
394
395             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
396             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
397             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
398             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
399             j_coord_offsetA  = DIM*jnrA;
400             j_coord_offsetB  = DIM*jnrB;
401             j_coord_offsetC  = DIM*jnrC;
402             j_coord_offsetD  = DIM*jnrD;
403
404             /* load j atom coordinates */
405             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
406                                                  x+j_coord_offsetC,x+j_coord_offsetD,
407                                                  &jx0,&jy0,&jz0);
408
409             /* Calculate displacement vector */
410             dx00             = _mm256_sub_pd(ix0,jx0);
411             dy00             = _mm256_sub_pd(iy0,jy0);
412             dz00             = _mm256_sub_pd(iz0,jz0);
413             dx10             = _mm256_sub_pd(ix1,jx0);
414             dy10             = _mm256_sub_pd(iy1,jy0);
415             dz10             = _mm256_sub_pd(iz1,jz0);
416             dx20             = _mm256_sub_pd(ix2,jx0);
417             dy20             = _mm256_sub_pd(iy2,jy0);
418             dz20             = _mm256_sub_pd(iz2,jz0);
419
420             /* Calculate squared distance and things based on it */
421             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
422             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
423             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
424
425             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
426             rinv10           = gmx_mm256_invsqrt_pd(rsq10);
427             rinv20           = gmx_mm256_invsqrt_pd(rsq20);
428
429             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
430             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
431             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
432
433             /* Load parameters for j particles */
434             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
435                                                                  charge+jnrC+0,charge+jnrD+0);
436
437             fjx0             = _mm256_setzero_pd();
438             fjy0             = _mm256_setzero_pd();
439             fjz0             = _mm256_setzero_pd();
440
441             /**************************
442              * CALCULATE INTERACTIONS *
443              **************************/
444
445             r00              = _mm256_mul_pd(rsq00,rinv00);
446             r00              = _mm256_andnot_pd(dummy_mask,r00);
447
448             /* Compute parameters for interactions between i and j atoms */
449             qq00             = _mm256_mul_pd(iq0,jq0);
450
451             /* EWALD ELECTROSTATICS */
452
453             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
454             ewrt             = _mm256_mul_pd(r00,ewtabscale);
455             ewitab           = _mm256_cvttpd_epi32(ewrt);
456             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
457             ewitab           = _mm_slli_epi32(ewitab,2);
458             ewtabF           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
459             ewtabD           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
460             ewtabV           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,2) );
461             ewtabFn          = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,3) );
462             GMX_MM256_FULLTRANSPOSE4_PD(ewtabF,ewtabD,ewtabV,ewtabFn);
463             felec            = _mm256_add_pd(ewtabF,_mm256_mul_pd(eweps,ewtabD));
464             velec            = _mm256_sub_pd(ewtabV,_mm256_mul_pd(_mm256_mul_pd(ewtabhalfspace,eweps),_mm256_add_pd(ewtabF,felec)));
465             velec            = _mm256_mul_pd(qq00,_mm256_sub_pd(rinv00,velec));
466             felec            = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
467
468             /* Update potential sum for this i atom from the interaction with this j atom. */
469             velec            = _mm256_andnot_pd(dummy_mask,velec);
470             velecsum         = _mm256_add_pd(velecsum,velec);
471
472             fscal            = felec;
473
474             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
475
476             /* Calculate temporary vectorial force */
477             tx               = _mm256_mul_pd(fscal,dx00);
478             ty               = _mm256_mul_pd(fscal,dy00);
479             tz               = _mm256_mul_pd(fscal,dz00);
480
481             /* Update vectorial force */
482             fix0             = _mm256_add_pd(fix0,tx);
483             fiy0             = _mm256_add_pd(fiy0,ty);
484             fiz0             = _mm256_add_pd(fiz0,tz);
485
486             fjx0             = _mm256_add_pd(fjx0,tx);
487             fjy0             = _mm256_add_pd(fjy0,ty);
488             fjz0             = _mm256_add_pd(fjz0,tz);
489
490             /**************************
491              * CALCULATE INTERACTIONS *
492              **************************/
493
494             r10              = _mm256_mul_pd(rsq10,rinv10);
495             r10              = _mm256_andnot_pd(dummy_mask,r10);
496
497             /* Compute parameters for interactions between i and j atoms */
498             qq10             = _mm256_mul_pd(iq1,jq0);
499
500             /* EWALD ELECTROSTATICS */
501
502             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
503             ewrt             = _mm256_mul_pd(r10,ewtabscale);
504             ewitab           = _mm256_cvttpd_epi32(ewrt);
505             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
506             ewitab           = _mm_slli_epi32(ewitab,2);
507             ewtabF           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
508             ewtabD           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
509             ewtabV           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,2) );
510             ewtabFn          = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,3) );
511             GMX_MM256_FULLTRANSPOSE4_PD(ewtabF,ewtabD,ewtabV,ewtabFn);
512             felec            = _mm256_add_pd(ewtabF,_mm256_mul_pd(eweps,ewtabD));
513             velec            = _mm256_sub_pd(ewtabV,_mm256_mul_pd(_mm256_mul_pd(ewtabhalfspace,eweps),_mm256_add_pd(ewtabF,felec)));
514             velec            = _mm256_mul_pd(qq10,_mm256_sub_pd(rinv10,velec));
515             felec            = _mm256_mul_pd(_mm256_mul_pd(qq10,rinv10),_mm256_sub_pd(rinvsq10,felec));
516
517             /* Update potential sum for this i atom from the interaction with this j atom. */
518             velec            = _mm256_andnot_pd(dummy_mask,velec);
519             velecsum         = _mm256_add_pd(velecsum,velec);
520
521             fscal            = felec;
522
523             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
524
525             /* Calculate temporary vectorial force */
526             tx               = _mm256_mul_pd(fscal,dx10);
527             ty               = _mm256_mul_pd(fscal,dy10);
528             tz               = _mm256_mul_pd(fscal,dz10);
529
530             /* Update vectorial force */
531             fix1             = _mm256_add_pd(fix1,tx);
532             fiy1             = _mm256_add_pd(fiy1,ty);
533             fiz1             = _mm256_add_pd(fiz1,tz);
534
535             fjx0             = _mm256_add_pd(fjx0,tx);
536             fjy0             = _mm256_add_pd(fjy0,ty);
537             fjz0             = _mm256_add_pd(fjz0,tz);
538
539             /**************************
540              * CALCULATE INTERACTIONS *
541              **************************/
542
543             r20              = _mm256_mul_pd(rsq20,rinv20);
544             r20              = _mm256_andnot_pd(dummy_mask,r20);
545
546             /* Compute parameters for interactions between i and j atoms */
547             qq20             = _mm256_mul_pd(iq2,jq0);
548
549             /* EWALD ELECTROSTATICS */
550
551             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
552             ewrt             = _mm256_mul_pd(r20,ewtabscale);
553             ewitab           = _mm256_cvttpd_epi32(ewrt);
554             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
555             ewitab           = _mm_slli_epi32(ewitab,2);
556             ewtabF           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
557             ewtabD           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
558             ewtabV           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,2) );
559             ewtabFn          = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,3) );
560             GMX_MM256_FULLTRANSPOSE4_PD(ewtabF,ewtabD,ewtabV,ewtabFn);
561             felec            = _mm256_add_pd(ewtabF,_mm256_mul_pd(eweps,ewtabD));
562             velec            = _mm256_sub_pd(ewtabV,_mm256_mul_pd(_mm256_mul_pd(ewtabhalfspace,eweps),_mm256_add_pd(ewtabF,felec)));
563             velec            = _mm256_mul_pd(qq20,_mm256_sub_pd(rinv20,velec));
564             felec            = _mm256_mul_pd(_mm256_mul_pd(qq20,rinv20),_mm256_sub_pd(rinvsq20,felec));
565
566             /* Update potential sum for this i atom from the interaction with this j atom. */
567             velec            = _mm256_andnot_pd(dummy_mask,velec);
568             velecsum         = _mm256_add_pd(velecsum,velec);
569
570             fscal            = felec;
571
572             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
573
574             /* Calculate temporary vectorial force */
575             tx               = _mm256_mul_pd(fscal,dx20);
576             ty               = _mm256_mul_pd(fscal,dy20);
577             tz               = _mm256_mul_pd(fscal,dz20);
578
579             /* Update vectorial force */
580             fix2             = _mm256_add_pd(fix2,tx);
581             fiy2             = _mm256_add_pd(fiy2,ty);
582             fiz2             = _mm256_add_pd(fiz2,tz);
583
584             fjx0             = _mm256_add_pd(fjx0,tx);
585             fjy0             = _mm256_add_pd(fjy0,ty);
586             fjz0             = _mm256_add_pd(fjz0,tz);
587
588             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
589             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
590             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
591             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
592
593             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
594
595             /* Inner loop uses 129 flops */
596         }
597
598         /* End of innermost loop */
599
600         gmx_mm256_update_iforce_3atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
601                                                  f+i_coord_offset,fshift+i_shift_offset);
602
603         ggid                        = gid[iidx];
604         /* Update potential energies */
605         gmx_mm256_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
606
607         /* Increment number of inner iterations */
608         inneriter                  += j_index_end - j_index_start;
609
610         /* Outer loop uses 19 flops */
611     }
612
613     /* Increment number of outer iterations */
614     outeriter        += nri;
615
616     /* Update outer/inner flops */
617
618     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3_VF,outeriter*19 + inneriter*129);
619 }
620 /*
621  * Gromacs nonbonded kernel:   nb_kernel_ElecEw_VdwNone_GeomW3P1_F_avx_256_double
622  * Electrostatics interaction: Ewald
623  * VdW interaction:            None
624  * Geometry:                   Water3-Particle
625  * Calculate force/pot:        Force
626  */
627 void
628 nb_kernel_ElecEw_VdwNone_GeomW3P1_F_avx_256_double
629                     (t_nblist                    * gmx_restrict       nlist,
630                      rvec                        * gmx_restrict          xx,
631                      rvec                        * gmx_restrict          ff,
632                      t_forcerec                  * gmx_restrict          fr,
633                      t_mdatoms                   * gmx_restrict     mdatoms,
634                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
635                      t_nrnb                      * gmx_restrict        nrnb)
636 {
637     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
638      * just 0 for non-waters.
639      * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
640      * jnr indices corresponding to data put in the four positions in the SIMD register.
641      */
642     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
643     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
644     int              jnrA,jnrB,jnrC,jnrD;
645     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
646     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
647     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
648     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
649     real             rcutoff_scalar;
650     real             *shiftvec,*fshift,*x,*f;
651     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
652     real             scratch[4*DIM];
653     __m256d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
654     real *           vdwioffsetptr0;
655     __m256d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
656     real *           vdwioffsetptr1;
657     __m256d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
658     real *           vdwioffsetptr2;
659     __m256d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
660     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
661     __m256d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
662     __m256d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
663     __m256d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
664     __m256d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
665     __m256d          velec,felec,velecsum,facel,crf,krf,krf2;
666     real             *charge;
667     __m128i          ewitab;
668     __m256d          ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
669     __m256d          beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
670     real             *ewtab;
671     __m256d          dummy_mask,cutoff_mask;
672     __m128           tmpmask0,tmpmask1;
673     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
674     __m256d          one     = _mm256_set1_pd(1.0);
675     __m256d          two     = _mm256_set1_pd(2.0);
676     x                = xx[0];
677     f                = ff[0];
678
679     nri              = nlist->nri;
680     iinr             = nlist->iinr;
681     jindex           = nlist->jindex;
682     jjnr             = nlist->jjnr;
683     shiftidx         = nlist->shift;
684     gid              = nlist->gid;
685     shiftvec         = fr->shift_vec[0];
686     fshift           = fr->fshift[0];
687     facel            = _mm256_set1_pd(fr->epsfac);
688     charge           = mdatoms->chargeA;
689
690     sh_ewald         = _mm256_set1_pd(fr->ic->sh_ewald);
691     beta             = _mm256_set1_pd(fr->ic->ewaldcoeff_q);
692     beta2            = _mm256_mul_pd(beta,beta);
693     beta3            = _mm256_mul_pd(beta,beta2);
694
695     ewtab            = fr->ic->tabq_coul_F;
696     ewtabscale       = _mm256_set1_pd(fr->ic->tabq_scale);
697     ewtabhalfspace   = _mm256_set1_pd(0.5/fr->ic->tabq_scale);
698
699     /* Setup water-specific parameters */
700     inr              = nlist->iinr[0];
701     iq0              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
702     iq1              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+1]));
703     iq2              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+2]));
704
705     /* Avoid stupid compiler warnings */
706     jnrA = jnrB = jnrC = jnrD = 0;
707     j_coord_offsetA = 0;
708     j_coord_offsetB = 0;
709     j_coord_offsetC = 0;
710     j_coord_offsetD = 0;
711
712     outeriter        = 0;
713     inneriter        = 0;
714
715     for(iidx=0;iidx<4*DIM;iidx++)
716     {
717         scratch[iidx] = 0.0;
718     }
719
720     /* Start outer loop over neighborlists */
721     for(iidx=0; iidx<nri; iidx++)
722     {
723         /* Load shift vector for this list */
724         i_shift_offset   = DIM*shiftidx[iidx];
725
726         /* Load limits for loop over neighbors */
727         j_index_start    = jindex[iidx];
728         j_index_end      = jindex[iidx+1];
729
730         /* Get outer coordinate index */
731         inr              = iinr[iidx];
732         i_coord_offset   = DIM*inr;
733
734         /* Load i particle coords and add shift vector */
735         gmx_mm256_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
736                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
737
738         fix0             = _mm256_setzero_pd();
739         fiy0             = _mm256_setzero_pd();
740         fiz0             = _mm256_setzero_pd();
741         fix1             = _mm256_setzero_pd();
742         fiy1             = _mm256_setzero_pd();
743         fiz1             = _mm256_setzero_pd();
744         fix2             = _mm256_setzero_pd();
745         fiy2             = _mm256_setzero_pd();
746         fiz2             = _mm256_setzero_pd();
747
748         /* Start inner kernel loop */
749         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
750         {
751
752             /* Get j neighbor index, and coordinate index */
753             jnrA             = jjnr[jidx];
754             jnrB             = jjnr[jidx+1];
755             jnrC             = jjnr[jidx+2];
756             jnrD             = jjnr[jidx+3];
757             j_coord_offsetA  = DIM*jnrA;
758             j_coord_offsetB  = DIM*jnrB;
759             j_coord_offsetC  = DIM*jnrC;
760             j_coord_offsetD  = DIM*jnrD;
761
762             /* load j atom coordinates */
763             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
764                                                  x+j_coord_offsetC,x+j_coord_offsetD,
765                                                  &jx0,&jy0,&jz0);
766
767             /* Calculate displacement vector */
768             dx00             = _mm256_sub_pd(ix0,jx0);
769             dy00             = _mm256_sub_pd(iy0,jy0);
770             dz00             = _mm256_sub_pd(iz0,jz0);
771             dx10             = _mm256_sub_pd(ix1,jx0);
772             dy10             = _mm256_sub_pd(iy1,jy0);
773             dz10             = _mm256_sub_pd(iz1,jz0);
774             dx20             = _mm256_sub_pd(ix2,jx0);
775             dy20             = _mm256_sub_pd(iy2,jy0);
776             dz20             = _mm256_sub_pd(iz2,jz0);
777
778             /* Calculate squared distance and things based on it */
779             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
780             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
781             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
782
783             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
784             rinv10           = gmx_mm256_invsqrt_pd(rsq10);
785             rinv20           = gmx_mm256_invsqrt_pd(rsq20);
786
787             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
788             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
789             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
790
791             /* Load parameters for j particles */
792             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
793                                                                  charge+jnrC+0,charge+jnrD+0);
794
795             fjx0             = _mm256_setzero_pd();
796             fjy0             = _mm256_setzero_pd();
797             fjz0             = _mm256_setzero_pd();
798
799             /**************************
800              * CALCULATE INTERACTIONS *
801              **************************/
802
803             r00              = _mm256_mul_pd(rsq00,rinv00);
804
805             /* Compute parameters for interactions between i and j atoms */
806             qq00             = _mm256_mul_pd(iq0,jq0);
807
808             /* EWALD ELECTROSTATICS */
809
810             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
811             ewrt             = _mm256_mul_pd(r00,ewtabscale);
812             ewitab           = _mm256_cvttpd_epi32(ewrt);
813             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
814             gmx_mm256_load_4pair_swizzle_pd(ewtab + _mm_extract_epi32(ewitab,0),ewtab + _mm_extract_epi32(ewitab,1),
815                                             ewtab + _mm_extract_epi32(ewitab,2),ewtab + _mm_extract_epi32(ewitab,3),
816                                             &ewtabF,&ewtabFn);
817             felec            = _mm256_add_pd(_mm256_mul_pd( _mm256_sub_pd(one,eweps),ewtabF),_mm256_mul_pd(eweps,ewtabFn));
818             felec            = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
819
820             fscal            = felec;
821
822             /* Calculate temporary vectorial force */
823             tx               = _mm256_mul_pd(fscal,dx00);
824             ty               = _mm256_mul_pd(fscal,dy00);
825             tz               = _mm256_mul_pd(fscal,dz00);
826
827             /* Update vectorial force */
828             fix0             = _mm256_add_pd(fix0,tx);
829             fiy0             = _mm256_add_pd(fiy0,ty);
830             fiz0             = _mm256_add_pd(fiz0,tz);
831
832             fjx0             = _mm256_add_pd(fjx0,tx);
833             fjy0             = _mm256_add_pd(fjy0,ty);
834             fjz0             = _mm256_add_pd(fjz0,tz);
835
836             /**************************
837              * CALCULATE INTERACTIONS *
838              **************************/
839
840             r10              = _mm256_mul_pd(rsq10,rinv10);
841
842             /* Compute parameters for interactions between i and j atoms */
843             qq10             = _mm256_mul_pd(iq1,jq0);
844
845             /* EWALD ELECTROSTATICS */
846
847             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
848             ewrt             = _mm256_mul_pd(r10,ewtabscale);
849             ewitab           = _mm256_cvttpd_epi32(ewrt);
850             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
851             gmx_mm256_load_4pair_swizzle_pd(ewtab + _mm_extract_epi32(ewitab,0),ewtab + _mm_extract_epi32(ewitab,1),
852                                             ewtab + _mm_extract_epi32(ewitab,2),ewtab + _mm_extract_epi32(ewitab,3),
853                                             &ewtabF,&ewtabFn);
854             felec            = _mm256_add_pd(_mm256_mul_pd( _mm256_sub_pd(one,eweps),ewtabF),_mm256_mul_pd(eweps,ewtabFn));
855             felec            = _mm256_mul_pd(_mm256_mul_pd(qq10,rinv10),_mm256_sub_pd(rinvsq10,felec));
856
857             fscal            = felec;
858
859             /* Calculate temporary vectorial force */
860             tx               = _mm256_mul_pd(fscal,dx10);
861             ty               = _mm256_mul_pd(fscal,dy10);
862             tz               = _mm256_mul_pd(fscal,dz10);
863
864             /* Update vectorial force */
865             fix1             = _mm256_add_pd(fix1,tx);
866             fiy1             = _mm256_add_pd(fiy1,ty);
867             fiz1             = _mm256_add_pd(fiz1,tz);
868
869             fjx0             = _mm256_add_pd(fjx0,tx);
870             fjy0             = _mm256_add_pd(fjy0,ty);
871             fjz0             = _mm256_add_pd(fjz0,tz);
872
873             /**************************
874              * CALCULATE INTERACTIONS *
875              **************************/
876
877             r20              = _mm256_mul_pd(rsq20,rinv20);
878
879             /* Compute parameters for interactions between i and j atoms */
880             qq20             = _mm256_mul_pd(iq2,jq0);
881
882             /* EWALD ELECTROSTATICS */
883
884             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
885             ewrt             = _mm256_mul_pd(r20,ewtabscale);
886             ewitab           = _mm256_cvttpd_epi32(ewrt);
887             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
888             gmx_mm256_load_4pair_swizzle_pd(ewtab + _mm_extract_epi32(ewitab,0),ewtab + _mm_extract_epi32(ewitab,1),
889                                             ewtab + _mm_extract_epi32(ewitab,2),ewtab + _mm_extract_epi32(ewitab,3),
890                                             &ewtabF,&ewtabFn);
891             felec            = _mm256_add_pd(_mm256_mul_pd( _mm256_sub_pd(one,eweps),ewtabF),_mm256_mul_pd(eweps,ewtabFn));
892             felec            = _mm256_mul_pd(_mm256_mul_pd(qq20,rinv20),_mm256_sub_pd(rinvsq20,felec));
893
894             fscal            = felec;
895
896             /* Calculate temporary vectorial force */
897             tx               = _mm256_mul_pd(fscal,dx20);
898             ty               = _mm256_mul_pd(fscal,dy20);
899             tz               = _mm256_mul_pd(fscal,dz20);
900
901             /* Update vectorial force */
902             fix2             = _mm256_add_pd(fix2,tx);
903             fiy2             = _mm256_add_pd(fiy2,ty);
904             fiz2             = _mm256_add_pd(fiz2,tz);
905
906             fjx0             = _mm256_add_pd(fjx0,tx);
907             fjy0             = _mm256_add_pd(fjy0,ty);
908             fjz0             = _mm256_add_pd(fjz0,tz);
909
910             fjptrA             = f+j_coord_offsetA;
911             fjptrB             = f+j_coord_offsetB;
912             fjptrC             = f+j_coord_offsetC;
913             fjptrD             = f+j_coord_offsetD;
914
915             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
916
917             /* Inner loop uses 111 flops */
918         }
919
920         if(jidx<j_index_end)
921         {
922
923             /* Get j neighbor index, and coordinate index */
924             jnrlistA         = jjnr[jidx];
925             jnrlistB         = jjnr[jidx+1];
926             jnrlistC         = jjnr[jidx+2];
927             jnrlistD         = jjnr[jidx+3];
928             /* Sign of each element will be negative for non-real atoms.
929              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
930              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
931              */
932             tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
933
934             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
935             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
936             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
937
938             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
939             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
940             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
941             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
942             j_coord_offsetA  = DIM*jnrA;
943             j_coord_offsetB  = DIM*jnrB;
944             j_coord_offsetC  = DIM*jnrC;
945             j_coord_offsetD  = DIM*jnrD;
946
947             /* load j atom coordinates */
948             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
949                                                  x+j_coord_offsetC,x+j_coord_offsetD,
950                                                  &jx0,&jy0,&jz0);
951
952             /* Calculate displacement vector */
953             dx00             = _mm256_sub_pd(ix0,jx0);
954             dy00             = _mm256_sub_pd(iy0,jy0);
955             dz00             = _mm256_sub_pd(iz0,jz0);
956             dx10             = _mm256_sub_pd(ix1,jx0);
957             dy10             = _mm256_sub_pd(iy1,jy0);
958             dz10             = _mm256_sub_pd(iz1,jz0);
959             dx20             = _mm256_sub_pd(ix2,jx0);
960             dy20             = _mm256_sub_pd(iy2,jy0);
961             dz20             = _mm256_sub_pd(iz2,jz0);
962
963             /* Calculate squared distance and things based on it */
964             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
965             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
966             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
967
968             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
969             rinv10           = gmx_mm256_invsqrt_pd(rsq10);
970             rinv20           = gmx_mm256_invsqrt_pd(rsq20);
971
972             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
973             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
974             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
975
976             /* Load parameters for j particles */
977             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
978                                                                  charge+jnrC+0,charge+jnrD+0);
979
980             fjx0             = _mm256_setzero_pd();
981             fjy0             = _mm256_setzero_pd();
982             fjz0             = _mm256_setzero_pd();
983
984             /**************************
985              * CALCULATE INTERACTIONS *
986              **************************/
987
988             r00              = _mm256_mul_pd(rsq00,rinv00);
989             r00              = _mm256_andnot_pd(dummy_mask,r00);
990
991             /* Compute parameters for interactions between i and j atoms */
992             qq00             = _mm256_mul_pd(iq0,jq0);
993
994             /* EWALD ELECTROSTATICS */
995
996             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
997             ewrt             = _mm256_mul_pd(r00,ewtabscale);
998             ewitab           = _mm256_cvttpd_epi32(ewrt);
999             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
1000             gmx_mm256_load_4pair_swizzle_pd(ewtab + _mm_extract_epi32(ewitab,0),ewtab + _mm_extract_epi32(ewitab,1),
1001                                             ewtab + _mm_extract_epi32(ewitab,2),ewtab + _mm_extract_epi32(ewitab,3),
1002                                             &ewtabF,&ewtabFn);
1003             felec            = _mm256_add_pd(_mm256_mul_pd( _mm256_sub_pd(one,eweps),ewtabF),_mm256_mul_pd(eweps,ewtabFn));
1004             felec            = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
1005
1006             fscal            = felec;
1007
1008             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1009
1010             /* Calculate temporary vectorial force */
1011             tx               = _mm256_mul_pd(fscal,dx00);
1012             ty               = _mm256_mul_pd(fscal,dy00);
1013             tz               = _mm256_mul_pd(fscal,dz00);
1014
1015             /* Update vectorial force */
1016             fix0             = _mm256_add_pd(fix0,tx);
1017             fiy0             = _mm256_add_pd(fiy0,ty);
1018             fiz0             = _mm256_add_pd(fiz0,tz);
1019
1020             fjx0             = _mm256_add_pd(fjx0,tx);
1021             fjy0             = _mm256_add_pd(fjy0,ty);
1022             fjz0             = _mm256_add_pd(fjz0,tz);
1023
1024             /**************************
1025              * CALCULATE INTERACTIONS *
1026              **************************/
1027
1028             r10              = _mm256_mul_pd(rsq10,rinv10);
1029             r10              = _mm256_andnot_pd(dummy_mask,r10);
1030
1031             /* Compute parameters for interactions between i and j atoms */
1032             qq10             = _mm256_mul_pd(iq1,jq0);
1033
1034             /* EWALD ELECTROSTATICS */
1035
1036             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1037             ewrt             = _mm256_mul_pd(r10,ewtabscale);
1038             ewitab           = _mm256_cvttpd_epi32(ewrt);
1039             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
1040             gmx_mm256_load_4pair_swizzle_pd(ewtab + _mm_extract_epi32(ewitab,0),ewtab + _mm_extract_epi32(ewitab,1),
1041                                             ewtab + _mm_extract_epi32(ewitab,2),ewtab + _mm_extract_epi32(ewitab,3),
1042                                             &ewtabF,&ewtabFn);
1043             felec            = _mm256_add_pd(_mm256_mul_pd( _mm256_sub_pd(one,eweps),ewtabF),_mm256_mul_pd(eweps,ewtabFn));
1044             felec            = _mm256_mul_pd(_mm256_mul_pd(qq10,rinv10),_mm256_sub_pd(rinvsq10,felec));
1045
1046             fscal            = felec;
1047
1048             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1049
1050             /* Calculate temporary vectorial force */
1051             tx               = _mm256_mul_pd(fscal,dx10);
1052             ty               = _mm256_mul_pd(fscal,dy10);
1053             tz               = _mm256_mul_pd(fscal,dz10);
1054
1055             /* Update vectorial force */
1056             fix1             = _mm256_add_pd(fix1,tx);
1057             fiy1             = _mm256_add_pd(fiy1,ty);
1058             fiz1             = _mm256_add_pd(fiz1,tz);
1059
1060             fjx0             = _mm256_add_pd(fjx0,tx);
1061             fjy0             = _mm256_add_pd(fjy0,ty);
1062             fjz0             = _mm256_add_pd(fjz0,tz);
1063
1064             /**************************
1065              * CALCULATE INTERACTIONS *
1066              **************************/
1067
1068             r20              = _mm256_mul_pd(rsq20,rinv20);
1069             r20              = _mm256_andnot_pd(dummy_mask,r20);
1070
1071             /* Compute parameters for interactions between i and j atoms */
1072             qq20             = _mm256_mul_pd(iq2,jq0);
1073
1074             /* EWALD ELECTROSTATICS */
1075
1076             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
1077             ewrt             = _mm256_mul_pd(r20,ewtabscale);
1078             ewitab           = _mm256_cvttpd_epi32(ewrt);
1079             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
1080             gmx_mm256_load_4pair_swizzle_pd(ewtab + _mm_extract_epi32(ewitab,0),ewtab + _mm_extract_epi32(ewitab,1),
1081                                             ewtab + _mm_extract_epi32(ewitab,2),ewtab + _mm_extract_epi32(ewitab,3),
1082                                             &ewtabF,&ewtabFn);
1083             felec            = _mm256_add_pd(_mm256_mul_pd( _mm256_sub_pd(one,eweps),ewtabF),_mm256_mul_pd(eweps,ewtabFn));
1084             felec            = _mm256_mul_pd(_mm256_mul_pd(qq20,rinv20),_mm256_sub_pd(rinvsq20,felec));
1085
1086             fscal            = felec;
1087
1088             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1089
1090             /* Calculate temporary vectorial force */
1091             tx               = _mm256_mul_pd(fscal,dx20);
1092             ty               = _mm256_mul_pd(fscal,dy20);
1093             tz               = _mm256_mul_pd(fscal,dz20);
1094
1095             /* Update vectorial force */
1096             fix2             = _mm256_add_pd(fix2,tx);
1097             fiy2             = _mm256_add_pd(fiy2,ty);
1098             fiz2             = _mm256_add_pd(fiz2,tz);
1099
1100             fjx0             = _mm256_add_pd(fjx0,tx);
1101             fjy0             = _mm256_add_pd(fjy0,ty);
1102             fjz0             = _mm256_add_pd(fjz0,tz);
1103
1104             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1105             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1106             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1107             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1108
1109             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1110
1111             /* Inner loop uses 114 flops */
1112         }
1113
1114         /* End of innermost loop */
1115
1116         gmx_mm256_update_iforce_3atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1117                                                  f+i_coord_offset,fshift+i_shift_offset);
1118
1119         /* Increment number of inner iterations */
1120         inneriter                  += j_index_end - j_index_start;
1121
1122         /* Outer loop uses 18 flops */
1123     }
1124
1125     /* Increment number of outer iterations */
1126     outeriter        += nri;
1127
1128     /* Update outer/inner flops */
1129
1130     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3_F,outeriter*18 + inneriter*114);
1131 }