Introduce gmxpre.h for truly global definitions
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_single / nb_kernel_ElecRF_VdwNone_GeomW3P1_avx_256_single.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
5  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6  * and including many others, as listed in the AUTHORS file in the
7  * top-level source directory and at http://www.gromacs.org.
8  *
9  * GROMACS is free software; you can redistribute it and/or
10  * modify it under the terms of the GNU Lesser General Public License
11  * as published by the Free Software Foundation; either version 2.1
12  * of the License, or (at your option) any later version.
13  *
14  * GROMACS is distributed in the hope that it will be useful,
15  * but WITHOUT ANY WARRANTY; without even the implied warranty of
16  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
17  * Lesser General Public License for more details.
18  *
19  * You should have received a copy of the GNU Lesser General Public
20  * License along with GROMACS; if not, see
21  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
23  *
24  * If you want to redistribute modifications to GROMACS, please
25  * consider that scientific software is very special. Version
26  * control is crucial - bugs must be traceable. We will be happy to
27  * consider code for inclusion in the official distribution, but
28  * derived work must not be called official GROMACS. Details are found
29  * in the README & COPYING files - if they are missing, get the
30  * official version at http://www.gromacs.org.
31  *
32  * To help us fund GROMACS development, we humbly ask that you cite
33  * the research papers on the package. Check out http://www.gromacs.org.
34  */
35 /*
36  * Note: this file was generated by the GROMACS avx_256_single kernel generator.
37  */
38 #include "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_ElecRF_VdwNone_GeomW3P1_VF_avx_256_single
54  * Electrostatics interaction: ReactionField
55  * VdW interaction:            None
56  * Geometry:                   Water3-Particle
57  * Calculate force/pot:        PotentialAndForce
58  */
59 void
60 nb_kernel_ElecRF_VdwNone_GeomW3P1_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     real *           vdwioffsetptr1;
91     __m256           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
92     real *           vdwioffsetptr2;
93     __m256           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
94     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
95     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
96     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
97     __m256           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
98     __m256           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
99     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
100     real             *charge;
101     __m256           dummy_mask,cutoff_mask;
102     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
103     __m256           one     = _mm256_set1_ps(1.0);
104     __m256           two     = _mm256_set1_ps(2.0);
105     x                = xx[0];
106     f                = ff[0];
107
108     nri              = nlist->nri;
109     iinr             = nlist->iinr;
110     jindex           = nlist->jindex;
111     jjnr             = nlist->jjnr;
112     shiftidx         = nlist->shift;
113     gid              = nlist->gid;
114     shiftvec         = fr->shift_vec[0];
115     fshift           = fr->fshift[0];
116     facel            = _mm256_set1_ps(fr->epsfac);
117     charge           = mdatoms->chargeA;
118     krf              = _mm256_set1_ps(fr->ic->k_rf);
119     krf2             = _mm256_set1_ps(fr->ic->k_rf*2.0);
120     crf              = _mm256_set1_ps(fr->ic->c_rf);
121
122     /* Setup water-specific parameters */
123     inr              = nlist->iinr[0];
124     iq0              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
125     iq1              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
126     iq2              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
127
128     /* Avoid stupid compiler warnings */
129     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
130     j_coord_offsetA = 0;
131     j_coord_offsetB = 0;
132     j_coord_offsetC = 0;
133     j_coord_offsetD = 0;
134     j_coord_offsetE = 0;
135     j_coord_offsetF = 0;
136     j_coord_offsetG = 0;
137     j_coord_offsetH = 0;
138
139     outeriter        = 0;
140     inneriter        = 0;
141
142     for(iidx=0;iidx<4*DIM;iidx++)
143     {
144         scratch[iidx] = 0.0;
145     }
146
147     /* Start outer loop over neighborlists */
148     for(iidx=0; iidx<nri; iidx++)
149     {
150         /* Load shift vector for this list */
151         i_shift_offset   = DIM*shiftidx[iidx];
152
153         /* Load limits for loop over neighbors */
154         j_index_start    = jindex[iidx];
155         j_index_end      = jindex[iidx+1];
156
157         /* Get outer coordinate index */
158         inr              = iinr[iidx];
159         i_coord_offset   = DIM*inr;
160
161         /* Load i particle coords and add shift vector */
162         gmx_mm256_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
163                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
164
165         fix0             = _mm256_setzero_ps();
166         fiy0             = _mm256_setzero_ps();
167         fiz0             = _mm256_setzero_ps();
168         fix1             = _mm256_setzero_ps();
169         fiy1             = _mm256_setzero_ps();
170         fiz1             = _mm256_setzero_ps();
171         fix2             = _mm256_setzero_ps();
172         fiy2             = _mm256_setzero_ps();
173         fiz2             = _mm256_setzero_ps();
174
175         /* Reset potential sums */
176         velecsum         = _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             dx10             = _mm256_sub_ps(ix1,jx0);
212             dy10             = _mm256_sub_ps(iy1,jy0);
213             dz10             = _mm256_sub_ps(iz1,jz0);
214             dx20             = _mm256_sub_ps(ix2,jx0);
215             dy20             = _mm256_sub_ps(iy2,jy0);
216             dz20             = _mm256_sub_ps(iz2,jz0);
217
218             /* Calculate squared distance and things based on it */
219             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
220             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
221             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
222
223             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
224             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
225             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
226
227             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
228             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
229             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
230
231             /* Load parameters for j particles */
232             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
233                                                                  charge+jnrC+0,charge+jnrD+0,
234                                                                  charge+jnrE+0,charge+jnrF+0,
235                                                                  charge+jnrG+0,charge+jnrH+0);
236
237             fjx0             = _mm256_setzero_ps();
238             fjy0             = _mm256_setzero_ps();
239             fjz0             = _mm256_setzero_ps();
240
241             /**************************
242              * CALCULATE INTERACTIONS *
243              **************************/
244
245             /* Compute parameters for interactions between i and j atoms */
246             qq00             = _mm256_mul_ps(iq0,jq0);
247
248             /* REACTION-FIELD ELECTROSTATICS */
249             velec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
250             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
251
252             /* Update potential sum for this i atom from the interaction with this j atom. */
253             velecsum         = _mm256_add_ps(velecsum,velec);
254
255             fscal            = felec;
256
257             /* Calculate temporary vectorial force */
258             tx               = _mm256_mul_ps(fscal,dx00);
259             ty               = _mm256_mul_ps(fscal,dy00);
260             tz               = _mm256_mul_ps(fscal,dz00);
261
262             /* Update vectorial force */
263             fix0             = _mm256_add_ps(fix0,tx);
264             fiy0             = _mm256_add_ps(fiy0,ty);
265             fiz0             = _mm256_add_ps(fiz0,tz);
266
267             fjx0             = _mm256_add_ps(fjx0,tx);
268             fjy0             = _mm256_add_ps(fjy0,ty);
269             fjz0             = _mm256_add_ps(fjz0,tz);
270
271             /**************************
272              * CALCULATE INTERACTIONS *
273              **************************/
274
275             /* Compute parameters for interactions between i and j atoms */
276             qq10             = _mm256_mul_ps(iq1,jq0);
277
278             /* REACTION-FIELD ELECTROSTATICS */
279             velec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_add_ps(rinv10,_mm256_mul_ps(krf,rsq10)),crf));
280             felec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_mul_ps(rinv10,rinvsq10),krf2));
281
282             /* Update potential sum for this i atom from the interaction with this j atom. */
283             velecsum         = _mm256_add_ps(velecsum,velec);
284
285             fscal            = felec;
286
287             /* Calculate temporary vectorial force */
288             tx               = _mm256_mul_ps(fscal,dx10);
289             ty               = _mm256_mul_ps(fscal,dy10);
290             tz               = _mm256_mul_ps(fscal,dz10);
291
292             /* Update vectorial force */
293             fix1             = _mm256_add_ps(fix1,tx);
294             fiy1             = _mm256_add_ps(fiy1,ty);
295             fiz1             = _mm256_add_ps(fiz1,tz);
296
297             fjx0             = _mm256_add_ps(fjx0,tx);
298             fjy0             = _mm256_add_ps(fjy0,ty);
299             fjz0             = _mm256_add_ps(fjz0,tz);
300
301             /**************************
302              * CALCULATE INTERACTIONS *
303              **************************/
304
305             /* Compute parameters for interactions between i and j atoms */
306             qq20             = _mm256_mul_ps(iq2,jq0);
307
308             /* REACTION-FIELD ELECTROSTATICS */
309             velec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_add_ps(rinv20,_mm256_mul_ps(krf,rsq20)),crf));
310             felec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_mul_ps(rinv20,rinvsq20),krf2));
311
312             /* Update potential sum for this i atom from the interaction with this j atom. */
313             velecsum         = _mm256_add_ps(velecsum,velec);
314
315             fscal            = felec;
316
317             /* Calculate temporary vectorial force */
318             tx               = _mm256_mul_ps(fscal,dx20);
319             ty               = _mm256_mul_ps(fscal,dy20);
320             tz               = _mm256_mul_ps(fscal,dz20);
321
322             /* Update vectorial force */
323             fix2             = _mm256_add_ps(fix2,tx);
324             fiy2             = _mm256_add_ps(fiy2,ty);
325             fiz2             = _mm256_add_ps(fiz2,tz);
326
327             fjx0             = _mm256_add_ps(fjx0,tx);
328             fjy0             = _mm256_add_ps(fjy0,ty);
329             fjz0             = _mm256_add_ps(fjz0,tz);
330
331             fjptrA             = f+j_coord_offsetA;
332             fjptrB             = f+j_coord_offsetB;
333             fjptrC             = f+j_coord_offsetC;
334             fjptrD             = f+j_coord_offsetD;
335             fjptrE             = f+j_coord_offsetE;
336             fjptrF             = f+j_coord_offsetF;
337             fjptrG             = f+j_coord_offsetG;
338             fjptrH             = f+j_coord_offsetH;
339
340             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
341
342             /* Inner loop uses 99 flops */
343         }
344
345         if(jidx<j_index_end)
346         {
347
348             /* Get j neighbor index, and coordinate index */
349             jnrlistA         = jjnr[jidx];
350             jnrlistB         = jjnr[jidx+1];
351             jnrlistC         = jjnr[jidx+2];
352             jnrlistD         = jjnr[jidx+3];
353             jnrlistE         = jjnr[jidx+4];
354             jnrlistF         = jjnr[jidx+5];
355             jnrlistG         = jjnr[jidx+6];
356             jnrlistH         = jjnr[jidx+7];
357             /* Sign of each element will be negative for non-real atoms.
358              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
359              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
360              */
361             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
362                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
363                                             
364             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
365             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
366             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
367             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
368             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
369             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
370             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
371             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
372             j_coord_offsetA  = DIM*jnrA;
373             j_coord_offsetB  = DIM*jnrB;
374             j_coord_offsetC  = DIM*jnrC;
375             j_coord_offsetD  = DIM*jnrD;
376             j_coord_offsetE  = DIM*jnrE;
377             j_coord_offsetF  = DIM*jnrF;
378             j_coord_offsetG  = DIM*jnrG;
379             j_coord_offsetH  = DIM*jnrH;
380
381             /* load j atom coordinates */
382             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
383                                                  x+j_coord_offsetC,x+j_coord_offsetD,
384                                                  x+j_coord_offsetE,x+j_coord_offsetF,
385                                                  x+j_coord_offsetG,x+j_coord_offsetH,
386                                                  &jx0,&jy0,&jz0);
387
388             /* Calculate displacement vector */
389             dx00             = _mm256_sub_ps(ix0,jx0);
390             dy00             = _mm256_sub_ps(iy0,jy0);
391             dz00             = _mm256_sub_ps(iz0,jz0);
392             dx10             = _mm256_sub_ps(ix1,jx0);
393             dy10             = _mm256_sub_ps(iy1,jy0);
394             dz10             = _mm256_sub_ps(iz1,jz0);
395             dx20             = _mm256_sub_ps(ix2,jx0);
396             dy20             = _mm256_sub_ps(iy2,jy0);
397             dz20             = _mm256_sub_ps(iz2,jz0);
398
399             /* Calculate squared distance and things based on it */
400             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
401             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
402             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
403
404             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
405             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
406             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
407
408             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
409             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
410             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
411
412             /* Load parameters for j particles */
413             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
414                                                                  charge+jnrC+0,charge+jnrD+0,
415                                                                  charge+jnrE+0,charge+jnrF+0,
416                                                                  charge+jnrG+0,charge+jnrH+0);
417
418             fjx0             = _mm256_setzero_ps();
419             fjy0             = _mm256_setzero_ps();
420             fjz0             = _mm256_setzero_ps();
421
422             /**************************
423              * CALCULATE INTERACTIONS *
424              **************************/
425
426             /* Compute parameters for interactions between i and j atoms */
427             qq00             = _mm256_mul_ps(iq0,jq0);
428
429             /* REACTION-FIELD ELECTROSTATICS */
430             velec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
431             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
432
433             /* Update potential sum for this i atom from the interaction with this j atom. */
434             velec            = _mm256_andnot_ps(dummy_mask,velec);
435             velecsum         = _mm256_add_ps(velecsum,velec);
436
437             fscal            = felec;
438
439             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
440
441             /* Calculate temporary vectorial force */
442             tx               = _mm256_mul_ps(fscal,dx00);
443             ty               = _mm256_mul_ps(fscal,dy00);
444             tz               = _mm256_mul_ps(fscal,dz00);
445
446             /* Update vectorial force */
447             fix0             = _mm256_add_ps(fix0,tx);
448             fiy0             = _mm256_add_ps(fiy0,ty);
449             fiz0             = _mm256_add_ps(fiz0,tz);
450
451             fjx0             = _mm256_add_ps(fjx0,tx);
452             fjy0             = _mm256_add_ps(fjy0,ty);
453             fjz0             = _mm256_add_ps(fjz0,tz);
454
455             /**************************
456              * CALCULATE INTERACTIONS *
457              **************************/
458
459             /* Compute parameters for interactions between i and j atoms */
460             qq10             = _mm256_mul_ps(iq1,jq0);
461
462             /* REACTION-FIELD ELECTROSTATICS */
463             velec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_add_ps(rinv10,_mm256_mul_ps(krf,rsq10)),crf));
464             felec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_mul_ps(rinv10,rinvsq10),krf2));
465
466             /* Update potential sum for this i atom from the interaction with this j atom. */
467             velec            = _mm256_andnot_ps(dummy_mask,velec);
468             velecsum         = _mm256_add_ps(velecsum,velec);
469
470             fscal            = felec;
471
472             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
473
474             /* Calculate temporary vectorial force */
475             tx               = _mm256_mul_ps(fscal,dx10);
476             ty               = _mm256_mul_ps(fscal,dy10);
477             tz               = _mm256_mul_ps(fscal,dz10);
478
479             /* Update vectorial force */
480             fix1             = _mm256_add_ps(fix1,tx);
481             fiy1             = _mm256_add_ps(fiy1,ty);
482             fiz1             = _mm256_add_ps(fiz1,tz);
483
484             fjx0             = _mm256_add_ps(fjx0,tx);
485             fjy0             = _mm256_add_ps(fjy0,ty);
486             fjz0             = _mm256_add_ps(fjz0,tz);
487
488             /**************************
489              * CALCULATE INTERACTIONS *
490              **************************/
491
492             /* Compute parameters for interactions between i and j atoms */
493             qq20             = _mm256_mul_ps(iq2,jq0);
494
495             /* REACTION-FIELD ELECTROSTATICS */
496             velec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_add_ps(rinv20,_mm256_mul_ps(krf,rsq20)),crf));
497             felec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_mul_ps(rinv20,rinvsq20),krf2));
498
499             /* Update potential sum for this i atom from the interaction with this j atom. */
500             velec            = _mm256_andnot_ps(dummy_mask,velec);
501             velecsum         = _mm256_add_ps(velecsum,velec);
502
503             fscal            = felec;
504
505             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
506
507             /* Calculate temporary vectorial force */
508             tx               = _mm256_mul_ps(fscal,dx20);
509             ty               = _mm256_mul_ps(fscal,dy20);
510             tz               = _mm256_mul_ps(fscal,dz20);
511
512             /* Update vectorial force */
513             fix2             = _mm256_add_ps(fix2,tx);
514             fiy2             = _mm256_add_ps(fiy2,ty);
515             fiz2             = _mm256_add_ps(fiz2,tz);
516
517             fjx0             = _mm256_add_ps(fjx0,tx);
518             fjy0             = _mm256_add_ps(fjy0,ty);
519             fjz0             = _mm256_add_ps(fjz0,tz);
520
521             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
522             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
523             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
524             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
525             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
526             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
527             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
528             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
529
530             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
531
532             /* Inner loop uses 99 flops */
533         }
534
535         /* End of innermost loop */
536
537         gmx_mm256_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
538                                                  f+i_coord_offset,fshift+i_shift_offset);
539
540         ggid                        = gid[iidx];
541         /* Update potential energies */
542         gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
543
544         /* Increment number of inner iterations */
545         inneriter                  += j_index_end - j_index_start;
546
547         /* Outer loop uses 19 flops */
548     }
549
550     /* Increment number of outer iterations */
551     outeriter        += nri;
552
553     /* Update outer/inner flops */
554
555     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3_VF,outeriter*19 + inneriter*99);
556 }
557 /*
558  * Gromacs nonbonded kernel:   nb_kernel_ElecRF_VdwNone_GeomW3P1_F_avx_256_single
559  * Electrostatics interaction: ReactionField
560  * VdW interaction:            None
561  * Geometry:                   Water3-Particle
562  * Calculate force/pot:        Force
563  */
564 void
565 nb_kernel_ElecRF_VdwNone_GeomW3P1_F_avx_256_single
566                     (t_nblist                    * gmx_restrict       nlist,
567                      rvec                        * gmx_restrict          xx,
568                      rvec                        * gmx_restrict          ff,
569                      t_forcerec                  * gmx_restrict          fr,
570                      t_mdatoms                   * gmx_restrict     mdatoms,
571                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
572                      t_nrnb                      * gmx_restrict        nrnb)
573 {
574     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
575      * just 0 for non-waters.
576      * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
577      * jnr indices corresponding to data put in the four positions in the SIMD register.
578      */
579     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
580     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
581     int              jnrA,jnrB,jnrC,jnrD;
582     int              jnrE,jnrF,jnrG,jnrH;
583     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
584     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
585     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
586     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
587     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
588     real             rcutoff_scalar;
589     real             *shiftvec,*fshift,*x,*f;
590     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
591     real             scratch[4*DIM];
592     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
593     real *           vdwioffsetptr0;
594     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
595     real *           vdwioffsetptr1;
596     __m256           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
597     real *           vdwioffsetptr2;
598     __m256           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
599     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
600     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
601     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
602     __m256           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
603     __m256           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
604     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
605     real             *charge;
606     __m256           dummy_mask,cutoff_mask;
607     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
608     __m256           one     = _mm256_set1_ps(1.0);
609     __m256           two     = _mm256_set1_ps(2.0);
610     x                = xx[0];
611     f                = ff[0];
612
613     nri              = nlist->nri;
614     iinr             = nlist->iinr;
615     jindex           = nlist->jindex;
616     jjnr             = nlist->jjnr;
617     shiftidx         = nlist->shift;
618     gid              = nlist->gid;
619     shiftvec         = fr->shift_vec[0];
620     fshift           = fr->fshift[0];
621     facel            = _mm256_set1_ps(fr->epsfac);
622     charge           = mdatoms->chargeA;
623     krf              = _mm256_set1_ps(fr->ic->k_rf);
624     krf2             = _mm256_set1_ps(fr->ic->k_rf*2.0);
625     crf              = _mm256_set1_ps(fr->ic->c_rf);
626
627     /* Setup water-specific parameters */
628     inr              = nlist->iinr[0];
629     iq0              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
630     iq1              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
631     iq2              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
632
633     /* Avoid stupid compiler warnings */
634     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
635     j_coord_offsetA = 0;
636     j_coord_offsetB = 0;
637     j_coord_offsetC = 0;
638     j_coord_offsetD = 0;
639     j_coord_offsetE = 0;
640     j_coord_offsetF = 0;
641     j_coord_offsetG = 0;
642     j_coord_offsetH = 0;
643
644     outeriter        = 0;
645     inneriter        = 0;
646
647     for(iidx=0;iidx<4*DIM;iidx++)
648     {
649         scratch[iidx] = 0.0;
650     }
651
652     /* Start outer loop over neighborlists */
653     for(iidx=0; iidx<nri; iidx++)
654     {
655         /* Load shift vector for this list */
656         i_shift_offset   = DIM*shiftidx[iidx];
657
658         /* Load limits for loop over neighbors */
659         j_index_start    = jindex[iidx];
660         j_index_end      = jindex[iidx+1];
661
662         /* Get outer coordinate index */
663         inr              = iinr[iidx];
664         i_coord_offset   = DIM*inr;
665
666         /* Load i particle coords and add shift vector */
667         gmx_mm256_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
668                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
669
670         fix0             = _mm256_setzero_ps();
671         fiy0             = _mm256_setzero_ps();
672         fiz0             = _mm256_setzero_ps();
673         fix1             = _mm256_setzero_ps();
674         fiy1             = _mm256_setzero_ps();
675         fiz1             = _mm256_setzero_ps();
676         fix2             = _mm256_setzero_ps();
677         fiy2             = _mm256_setzero_ps();
678         fiz2             = _mm256_setzero_ps();
679
680         /* Start inner kernel loop */
681         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
682         {
683
684             /* Get j neighbor index, and coordinate index */
685             jnrA             = jjnr[jidx];
686             jnrB             = jjnr[jidx+1];
687             jnrC             = jjnr[jidx+2];
688             jnrD             = jjnr[jidx+3];
689             jnrE             = jjnr[jidx+4];
690             jnrF             = jjnr[jidx+5];
691             jnrG             = jjnr[jidx+6];
692             jnrH             = jjnr[jidx+7];
693             j_coord_offsetA  = DIM*jnrA;
694             j_coord_offsetB  = DIM*jnrB;
695             j_coord_offsetC  = DIM*jnrC;
696             j_coord_offsetD  = DIM*jnrD;
697             j_coord_offsetE  = DIM*jnrE;
698             j_coord_offsetF  = DIM*jnrF;
699             j_coord_offsetG  = DIM*jnrG;
700             j_coord_offsetH  = DIM*jnrH;
701
702             /* load j atom coordinates */
703             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
704                                                  x+j_coord_offsetC,x+j_coord_offsetD,
705                                                  x+j_coord_offsetE,x+j_coord_offsetF,
706                                                  x+j_coord_offsetG,x+j_coord_offsetH,
707                                                  &jx0,&jy0,&jz0);
708
709             /* Calculate displacement vector */
710             dx00             = _mm256_sub_ps(ix0,jx0);
711             dy00             = _mm256_sub_ps(iy0,jy0);
712             dz00             = _mm256_sub_ps(iz0,jz0);
713             dx10             = _mm256_sub_ps(ix1,jx0);
714             dy10             = _mm256_sub_ps(iy1,jy0);
715             dz10             = _mm256_sub_ps(iz1,jz0);
716             dx20             = _mm256_sub_ps(ix2,jx0);
717             dy20             = _mm256_sub_ps(iy2,jy0);
718             dz20             = _mm256_sub_ps(iz2,jz0);
719
720             /* Calculate squared distance and things based on it */
721             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
722             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
723             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
724
725             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
726             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
727             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
728
729             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
730             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
731             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
732
733             /* Load parameters for j particles */
734             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
735                                                                  charge+jnrC+0,charge+jnrD+0,
736                                                                  charge+jnrE+0,charge+jnrF+0,
737                                                                  charge+jnrG+0,charge+jnrH+0);
738
739             fjx0             = _mm256_setzero_ps();
740             fjy0             = _mm256_setzero_ps();
741             fjz0             = _mm256_setzero_ps();
742
743             /**************************
744              * CALCULATE INTERACTIONS *
745              **************************/
746
747             /* Compute parameters for interactions between i and j atoms */
748             qq00             = _mm256_mul_ps(iq0,jq0);
749
750             /* REACTION-FIELD ELECTROSTATICS */
751             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
752
753             fscal            = felec;
754
755             /* Calculate temporary vectorial force */
756             tx               = _mm256_mul_ps(fscal,dx00);
757             ty               = _mm256_mul_ps(fscal,dy00);
758             tz               = _mm256_mul_ps(fscal,dz00);
759
760             /* Update vectorial force */
761             fix0             = _mm256_add_ps(fix0,tx);
762             fiy0             = _mm256_add_ps(fiy0,ty);
763             fiz0             = _mm256_add_ps(fiz0,tz);
764
765             fjx0             = _mm256_add_ps(fjx0,tx);
766             fjy0             = _mm256_add_ps(fjy0,ty);
767             fjz0             = _mm256_add_ps(fjz0,tz);
768
769             /**************************
770              * CALCULATE INTERACTIONS *
771              **************************/
772
773             /* Compute parameters for interactions between i and j atoms */
774             qq10             = _mm256_mul_ps(iq1,jq0);
775
776             /* REACTION-FIELD ELECTROSTATICS */
777             felec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_mul_ps(rinv10,rinvsq10),krf2));
778
779             fscal            = felec;
780
781             /* Calculate temporary vectorial force */
782             tx               = _mm256_mul_ps(fscal,dx10);
783             ty               = _mm256_mul_ps(fscal,dy10);
784             tz               = _mm256_mul_ps(fscal,dz10);
785
786             /* Update vectorial force */
787             fix1             = _mm256_add_ps(fix1,tx);
788             fiy1             = _mm256_add_ps(fiy1,ty);
789             fiz1             = _mm256_add_ps(fiz1,tz);
790
791             fjx0             = _mm256_add_ps(fjx0,tx);
792             fjy0             = _mm256_add_ps(fjy0,ty);
793             fjz0             = _mm256_add_ps(fjz0,tz);
794
795             /**************************
796              * CALCULATE INTERACTIONS *
797              **************************/
798
799             /* Compute parameters for interactions between i and j atoms */
800             qq20             = _mm256_mul_ps(iq2,jq0);
801
802             /* REACTION-FIELD ELECTROSTATICS */
803             felec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_mul_ps(rinv20,rinvsq20),krf2));
804
805             fscal            = felec;
806
807             /* Calculate temporary vectorial force */
808             tx               = _mm256_mul_ps(fscal,dx20);
809             ty               = _mm256_mul_ps(fscal,dy20);
810             tz               = _mm256_mul_ps(fscal,dz20);
811
812             /* Update vectorial force */
813             fix2             = _mm256_add_ps(fix2,tx);
814             fiy2             = _mm256_add_ps(fiy2,ty);
815             fiz2             = _mm256_add_ps(fiz2,tz);
816
817             fjx0             = _mm256_add_ps(fjx0,tx);
818             fjy0             = _mm256_add_ps(fjy0,ty);
819             fjz0             = _mm256_add_ps(fjz0,tz);
820
821             fjptrA             = f+j_coord_offsetA;
822             fjptrB             = f+j_coord_offsetB;
823             fjptrC             = f+j_coord_offsetC;
824             fjptrD             = f+j_coord_offsetD;
825             fjptrE             = f+j_coord_offsetE;
826             fjptrF             = f+j_coord_offsetF;
827             fjptrG             = f+j_coord_offsetG;
828             fjptrH             = f+j_coord_offsetH;
829
830             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
831
832             /* Inner loop uses 84 flops */
833         }
834
835         if(jidx<j_index_end)
836         {
837
838             /* Get j neighbor index, and coordinate index */
839             jnrlistA         = jjnr[jidx];
840             jnrlistB         = jjnr[jidx+1];
841             jnrlistC         = jjnr[jidx+2];
842             jnrlistD         = jjnr[jidx+3];
843             jnrlistE         = jjnr[jidx+4];
844             jnrlistF         = jjnr[jidx+5];
845             jnrlistG         = jjnr[jidx+6];
846             jnrlistH         = jjnr[jidx+7];
847             /* Sign of each element will be negative for non-real atoms.
848              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
849              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
850              */
851             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
852                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
853                                             
854             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
855             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
856             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
857             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
858             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
859             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
860             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
861             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
862             j_coord_offsetA  = DIM*jnrA;
863             j_coord_offsetB  = DIM*jnrB;
864             j_coord_offsetC  = DIM*jnrC;
865             j_coord_offsetD  = DIM*jnrD;
866             j_coord_offsetE  = DIM*jnrE;
867             j_coord_offsetF  = DIM*jnrF;
868             j_coord_offsetG  = DIM*jnrG;
869             j_coord_offsetH  = DIM*jnrH;
870
871             /* load j atom coordinates */
872             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
873                                                  x+j_coord_offsetC,x+j_coord_offsetD,
874                                                  x+j_coord_offsetE,x+j_coord_offsetF,
875                                                  x+j_coord_offsetG,x+j_coord_offsetH,
876                                                  &jx0,&jy0,&jz0);
877
878             /* Calculate displacement vector */
879             dx00             = _mm256_sub_ps(ix0,jx0);
880             dy00             = _mm256_sub_ps(iy0,jy0);
881             dz00             = _mm256_sub_ps(iz0,jz0);
882             dx10             = _mm256_sub_ps(ix1,jx0);
883             dy10             = _mm256_sub_ps(iy1,jy0);
884             dz10             = _mm256_sub_ps(iz1,jz0);
885             dx20             = _mm256_sub_ps(ix2,jx0);
886             dy20             = _mm256_sub_ps(iy2,jy0);
887             dz20             = _mm256_sub_ps(iz2,jz0);
888
889             /* Calculate squared distance and things based on it */
890             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
891             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
892             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
893
894             rinv00           = gmx_mm256_invsqrt_ps(rsq00);
895             rinv10           = gmx_mm256_invsqrt_ps(rsq10);
896             rinv20           = gmx_mm256_invsqrt_ps(rsq20);
897
898             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
899             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
900             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
901
902             /* Load parameters for j particles */
903             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
904                                                                  charge+jnrC+0,charge+jnrD+0,
905                                                                  charge+jnrE+0,charge+jnrF+0,
906                                                                  charge+jnrG+0,charge+jnrH+0);
907
908             fjx0             = _mm256_setzero_ps();
909             fjy0             = _mm256_setzero_ps();
910             fjz0             = _mm256_setzero_ps();
911
912             /**************************
913              * CALCULATE INTERACTIONS *
914              **************************/
915
916             /* Compute parameters for interactions between i and j atoms */
917             qq00             = _mm256_mul_ps(iq0,jq0);
918
919             /* REACTION-FIELD ELECTROSTATICS */
920             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
921
922             fscal            = felec;
923
924             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
925
926             /* Calculate temporary vectorial force */
927             tx               = _mm256_mul_ps(fscal,dx00);
928             ty               = _mm256_mul_ps(fscal,dy00);
929             tz               = _mm256_mul_ps(fscal,dz00);
930
931             /* Update vectorial force */
932             fix0             = _mm256_add_ps(fix0,tx);
933             fiy0             = _mm256_add_ps(fiy0,ty);
934             fiz0             = _mm256_add_ps(fiz0,tz);
935
936             fjx0             = _mm256_add_ps(fjx0,tx);
937             fjy0             = _mm256_add_ps(fjy0,ty);
938             fjz0             = _mm256_add_ps(fjz0,tz);
939
940             /**************************
941              * CALCULATE INTERACTIONS *
942              **************************/
943
944             /* Compute parameters for interactions between i and j atoms */
945             qq10             = _mm256_mul_ps(iq1,jq0);
946
947             /* REACTION-FIELD ELECTROSTATICS */
948             felec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_mul_ps(rinv10,rinvsq10),krf2));
949
950             fscal            = felec;
951
952             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
953
954             /* Calculate temporary vectorial force */
955             tx               = _mm256_mul_ps(fscal,dx10);
956             ty               = _mm256_mul_ps(fscal,dy10);
957             tz               = _mm256_mul_ps(fscal,dz10);
958
959             /* Update vectorial force */
960             fix1             = _mm256_add_ps(fix1,tx);
961             fiy1             = _mm256_add_ps(fiy1,ty);
962             fiz1             = _mm256_add_ps(fiz1,tz);
963
964             fjx0             = _mm256_add_ps(fjx0,tx);
965             fjy0             = _mm256_add_ps(fjy0,ty);
966             fjz0             = _mm256_add_ps(fjz0,tz);
967
968             /**************************
969              * CALCULATE INTERACTIONS *
970              **************************/
971
972             /* Compute parameters for interactions between i and j atoms */
973             qq20             = _mm256_mul_ps(iq2,jq0);
974
975             /* REACTION-FIELD ELECTROSTATICS */
976             felec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_mul_ps(rinv20,rinvsq20),krf2));
977
978             fscal            = felec;
979
980             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
981
982             /* Calculate temporary vectorial force */
983             tx               = _mm256_mul_ps(fscal,dx20);
984             ty               = _mm256_mul_ps(fscal,dy20);
985             tz               = _mm256_mul_ps(fscal,dz20);
986
987             /* Update vectorial force */
988             fix2             = _mm256_add_ps(fix2,tx);
989             fiy2             = _mm256_add_ps(fiy2,ty);
990             fiz2             = _mm256_add_ps(fiz2,tz);
991
992             fjx0             = _mm256_add_ps(fjx0,tx);
993             fjy0             = _mm256_add_ps(fjy0,ty);
994             fjz0             = _mm256_add_ps(fjz0,tz);
995
996             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
997             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
998             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
999             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1000             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1001             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1002             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1003             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1004
1005             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
1006
1007             /* Inner loop uses 84 flops */
1008         }
1009
1010         /* End of innermost loop */
1011
1012         gmx_mm256_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1013                                                  f+i_coord_offset,fshift+i_shift_offset);
1014
1015         /* Increment number of inner iterations */
1016         inneriter                  += j_index_end - j_index_start;
1017
1018         /* Outer loop uses 18 flops */
1019     }
1020
1021     /* Increment number of outer iterations */
1022     outeriter        += nri;
1023
1024     /* Update outer/inner flops */
1025
1026     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3_F,outeriter*18 + inneriter*84);
1027 }