Merge release-5-0 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_single / nb_kernel_ElecGB_VdwNone_GeomP1P1_avx_256_single.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS avx_256_single kernel generator.
37  */
38 #include "gmxpre.h"
39
40 #include "config.h"
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "gromacs/legacyheaders/types/simple.h"
46 #include "gromacs/math/vec.h"
47 #include "gromacs/legacyheaders/nrnb.h"
48
49 #include "gromacs/simd/math_x86_avx_256_single.h"
50 #include "kernelutil_x86_avx_256_single.h"
51
52 /*
53  * Gromacs nonbonded kernel:   nb_kernel_ElecGB_VdwNone_GeomP1P1_VF_avx_256_single
54  * Electrostatics interaction: GeneralizedBorn
55  * VdW interaction:            None
56  * Geometry:                   Particle-Particle
57  * Calculate force/pot:        PotentialAndForce
58  */
59 void
60 nb_kernel_ElecGB_VdwNone_GeomP1P1_VF_avx_256_single
61                     (t_nblist                    * gmx_restrict       nlist,
62                      rvec                        * gmx_restrict          xx,
63                      rvec                        * gmx_restrict          ff,
64                      t_forcerec                  * gmx_restrict          fr,
65                      t_mdatoms                   * gmx_restrict     mdatoms,
66                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
67                      t_nrnb                      * gmx_restrict        nrnb)
68 {
69     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
70      * just 0 for non-waters.
71      * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
72      * jnr indices corresponding to data put in the four positions in the SIMD register.
73      */
74     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
75     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
76     int              jnrA,jnrB,jnrC,jnrD;
77     int              jnrE,jnrF,jnrG,jnrH;
78     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
79     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
80     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
81     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
82     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
83     real             rcutoff_scalar;
84     real             *shiftvec,*fshift,*x,*f;
85     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
86     real             scratch[4*DIM];
87     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
88     real *           vdwioffsetptr0;
89     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
90     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
91     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
92     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
93     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
94     real             *charge;
95     __m256i          gbitab;
96     __m128i          gbitab_lo,gbitab_hi;
97     __m256           vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
98     __m256           minushalf = _mm256_set1_ps(-0.5);
99     real             *invsqrta,*dvda,*gbtab;
100     __m256i          vfitab;
101     __m128i          vfitab_lo,vfitab_hi;
102     __m128i          ifour       = _mm_set1_epi32(4);
103     __m256           rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
104     real             *vftab;
105     __m256           dummy_mask,cutoff_mask;
106     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
107     __m256           one     = _mm256_set1_ps(1.0);
108     __m256           two     = _mm256_set1_ps(2.0);
109     x                = xx[0];
110     f                = ff[0];
111
112     nri              = nlist->nri;
113     iinr             = nlist->iinr;
114     jindex           = nlist->jindex;
115     jjnr             = nlist->jjnr;
116     shiftidx         = nlist->shift;
117     gid              = nlist->gid;
118     shiftvec         = fr->shift_vec[0];
119     fshift           = fr->fshift[0];
120     facel            = _mm256_set1_ps(fr->epsfac);
121     charge           = mdatoms->chargeA;
122
123     invsqrta         = fr->invsqrta;
124     dvda             = fr->dvda;
125     gbtabscale       = _mm256_set1_ps(fr->gbtab.scale);
126     gbtab            = fr->gbtab.data;
127     gbinvepsdiff     = _mm256_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
128
129     /* Avoid stupid compiler warnings */
130     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
131     j_coord_offsetA = 0;
132     j_coord_offsetB = 0;
133     j_coord_offsetC = 0;
134     j_coord_offsetD = 0;
135     j_coord_offsetE = 0;
136     j_coord_offsetF = 0;
137     j_coord_offsetG = 0;
138     j_coord_offsetH = 0;
139
140     outeriter        = 0;
141     inneriter        = 0;
142
143     for(iidx=0;iidx<4*DIM;iidx++)
144     {
145         scratch[iidx] = 0.0;
146     }
147
148     /* Start outer loop over neighborlists */
149     for(iidx=0; iidx<nri; iidx++)
150     {
151         /* Load shift vector for this list */
152         i_shift_offset   = DIM*shiftidx[iidx];
153
154         /* Load limits for loop over neighbors */
155         j_index_start    = jindex[iidx];
156         j_index_end      = jindex[iidx+1];
157
158         /* Get outer coordinate index */
159         inr              = iinr[iidx];
160         i_coord_offset   = DIM*inr;
161
162         /* Load i particle coords and add shift vector */
163         gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
164
165         fix0             = _mm256_setzero_ps();
166         fiy0             = _mm256_setzero_ps();
167         fiz0             = _mm256_setzero_ps();
168
169         /* Load parameters for i particles */
170         iq0              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
171         isai0            = _mm256_set1_ps(invsqrta[inr+0]);
172
173         /* Reset potential sums */
174         velecsum         = _mm256_setzero_ps();
175         vgbsum           = _mm256_setzero_ps();
176         dvdasum          = _mm256_setzero_ps();
177
178         /* Start inner kernel loop */
179         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
180         {
181
182             /* Get j neighbor index, and coordinate index */
183             jnrA             = jjnr[jidx];
184             jnrB             = jjnr[jidx+1];
185             jnrC             = jjnr[jidx+2];
186             jnrD             = jjnr[jidx+3];
187             jnrE             = jjnr[jidx+4];
188             jnrF             = jjnr[jidx+5];
189             jnrG             = jjnr[jidx+6];
190             jnrH             = jjnr[jidx+7];
191             j_coord_offsetA  = DIM*jnrA;
192             j_coord_offsetB  = DIM*jnrB;
193             j_coord_offsetC  = DIM*jnrC;
194             j_coord_offsetD  = DIM*jnrD;
195             j_coord_offsetE  = DIM*jnrE;
196             j_coord_offsetF  = DIM*jnrF;
197             j_coord_offsetG  = DIM*jnrG;
198             j_coord_offsetH  = DIM*jnrH;
199
200             /* load j atom coordinates */
201             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
202                                                  x+j_coord_offsetC,x+j_coord_offsetD,
203                                                  x+j_coord_offsetE,x+j_coord_offsetF,
204                                                  x+j_coord_offsetG,x+j_coord_offsetH,
205                                                  &jx0,&jy0,&jz0);
206
207             /* Calculate displacement vector */
208             dx00             = _mm256_sub_ps(ix0,jx0);
209             dy00             = _mm256_sub_ps(iy0,jy0);
210             dz00             = _mm256_sub_ps(iz0,jz0);
211
212             /* Calculate squared distance and things based on it */
213             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
214
215             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
216
217             /* Load parameters for j particles */
218             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
219                                                                  charge+jnrC+0,charge+jnrD+0,
220                                                                  charge+jnrE+0,charge+jnrF+0,
221                                                                  charge+jnrG+0,charge+jnrH+0);
222             isaj0            = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
223                                                                  invsqrta+jnrC+0,invsqrta+jnrD+0,
224                                                                  invsqrta+jnrE+0,invsqrta+jnrF+0,
225                                                                  invsqrta+jnrG+0,invsqrta+jnrH+0);
226
227             /**************************
228              * CALCULATE INTERACTIONS *
229              **************************/
230
231             r00              = _mm256_mul_ps(rsq00,rinv00);
232
233             /* Compute parameters for interactions between i and j atoms */
234             qq00             = _mm256_mul_ps(iq0,jq0);
235
236             /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
237             isaprod          = _mm256_mul_ps(isai0,isaj0);
238             gbqqfactor       = _mm256_xor_ps(signbit,_mm256_mul_ps(qq00,_mm256_mul_ps(isaprod,gbinvepsdiff)));
239             gbscale          = _mm256_mul_ps(isaprod,gbtabscale);
240
241             /* Calculate generalized born table index - this is a separate table from the normal one,
242              * but we use the same procedure by multiplying r with scale and truncating to integer.
243              */
244             rt               = _mm256_mul_ps(r00,gbscale);
245             gbitab           = _mm256_cvttps_epi32(rt);
246             gbeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
247             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
248             gbitab_lo        = _mm256_extractf128_si256(gbitab,0x0);
249             gbitab_hi        = _mm256_extractf128_si256(gbitab,0x1);
250             gbitab_lo        = _mm_slli_epi32(gbitab_lo,2);
251             gbitab_hi        = _mm_slli_epi32(gbitab_hi,2);
252             Y                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
253                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
254             F                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
255                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
256             G                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
257                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
258             H                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
259                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
260             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
261             Heps             = _mm256_mul_ps(gbeps,H);
262             Fp               = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
263             VV               = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
264             vgb              = _mm256_mul_ps(gbqqfactor,VV);
265
266             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
267             fgb              = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
268             dvdatmp          = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r00)));
269             dvdasum          = _mm256_add_ps(dvdasum,dvdatmp);
270             fjptrA           = dvda+jnrA;
271             fjptrB           = dvda+jnrB;
272             fjptrC           = dvda+jnrC;
273             fjptrD           = dvda+jnrD;
274             fjptrE           = dvda+jnrE;
275             fjptrF           = dvda+jnrF;
276             fjptrG           = dvda+jnrG;
277             fjptrH           = dvda+jnrH;
278             gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
279                                                  _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj0,isaj0)));
280             velec            = _mm256_mul_ps(qq00,rinv00);
281             felec            = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv00),fgb),rinv00);
282
283             /* Update potential sum for this i atom from the interaction with this j atom. */
284             velecsum         = _mm256_add_ps(velecsum,velec);
285             vgbsum           = _mm256_add_ps(vgbsum,vgb);
286
287             fscal            = felec;
288
289             /* Calculate temporary vectorial force */
290             tx               = _mm256_mul_ps(fscal,dx00);
291             ty               = _mm256_mul_ps(fscal,dy00);
292             tz               = _mm256_mul_ps(fscal,dz00);
293
294             /* Update vectorial force */
295             fix0             = _mm256_add_ps(fix0,tx);
296             fiy0             = _mm256_add_ps(fiy0,ty);
297             fiz0             = _mm256_add_ps(fiz0,tz);
298
299             fjptrA             = f+j_coord_offsetA;
300             fjptrB             = f+j_coord_offsetB;
301             fjptrC             = f+j_coord_offsetC;
302             fjptrD             = f+j_coord_offsetD;
303             fjptrE             = f+j_coord_offsetE;
304             fjptrF             = f+j_coord_offsetF;
305             fjptrG             = f+j_coord_offsetG;
306             fjptrH             = f+j_coord_offsetH;
307             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
308
309             /* Inner loop uses 57 flops */
310         }
311
312         if(jidx<j_index_end)
313         {
314
315             /* Get j neighbor index, and coordinate index */
316             jnrlistA         = jjnr[jidx];
317             jnrlistB         = jjnr[jidx+1];
318             jnrlistC         = jjnr[jidx+2];
319             jnrlistD         = jjnr[jidx+3];
320             jnrlistE         = jjnr[jidx+4];
321             jnrlistF         = jjnr[jidx+5];
322             jnrlistG         = jjnr[jidx+6];
323             jnrlistH         = jjnr[jidx+7];
324             /* Sign of each element will be negative for non-real atoms.
325              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
326              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
327              */
328             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
329                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
330                                             
331             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
332             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
333             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
334             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
335             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
336             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
337             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
338             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
339             j_coord_offsetA  = DIM*jnrA;
340             j_coord_offsetB  = DIM*jnrB;
341             j_coord_offsetC  = DIM*jnrC;
342             j_coord_offsetD  = DIM*jnrD;
343             j_coord_offsetE  = DIM*jnrE;
344             j_coord_offsetF  = DIM*jnrF;
345             j_coord_offsetG  = DIM*jnrG;
346             j_coord_offsetH  = DIM*jnrH;
347
348             /* load j atom coordinates */
349             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
350                                                  x+j_coord_offsetC,x+j_coord_offsetD,
351                                                  x+j_coord_offsetE,x+j_coord_offsetF,
352                                                  x+j_coord_offsetG,x+j_coord_offsetH,
353                                                  &jx0,&jy0,&jz0);
354
355             /* Calculate displacement vector */
356             dx00             = _mm256_sub_ps(ix0,jx0);
357             dy00             = _mm256_sub_ps(iy0,jy0);
358             dz00             = _mm256_sub_ps(iz0,jz0);
359
360             /* Calculate squared distance and things based on it */
361             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
362
363             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
364
365             /* Load parameters for j particles */
366             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
367                                                                  charge+jnrC+0,charge+jnrD+0,
368                                                                  charge+jnrE+0,charge+jnrF+0,
369                                                                  charge+jnrG+0,charge+jnrH+0);
370             isaj0            = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
371                                                                  invsqrta+jnrC+0,invsqrta+jnrD+0,
372                                                                  invsqrta+jnrE+0,invsqrta+jnrF+0,
373                                                                  invsqrta+jnrG+0,invsqrta+jnrH+0);
374
375             /**************************
376              * CALCULATE INTERACTIONS *
377              **************************/
378
379             r00              = _mm256_mul_ps(rsq00,rinv00);
380             r00              = _mm256_andnot_ps(dummy_mask,r00);
381
382             /* Compute parameters for interactions between i and j atoms */
383             qq00             = _mm256_mul_ps(iq0,jq0);
384
385             /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
386             isaprod          = _mm256_mul_ps(isai0,isaj0);
387             gbqqfactor       = _mm256_xor_ps(signbit,_mm256_mul_ps(qq00,_mm256_mul_ps(isaprod,gbinvepsdiff)));
388             gbscale          = _mm256_mul_ps(isaprod,gbtabscale);
389
390             /* Calculate generalized born table index - this is a separate table from the normal one,
391              * but we use the same procedure by multiplying r with scale and truncating to integer.
392              */
393             rt               = _mm256_mul_ps(r00,gbscale);
394             gbitab           = _mm256_cvttps_epi32(rt);
395             gbeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
396             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
397             gbitab_lo        = _mm256_extractf128_si256(gbitab,0x0);
398             gbitab_hi        = _mm256_extractf128_si256(gbitab,0x1);
399             gbitab_lo        = _mm_slli_epi32(gbitab_lo,2);
400             gbitab_hi        = _mm_slli_epi32(gbitab_hi,2);
401             Y                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
402                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
403             F                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
404                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
405             G                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
406                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
407             H                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
408                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
409             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
410             Heps             = _mm256_mul_ps(gbeps,H);
411             Fp               = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
412             VV               = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
413             vgb              = _mm256_mul_ps(gbqqfactor,VV);
414
415             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
416             fgb              = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
417             dvdatmp          = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r00)));
418             dvdatmp          = _mm256_andnot_ps(dummy_mask,dvdatmp);
419             dvdasum          = _mm256_add_ps(dvdasum,dvdatmp);
420             /* The pointers to scratch make sure that this code with compilers that take gmx_restrict seriously (e.g. icc 13) really can't screw things up. */
421             fjptrA             = (jnrlistA>=0) ? dvda+jnrA : scratch;
422             fjptrB             = (jnrlistB>=0) ? dvda+jnrB : scratch;
423             fjptrC             = (jnrlistC>=0) ? dvda+jnrC : scratch;
424             fjptrD             = (jnrlistD>=0) ? dvda+jnrD : scratch;
425             fjptrE             = (jnrlistE>=0) ? dvda+jnrE : scratch;
426             fjptrF             = (jnrlistF>=0) ? dvda+jnrF : scratch;
427             fjptrG             = (jnrlistG>=0) ? dvda+jnrG : scratch;
428             fjptrH             = (jnrlistH>=0) ? dvda+jnrH : scratch;
429             gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
430                                                  _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj0,isaj0)));
431             velec            = _mm256_mul_ps(qq00,rinv00);
432             felec            = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv00),fgb),rinv00);
433
434             /* Update potential sum for this i atom from the interaction with this j atom. */
435             velec            = _mm256_andnot_ps(dummy_mask,velec);
436             velecsum         = _mm256_add_ps(velecsum,velec);
437             vgb              = _mm256_andnot_ps(dummy_mask,vgb);
438             vgbsum           = _mm256_add_ps(vgbsum,vgb);
439
440             fscal            = felec;
441
442             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
443
444             /* Calculate temporary vectorial force */
445             tx               = _mm256_mul_ps(fscal,dx00);
446             ty               = _mm256_mul_ps(fscal,dy00);
447             tz               = _mm256_mul_ps(fscal,dz00);
448
449             /* Update vectorial force */
450             fix0             = _mm256_add_ps(fix0,tx);
451             fiy0             = _mm256_add_ps(fiy0,ty);
452             fiz0             = _mm256_add_ps(fiz0,tz);
453
454             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
455             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
456             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
457             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
458             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
459             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
460             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
461             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
462             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
463
464             /* Inner loop uses 58 flops */
465         }
466
467         /* End of innermost loop */
468
469         gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
470                                                  f+i_coord_offset,fshift+i_shift_offset);
471
472         ggid                        = gid[iidx];
473         /* Update potential energies */
474         gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
475         gmx_mm256_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
476         dvdasum = _mm256_mul_ps(dvdasum, _mm256_mul_ps(isai0,isai0));
477         gmx_mm256_update_1pot_ps(dvdasum,dvda+inr);
478
479         /* Increment number of inner iterations */
480         inneriter                  += j_index_end - j_index_start;
481
482         /* Outer loop uses 9 flops */
483     }
484
485     /* Increment number of outer iterations */
486     outeriter        += nri;
487
488     /* Update outer/inner flops */
489
490     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VF,outeriter*9 + inneriter*58);
491 }
492 /*
493  * Gromacs nonbonded kernel:   nb_kernel_ElecGB_VdwNone_GeomP1P1_F_avx_256_single
494  * Electrostatics interaction: GeneralizedBorn
495  * VdW interaction:            None
496  * Geometry:                   Particle-Particle
497  * Calculate force/pot:        Force
498  */
499 void
500 nb_kernel_ElecGB_VdwNone_GeomP1P1_F_avx_256_single
501                     (t_nblist                    * gmx_restrict       nlist,
502                      rvec                        * gmx_restrict          xx,
503                      rvec                        * gmx_restrict          ff,
504                      t_forcerec                  * gmx_restrict          fr,
505                      t_mdatoms                   * gmx_restrict     mdatoms,
506                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
507                      t_nrnb                      * gmx_restrict        nrnb)
508 {
509     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
510      * just 0 for non-waters.
511      * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
512      * jnr indices corresponding to data put in the four positions in the SIMD register.
513      */
514     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
515     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
516     int              jnrA,jnrB,jnrC,jnrD;
517     int              jnrE,jnrF,jnrG,jnrH;
518     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
519     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
520     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
521     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
522     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
523     real             rcutoff_scalar;
524     real             *shiftvec,*fshift,*x,*f;
525     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
526     real             scratch[4*DIM];
527     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
528     real *           vdwioffsetptr0;
529     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
530     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
531     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
532     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
533     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
534     real             *charge;
535     __m256i          gbitab;
536     __m128i          gbitab_lo,gbitab_hi;
537     __m256           vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
538     __m256           minushalf = _mm256_set1_ps(-0.5);
539     real             *invsqrta,*dvda,*gbtab;
540     __m256i          vfitab;
541     __m128i          vfitab_lo,vfitab_hi;
542     __m128i          ifour       = _mm_set1_epi32(4);
543     __m256           rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
544     real             *vftab;
545     __m256           dummy_mask,cutoff_mask;
546     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
547     __m256           one     = _mm256_set1_ps(1.0);
548     __m256           two     = _mm256_set1_ps(2.0);
549     x                = xx[0];
550     f                = ff[0];
551
552     nri              = nlist->nri;
553     iinr             = nlist->iinr;
554     jindex           = nlist->jindex;
555     jjnr             = nlist->jjnr;
556     shiftidx         = nlist->shift;
557     gid              = nlist->gid;
558     shiftvec         = fr->shift_vec[0];
559     fshift           = fr->fshift[0];
560     facel            = _mm256_set1_ps(fr->epsfac);
561     charge           = mdatoms->chargeA;
562
563     invsqrta         = fr->invsqrta;
564     dvda             = fr->dvda;
565     gbtabscale       = _mm256_set1_ps(fr->gbtab.scale);
566     gbtab            = fr->gbtab.data;
567     gbinvepsdiff     = _mm256_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
568
569     /* Avoid stupid compiler warnings */
570     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
571     j_coord_offsetA = 0;
572     j_coord_offsetB = 0;
573     j_coord_offsetC = 0;
574     j_coord_offsetD = 0;
575     j_coord_offsetE = 0;
576     j_coord_offsetF = 0;
577     j_coord_offsetG = 0;
578     j_coord_offsetH = 0;
579
580     outeriter        = 0;
581     inneriter        = 0;
582
583     for(iidx=0;iidx<4*DIM;iidx++)
584     {
585         scratch[iidx] = 0.0;
586     }
587
588     /* Start outer loop over neighborlists */
589     for(iidx=0; iidx<nri; iidx++)
590     {
591         /* Load shift vector for this list */
592         i_shift_offset   = DIM*shiftidx[iidx];
593
594         /* Load limits for loop over neighbors */
595         j_index_start    = jindex[iidx];
596         j_index_end      = jindex[iidx+1];
597
598         /* Get outer coordinate index */
599         inr              = iinr[iidx];
600         i_coord_offset   = DIM*inr;
601
602         /* Load i particle coords and add shift vector */
603         gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
604
605         fix0             = _mm256_setzero_ps();
606         fiy0             = _mm256_setzero_ps();
607         fiz0             = _mm256_setzero_ps();
608
609         /* Load parameters for i particles */
610         iq0              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
611         isai0            = _mm256_set1_ps(invsqrta[inr+0]);
612
613         dvdasum          = _mm256_setzero_ps();
614
615         /* Start inner kernel loop */
616         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
617         {
618
619             /* Get j neighbor index, and coordinate index */
620             jnrA             = jjnr[jidx];
621             jnrB             = jjnr[jidx+1];
622             jnrC             = jjnr[jidx+2];
623             jnrD             = jjnr[jidx+3];
624             jnrE             = jjnr[jidx+4];
625             jnrF             = jjnr[jidx+5];
626             jnrG             = jjnr[jidx+6];
627             jnrH             = jjnr[jidx+7];
628             j_coord_offsetA  = DIM*jnrA;
629             j_coord_offsetB  = DIM*jnrB;
630             j_coord_offsetC  = DIM*jnrC;
631             j_coord_offsetD  = DIM*jnrD;
632             j_coord_offsetE  = DIM*jnrE;
633             j_coord_offsetF  = DIM*jnrF;
634             j_coord_offsetG  = DIM*jnrG;
635             j_coord_offsetH  = DIM*jnrH;
636
637             /* load j atom coordinates */
638             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
639                                                  x+j_coord_offsetC,x+j_coord_offsetD,
640                                                  x+j_coord_offsetE,x+j_coord_offsetF,
641                                                  x+j_coord_offsetG,x+j_coord_offsetH,
642                                                  &jx0,&jy0,&jz0);
643
644             /* Calculate displacement vector */
645             dx00             = _mm256_sub_ps(ix0,jx0);
646             dy00             = _mm256_sub_ps(iy0,jy0);
647             dz00             = _mm256_sub_ps(iz0,jz0);
648
649             /* Calculate squared distance and things based on it */
650             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
651
652             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
653
654             /* Load parameters for j particles */
655             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
656                                                                  charge+jnrC+0,charge+jnrD+0,
657                                                                  charge+jnrE+0,charge+jnrF+0,
658                                                                  charge+jnrG+0,charge+jnrH+0);
659             isaj0            = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
660                                                                  invsqrta+jnrC+0,invsqrta+jnrD+0,
661                                                                  invsqrta+jnrE+0,invsqrta+jnrF+0,
662                                                                  invsqrta+jnrG+0,invsqrta+jnrH+0);
663
664             /**************************
665              * CALCULATE INTERACTIONS *
666              **************************/
667
668             r00              = _mm256_mul_ps(rsq00,rinv00);
669
670             /* Compute parameters for interactions between i and j atoms */
671             qq00             = _mm256_mul_ps(iq0,jq0);
672
673             /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
674             isaprod          = _mm256_mul_ps(isai0,isaj0);
675             gbqqfactor       = _mm256_xor_ps(signbit,_mm256_mul_ps(qq00,_mm256_mul_ps(isaprod,gbinvepsdiff)));
676             gbscale          = _mm256_mul_ps(isaprod,gbtabscale);
677
678             /* Calculate generalized born table index - this is a separate table from the normal one,
679              * but we use the same procedure by multiplying r with scale and truncating to integer.
680              */
681             rt               = _mm256_mul_ps(r00,gbscale);
682             gbitab           = _mm256_cvttps_epi32(rt);
683             gbeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
684             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
685             gbitab_lo        = _mm256_extractf128_si256(gbitab,0x0);
686             gbitab_hi        = _mm256_extractf128_si256(gbitab,0x1);
687             gbitab_lo        = _mm_slli_epi32(gbitab_lo,2);
688             gbitab_hi        = _mm_slli_epi32(gbitab_hi,2);
689             Y                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
690                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
691             F                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
692                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
693             G                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
694                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
695             H                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
696                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
697             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
698             Heps             = _mm256_mul_ps(gbeps,H);
699             Fp               = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
700             VV               = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
701             vgb              = _mm256_mul_ps(gbqqfactor,VV);
702
703             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
704             fgb              = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
705             dvdatmp          = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r00)));
706             dvdasum          = _mm256_add_ps(dvdasum,dvdatmp);
707             fjptrA           = dvda+jnrA;
708             fjptrB           = dvda+jnrB;
709             fjptrC           = dvda+jnrC;
710             fjptrD           = dvda+jnrD;
711             fjptrE           = dvda+jnrE;
712             fjptrF           = dvda+jnrF;
713             fjptrG           = dvda+jnrG;
714             fjptrH           = dvda+jnrH;
715             gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
716                                                  _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj0,isaj0)));
717             velec            = _mm256_mul_ps(qq00,rinv00);
718             felec            = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv00),fgb),rinv00);
719
720             fscal            = felec;
721
722             /* Calculate temporary vectorial force */
723             tx               = _mm256_mul_ps(fscal,dx00);
724             ty               = _mm256_mul_ps(fscal,dy00);
725             tz               = _mm256_mul_ps(fscal,dz00);
726
727             /* Update vectorial force */
728             fix0             = _mm256_add_ps(fix0,tx);
729             fiy0             = _mm256_add_ps(fiy0,ty);
730             fiz0             = _mm256_add_ps(fiz0,tz);
731
732             fjptrA             = f+j_coord_offsetA;
733             fjptrB             = f+j_coord_offsetB;
734             fjptrC             = f+j_coord_offsetC;
735             fjptrD             = f+j_coord_offsetD;
736             fjptrE             = f+j_coord_offsetE;
737             fjptrF             = f+j_coord_offsetF;
738             fjptrG             = f+j_coord_offsetG;
739             fjptrH             = f+j_coord_offsetH;
740             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
741
742             /* Inner loop uses 55 flops */
743         }
744
745         if(jidx<j_index_end)
746         {
747
748             /* Get j neighbor index, and coordinate index */
749             jnrlistA         = jjnr[jidx];
750             jnrlistB         = jjnr[jidx+1];
751             jnrlistC         = jjnr[jidx+2];
752             jnrlistD         = jjnr[jidx+3];
753             jnrlistE         = jjnr[jidx+4];
754             jnrlistF         = jjnr[jidx+5];
755             jnrlistG         = jjnr[jidx+6];
756             jnrlistH         = jjnr[jidx+7];
757             /* Sign of each element will be negative for non-real atoms.
758              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
759              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
760              */
761             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
762                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
763                                             
764             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
765             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
766             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
767             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
768             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
769             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
770             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
771             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
772             j_coord_offsetA  = DIM*jnrA;
773             j_coord_offsetB  = DIM*jnrB;
774             j_coord_offsetC  = DIM*jnrC;
775             j_coord_offsetD  = DIM*jnrD;
776             j_coord_offsetE  = DIM*jnrE;
777             j_coord_offsetF  = DIM*jnrF;
778             j_coord_offsetG  = DIM*jnrG;
779             j_coord_offsetH  = DIM*jnrH;
780
781             /* load j atom coordinates */
782             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
783                                                  x+j_coord_offsetC,x+j_coord_offsetD,
784                                                  x+j_coord_offsetE,x+j_coord_offsetF,
785                                                  x+j_coord_offsetG,x+j_coord_offsetH,
786                                                  &jx0,&jy0,&jz0);
787
788             /* Calculate displacement vector */
789             dx00             = _mm256_sub_ps(ix0,jx0);
790             dy00             = _mm256_sub_ps(iy0,jy0);
791             dz00             = _mm256_sub_ps(iz0,jz0);
792
793             /* Calculate squared distance and things based on it */
794             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
795
796             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
797
798             /* Load parameters for j particles */
799             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
800                                                                  charge+jnrC+0,charge+jnrD+0,
801                                                                  charge+jnrE+0,charge+jnrF+0,
802                                                                  charge+jnrG+0,charge+jnrH+0);
803             isaj0            = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+0,invsqrta+jnrB+0,
804                                                                  invsqrta+jnrC+0,invsqrta+jnrD+0,
805                                                                  invsqrta+jnrE+0,invsqrta+jnrF+0,
806                                                                  invsqrta+jnrG+0,invsqrta+jnrH+0);
807
808             /**************************
809              * CALCULATE INTERACTIONS *
810              **************************/
811
812             r00              = _mm256_mul_ps(rsq00,rinv00);
813             r00              = _mm256_andnot_ps(dummy_mask,r00);
814
815             /* Compute parameters for interactions between i and j atoms */
816             qq00             = _mm256_mul_ps(iq0,jq0);
817
818             /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
819             isaprod          = _mm256_mul_ps(isai0,isaj0);
820             gbqqfactor       = _mm256_xor_ps(signbit,_mm256_mul_ps(qq00,_mm256_mul_ps(isaprod,gbinvepsdiff)));
821             gbscale          = _mm256_mul_ps(isaprod,gbtabscale);
822
823             /* Calculate generalized born table index - this is a separate table from the normal one,
824              * but we use the same procedure by multiplying r with scale and truncating to integer.
825              */
826             rt               = _mm256_mul_ps(r00,gbscale);
827             gbitab           = _mm256_cvttps_epi32(rt);
828             gbeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
829             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
830             gbitab_lo        = _mm256_extractf128_si256(gbitab,0x0);
831             gbitab_hi        = _mm256_extractf128_si256(gbitab,0x1);
832             gbitab_lo        = _mm_slli_epi32(gbitab_lo,2);
833             gbitab_hi        = _mm_slli_epi32(gbitab_hi,2);
834             Y                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
835                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
836             F                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
837                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
838             G                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
839                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
840             H                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
841                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
842             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
843             Heps             = _mm256_mul_ps(gbeps,H);
844             Fp               = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
845             VV               = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
846             vgb              = _mm256_mul_ps(gbqqfactor,VV);
847
848             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
849             fgb              = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
850             dvdatmp          = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r00)));
851             dvdatmp          = _mm256_andnot_ps(dummy_mask,dvdatmp);
852             dvdasum          = _mm256_add_ps(dvdasum,dvdatmp);
853             /* The pointers to scratch make sure that this code with compilers that take gmx_restrict seriously (e.g. icc 13) really can't screw things up. */
854             fjptrA             = (jnrlistA>=0) ? dvda+jnrA : scratch;
855             fjptrB             = (jnrlistB>=0) ? dvda+jnrB : scratch;
856             fjptrC             = (jnrlistC>=0) ? dvda+jnrC : scratch;
857             fjptrD             = (jnrlistD>=0) ? dvda+jnrD : scratch;
858             fjptrE             = (jnrlistE>=0) ? dvda+jnrE : scratch;
859             fjptrF             = (jnrlistF>=0) ? dvda+jnrF : scratch;
860             fjptrG             = (jnrlistG>=0) ? dvda+jnrG : scratch;
861             fjptrH             = (jnrlistH>=0) ? dvda+jnrH : scratch;
862             gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
863                                                  _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj0,isaj0)));
864             velec            = _mm256_mul_ps(qq00,rinv00);
865             felec            = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv00),fgb),rinv00);
866
867             fscal            = felec;
868
869             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
870
871             /* Calculate temporary vectorial force */
872             tx               = _mm256_mul_ps(fscal,dx00);
873             ty               = _mm256_mul_ps(fscal,dy00);
874             tz               = _mm256_mul_ps(fscal,dz00);
875
876             /* Update vectorial force */
877             fix0             = _mm256_add_ps(fix0,tx);
878             fiy0             = _mm256_add_ps(fiy0,ty);
879             fiz0             = _mm256_add_ps(fiz0,tz);
880
881             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
882             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
883             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
884             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
885             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
886             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
887             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
888             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
889             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
890
891             /* Inner loop uses 56 flops */
892         }
893
894         /* End of innermost loop */
895
896         gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
897                                                  f+i_coord_offset,fshift+i_shift_offset);
898
899         dvdasum = _mm256_mul_ps(dvdasum, _mm256_mul_ps(isai0,isai0));
900         gmx_mm256_update_1pot_ps(dvdasum,dvda+inr);
901
902         /* Increment number of inner iterations */
903         inneriter                  += j_index_end - j_index_start;
904
905         /* Outer loop uses 7 flops */
906     }
907
908     /* Increment number of outer iterations */
909     outeriter        += nri;
910
911     /* Update outer/inner flops */
912
913     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_F,outeriter*7 + inneriter*56);
914 }