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