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