Compile nonbonded kernels as C++
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_double / nb_kernel_ElecEwSh_VdwNone_GeomP1P1_avx_256_double.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014,2015,2017,2018, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS avx_256_double kernel generator.
37  */
38 #include "gmxpre.h"
39
40 #include "config.h"
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "gromacs/gmxlib/nrnb.h"
46
47 #include "kernelutil_x86_avx_256_double.h"
48
49 /*
50  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwNone_GeomP1P1_VF_avx_256_double
51  * Electrostatics interaction: Ewald
52  * VdW interaction:            None
53  * Geometry:                   Particle-Particle
54  * Calculate force/pot:        PotentialAndForce
55  */
56 void
57 nb_kernel_ElecEwSh_VdwNone_GeomP1P1_VF_avx_256_double
58                     (t_nblist                    * gmx_restrict       nlist,
59                      rvec                        * gmx_restrict          xx,
60                      rvec                        * gmx_restrict          ff,
61                      struct t_forcerec           * gmx_restrict          fr,
62                      t_mdatoms                   * gmx_restrict     mdatoms,
63                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
64                      t_nrnb                      * gmx_restrict        nrnb)
65 {
66     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
67      * just 0 for non-waters.
68      * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
69      * jnr indices corresponding to data put in the four positions in the SIMD register.
70      */
71     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
72     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
73     int              jnrA,jnrB,jnrC,jnrD;
74     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
75     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
76     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
77     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
78     real             rcutoff_scalar;
79     real             *shiftvec,*fshift,*x,*f;
80     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
81     real             scratch[4*DIM];
82     __m256d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
83     real *           vdwioffsetptr0;
84     __m256d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
85     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
86     __m256d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
87     __m256d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
88     __m256d          velec,felec,velecsum,facel,crf,krf,krf2;
89     real             *charge;
90     __m128i          ewitab;
91     __m256d          ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
92     __m256d          beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
93     real             *ewtab;
94     __m256d          dummy_mask,cutoff_mask;
95     __m128           tmpmask0,tmpmask1;
96     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
97     __m256d          one     = _mm256_set1_pd(1.0);
98     __m256d          two     = _mm256_set1_pd(2.0);
99     x                = xx[0];
100     f                = ff[0];
101
102     nri              = nlist->nri;
103     iinr             = nlist->iinr;
104     jindex           = nlist->jindex;
105     jjnr             = nlist->jjnr;
106     shiftidx         = nlist->shift;
107     gid              = nlist->gid;
108     shiftvec         = fr->shift_vec[0];
109     fshift           = fr->fshift[0];
110     facel            = _mm256_set1_pd(fr->ic->epsfac);
111     charge           = mdatoms->chargeA;
112
113     sh_ewald         = _mm256_set1_pd(fr->ic->sh_ewald);
114     beta             = _mm256_set1_pd(fr->ic->ewaldcoeff_q);
115     beta2            = _mm256_mul_pd(beta,beta);
116     beta3            = _mm256_mul_pd(beta,beta2);
117
118     ewtab            = fr->ic->tabq_coul_FDV0;
119     ewtabscale       = _mm256_set1_pd(fr->ic->tabq_scale);
120     ewtabhalfspace   = _mm256_set1_pd(0.5/fr->ic->tabq_scale);
121
122     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
123     rcutoff_scalar   = fr->ic->rcoulomb;
124     rcutoff          = _mm256_set1_pd(rcutoff_scalar);
125     rcutoff2         = _mm256_mul_pd(rcutoff,rcutoff);
126
127     /* Avoid stupid compiler warnings */
128     jnrA = jnrB = jnrC = jnrD = 0;
129     j_coord_offsetA = 0;
130     j_coord_offsetB = 0;
131     j_coord_offsetC = 0;
132     j_coord_offsetD = 0;
133
134     outeriter        = 0;
135     inneriter        = 0;
136
137     for(iidx=0;iidx<4*DIM;iidx++)
138     {
139         scratch[iidx] = 0.0;
140     }
141
142     /* Start outer loop over neighborlists */
143     for(iidx=0; iidx<nri; iidx++)
144     {
145         /* Load shift vector for this list */
146         i_shift_offset   = DIM*shiftidx[iidx];
147
148         /* Load limits for loop over neighbors */
149         j_index_start    = jindex[iidx];
150         j_index_end      = jindex[iidx+1];
151
152         /* Get outer coordinate index */
153         inr              = iinr[iidx];
154         i_coord_offset   = DIM*inr;
155
156         /* Load i particle coords and add shift vector */
157         gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
158
159         fix0             = _mm256_setzero_pd();
160         fiy0             = _mm256_setzero_pd();
161         fiz0             = _mm256_setzero_pd();
162
163         /* Load parameters for i particles */
164         iq0              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
165
166         /* Reset potential sums */
167         velecsum         = _mm256_setzero_pd();
168
169         /* Start inner kernel loop */
170         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
171         {
172
173             /* Get j neighbor index, and coordinate index */
174             jnrA             = jjnr[jidx];
175             jnrB             = jjnr[jidx+1];
176             jnrC             = jjnr[jidx+2];
177             jnrD             = jjnr[jidx+3];
178             j_coord_offsetA  = DIM*jnrA;
179             j_coord_offsetB  = DIM*jnrB;
180             j_coord_offsetC  = DIM*jnrC;
181             j_coord_offsetD  = DIM*jnrD;
182
183             /* load j atom coordinates */
184             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
185                                                  x+j_coord_offsetC,x+j_coord_offsetD,
186                                                  &jx0,&jy0,&jz0);
187
188             /* Calculate displacement vector */
189             dx00             = _mm256_sub_pd(ix0,jx0);
190             dy00             = _mm256_sub_pd(iy0,jy0);
191             dz00             = _mm256_sub_pd(iz0,jz0);
192
193             /* Calculate squared distance and things based on it */
194             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
195
196             rinv00           = avx256_invsqrt_d(rsq00);
197
198             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
199
200             /* Load parameters for j particles */
201             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
202                                                                  charge+jnrC+0,charge+jnrD+0);
203
204             /**************************
205              * CALCULATE INTERACTIONS *
206              **************************/
207
208             if (gmx_mm256_any_lt(rsq00,rcutoff2))
209             {
210
211             r00              = _mm256_mul_pd(rsq00,rinv00);
212
213             /* Compute parameters for interactions between i and j atoms */
214             qq00             = _mm256_mul_pd(iq0,jq0);
215
216             /* EWALD ELECTROSTATICS */
217
218             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
219             ewrt             = _mm256_mul_pd(r00,ewtabscale);
220             ewitab           = _mm256_cvttpd_epi32(ewrt);
221             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
222             ewitab           = _mm_slli_epi32(ewitab,2);
223             ewtabF           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
224             ewtabD           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
225             ewtabV           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,2) );
226             ewtabFn          = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,3) );
227             GMX_MM256_FULLTRANSPOSE4_PD(ewtabF,ewtabD,ewtabV,ewtabFn);
228             felec            = _mm256_add_pd(ewtabF,_mm256_mul_pd(eweps,ewtabD));
229             velec            = _mm256_sub_pd(ewtabV,_mm256_mul_pd(_mm256_mul_pd(ewtabhalfspace,eweps),_mm256_add_pd(ewtabF,felec)));
230             velec            = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_sub_pd(rinv00,sh_ewald),velec));
231             felec            = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
232
233             cutoff_mask      = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
234
235             /* Update potential sum for this i atom from the interaction with this j atom. */
236             velec            = _mm256_and_pd(velec,cutoff_mask);
237             velecsum         = _mm256_add_pd(velecsum,velec);
238
239             fscal            = felec;
240
241             fscal            = _mm256_and_pd(fscal,cutoff_mask);
242
243             /* Calculate temporary vectorial force */
244             tx               = _mm256_mul_pd(fscal,dx00);
245             ty               = _mm256_mul_pd(fscal,dy00);
246             tz               = _mm256_mul_pd(fscal,dz00);
247
248             /* Update vectorial force */
249             fix0             = _mm256_add_pd(fix0,tx);
250             fiy0             = _mm256_add_pd(fiy0,ty);
251             fiz0             = _mm256_add_pd(fiz0,tz);
252
253             fjptrA             = f+j_coord_offsetA;
254             fjptrB             = f+j_coord_offsetB;
255             fjptrC             = f+j_coord_offsetC;
256             fjptrD             = f+j_coord_offsetD;
257             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
258
259             }
260
261             /* Inner loop uses 46 flops */
262         }
263
264         if(jidx<j_index_end)
265         {
266
267             /* Get j neighbor index, and coordinate index */
268             jnrlistA         = jjnr[jidx];
269             jnrlistB         = jjnr[jidx+1];
270             jnrlistC         = jjnr[jidx+2];
271             jnrlistD         = jjnr[jidx+3];
272             /* Sign of each element will be negative for non-real atoms.
273              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
274              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
275              */
276             tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
277
278             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
279             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
280             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
281
282             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
283             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
284             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
285             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
286             j_coord_offsetA  = DIM*jnrA;
287             j_coord_offsetB  = DIM*jnrB;
288             j_coord_offsetC  = DIM*jnrC;
289             j_coord_offsetD  = DIM*jnrD;
290
291             /* load j atom coordinates */
292             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
293                                                  x+j_coord_offsetC,x+j_coord_offsetD,
294                                                  &jx0,&jy0,&jz0);
295
296             /* Calculate displacement vector */
297             dx00             = _mm256_sub_pd(ix0,jx0);
298             dy00             = _mm256_sub_pd(iy0,jy0);
299             dz00             = _mm256_sub_pd(iz0,jz0);
300
301             /* Calculate squared distance and things based on it */
302             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
303
304             rinv00           = avx256_invsqrt_d(rsq00);
305
306             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
307
308             /* Load parameters for j particles */
309             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
310                                                                  charge+jnrC+0,charge+jnrD+0);
311
312             /**************************
313              * CALCULATE INTERACTIONS *
314              **************************/
315
316             if (gmx_mm256_any_lt(rsq00,rcutoff2))
317             {
318
319             r00              = _mm256_mul_pd(rsq00,rinv00);
320             r00              = _mm256_andnot_pd(dummy_mask,r00);
321
322             /* Compute parameters for interactions between i and j atoms */
323             qq00             = _mm256_mul_pd(iq0,jq0);
324
325             /* EWALD ELECTROSTATICS */
326
327             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
328             ewrt             = _mm256_mul_pd(r00,ewtabscale);
329             ewitab           = _mm256_cvttpd_epi32(ewrt);
330             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
331             ewitab           = _mm_slli_epi32(ewitab,2);
332             ewtabF           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
333             ewtabD           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
334             ewtabV           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,2) );
335             ewtabFn          = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,3) );
336             GMX_MM256_FULLTRANSPOSE4_PD(ewtabF,ewtabD,ewtabV,ewtabFn);
337             felec            = _mm256_add_pd(ewtabF,_mm256_mul_pd(eweps,ewtabD));
338             velec            = _mm256_sub_pd(ewtabV,_mm256_mul_pd(_mm256_mul_pd(ewtabhalfspace,eweps),_mm256_add_pd(ewtabF,felec)));
339             velec            = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_sub_pd(rinv00,sh_ewald),velec));
340             felec            = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
341
342             cutoff_mask      = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
343
344             /* Update potential sum for this i atom from the interaction with this j atom. */
345             velec            = _mm256_and_pd(velec,cutoff_mask);
346             velec            = _mm256_andnot_pd(dummy_mask,velec);
347             velecsum         = _mm256_add_pd(velecsum,velec);
348
349             fscal            = felec;
350
351             fscal            = _mm256_and_pd(fscal,cutoff_mask);
352
353             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
354
355             /* Calculate temporary vectorial force */
356             tx               = _mm256_mul_pd(fscal,dx00);
357             ty               = _mm256_mul_pd(fscal,dy00);
358             tz               = _mm256_mul_pd(fscal,dz00);
359
360             /* Update vectorial force */
361             fix0             = _mm256_add_pd(fix0,tx);
362             fiy0             = _mm256_add_pd(fiy0,ty);
363             fiz0             = _mm256_add_pd(fiz0,tz);
364
365             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
366             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
367             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
368             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
369             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
370
371             }
372
373             /* Inner loop uses 47 flops */
374         }
375
376         /* End of innermost loop */
377
378         gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
379                                                  f+i_coord_offset,fshift+i_shift_offset);
380
381         ggid                        = gid[iidx];
382         /* Update potential energies */
383         gmx_mm256_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
384
385         /* Increment number of inner iterations */
386         inneriter                  += j_index_end - j_index_start;
387
388         /* Outer loop uses 8 flops */
389     }
390
391     /* Increment number of outer iterations */
392     outeriter        += nri;
393
394     /* Update outer/inner flops */
395
396     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*8 + inneriter*47);
397 }
398 /*
399  * Gromacs nonbonded kernel:   nb_kernel_ElecEwSh_VdwNone_GeomP1P1_F_avx_256_double
400  * Electrostatics interaction: Ewald
401  * VdW interaction:            None
402  * Geometry:                   Particle-Particle
403  * Calculate force/pot:        Force
404  */
405 void
406 nb_kernel_ElecEwSh_VdwNone_GeomP1P1_F_avx_256_double
407                     (t_nblist                    * gmx_restrict       nlist,
408                      rvec                        * gmx_restrict          xx,
409                      rvec                        * gmx_restrict          ff,
410                      struct t_forcerec           * gmx_restrict          fr,
411                      t_mdatoms                   * gmx_restrict     mdatoms,
412                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
413                      t_nrnb                      * gmx_restrict        nrnb)
414 {
415     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
416      * just 0 for non-waters.
417      * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
418      * jnr indices corresponding to data put in the four positions in the SIMD register.
419      */
420     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
421     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
422     int              jnrA,jnrB,jnrC,jnrD;
423     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
424     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
425     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
426     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
427     real             rcutoff_scalar;
428     real             *shiftvec,*fshift,*x,*f;
429     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
430     real             scratch[4*DIM];
431     __m256d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
432     real *           vdwioffsetptr0;
433     __m256d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
434     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
435     __m256d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
436     __m256d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
437     __m256d          velec,felec,velecsum,facel,crf,krf,krf2;
438     real             *charge;
439     __m128i          ewitab;
440     __m256d          ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
441     __m256d          beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
442     real             *ewtab;
443     __m256d          dummy_mask,cutoff_mask;
444     __m128           tmpmask0,tmpmask1;
445     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
446     __m256d          one     = _mm256_set1_pd(1.0);
447     __m256d          two     = _mm256_set1_pd(2.0);
448     x                = xx[0];
449     f                = ff[0];
450
451     nri              = nlist->nri;
452     iinr             = nlist->iinr;
453     jindex           = nlist->jindex;
454     jjnr             = nlist->jjnr;
455     shiftidx         = nlist->shift;
456     gid              = nlist->gid;
457     shiftvec         = fr->shift_vec[0];
458     fshift           = fr->fshift[0];
459     facel            = _mm256_set1_pd(fr->ic->epsfac);
460     charge           = mdatoms->chargeA;
461
462     sh_ewald         = _mm256_set1_pd(fr->ic->sh_ewald);
463     beta             = _mm256_set1_pd(fr->ic->ewaldcoeff_q);
464     beta2            = _mm256_mul_pd(beta,beta);
465     beta3            = _mm256_mul_pd(beta,beta2);
466
467     ewtab            = fr->ic->tabq_coul_F;
468     ewtabscale       = _mm256_set1_pd(fr->ic->tabq_scale);
469     ewtabhalfspace   = _mm256_set1_pd(0.5/fr->ic->tabq_scale);
470
471     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
472     rcutoff_scalar   = fr->ic->rcoulomb;
473     rcutoff          = _mm256_set1_pd(rcutoff_scalar);
474     rcutoff2         = _mm256_mul_pd(rcutoff,rcutoff);
475
476     /* Avoid stupid compiler warnings */
477     jnrA = jnrB = jnrC = jnrD = 0;
478     j_coord_offsetA = 0;
479     j_coord_offsetB = 0;
480     j_coord_offsetC = 0;
481     j_coord_offsetD = 0;
482
483     outeriter        = 0;
484     inneriter        = 0;
485
486     for(iidx=0;iidx<4*DIM;iidx++)
487     {
488         scratch[iidx] = 0.0;
489     }
490
491     /* Start outer loop over neighborlists */
492     for(iidx=0; iidx<nri; iidx++)
493     {
494         /* Load shift vector for this list */
495         i_shift_offset   = DIM*shiftidx[iidx];
496
497         /* Load limits for loop over neighbors */
498         j_index_start    = jindex[iidx];
499         j_index_end      = jindex[iidx+1];
500
501         /* Get outer coordinate index */
502         inr              = iinr[iidx];
503         i_coord_offset   = DIM*inr;
504
505         /* Load i particle coords and add shift vector */
506         gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
507
508         fix0             = _mm256_setzero_pd();
509         fiy0             = _mm256_setzero_pd();
510         fiz0             = _mm256_setzero_pd();
511
512         /* Load parameters for i particles */
513         iq0              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
514
515         /* Start inner kernel loop */
516         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
517         {
518
519             /* Get j neighbor index, and coordinate index */
520             jnrA             = jjnr[jidx];
521             jnrB             = jjnr[jidx+1];
522             jnrC             = jjnr[jidx+2];
523             jnrD             = jjnr[jidx+3];
524             j_coord_offsetA  = DIM*jnrA;
525             j_coord_offsetB  = DIM*jnrB;
526             j_coord_offsetC  = DIM*jnrC;
527             j_coord_offsetD  = DIM*jnrD;
528
529             /* load j atom coordinates */
530             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
531                                                  x+j_coord_offsetC,x+j_coord_offsetD,
532                                                  &jx0,&jy0,&jz0);
533
534             /* Calculate displacement vector */
535             dx00             = _mm256_sub_pd(ix0,jx0);
536             dy00             = _mm256_sub_pd(iy0,jy0);
537             dz00             = _mm256_sub_pd(iz0,jz0);
538
539             /* Calculate squared distance and things based on it */
540             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
541
542             rinv00           = avx256_invsqrt_d(rsq00);
543
544             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
545
546             /* Load parameters for j particles */
547             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
548                                                                  charge+jnrC+0,charge+jnrD+0);
549
550             /**************************
551              * CALCULATE INTERACTIONS *
552              **************************/
553
554             if (gmx_mm256_any_lt(rsq00,rcutoff2))
555             {
556
557             r00              = _mm256_mul_pd(rsq00,rinv00);
558
559             /* Compute parameters for interactions between i and j atoms */
560             qq00             = _mm256_mul_pd(iq0,jq0);
561
562             /* EWALD ELECTROSTATICS */
563
564             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
565             ewrt             = _mm256_mul_pd(r00,ewtabscale);
566             ewitab           = _mm256_cvttpd_epi32(ewrt);
567             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
568             gmx_mm256_load_4pair_swizzle_pd(ewtab + _mm_extract_epi32(ewitab,0),ewtab + _mm_extract_epi32(ewitab,1),
569                                             ewtab + _mm_extract_epi32(ewitab,2),ewtab + _mm_extract_epi32(ewitab,3),
570                                             &ewtabF,&ewtabFn);
571             felec            = _mm256_add_pd(_mm256_mul_pd( _mm256_sub_pd(one,eweps),ewtabF),_mm256_mul_pd(eweps,ewtabFn));
572             felec            = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
573
574             cutoff_mask      = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
575
576             fscal            = felec;
577
578             fscal            = _mm256_and_pd(fscal,cutoff_mask);
579
580             /* Calculate temporary vectorial force */
581             tx               = _mm256_mul_pd(fscal,dx00);
582             ty               = _mm256_mul_pd(fscal,dy00);
583             tz               = _mm256_mul_pd(fscal,dz00);
584
585             /* Update vectorial force */
586             fix0             = _mm256_add_pd(fix0,tx);
587             fiy0             = _mm256_add_pd(fiy0,ty);
588             fiz0             = _mm256_add_pd(fiz0,tz);
589
590             fjptrA             = f+j_coord_offsetA;
591             fjptrB             = f+j_coord_offsetB;
592             fjptrC             = f+j_coord_offsetC;
593             fjptrD             = f+j_coord_offsetD;
594             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
595
596             }
597
598             /* Inner loop uses 39 flops */
599         }
600
601         if(jidx<j_index_end)
602         {
603
604             /* Get j neighbor index, and coordinate index */
605             jnrlistA         = jjnr[jidx];
606             jnrlistB         = jjnr[jidx+1];
607             jnrlistC         = jjnr[jidx+2];
608             jnrlistD         = jjnr[jidx+3];
609             /* Sign of each element will be negative for non-real atoms.
610              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
611              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
612              */
613             tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
614
615             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
616             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
617             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
618
619             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
620             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
621             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
622             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
623             j_coord_offsetA  = DIM*jnrA;
624             j_coord_offsetB  = DIM*jnrB;
625             j_coord_offsetC  = DIM*jnrC;
626             j_coord_offsetD  = DIM*jnrD;
627
628             /* load j atom coordinates */
629             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
630                                                  x+j_coord_offsetC,x+j_coord_offsetD,
631                                                  &jx0,&jy0,&jz0);
632
633             /* Calculate displacement vector */
634             dx00             = _mm256_sub_pd(ix0,jx0);
635             dy00             = _mm256_sub_pd(iy0,jy0);
636             dz00             = _mm256_sub_pd(iz0,jz0);
637
638             /* Calculate squared distance and things based on it */
639             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
640
641             rinv00           = avx256_invsqrt_d(rsq00);
642
643             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
644
645             /* Load parameters for j particles */
646             jq0              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+0,charge+jnrB+0,
647                                                                  charge+jnrC+0,charge+jnrD+0);
648
649             /**************************
650              * CALCULATE INTERACTIONS *
651              **************************/
652
653             if (gmx_mm256_any_lt(rsq00,rcutoff2))
654             {
655
656             r00              = _mm256_mul_pd(rsq00,rinv00);
657             r00              = _mm256_andnot_pd(dummy_mask,r00);
658
659             /* Compute parameters for interactions between i and j atoms */
660             qq00             = _mm256_mul_pd(iq0,jq0);
661
662             /* EWALD ELECTROSTATICS */
663
664             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
665             ewrt             = _mm256_mul_pd(r00,ewtabscale);
666             ewitab           = _mm256_cvttpd_epi32(ewrt);
667             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
668             gmx_mm256_load_4pair_swizzle_pd(ewtab + _mm_extract_epi32(ewitab,0),ewtab + _mm_extract_epi32(ewitab,1),
669                                             ewtab + _mm_extract_epi32(ewitab,2),ewtab + _mm_extract_epi32(ewitab,3),
670                                             &ewtabF,&ewtabFn);
671             felec            = _mm256_add_pd(_mm256_mul_pd( _mm256_sub_pd(one,eweps),ewtabF),_mm256_mul_pd(eweps,ewtabFn));
672             felec            = _mm256_mul_pd(_mm256_mul_pd(qq00,rinv00),_mm256_sub_pd(rinvsq00,felec));
673
674             cutoff_mask      = _mm256_cmp_pd(rsq00,rcutoff2,_CMP_LT_OQ);
675
676             fscal            = felec;
677
678             fscal            = _mm256_and_pd(fscal,cutoff_mask);
679
680             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
681
682             /* Calculate temporary vectorial force */
683             tx               = _mm256_mul_pd(fscal,dx00);
684             ty               = _mm256_mul_pd(fscal,dy00);
685             tz               = _mm256_mul_pd(fscal,dz00);
686
687             /* Update vectorial force */
688             fix0             = _mm256_add_pd(fix0,tx);
689             fiy0             = _mm256_add_pd(fiy0,ty);
690             fiz0             = _mm256_add_pd(fiz0,tz);
691
692             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
693             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
694             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
695             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
696             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
697
698             }
699
700             /* Inner loop uses 40 flops */
701         }
702
703         /* End of innermost loop */
704
705         gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
706                                                  f+i_coord_offset,fshift+i_shift_offset);
707
708         /* Increment number of inner iterations */
709         inneriter                  += j_index_end - j_index_start;
710
711         /* Outer loop uses 7 flops */
712     }
713
714     /* Increment number of outer iterations */
715     outeriter        += nri;
716
717     /* Update outer/inner flops */
718
719     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*40);
720 }