Merge release-4-6 into master
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_double / nb_kernel_ElecRF_VdwNone_GeomW3W3_avx_256_double.c
1 /*
2  * Note: this file was generated by the Gromacs avx_256_double kernel generator.
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  * Copyright (c) 2001-2012, The GROMACS Development Team
9  *
10  * Gromacs is a library for molecular simulation and trajectory analysis,
11  * written by Erik Lindahl, David van der Spoel, Berk Hess, and others - for
12  * a full list of developers and information, check out http://www.gromacs.org
13  *
14  * This program is free software; you can redistribute it and/or modify it under
15  * the terms of the GNU Lesser General Public License as published by the Free
16  * Software Foundation; either version 2 of the License, or (at your option) any
17  * later version.
18  *
19  * To help fund GROMACS development, we humbly ask that you cite
20  * the papers people have written on it - you can find them on the website.
21  */
22 #ifdef HAVE_CONFIG_H
23 #include <config.h>
24 #endif
25
26 #include <math.h>
27
28 #include "../nb_kernel.h"
29 #include "types/simple.h"
30 #include "vec.h"
31 #include "nrnb.h"
32
33 #include "gmx_math_x86_avx_256_double.h"
34 #include "kernelutil_x86_avx_256_double.h"
35
36 /*
37  * Gromacs nonbonded kernel:   nb_kernel_ElecRF_VdwNone_GeomW3W3_VF_avx_256_double
38  * Electrostatics interaction: ReactionField
39  * VdW interaction:            None
40  * Geometry:                   Water3-Water3
41  * Calculate force/pot:        PotentialAndForce
42  */
43 void
44 nb_kernel_ElecRF_VdwNone_GeomW3W3_VF_avx_256_double
45                     (t_nblist * gmx_restrict                nlist,
46                      rvec * gmx_restrict                    xx,
47                      rvec * gmx_restrict                    ff,
48                      t_forcerec * gmx_restrict              fr,
49                      t_mdatoms * gmx_restrict               mdatoms,
50                      nb_kernel_data_t * gmx_restrict        kernel_data,
51                      t_nrnb * gmx_restrict                  nrnb)
52 {
53     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
54      * just 0 for non-waters.
55      * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
56      * jnr indices corresponding to data put in the four positions in the SIMD register.
57      */
58     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
59     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
60     int              jnrA,jnrB,jnrC,jnrD;
61     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
62     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
63     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
64     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
65     real             rcutoff_scalar;
66     real             *shiftvec,*fshift,*x,*f;
67     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
68     real             scratch[4*DIM];
69     __m256d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
70     real *           vdwioffsetptr0;
71     __m256d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
72     real *           vdwioffsetptr1;
73     __m256d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
74     real *           vdwioffsetptr2;
75     __m256d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
76     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
77     __m256d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
78     int              vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D;
79     __m256d          jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
80     int              vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D;
81     __m256d          jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
82     __m256d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
83     __m256d          dx01,dy01,dz01,rsq01,rinv01,rinvsq01,r01,qq01,c6_01,c12_01;
84     __m256d          dx02,dy02,dz02,rsq02,rinv02,rinvsq02,r02,qq02,c6_02,c12_02;
85     __m256d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
86     __m256d          dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
87     __m256d          dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
88     __m256d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
89     __m256d          dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
90     __m256d          dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
91     __m256d          velec,felec,velecsum,facel,crf,krf,krf2;
92     real             *charge;
93     __m256d          dummy_mask,cutoff_mask;
94     __m128           tmpmask0,tmpmask1;
95     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
96     __m256d          one     = _mm256_set1_pd(1.0);
97     __m256d          two     = _mm256_set1_pd(2.0);
98     x                = xx[0];
99     f                = ff[0];
100
101     nri              = nlist->nri;
102     iinr             = nlist->iinr;
103     jindex           = nlist->jindex;
104     jjnr             = nlist->jjnr;
105     shiftidx         = nlist->shift;
106     gid              = nlist->gid;
107     shiftvec         = fr->shift_vec[0];
108     fshift           = fr->fshift[0];
109     facel            = _mm256_set1_pd(fr->epsfac);
110     charge           = mdatoms->chargeA;
111     krf              = _mm256_set1_pd(fr->ic->k_rf);
112     krf2             = _mm256_set1_pd(fr->ic->k_rf*2.0);
113     crf              = _mm256_set1_pd(fr->ic->c_rf);
114
115     /* Setup water-specific parameters */
116     inr              = nlist->iinr[0];
117     iq0              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
118     iq1              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+1]));
119     iq2              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+2]));
120
121     jq0              = _mm256_set1_pd(charge[inr+0]);
122     jq1              = _mm256_set1_pd(charge[inr+1]);
123     jq2              = _mm256_set1_pd(charge[inr+2]);
124     qq00             = _mm256_mul_pd(iq0,jq0);
125     qq01             = _mm256_mul_pd(iq0,jq1);
126     qq02             = _mm256_mul_pd(iq0,jq2);
127     qq10             = _mm256_mul_pd(iq1,jq0);
128     qq11             = _mm256_mul_pd(iq1,jq1);
129     qq12             = _mm256_mul_pd(iq1,jq2);
130     qq20             = _mm256_mul_pd(iq2,jq0);
131     qq21             = _mm256_mul_pd(iq2,jq1);
132     qq22             = _mm256_mul_pd(iq2,jq2);
133
134     /* Avoid stupid compiler warnings */
135     jnrA = jnrB = jnrC = jnrD = 0;
136     j_coord_offsetA = 0;
137     j_coord_offsetB = 0;
138     j_coord_offsetC = 0;
139     j_coord_offsetD = 0;
140
141     outeriter        = 0;
142     inneriter        = 0;
143
144     for(iidx=0;iidx<4*DIM;iidx++)
145     {
146         scratch[iidx] = 0.0;
147     }
148
149     /* Start outer loop over neighborlists */
150     for(iidx=0; iidx<nri; iidx++)
151     {
152         /* Load shift vector for this list */
153         i_shift_offset   = DIM*shiftidx[iidx];
154
155         /* Load limits for loop over neighbors */
156         j_index_start    = jindex[iidx];
157         j_index_end      = jindex[iidx+1];
158
159         /* Get outer coordinate index */
160         inr              = iinr[iidx];
161         i_coord_offset   = DIM*inr;
162
163         /* Load i particle coords and add shift vector */
164         gmx_mm256_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
165                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
166
167         fix0             = _mm256_setzero_pd();
168         fiy0             = _mm256_setzero_pd();
169         fiz0             = _mm256_setzero_pd();
170         fix1             = _mm256_setzero_pd();
171         fiy1             = _mm256_setzero_pd();
172         fiz1             = _mm256_setzero_pd();
173         fix2             = _mm256_setzero_pd();
174         fiy2             = _mm256_setzero_pd();
175         fiz2             = _mm256_setzero_pd();
176
177         /* Reset potential sums */
178         velecsum         = _mm256_setzero_pd();
179
180         /* Start inner kernel loop */
181         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
182         {
183
184             /* Get j neighbor index, and coordinate index */
185             jnrA             = jjnr[jidx];
186             jnrB             = jjnr[jidx+1];
187             jnrC             = jjnr[jidx+2];
188             jnrD             = jjnr[jidx+3];
189             j_coord_offsetA  = DIM*jnrA;
190             j_coord_offsetB  = DIM*jnrB;
191             j_coord_offsetC  = DIM*jnrC;
192             j_coord_offsetD  = DIM*jnrD;
193
194             /* load j atom coordinates */
195             gmx_mm256_load_3rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
196                                                  x+j_coord_offsetC,x+j_coord_offsetD,
197                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
198
199             /* Calculate displacement vector */
200             dx00             = _mm256_sub_pd(ix0,jx0);
201             dy00             = _mm256_sub_pd(iy0,jy0);
202             dz00             = _mm256_sub_pd(iz0,jz0);
203             dx01             = _mm256_sub_pd(ix0,jx1);
204             dy01             = _mm256_sub_pd(iy0,jy1);
205             dz01             = _mm256_sub_pd(iz0,jz1);
206             dx02             = _mm256_sub_pd(ix0,jx2);
207             dy02             = _mm256_sub_pd(iy0,jy2);
208             dz02             = _mm256_sub_pd(iz0,jz2);
209             dx10             = _mm256_sub_pd(ix1,jx0);
210             dy10             = _mm256_sub_pd(iy1,jy0);
211             dz10             = _mm256_sub_pd(iz1,jz0);
212             dx11             = _mm256_sub_pd(ix1,jx1);
213             dy11             = _mm256_sub_pd(iy1,jy1);
214             dz11             = _mm256_sub_pd(iz1,jz1);
215             dx12             = _mm256_sub_pd(ix1,jx2);
216             dy12             = _mm256_sub_pd(iy1,jy2);
217             dz12             = _mm256_sub_pd(iz1,jz2);
218             dx20             = _mm256_sub_pd(ix2,jx0);
219             dy20             = _mm256_sub_pd(iy2,jy0);
220             dz20             = _mm256_sub_pd(iz2,jz0);
221             dx21             = _mm256_sub_pd(ix2,jx1);
222             dy21             = _mm256_sub_pd(iy2,jy1);
223             dz21             = _mm256_sub_pd(iz2,jz1);
224             dx22             = _mm256_sub_pd(ix2,jx2);
225             dy22             = _mm256_sub_pd(iy2,jy2);
226             dz22             = _mm256_sub_pd(iz2,jz2);
227
228             /* Calculate squared distance and things based on it */
229             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
230             rsq01            = gmx_mm256_calc_rsq_pd(dx01,dy01,dz01);
231             rsq02            = gmx_mm256_calc_rsq_pd(dx02,dy02,dz02);
232             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
233             rsq11            = gmx_mm256_calc_rsq_pd(dx11,dy11,dz11);
234             rsq12            = gmx_mm256_calc_rsq_pd(dx12,dy12,dz12);
235             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
236             rsq21            = gmx_mm256_calc_rsq_pd(dx21,dy21,dz21);
237             rsq22            = gmx_mm256_calc_rsq_pd(dx22,dy22,dz22);
238
239             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
240             rinv01           = gmx_mm256_invsqrt_pd(rsq01);
241             rinv02           = gmx_mm256_invsqrt_pd(rsq02);
242             rinv10           = gmx_mm256_invsqrt_pd(rsq10);
243             rinv11           = gmx_mm256_invsqrt_pd(rsq11);
244             rinv12           = gmx_mm256_invsqrt_pd(rsq12);
245             rinv20           = gmx_mm256_invsqrt_pd(rsq20);
246             rinv21           = gmx_mm256_invsqrt_pd(rsq21);
247             rinv22           = gmx_mm256_invsqrt_pd(rsq22);
248
249             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
250             rinvsq01         = _mm256_mul_pd(rinv01,rinv01);
251             rinvsq02         = _mm256_mul_pd(rinv02,rinv02);
252             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
253             rinvsq11         = _mm256_mul_pd(rinv11,rinv11);
254             rinvsq12         = _mm256_mul_pd(rinv12,rinv12);
255             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
256             rinvsq21         = _mm256_mul_pd(rinv21,rinv21);
257             rinvsq22         = _mm256_mul_pd(rinv22,rinv22);
258
259             fjx0             = _mm256_setzero_pd();
260             fjy0             = _mm256_setzero_pd();
261             fjz0             = _mm256_setzero_pd();
262             fjx1             = _mm256_setzero_pd();
263             fjy1             = _mm256_setzero_pd();
264             fjz1             = _mm256_setzero_pd();
265             fjx2             = _mm256_setzero_pd();
266             fjy2             = _mm256_setzero_pd();
267             fjz2             = _mm256_setzero_pd();
268
269             /**************************
270              * CALCULATE INTERACTIONS *
271              **************************/
272
273             /* REACTION-FIELD ELECTROSTATICS */
274             velec            = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_add_pd(rinv00,_mm256_mul_pd(krf,rsq00)),crf));
275             felec            = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_mul_pd(rinv00,rinvsq00),krf2));
276
277             /* Update potential sum for this i atom from the interaction with this j atom. */
278             velecsum         = _mm256_add_pd(velecsum,velec);
279
280             fscal            = felec;
281
282             /* Calculate temporary vectorial force */
283             tx               = _mm256_mul_pd(fscal,dx00);
284             ty               = _mm256_mul_pd(fscal,dy00);
285             tz               = _mm256_mul_pd(fscal,dz00);
286
287             /* Update vectorial force */
288             fix0             = _mm256_add_pd(fix0,tx);
289             fiy0             = _mm256_add_pd(fiy0,ty);
290             fiz0             = _mm256_add_pd(fiz0,tz);
291
292             fjx0             = _mm256_add_pd(fjx0,tx);
293             fjy0             = _mm256_add_pd(fjy0,ty);
294             fjz0             = _mm256_add_pd(fjz0,tz);
295
296             /**************************
297              * CALCULATE INTERACTIONS *
298              **************************/
299
300             /* REACTION-FIELD ELECTROSTATICS */
301             velec            = _mm256_mul_pd(qq01,_mm256_sub_pd(_mm256_add_pd(rinv01,_mm256_mul_pd(krf,rsq01)),crf));
302             felec            = _mm256_mul_pd(qq01,_mm256_sub_pd(_mm256_mul_pd(rinv01,rinvsq01),krf2));
303
304             /* Update potential sum for this i atom from the interaction with this j atom. */
305             velecsum         = _mm256_add_pd(velecsum,velec);
306
307             fscal            = felec;
308
309             /* Calculate temporary vectorial force */
310             tx               = _mm256_mul_pd(fscal,dx01);
311             ty               = _mm256_mul_pd(fscal,dy01);
312             tz               = _mm256_mul_pd(fscal,dz01);
313
314             /* Update vectorial force */
315             fix0             = _mm256_add_pd(fix0,tx);
316             fiy0             = _mm256_add_pd(fiy0,ty);
317             fiz0             = _mm256_add_pd(fiz0,tz);
318
319             fjx1             = _mm256_add_pd(fjx1,tx);
320             fjy1             = _mm256_add_pd(fjy1,ty);
321             fjz1             = _mm256_add_pd(fjz1,tz);
322
323             /**************************
324              * CALCULATE INTERACTIONS *
325              **************************/
326
327             /* REACTION-FIELD ELECTROSTATICS */
328             velec            = _mm256_mul_pd(qq02,_mm256_sub_pd(_mm256_add_pd(rinv02,_mm256_mul_pd(krf,rsq02)),crf));
329             felec            = _mm256_mul_pd(qq02,_mm256_sub_pd(_mm256_mul_pd(rinv02,rinvsq02),krf2));
330
331             /* Update potential sum for this i atom from the interaction with this j atom. */
332             velecsum         = _mm256_add_pd(velecsum,velec);
333
334             fscal            = felec;
335
336             /* Calculate temporary vectorial force */
337             tx               = _mm256_mul_pd(fscal,dx02);
338             ty               = _mm256_mul_pd(fscal,dy02);
339             tz               = _mm256_mul_pd(fscal,dz02);
340
341             /* Update vectorial force */
342             fix0             = _mm256_add_pd(fix0,tx);
343             fiy0             = _mm256_add_pd(fiy0,ty);
344             fiz0             = _mm256_add_pd(fiz0,tz);
345
346             fjx2             = _mm256_add_pd(fjx2,tx);
347             fjy2             = _mm256_add_pd(fjy2,ty);
348             fjz2             = _mm256_add_pd(fjz2,tz);
349
350             /**************************
351              * CALCULATE INTERACTIONS *
352              **************************/
353
354             /* REACTION-FIELD ELECTROSTATICS */
355             velec            = _mm256_mul_pd(qq10,_mm256_sub_pd(_mm256_add_pd(rinv10,_mm256_mul_pd(krf,rsq10)),crf));
356             felec            = _mm256_mul_pd(qq10,_mm256_sub_pd(_mm256_mul_pd(rinv10,rinvsq10),krf2));
357
358             /* Update potential sum for this i atom from the interaction with this j atom. */
359             velecsum         = _mm256_add_pd(velecsum,velec);
360
361             fscal            = felec;
362
363             /* Calculate temporary vectorial force */
364             tx               = _mm256_mul_pd(fscal,dx10);
365             ty               = _mm256_mul_pd(fscal,dy10);
366             tz               = _mm256_mul_pd(fscal,dz10);
367
368             /* Update vectorial force */
369             fix1             = _mm256_add_pd(fix1,tx);
370             fiy1             = _mm256_add_pd(fiy1,ty);
371             fiz1             = _mm256_add_pd(fiz1,tz);
372
373             fjx0             = _mm256_add_pd(fjx0,tx);
374             fjy0             = _mm256_add_pd(fjy0,ty);
375             fjz0             = _mm256_add_pd(fjz0,tz);
376
377             /**************************
378              * CALCULATE INTERACTIONS *
379              **************************/
380
381             /* REACTION-FIELD ELECTROSTATICS */
382             velec            = _mm256_mul_pd(qq11,_mm256_sub_pd(_mm256_add_pd(rinv11,_mm256_mul_pd(krf,rsq11)),crf));
383             felec            = _mm256_mul_pd(qq11,_mm256_sub_pd(_mm256_mul_pd(rinv11,rinvsq11),krf2));
384
385             /* Update potential sum for this i atom from the interaction with this j atom. */
386             velecsum         = _mm256_add_pd(velecsum,velec);
387
388             fscal            = felec;
389
390             /* Calculate temporary vectorial force */
391             tx               = _mm256_mul_pd(fscal,dx11);
392             ty               = _mm256_mul_pd(fscal,dy11);
393             tz               = _mm256_mul_pd(fscal,dz11);
394
395             /* Update vectorial force */
396             fix1             = _mm256_add_pd(fix1,tx);
397             fiy1             = _mm256_add_pd(fiy1,ty);
398             fiz1             = _mm256_add_pd(fiz1,tz);
399
400             fjx1             = _mm256_add_pd(fjx1,tx);
401             fjy1             = _mm256_add_pd(fjy1,ty);
402             fjz1             = _mm256_add_pd(fjz1,tz);
403
404             /**************************
405              * CALCULATE INTERACTIONS *
406              **************************/
407
408             /* REACTION-FIELD ELECTROSTATICS */
409             velec            = _mm256_mul_pd(qq12,_mm256_sub_pd(_mm256_add_pd(rinv12,_mm256_mul_pd(krf,rsq12)),crf));
410             felec            = _mm256_mul_pd(qq12,_mm256_sub_pd(_mm256_mul_pd(rinv12,rinvsq12),krf2));
411
412             /* Update potential sum for this i atom from the interaction with this j atom. */
413             velecsum         = _mm256_add_pd(velecsum,velec);
414
415             fscal            = felec;
416
417             /* Calculate temporary vectorial force */
418             tx               = _mm256_mul_pd(fscal,dx12);
419             ty               = _mm256_mul_pd(fscal,dy12);
420             tz               = _mm256_mul_pd(fscal,dz12);
421
422             /* Update vectorial force */
423             fix1             = _mm256_add_pd(fix1,tx);
424             fiy1             = _mm256_add_pd(fiy1,ty);
425             fiz1             = _mm256_add_pd(fiz1,tz);
426
427             fjx2             = _mm256_add_pd(fjx2,tx);
428             fjy2             = _mm256_add_pd(fjy2,ty);
429             fjz2             = _mm256_add_pd(fjz2,tz);
430
431             /**************************
432              * CALCULATE INTERACTIONS *
433              **************************/
434
435             /* REACTION-FIELD ELECTROSTATICS */
436             velec            = _mm256_mul_pd(qq20,_mm256_sub_pd(_mm256_add_pd(rinv20,_mm256_mul_pd(krf,rsq20)),crf));
437             felec            = _mm256_mul_pd(qq20,_mm256_sub_pd(_mm256_mul_pd(rinv20,rinvsq20),krf2));
438
439             /* Update potential sum for this i atom from the interaction with this j atom. */
440             velecsum         = _mm256_add_pd(velecsum,velec);
441
442             fscal            = felec;
443
444             /* Calculate temporary vectorial force */
445             tx               = _mm256_mul_pd(fscal,dx20);
446             ty               = _mm256_mul_pd(fscal,dy20);
447             tz               = _mm256_mul_pd(fscal,dz20);
448
449             /* Update vectorial force */
450             fix2             = _mm256_add_pd(fix2,tx);
451             fiy2             = _mm256_add_pd(fiy2,ty);
452             fiz2             = _mm256_add_pd(fiz2,tz);
453
454             fjx0             = _mm256_add_pd(fjx0,tx);
455             fjy0             = _mm256_add_pd(fjy0,ty);
456             fjz0             = _mm256_add_pd(fjz0,tz);
457
458             /**************************
459              * CALCULATE INTERACTIONS *
460              **************************/
461
462             /* REACTION-FIELD ELECTROSTATICS */
463             velec            = _mm256_mul_pd(qq21,_mm256_sub_pd(_mm256_add_pd(rinv21,_mm256_mul_pd(krf,rsq21)),crf));
464             felec            = _mm256_mul_pd(qq21,_mm256_sub_pd(_mm256_mul_pd(rinv21,rinvsq21),krf2));
465
466             /* Update potential sum for this i atom from the interaction with this j atom. */
467             velecsum         = _mm256_add_pd(velecsum,velec);
468
469             fscal            = felec;
470
471             /* Calculate temporary vectorial force */
472             tx               = _mm256_mul_pd(fscal,dx21);
473             ty               = _mm256_mul_pd(fscal,dy21);
474             tz               = _mm256_mul_pd(fscal,dz21);
475
476             /* Update vectorial force */
477             fix2             = _mm256_add_pd(fix2,tx);
478             fiy2             = _mm256_add_pd(fiy2,ty);
479             fiz2             = _mm256_add_pd(fiz2,tz);
480
481             fjx1             = _mm256_add_pd(fjx1,tx);
482             fjy1             = _mm256_add_pd(fjy1,ty);
483             fjz1             = _mm256_add_pd(fjz1,tz);
484
485             /**************************
486              * CALCULATE INTERACTIONS *
487              **************************/
488
489             /* REACTION-FIELD ELECTROSTATICS */
490             velec            = _mm256_mul_pd(qq22,_mm256_sub_pd(_mm256_add_pd(rinv22,_mm256_mul_pd(krf,rsq22)),crf));
491             felec            = _mm256_mul_pd(qq22,_mm256_sub_pd(_mm256_mul_pd(rinv22,rinvsq22),krf2));
492
493             /* Update potential sum for this i atom from the interaction with this j atom. */
494             velecsum         = _mm256_add_pd(velecsum,velec);
495
496             fscal            = felec;
497
498             /* Calculate temporary vectorial force */
499             tx               = _mm256_mul_pd(fscal,dx22);
500             ty               = _mm256_mul_pd(fscal,dy22);
501             tz               = _mm256_mul_pd(fscal,dz22);
502
503             /* Update vectorial force */
504             fix2             = _mm256_add_pd(fix2,tx);
505             fiy2             = _mm256_add_pd(fiy2,ty);
506             fiz2             = _mm256_add_pd(fiz2,tz);
507
508             fjx2             = _mm256_add_pd(fjx2,tx);
509             fjy2             = _mm256_add_pd(fjy2,ty);
510             fjz2             = _mm256_add_pd(fjz2,tz);
511
512             fjptrA             = f+j_coord_offsetA;
513             fjptrB             = f+j_coord_offsetB;
514             fjptrC             = f+j_coord_offsetC;
515             fjptrD             = f+j_coord_offsetD;
516
517             gmx_mm256_decrement_3rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
518                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
519
520             /* Inner loop uses 288 flops */
521         }
522
523         if(jidx<j_index_end)
524         {
525
526             /* Get j neighbor index, and coordinate index */
527             jnrlistA         = jjnr[jidx];
528             jnrlistB         = jjnr[jidx+1];
529             jnrlistC         = jjnr[jidx+2];
530             jnrlistD         = jjnr[jidx+3];
531             /* Sign of each element will be negative for non-real atoms.
532              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
533              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
534              */
535             tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
536
537             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
538             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
539             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
540
541             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
542             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
543             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
544             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
545             j_coord_offsetA  = DIM*jnrA;
546             j_coord_offsetB  = DIM*jnrB;
547             j_coord_offsetC  = DIM*jnrC;
548             j_coord_offsetD  = DIM*jnrD;
549
550             /* load j atom coordinates */
551             gmx_mm256_load_3rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
552                                                  x+j_coord_offsetC,x+j_coord_offsetD,
553                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
554
555             /* Calculate displacement vector */
556             dx00             = _mm256_sub_pd(ix0,jx0);
557             dy00             = _mm256_sub_pd(iy0,jy0);
558             dz00             = _mm256_sub_pd(iz0,jz0);
559             dx01             = _mm256_sub_pd(ix0,jx1);
560             dy01             = _mm256_sub_pd(iy0,jy1);
561             dz01             = _mm256_sub_pd(iz0,jz1);
562             dx02             = _mm256_sub_pd(ix0,jx2);
563             dy02             = _mm256_sub_pd(iy0,jy2);
564             dz02             = _mm256_sub_pd(iz0,jz2);
565             dx10             = _mm256_sub_pd(ix1,jx0);
566             dy10             = _mm256_sub_pd(iy1,jy0);
567             dz10             = _mm256_sub_pd(iz1,jz0);
568             dx11             = _mm256_sub_pd(ix1,jx1);
569             dy11             = _mm256_sub_pd(iy1,jy1);
570             dz11             = _mm256_sub_pd(iz1,jz1);
571             dx12             = _mm256_sub_pd(ix1,jx2);
572             dy12             = _mm256_sub_pd(iy1,jy2);
573             dz12             = _mm256_sub_pd(iz1,jz2);
574             dx20             = _mm256_sub_pd(ix2,jx0);
575             dy20             = _mm256_sub_pd(iy2,jy0);
576             dz20             = _mm256_sub_pd(iz2,jz0);
577             dx21             = _mm256_sub_pd(ix2,jx1);
578             dy21             = _mm256_sub_pd(iy2,jy1);
579             dz21             = _mm256_sub_pd(iz2,jz1);
580             dx22             = _mm256_sub_pd(ix2,jx2);
581             dy22             = _mm256_sub_pd(iy2,jy2);
582             dz22             = _mm256_sub_pd(iz2,jz2);
583
584             /* Calculate squared distance and things based on it */
585             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
586             rsq01            = gmx_mm256_calc_rsq_pd(dx01,dy01,dz01);
587             rsq02            = gmx_mm256_calc_rsq_pd(dx02,dy02,dz02);
588             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
589             rsq11            = gmx_mm256_calc_rsq_pd(dx11,dy11,dz11);
590             rsq12            = gmx_mm256_calc_rsq_pd(dx12,dy12,dz12);
591             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
592             rsq21            = gmx_mm256_calc_rsq_pd(dx21,dy21,dz21);
593             rsq22            = gmx_mm256_calc_rsq_pd(dx22,dy22,dz22);
594
595             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
596             rinv01           = gmx_mm256_invsqrt_pd(rsq01);
597             rinv02           = gmx_mm256_invsqrt_pd(rsq02);
598             rinv10           = gmx_mm256_invsqrt_pd(rsq10);
599             rinv11           = gmx_mm256_invsqrt_pd(rsq11);
600             rinv12           = gmx_mm256_invsqrt_pd(rsq12);
601             rinv20           = gmx_mm256_invsqrt_pd(rsq20);
602             rinv21           = gmx_mm256_invsqrt_pd(rsq21);
603             rinv22           = gmx_mm256_invsqrt_pd(rsq22);
604
605             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
606             rinvsq01         = _mm256_mul_pd(rinv01,rinv01);
607             rinvsq02         = _mm256_mul_pd(rinv02,rinv02);
608             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
609             rinvsq11         = _mm256_mul_pd(rinv11,rinv11);
610             rinvsq12         = _mm256_mul_pd(rinv12,rinv12);
611             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
612             rinvsq21         = _mm256_mul_pd(rinv21,rinv21);
613             rinvsq22         = _mm256_mul_pd(rinv22,rinv22);
614
615             fjx0             = _mm256_setzero_pd();
616             fjy0             = _mm256_setzero_pd();
617             fjz0             = _mm256_setzero_pd();
618             fjx1             = _mm256_setzero_pd();
619             fjy1             = _mm256_setzero_pd();
620             fjz1             = _mm256_setzero_pd();
621             fjx2             = _mm256_setzero_pd();
622             fjy2             = _mm256_setzero_pd();
623             fjz2             = _mm256_setzero_pd();
624
625             /**************************
626              * CALCULATE INTERACTIONS *
627              **************************/
628
629             /* REACTION-FIELD ELECTROSTATICS */
630             velec            = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_add_pd(rinv00,_mm256_mul_pd(krf,rsq00)),crf));
631             felec            = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_mul_pd(rinv00,rinvsq00),krf2));
632
633             /* Update potential sum for this i atom from the interaction with this j atom. */
634             velec            = _mm256_andnot_pd(dummy_mask,velec);
635             velecsum         = _mm256_add_pd(velecsum,velec);
636
637             fscal            = felec;
638
639             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
640
641             /* Calculate temporary vectorial force */
642             tx               = _mm256_mul_pd(fscal,dx00);
643             ty               = _mm256_mul_pd(fscal,dy00);
644             tz               = _mm256_mul_pd(fscal,dz00);
645
646             /* Update vectorial force */
647             fix0             = _mm256_add_pd(fix0,tx);
648             fiy0             = _mm256_add_pd(fiy0,ty);
649             fiz0             = _mm256_add_pd(fiz0,tz);
650
651             fjx0             = _mm256_add_pd(fjx0,tx);
652             fjy0             = _mm256_add_pd(fjy0,ty);
653             fjz0             = _mm256_add_pd(fjz0,tz);
654
655             /**************************
656              * CALCULATE INTERACTIONS *
657              **************************/
658
659             /* REACTION-FIELD ELECTROSTATICS */
660             velec            = _mm256_mul_pd(qq01,_mm256_sub_pd(_mm256_add_pd(rinv01,_mm256_mul_pd(krf,rsq01)),crf));
661             felec            = _mm256_mul_pd(qq01,_mm256_sub_pd(_mm256_mul_pd(rinv01,rinvsq01),krf2));
662
663             /* Update potential sum for this i atom from the interaction with this j atom. */
664             velec            = _mm256_andnot_pd(dummy_mask,velec);
665             velecsum         = _mm256_add_pd(velecsum,velec);
666
667             fscal            = felec;
668
669             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
670
671             /* Calculate temporary vectorial force */
672             tx               = _mm256_mul_pd(fscal,dx01);
673             ty               = _mm256_mul_pd(fscal,dy01);
674             tz               = _mm256_mul_pd(fscal,dz01);
675
676             /* Update vectorial force */
677             fix0             = _mm256_add_pd(fix0,tx);
678             fiy0             = _mm256_add_pd(fiy0,ty);
679             fiz0             = _mm256_add_pd(fiz0,tz);
680
681             fjx1             = _mm256_add_pd(fjx1,tx);
682             fjy1             = _mm256_add_pd(fjy1,ty);
683             fjz1             = _mm256_add_pd(fjz1,tz);
684
685             /**************************
686              * CALCULATE INTERACTIONS *
687              **************************/
688
689             /* REACTION-FIELD ELECTROSTATICS */
690             velec            = _mm256_mul_pd(qq02,_mm256_sub_pd(_mm256_add_pd(rinv02,_mm256_mul_pd(krf,rsq02)),crf));
691             felec            = _mm256_mul_pd(qq02,_mm256_sub_pd(_mm256_mul_pd(rinv02,rinvsq02),krf2));
692
693             /* Update potential sum for this i atom from the interaction with this j atom. */
694             velec            = _mm256_andnot_pd(dummy_mask,velec);
695             velecsum         = _mm256_add_pd(velecsum,velec);
696
697             fscal            = felec;
698
699             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
700
701             /* Calculate temporary vectorial force */
702             tx               = _mm256_mul_pd(fscal,dx02);
703             ty               = _mm256_mul_pd(fscal,dy02);
704             tz               = _mm256_mul_pd(fscal,dz02);
705
706             /* Update vectorial force */
707             fix0             = _mm256_add_pd(fix0,tx);
708             fiy0             = _mm256_add_pd(fiy0,ty);
709             fiz0             = _mm256_add_pd(fiz0,tz);
710
711             fjx2             = _mm256_add_pd(fjx2,tx);
712             fjy2             = _mm256_add_pd(fjy2,ty);
713             fjz2             = _mm256_add_pd(fjz2,tz);
714
715             /**************************
716              * CALCULATE INTERACTIONS *
717              **************************/
718
719             /* REACTION-FIELD ELECTROSTATICS */
720             velec            = _mm256_mul_pd(qq10,_mm256_sub_pd(_mm256_add_pd(rinv10,_mm256_mul_pd(krf,rsq10)),crf));
721             felec            = _mm256_mul_pd(qq10,_mm256_sub_pd(_mm256_mul_pd(rinv10,rinvsq10),krf2));
722
723             /* Update potential sum for this i atom from the interaction with this j atom. */
724             velec            = _mm256_andnot_pd(dummy_mask,velec);
725             velecsum         = _mm256_add_pd(velecsum,velec);
726
727             fscal            = felec;
728
729             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
730
731             /* Calculate temporary vectorial force */
732             tx               = _mm256_mul_pd(fscal,dx10);
733             ty               = _mm256_mul_pd(fscal,dy10);
734             tz               = _mm256_mul_pd(fscal,dz10);
735
736             /* Update vectorial force */
737             fix1             = _mm256_add_pd(fix1,tx);
738             fiy1             = _mm256_add_pd(fiy1,ty);
739             fiz1             = _mm256_add_pd(fiz1,tz);
740
741             fjx0             = _mm256_add_pd(fjx0,tx);
742             fjy0             = _mm256_add_pd(fjy0,ty);
743             fjz0             = _mm256_add_pd(fjz0,tz);
744
745             /**************************
746              * CALCULATE INTERACTIONS *
747              **************************/
748
749             /* REACTION-FIELD ELECTROSTATICS */
750             velec            = _mm256_mul_pd(qq11,_mm256_sub_pd(_mm256_add_pd(rinv11,_mm256_mul_pd(krf,rsq11)),crf));
751             felec            = _mm256_mul_pd(qq11,_mm256_sub_pd(_mm256_mul_pd(rinv11,rinvsq11),krf2));
752
753             /* Update potential sum for this i atom from the interaction with this j atom. */
754             velec            = _mm256_andnot_pd(dummy_mask,velec);
755             velecsum         = _mm256_add_pd(velecsum,velec);
756
757             fscal            = felec;
758
759             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
760
761             /* Calculate temporary vectorial force */
762             tx               = _mm256_mul_pd(fscal,dx11);
763             ty               = _mm256_mul_pd(fscal,dy11);
764             tz               = _mm256_mul_pd(fscal,dz11);
765
766             /* Update vectorial force */
767             fix1             = _mm256_add_pd(fix1,tx);
768             fiy1             = _mm256_add_pd(fiy1,ty);
769             fiz1             = _mm256_add_pd(fiz1,tz);
770
771             fjx1             = _mm256_add_pd(fjx1,tx);
772             fjy1             = _mm256_add_pd(fjy1,ty);
773             fjz1             = _mm256_add_pd(fjz1,tz);
774
775             /**************************
776              * CALCULATE INTERACTIONS *
777              **************************/
778
779             /* REACTION-FIELD ELECTROSTATICS */
780             velec            = _mm256_mul_pd(qq12,_mm256_sub_pd(_mm256_add_pd(rinv12,_mm256_mul_pd(krf,rsq12)),crf));
781             felec            = _mm256_mul_pd(qq12,_mm256_sub_pd(_mm256_mul_pd(rinv12,rinvsq12),krf2));
782
783             /* Update potential sum for this i atom from the interaction with this j atom. */
784             velec            = _mm256_andnot_pd(dummy_mask,velec);
785             velecsum         = _mm256_add_pd(velecsum,velec);
786
787             fscal            = felec;
788
789             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
790
791             /* Calculate temporary vectorial force */
792             tx               = _mm256_mul_pd(fscal,dx12);
793             ty               = _mm256_mul_pd(fscal,dy12);
794             tz               = _mm256_mul_pd(fscal,dz12);
795
796             /* Update vectorial force */
797             fix1             = _mm256_add_pd(fix1,tx);
798             fiy1             = _mm256_add_pd(fiy1,ty);
799             fiz1             = _mm256_add_pd(fiz1,tz);
800
801             fjx2             = _mm256_add_pd(fjx2,tx);
802             fjy2             = _mm256_add_pd(fjy2,ty);
803             fjz2             = _mm256_add_pd(fjz2,tz);
804
805             /**************************
806              * CALCULATE INTERACTIONS *
807              **************************/
808
809             /* REACTION-FIELD ELECTROSTATICS */
810             velec            = _mm256_mul_pd(qq20,_mm256_sub_pd(_mm256_add_pd(rinv20,_mm256_mul_pd(krf,rsq20)),crf));
811             felec            = _mm256_mul_pd(qq20,_mm256_sub_pd(_mm256_mul_pd(rinv20,rinvsq20),krf2));
812
813             /* Update potential sum for this i atom from the interaction with this j atom. */
814             velec            = _mm256_andnot_pd(dummy_mask,velec);
815             velecsum         = _mm256_add_pd(velecsum,velec);
816
817             fscal            = felec;
818
819             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
820
821             /* Calculate temporary vectorial force */
822             tx               = _mm256_mul_pd(fscal,dx20);
823             ty               = _mm256_mul_pd(fscal,dy20);
824             tz               = _mm256_mul_pd(fscal,dz20);
825
826             /* Update vectorial force */
827             fix2             = _mm256_add_pd(fix2,tx);
828             fiy2             = _mm256_add_pd(fiy2,ty);
829             fiz2             = _mm256_add_pd(fiz2,tz);
830
831             fjx0             = _mm256_add_pd(fjx0,tx);
832             fjy0             = _mm256_add_pd(fjy0,ty);
833             fjz0             = _mm256_add_pd(fjz0,tz);
834
835             /**************************
836              * CALCULATE INTERACTIONS *
837              **************************/
838
839             /* REACTION-FIELD ELECTROSTATICS */
840             velec            = _mm256_mul_pd(qq21,_mm256_sub_pd(_mm256_add_pd(rinv21,_mm256_mul_pd(krf,rsq21)),crf));
841             felec            = _mm256_mul_pd(qq21,_mm256_sub_pd(_mm256_mul_pd(rinv21,rinvsq21),krf2));
842
843             /* Update potential sum for this i atom from the interaction with this j atom. */
844             velec            = _mm256_andnot_pd(dummy_mask,velec);
845             velecsum         = _mm256_add_pd(velecsum,velec);
846
847             fscal            = felec;
848
849             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
850
851             /* Calculate temporary vectorial force */
852             tx               = _mm256_mul_pd(fscal,dx21);
853             ty               = _mm256_mul_pd(fscal,dy21);
854             tz               = _mm256_mul_pd(fscal,dz21);
855
856             /* Update vectorial force */
857             fix2             = _mm256_add_pd(fix2,tx);
858             fiy2             = _mm256_add_pd(fiy2,ty);
859             fiz2             = _mm256_add_pd(fiz2,tz);
860
861             fjx1             = _mm256_add_pd(fjx1,tx);
862             fjy1             = _mm256_add_pd(fjy1,ty);
863             fjz1             = _mm256_add_pd(fjz1,tz);
864
865             /**************************
866              * CALCULATE INTERACTIONS *
867              **************************/
868
869             /* REACTION-FIELD ELECTROSTATICS */
870             velec            = _mm256_mul_pd(qq22,_mm256_sub_pd(_mm256_add_pd(rinv22,_mm256_mul_pd(krf,rsq22)),crf));
871             felec            = _mm256_mul_pd(qq22,_mm256_sub_pd(_mm256_mul_pd(rinv22,rinvsq22),krf2));
872
873             /* Update potential sum for this i atom from the interaction with this j atom. */
874             velec            = _mm256_andnot_pd(dummy_mask,velec);
875             velecsum         = _mm256_add_pd(velecsum,velec);
876
877             fscal            = felec;
878
879             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
880
881             /* Calculate temporary vectorial force */
882             tx               = _mm256_mul_pd(fscal,dx22);
883             ty               = _mm256_mul_pd(fscal,dy22);
884             tz               = _mm256_mul_pd(fscal,dz22);
885
886             /* Update vectorial force */
887             fix2             = _mm256_add_pd(fix2,tx);
888             fiy2             = _mm256_add_pd(fiy2,ty);
889             fiz2             = _mm256_add_pd(fiz2,tz);
890
891             fjx2             = _mm256_add_pd(fjx2,tx);
892             fjy2             = _mm256_add_pd(fjy2,ty);
893             fjz2             = _mm256_add_pd(fjz2,tz);
894
895             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
896             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
897             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
898             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
899
900             gmx_mm256_decrement_3rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
901                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
902
903             /* Inner loop uses 288 flops */
904         }
905
906         /* End of innermost loop */
907
908         gmx_mm256_update_iforce_3atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
909                                                  f+i_coord_offset,fshift+i_shift_offset);
910
911         ggid                        = gid[iidx];
912         /* Update potential energies */
913         gmx_mm256_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
914
915         /* Increment number of inner iterations */
916         inneriter                  += j_index_end - j_index_start;
917
918         /* Outer loop uses 19 flops */
919     }
920
921     /* Increment number of outer iterations */
922     outeriter        += nri;
923
924     /* Update outer/inner flops */
925
926     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3W3_VF,outeriter*19 + inneriter*288);
927 }
928 /*
929  * Gromacs nonbonded kernel:   nb_kernel_ElecRF_VdwNone_GeomW3W3_F_avx_256_double
930  * Electrostatics interaction: ReactionField
931  * VdW interaction:            None
932  * Geometry:                   Water3-Water3
933  * Calculate force/pot:        Force
934  */
935 void
936 nb_kernel_ElecRF_VdwNone_GeomW3W3_F_avx_256_double
937                     (t_nblist * gmx_restrict                nlist,
938                      rvec * gmx_restrict                    xx,
939                      rvec * gmx_restrict                    ff,
940                      t_forcerec * gmx_restrict              fr,
941                      t_mdatoms * gmx_restrict               mdatoms,
942                      nb_kernel_data_t * gmx_restrict        kernel_data,
943                      t_nrnb * gmx_restrict                  nrnb)
944 {
945     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
946      * just 0 for non-waters.
947      * Suffixes A,B,C,D refer to j loop unrolling done with AVX, e.g. for the four different
948      * jnr indices corresponding to data put in the four positions in the SIMD register.
949      */
950     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
951     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
952     int              jnrA,jnrB,jnrC,jnrD;
953     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
954     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
955     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
956     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
957     real             rcutoff_scalar;
958     real             *shiftvec,*fshift,*x,*f;
959     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
960     real             scratch[4*DIM];
961     __m256d          tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
962     real *           vdwioffsetptr0;
963     __m256d          ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
964     real *           vdwioffsetptr1;
965     __m256d          ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
966     real *           vdwioffsetptr2;
967     __m256d          ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
968     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D;
969     __m256d          jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
970     int              vdwjidx1A,vdwjidx1B,vdwjidx1C,vdwjidx1D;
971     __m256d          jx1,jy1,jz1,fjx1,fjy1,fjz1,jq1,isaj1;
972     int              vdwjidx2A,vdwjidx2B,vdwjidx2C,vdwjidx2D;
973     __m256d          jx2,jy2,jz2,fjx2,fjy2,fjz2,jq2,isaj2;
974     __m256d          dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
975     __m256d          dx01,dy01,dz01,rsq01,rinv01,rinvsq01,r01,qq01,c6_01,c12_01;
976     __m256d          dx02,dy02,dz02,rsq02,rinv02,rinvsq02,r02,qq02,c6_02,c12_02;
977     __m256d          dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
978     __m256d          dx11,dy11,dz11,rsq11,rinv11,rinvsq11,r11,qq11,c6_11,c12_11;
979     __m256d          dx12,dy12,dz12,rsq12,rinv12,rinvsq12,r12,qq12,c6_12,c12_12;
980     __m256d          dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
981     __m256d          dx21,dy21,dz21,rsq21,rinv21,rinvsq21,r21,qq21,c6_21,c12_21;
982     __m256d          dx22,dy22,dz22,rsq22,rinv22,rinvsq22,r22,qq22,c6_22,c12_22;
983     __m256d          velec,felec,velecsum,facel,crf,krf,krf2;
984     real             *charge;
985     __m256d          dummy_mask,cutoff_mask;
986     __m128           tmpmask0,tmpmask1;
987     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
988     __m256d          one     = _mm256_set1_pd(1.0);
989     __m256d          two     = _mm256_set1_pd(2.0);
990     x                = xx[0];
991     f                = ff[0];
992
993     nri              = nlist->nri;
994     iinr             = nlist->iinr;
995     jindex           = nlist->jindex;
996     jjnr             = nlist->jjnr;
997     shiftidx         = nlist->shift;
998     gid              = nlist->gid;
999     shiftvec         = fr->shift_vec[0];
1000     fshift           = fr->fshift[0];
1001     facel            = _mm256_set1_pd(fr->epsfac);
1002     charge           = mdatoms->chargeA;
1003     krf              = _mm256_set1_pd(fr->ic->k_rf);
1004     krf2             = _mm256_set1_pd(fr->ic->k_rf*2.0);
1005     crf              = _mm256_set1_pd(fr->ic->c_rf);
1006
1007     /* Setup water-specific parameters */
1008     inr              = nlist->iinr[0];
1009     iq0              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+0]));
1010     iq1              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+1]));
1011     iq2              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+2]));
1012
1013     jq0              = _mm256_set1_pd(charge[inr+0]);
1014     jq1              = _mm256_set1_pd(charge[inr+1]);
1015     jq2              = _mm256_set1_pd(charge[inr+2]);
1016     qq00             = _mm256_mul_pd(iq0,jq0);
1017     qq01             = _mm256_mul_pd(iq0,jq1);
1018     qq02             = _mm256_mul_pd(iq0,jq2);
1019     qq10             = _mm256_mul_pd(iq1,jq0);
1020     qq11             = _mm256_mul_pd(iq1,jq1);
1021     qq12             = _mm256_mul_pd(iq1,jq2);
1022     qq20             = _mm256_mul_pd(iq2,jq0);
1023     qq21             = _mm256_mul_pd(iq2,jq1);
1024     qq22             = _mm256_mul_pd(iq2,jq2);
1025
1026     /* Avoid stupid compiler warnings */
1027     jnrA = jnrB = jnrC = jnrD = 0;
1028     j_coord_offsetA = 0;
1029     j_coord_offsetB = 0;
1030     j_coord_offsetC = 0;
1031     j_coord_offsetD = 0;
1032
1033     outeriter        = 0;
1034     inneriter        = 0;
1035
1036     for(iidx=0;iidx<4*DIM;iidx++)
1037     {
1038         scratch[iidx] = 0.0;
1039     }
1040
1041     /* Start outer loop over neighborlists */
1042     for(iidx=0; iidx<nri; iidx++)
1043     {
1044         /* Load shift vector for this list */
1045         i_shift_offset   = DIM*shiftidx[iidx];
1046
1047         /* Load limits for loop over neighbors */
1048         j_index_start    = jindex[iidx];
1049         j_index_end      = jindex[iidx+1];
1050
1051         /* Get outer coordinate index */
1052         inr              = iinr[iidx];
1053         i_coord_offset   = DIM*inr;
1054
1055         /* Load i particle coords and add shift vector */
1056         gmx_mm256_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
1057                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
1058
1059         fix0             = _mm256_setzero_pd();
1060         fiy0             = _mm256_setzero_pd();
1061         fiz0             = _mm256_setzero_pd();
1062         fix1             = _mm256_setzero_pd();
1063         fiy1             = _mm256_setzero_pd();
1064         fiz1             = _mm256_setzero_pd();
1065         fix2             = _mm256_setzero_pd();
1066         fiy2             = _mm256_setzero_pd();
1067         fiz2             = _mm256_setzero_pd();
1068
1069         /* Start inner kernel loop */
1070         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
1071         {
1072
1073             /* Get j neighbor index, and coordinate index */
1074             jnrA             = jjnr[jidx];
1075             jnrB             = jjnr[jidx+1];
1076             jnrC             = jjnr[jidx+2];
1077             jnrD             = jjnr[jidx+3];
1078             j_coord_offsetA  = DIM*jnrA;
1079             j_coord_offsetB  = DIM*jnrB;
1080             j_coord_offsetC  = DIM*jnrC;
1081             j_coord_offsetD  = DIM*jnrD;
1082
1083             /* load j atom coordinates */
1084             gmx_mm256_load_3rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
1085                                                  x+j_coord_offsetC,x+j_coord_offsetD,
1086                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
1087
1088             /* Calculate displacement vector */
1089             dx00             = _mm256_sub_pd(ix0,jx0);
1090             dy00             = _mm256_sub_pd(iy0,jy0);
1091             dz00             = _mm256_sub_pd(iz0,jz0);
1092             dx01             = _mm256_sub_pd(ix0,jx1);
1093             dy01             = _mm256_sub_pd(iy0,jy1);
1094             dz01             = _mm256_sub_pd(iz0,jz1);
1095             dx02             = _mm256_sub_pd(ix0,jx2);
1096             dy02             = _mm256_sub_pd(iy0,jy2);
1097             dz02             = _mm256_sub_pd(iz0,jz2);
1098             dx10             = _mm256_sub_pd(ix1,jx0);
1099             dy10             = _mm256_sub_pd(iy1,jy0);
1100             dz10             = _mm256_sub_pd(iz1,jz0);
1101             dx11             = _mm256_sub_pd(ix1,jx1);
1102             dy11             = _mm256_sub_pd(iy1,jy1);
1103             dz11             = _mm256_sub_pd(iz1,jz1);
1104             dx12             = _mm256_sub_pd(ix1,jx2);
1105             dy12             = _mm256_sub_pd(iy1,jy2);
1106             dz12             = _mm256_sub_pd(iz1,jz2);
1107             dx20             = _mm256_sub_pd(ix2,jx0);
1108             dy20             = _mm256_sub_pd(iy2,jy0);
1109             dz20             = _mm256_sub_pd(iz2,jz0);
1110             dx21             = _mm256_sub_pd(ix2,jx1);
1111             dy21             = _mm256_sub_pd(iy2,jy1);
1112             dz21             = _mm256_sub_pd(iz2,jz1);
1113             dx22             = _mm256_sub_pd(ix2,jx2);
1114             dy22             = _mm256_sub_pd(iy2,jy2);
1115             dz22             = _mm256_sub_pd(iz2,jz2);
1116
1117             /* Calculate squared distance and things based on it */
1118             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
1119             rsq01            = gmx_mm256_calc_rsq_pd(dx01,dy01,dz01);
1120             rsq02            = gmx_mm256_calc_rsq_pd(dx02,dy02,dz02);
1121             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
1122             rsq11            = gmx_mm256_calc_rsq_pd(dx11,dy11,dz11);
1123             rsq12            = gmx_mm256_calc_rsq_pd(dx12,dy12,dz12);
1124             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
1125             rsq21            = gmx_mm256_calc_rsq_pd(dx21,dy21,dz21);
1126             rsq22            = gmx_mm256_calc_rsq_pd(dx22,dy22,dz22);
1127
1128             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
1129             rinv01           = gmx_mm256_invsqrt_pd(rsq01);
1130             rinv02           = gmx_mm256_invsqrt_pd(rsq02);
1131             rinv10           = gmx_mm256_invsqrt_pd(rsq10);
1132             rinv11           = gmx_mm256_invsqrt_pd(rsq11);
1133             rinv12           = gmx_mm256_invsqrt_pd(rsq12);
1134             rinv20           = gmx_mm256_invsqrt_pd(rsq20);
1135             rinv21           = gmx_mm256_invsqrt_pd(rsq21);
1136             rinv22           = gmx_mm256_invsqrt_pd(rsq22);
1137
1138             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
1139             rinvsq01         = _mm256_mul_pd(rinv01,rinv01);
1140             rinvsq02         = _mm256_mul_pd(rinv02,rinv02);
1141             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
1142             rinvsq11         = _mm256_mul_pd(rinv11,rinv11);
1143             rinvsq12         = _mm256_mul_pd(rinv12,rinv12);
1144             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
1145             rinvsq21         = _mm256_mul_pd(rinv21,rinv21);
1146             rinvsq22         = _mm256_mul_pd(rinv22,rinv22);
1147
1148             fjx0             = _mm256_setzero_pd();
1149             fjy0             = _mm256_setzero_pd();
1150             fjz0             = _mm256_setzero_pd();
1151             fjx1             = _mm256_setzero_pd();
1152             fjy1             = _mm256_setzero_pd();
1153             fjz1             = _mm256_setzero_pd();
1154             fjx2             = _mm256_setzero_pd();
1155             fjy2             = _mm256_setzero_pd();
1156             fjz2             = _mm256_setzero_pd();
1157
1158             /**************************
1159              * CALCULATE INTERACTIONS *
1160              **************************/
1161
1162             /* REACTION-FIELD ELECTROSTATICS */
1163             felec            = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_mul_pd(rinv00,rinvsq00),krf2));
1164
1165             fscal            = felec;
1166
1167             /* Calculate temporary vectorial force */
1168             tx               = _mm256_mul_pd(fscal,dx00);
1169             ty               = _mm256_mul_pd(fscal,dy00);
1170             tz               = _mm256_mul_pd(fscal,dz00);
1171
1172             /* Update vectorial force */
1173             fix0             = _mm256_add_pd(fix0,tx);
1174             fiy0             = _mm256_add_pd(fiy0,ty);
1175             fiz0             = _mm256_add_pd(fiz0,tz);
1176
1177             fjx0             = _mm256_add_pd(fjx0,tx);
1178             fjy0             = _mm256_add_pd(fjy0,ty);
1179             fjz0             = _mm256_add_pd(fjz0,tz);
1180
1181             /**************************
1182              * CALCULATE INTERACTIONS *
1183              **************************/
1184
1185             /* REACTION-FIELD ELECTROSTATICS */
1186             felec            = _mm256_mul_pd(qq01,_mm256_sub_pd(_mm256_mul_pd(rinv01,rinvsq01),krf2));
1187
1188             fscal            = felec;
1189
1190             /* Calculate temporary vectorial force */
1191             tx               = _mm256_mul_pd(fscal,dx01);
1192             ty               = _mm256_mul_pd(fscal,dy01);
1193             tz               = _mm256_mul_pd(fscal,dz01);
1194
1195             /* Update vectorial force */
1196             fix0             = _mm256_add_pd(fix0,tx);
1197             fiy0             = _mm256_add_pd(fiy0,ty);
1198             fiz0             = _mm256_add_pd(fiz0,tz);
1199
1200             fjx1             = _mm256_add_pd(fjx1,tx);
1201             fjy1             = _mm256_add_pd(fjy1,ty);
1202             fjz1             = _mm256_add_pd(fjz1,tz);
1203
1204             /**************************
1205              * CALCULATE INTERACTIONS *
1206              **************************/
1207
1208             /* REACTION-FIELD ELECTROSTATICS */
1209             felec            = _mm256_mul_pd(qq02,_mm256_sub_pd(_mm256_mul_pd(rinv02,rinvsq02),krf2));
1210
1211             fscal            = felec;
1212
1213             /* Calculate temporary vectorial force */
1214             tx               = _mm256_mul_pd(fscal,dx02);
1215             ty               = _mm256_mul_pd(fscal,dy02);
1216             tz               = _mm256_mul_pd(fscal,dz02);
1217
1218             /* Update vectorial force */
1219             fix0             = _mm256_add_pd(fix0,tx);
1220             fiy0             = _mm256_add_pd(fiy0,ty);
1221             fiz0             = _mm256_add_pd(fiz0,tz);
1222
1223             fjx2             = _mm256_add_pd(fjx2,tx);
1224             fjy2             = _mm256_add_pd(fjy2,ty);
1225             fjz2             = _mm256_add_pd(fjz2,tz);
1226
1227             /**************************
1228              * CALCULATE INTERACTIONS *
1229              **************************/
1230
1231             /* REACTION-FIELD ELECTROSTATICS */
1232             felec            = _mm256_mul_pd(qq10,_mm256_sub_pd(_mm256_mul_pd(rinv10,rinvsq10),krf2));
1233
1234             fscal            = felec;
1235
1236             /* Calculate temporary vectorial force */
1237             tx               = _mm256_mul_pd(fscal,dx10);
1238             ty               = _mm256_mul_pd(fscal,dy10);
1239             tz               = _mm256_mul_pd(fscal,dz10);
1240
1241             /* Update vectorial force */
1242             fix1             = _mm256_add_pd(fix1,tx);
1243             fiy1             = _mm256_add_pd(fiy1,ty);
1244             fiz1             = _mm256_add_pd(fiz1,tz);
1245
1246             fjx0             = _mm256_add_pd(fjx0,tx);
1247             fjy0             = _mm256_add_pd(fjy0,ty);
1248             fjz0             = _mm256_add_pd(fjz0,tz);
1249
1250             /**************************
1251              * CALCULATE INTERACTIONS *
1252              **************************/
1253
1254             /* REACTION-FIELD ELECTROSTATICS */
1255             felec            = _mm256_mul_pd(qq11,_mm256_sub_pd(_mm256_mul_pd(rinv11,rinvsq11),krf2));
1256
1257             fscal            = felec;
1258
1259             /* Calculate temporary vectorial force */
1260             tx               = _mm256_mul_pd(fscal,dx11);
1261             ty               = _mm256_mul_pd(fscal,dy11);
1262             tz               = _mm256_mul_pd(fscal,dz11);
1263
1264             /* Update vectorial force */
1265             fix1             = _mm256_add_pd(fix1,tx);
1266             fiy1             = _mm256_add_pd(fiy1,ty);
1267             fiz1             = _mm256_add_pd(fiz1,tz);
1268
1269             fjx1             = _mm256_add_pd(fjx1,tx);
1270             fjy1             = _mm256_add_pd(fjy1,ty);
1271             fjz1             = _mm256_add_pd(fjz1,tz);
1272
1273             /**************************
1274              * CALCULATE INTERACTIONS *
1275              **************************/
1276
1277             /* REACTION-FIELD ELECTROSTATICS */
1278             felec            = _mm256_mul_pd(qq12,_mm256_sub_pd(_mm256_mul_pd(rinv12,rinvsq12),krf2));
1279
1280             fscal            = felec;
1281
1282             /* Calculate temporary vectorial force */
1283             tx               = _mm256_mul_pd(fscal,dx12);
1284             ty               = _mm256_mul_pd(fscal,dy12);
1285             tz               = _mm256_mul_pd(fscal,dz12);
1286
1287             /* Update vectorial force */
1288             fix1             = _mm256_add_pd(fix1,tx);
1289             fiy1             = _mm256_add_pd(fiy1,ty);
1290             fiz1             = _mm256_add_pd(fiz1,tz);
1291
1292             fjx2             = _mm256_add_pd(fjx2,tx);
1293             fjy2             = _mm256_add_pd(fjy2,ty);
1294             fjz2             = _mm256_add_pd(fjz2,tz);
1295
1296             /**************************
1297              * CALCULATE INTERACTIONS *
1298              **************************/
1299
1300             /* REACTION-FIELD ELECTROSTATICS */
1301             felec            = _mm256_mul_pd(qq20,_mm256_sub_pd(_mm256_mul_pd(rinv20,rinvsq20),krf2));
1302
1303             fscal            = felec;
1304
1305             /* Calculate temporary vectorial force */
1306             tx               = _mm256_mul_pd(fscal,dx20);
1307             ty               = _mm256_mul_pd(fscal,dy20);
1308             tz               = _mm256_mul_pd(fscal,dz20);
1309
1310             /* Update vectorial force */
1311             fix2             = _mm256_add_pd(fix2,tx);
1312             fiy2             = _mm256_add_pd(fiy2,ty);
1313             fiz2             = _mm256_add_pd(fiz2,tz);
1314
1315             fjx0             = _mm256_add_pd(fjx0,tx);
1316             fjy0             = _mm256_add_pd(fjy0,ty);
1317             fjz0             = _mm256_add_pd(fjz0,tz);
1318
1319             /**************************
1320              * CALCULATE INTERACTIONS *
1321              **************************/
1322
1323             /* REACTION-FIELD ELECTROSTATICS */
1324             felec            = _mm256_mul_pd(qq21,_mm256_sub_pd(_mm256_mul_pd(rinv21,rinvsq21),krf2));
1325
1326             fscal            = felec;
1327
1328             /* Calculate temporary vectorial force */
1329             tx               = _mm256_mul_pd(fscal,dx21);
1330             ty               = _mm256_mul_pd(fscal,dy21);
1331             tz               = _mm256_mul_pd(fscal,dz21);
1332
1333             /* Update vectorial force */
1334             fix2             = _mm256_add_pd(fix2,tx);
1335             fiy2             = _mm256_add_pd(fiy2,ty);
1336             fiz2             = _mm256_add_pd(fiz2,tz);
1337
1338             fjx1             = _mm256_add_pd(fjx1,tx);
1339             fjy1             = _mm256_add_pd(fjy1,ty);
1340             fjz1             = _mm256_add_pd(fjz1,tz);
1341
1342             /**************************
1343              * CALCULATE INTERACTIONS *
1344              **************************/
1345
1346             /* REACTION-FIELD ELECTROSTATICS */
1347             felec            = _mm256_mul_pd(qq22,_mm256_sub_pd(_mm256_mul_pd(rinv22,rinvsq22),krf2));
1348
1349             fscal            = felec;
1350
1351             /* Calculate temporary vectorial force */
1352             tx               = _mm256_mul_pd(fscal,dx22);
1353             ty               = _mm256_mul_pd(fscal,dy22);
1354             tz               = _mm256_mul_pd(fscal,dz22);
1355
1356             /* Update vectorial force */
1357             fix2             = _mm256_add_pd(fix2,tx);
1358             fiy2             = _mm256_add_pd(fiy2,ty);
1359             fiz2             = _mm256_add_pd(fiz2,tz);
1360
1361             fjx2             = _mm256_add_pd(fjx2,tx);
1362             fjy2             = _mm256_add_pd(fjy2,ty);
1363             fjz2             = _mm256_add_pd(fjz2,tz);
1364
1365             fjptrA             = f+j_coord_offsetA;
1366             fjptrB             = f+j_coord_offsetB;
1367             fjptrC             = f+j_coord_offsetC;
1368             fjptrD             = f+j_coord_offsetD;
1369
1370             gmx_mm256_decrement_3rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
1371                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1372
1373             /* Inner loop uses 243 flops */
1374         }
1375
1376         if(jidx<j_index_end)
1377         {
1378
1379             /* Get j neighbor index, and coordinate index */
1380             jnrlistA         = jjnr[jidx];
1381             jnrlistB         = jjnr[jidx+1];
1382             jnrlistC         = jjnr[jidx+2];
1383             jnrlistD         = jjnr[jidx+3];
1384             /* Sign of each element will be negative for non-real atoms.
1385              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
1386              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
1387              */
1388             tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
1389
1390             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
1391             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
1392             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
1393
1394             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
1395             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
1396             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
1397             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
1398             j_coord_offsetA  = DIM*jnrA;
1399             j_coord_offsetB  = DIM*jnrB;
1400             j_coord_offsetC  = DIM*jnrC;
1401             j_coord_offsetD  = DIM*jnrD;
1402
1403             /* load j atom coordinates */
1404             gmx_mm256_load_3rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
1405                                                  x+j_coord_offsetC,x+j_coord_offsetD,
1406                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
1407
1408             /* Calculate displacement vector */
1409             dx00             = _mm256_sub_pd(ix0,jx0);
1410             dy00             = _mm256_sub_pd(iy0,jy0);
1411             dz00             = _mm256_sub_pd(iz0,jz0);
1412             dx01             = _mm256_sub_pd(ix0,jx1);
1413             dy01             = _mm256_sub_pd(iy0,jy1);
1414             dz01             = _mm256_sub_pd(iz0,jz1);
1415             dx02             = _mm256_sub_pd(ix0,jx2);
1416             dy02             = _mm256_sub_pd(iy0,jy2);
1417             dz02             = _mm256_sub_pd(iz0,jz2);
1418             dx10             = _mm256_sub_pd(ix1,jx0);
1419             dy10             = _mm256_sub_pd(iy1,jy0);
1420             dz10             = _mm256_sub_pd(iz1,jz0);
1421             dx11             = _mm256_sub_pd(ix1,jx1);
1422             dy11             = _mm256_sub_pd(iy1,jy1);
1423             dz11             = _mm256_sub_pd(iz1,jz1);
1424             dx12             = _mm256_sub_pd(ix1,jx2);
1425             dy12             = _mm256_sub_pd(iy1,jy2);
1426             dz12             = _mm256_sub_pd(iz1,jz2);
1427             dx20             = _mm256_sub_pd(ix2,jx0);
1428             dy20             = _mm256_sub_pd(iy2,jy0);
1429             dz20             = _mm256_sub_pd(iz2,jz0);
1430             dx21             = _mm256_sub_pd(ix2,jx1);
1431             dy21             = _mm256_sub_pd(iy2,jy1);
1432             dz21             = _mm256_sub_pd(iz2,jz1);
1433             dx22             = _mm256_sub_pd(ix2,jx2);
1434             dy22             = _mm256_sub_pd(iy2,jy2);
1435             dz22             = _mm256_sub_pd(iz2,jz2);
1436
1437             /* Calculate squared distance and things based on it */
1438             rsq00            = gmx_mm256_calc_rsq_pd(dx00,dy00,dz00);
1439             rsq01            = gmx_mm256_calc_rsq_pd(dx01,dy01,dz01);
1440             rsq02            = gmx_mm256_calc_rsq_pd(dx02,dy02,dz02);
1441             rsq10            = gmx_mm256_calc_rsq_pd(dx10,dy10,dz10);
1442             rsq11            = gmx_mm256_calc_rsq_pd(dx11,dy11,dz11);
1443             rsq12            = gmx_mm256_calc_rsq_pd(dx12,dy12,dz12);
1444             rsq20            = gmx_mm256_calc_rsq_pd(dx20,dy20,dz20);
1445             rsq21            = gmx_mm256_calc_rsq_pd(dx21,dy21,dz21);
1446             rsq22            = gmx_mm256_calc_rsq_pd(dx22,dy22,dz22);
1447
1448             rinv00           = gmx_mm256_invsqrt_pd(rsq00);
1449             rinv01           = gmx_mm256_invsqrt_pd(rsq01);
1450             rinv02           = gmx_mm256_invsqrt_pd(rsq02);
1451             rinv10           = gmx_mm256_invsqrt_pd(rsq10);
1452             rinv11           = gmx_mm256_invsqrt_pd(rsq11);
1453             rinv12           = gmx_mm256_invsqrt_pd(rsq12);
1454             rinv20           = gmx_mm256_invsqrt_pd(rsq20);
1455             rinv21           = gmx_mm256_invsqrt_pd(rsq21);
1456             rinv22           = gmx_mm256_invsqrt_pd(rsq22);
1457
1458             rinvsq00         = _mm256_mul_pd(rinv00,rinv00);
1459             rinvsq01         = _mm256_mul_pd(rinv01,rinv01);
1460             rinvsq02         = _mm256_mul_pd(rinv02,rinv02);
1461             rinvsq10         = _mm256_mul_pd(rinv10,rinv10);
1462             rinvsq11         = _mm256_mul_pd(rinv11,rinv11);
1463             rinvsq12         = _mm256_mul_pd(rinv12,rinv12);
1464             rinvsq20         = _mm256_mul_pd(rinv20,rinv20);
1465             rinvsq21         = _mm256_mul_pd(rinv21,rinv21);
1466             rinvsq22         = _mm256_mul_pd(rinv22,rinv22);
1467
1468             fjx0             = _mm256_setzero_pd();
1469             fjy0             = _mm256_setzero_pd();
1470             fjz0             = _mm256_setzero_pd();
1471             fjx1             = _mm256_setzero_pd();
1472             fjy1             = _mm256_setzero_pd();
1473             fjz1             = _mm256_setzero_pd();
1474             fjx2             = _mm256_setzero_pd();
1475             fjy2             = _mm256_setzero_pd();
1476             fjz2             = _mm256_setzero_pd();
1477
1478             /**************************
1479              * CALCULATE INTERACTIONS *
1480              **************************/
1481
1482             /* REACTION-FIELD ELECTROSTATICS */
1483             felec            = _mm256_mul_pd(qq00,_mm256_sub_pd(_mm256_mul_pd(rinv00,rinvsq00),krf2));
1484
1485             fscal            = felec;
1486
1487             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1488
1489             /* Calculate temporary vectorial force */
1490             tx               = _mm256_mul_pd(fscal,dx00);
1491             ty               = _mm256_mul_pd(fscal,dy00);
1492             tz               = _mm256_mul_pd(fscal,dz00);
1493
1494             /* Update vectorial force */
1495             fix0             = _mm256_add_pd(fix0,tx);
1496             fiy0             = _mm256_add_pd(fiy0,ty);
1497             fiz0             = _mm256_add_pd(fiz0,tz);
1498
1499             fjx0             = _mm256_add_pd(fjx0,tx);
1500             fjy0             = _mm256_add_pd(fjy0,ty);
1501             fjz0             = _mm256_add_pd(fjz0,tz);
1502
1503             /**************************
1504              * CALCULATE INTERACTIONS *
1505              **************************/
1506
1507             /* REACTION-FIELD ELECTROSTATICS */
1508             felec            = _mm256_mul_pd(qq01,_mm256_sub_pd(_mm256_mul_pd(rinv01,rinvsq01),krf2));
1509
1510             fscal            = felec;
1511
1512             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1513
1514             /* Calculate temporary vectorial force */
1515             tx               = _mm256_mul_pd(fscal,dx01);
1516             ty               = _mm256_mul_pd(fscal,dy01);
1517             tz               = _mm256_mul_pd(fscal,dz01);
1518
1519             /* Update vectorial force */
1520             fix0             = _mm256_add_pd(fix0,tx);
1521             fiy0             = _mm256_add_pd(fiy0,ty);
1522             fiz0             = _mm256_add_pd(fiz0,tz);
1523
1524             fjx1             = _mm256_add_pd(fjx1,tx);
1525             fjy1             = _mm256_add_pd(fjy1,ty);
1526             fjz1             = _mm256_add_pd(fjz1,tz);
1527
1528             /**************************
1529              * CALCULATE INTERACTIONS *
1530              **************************/
1531
1532             /* REACTION-FIELD ELECTROSTATICS */
1533             felec            = _mm256_mul_pd(qq02,_mm256_sub_pd(_mm256_mul_pd(rinv02,rinvsq02),krf2));
1534
1535             fscal            = felec;
1536
1537             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1538
1539             /* Calculate temporary vectorial force */
1540             tx               = _mm256_mul_pd(fscal,dx02);
1541             ty               = _mm256_mul_pd(fscal,dy02);
1542             tz               = _mm256_mul_pd(fscal,dz02);
1543
1544             /* Update vectorial force */
1545             fix0             = _mm256_add_pd(fix0,tx);
1546             fiy0             = _mm256_add_pd(fiy0,ty);
1547             fiz0             = _mm256_add_pd(fiz0,tz);
1548
1549             fjx2             = _mm256_add_pd(fjx2,tx);
1550             fjy2             = _mm256_add_pd(fjy2,ty);
1551             fjz2             = _mm256_add_pd(fjz2,tz);
1552
1553             /**************************
1554              * CALCULATE INTERACTIONS *
1555              **************************/
1556
1557             /* REACTION-FIELD ELECTROSTATICS */
1558             felec            = _mm256_mul_pd(qq10,_mm256_sub_pd(_mm256_mul_pd(rinv10,rinvsq10),krf2));
1559
1560             fscal            = felec;
1561
1562             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1563
1564             /* Calculate temporary vectorial force */
1565             tx               = _mm256_mul_pd(fscal,dx10);
1566             ty               = _mm256_mul_pd(fscal,dy10);
1567             tz               = _mm256_mul_pd(fscal,dz10);
1568
1569             /* Update vectorial force */
1570             fix1             = _mm256_add_pd(fix1,tx);
1571             fiy1             = _mm256_add_pd(fiy1,ty);
1572             fiz1             = _mm256_add_pd(fiz1,tz);
1573
1574             fjx0             = _mm256_add_pd(fjx0,tx);
1575             fjy0             = _mm256_add_pd(fjy0,ty);
1576             fjz0             = _mm256_add_pd(fjz0,tz);
1577
1578             /**************************
1579              * CALCULATE INTERACTIONS *
1580              **************************/
1581
1582             /* REACTION-FIELD ELECTROSTATICS */
1583             felec            = _mm256_mul_pd(qq11,_mm256_sub_pd(_mm256_mul_pd(rinv11,rinvsq11),krf2));
1584
1585             fscal            = felec;
1586
1587             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1588
1589             /* Calculate temporary vectorial force */
1590             tx               = _mm256_mul_pd(fscal,dx11);
1591             ty               = _mm256_mul_pd(fscal,dy11);
1592             tz               = _mm256_mul_pd(fscal,dz11);
1593
1594             /* Update vectorial force */
1595             fix1             = _mm256_add_pd(fix1,tx);
1596             fiy1             = _mm256_add_pd(fiy1,ty);
1597             fiz1             = _mm256_add_pd(fiz1,tz);
1598
1599             fjx1             = _mm256_add_pd(fjx1,tx);
1600             fjy1             = _mm256_add_pd(fjy1,ty);
1601             fjz1             = _mm256_add_pd(fjz1,tz);
1602
1603             /**************************
1604              * CALCULATE INTERACTIONS *
1605              **************************/
1606
1607             /* REACTION-FIELD ELECTROSTATICS */
1608             felec            = _mm256_mul_pd(qq12,_mm256_sub_pd(_mm256_mul_pd(rinv12,rinvsq12),krf2));
1609
1610             fscal            = felec;
1611
1612             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1613
1614             /* Calculate temporary vectorial force */
1615             tx               = _mm256_mul_pd(fscal,dx12);
1616             ty               = _mm256_mul_pd(fscal,dy12);
1617             tz               = _mm256_mul_pd(fscal,dz12);
1618
1619             /* Update vectorial force */
1620             fix1             = _mm256_add_pd(fix1,tx);
1621             fiy1             = _mm256_add_pd(fiy1,ty);
1622             fiz1             = _mm256_add_pd(fiz1,tz);
1623
1624             fjx2             = _mm256_add_pd(fjx2,tx);
1625             fjy2             = _mm256_add_pd(fjy2,ty);
1626             fjz2             = _mm256_add_pd(fjz2,tz);
1627
1628             /**************************
1629              * CALCULATE INTERACTIONS *
1630              **************************/
1631
1632             /* REACTION-FIELD ELECTROSTATICS */
1633             felec            = _mm256_mul_pd(qq20,_mm256_sub_pd(_mm256_mul_pd(rinv20,rinvsq20),krf2));
1634
1635             fscal            = felec;
1636
1637             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1638
1639             /* Calculate temporary vectorial force */
1640             tx               = _mm256_mul_pd(fscal,dx20);
1641             ty               = _mm256_mul_pd(fscal,dy20);
1642             tz               = _mm256_mul_pd(fscal,dz20);
1643
1644             /* Update vectorial force */
1645             fix2             = _mm256_add_pd(fix2,tx);
1646             fiy2             = _mm256_add_pd(fiy2,ty);
1647             fiz2             = _mm256_add_pd(fiz2,tz);
1648
1649             fjx0             = _mm256_add_pd(fjx0,tx);
1650             fjy0             = _mm256_add_pd(fjy0,ty);
1651             fjz0             = _mm256_add_pd(fjz0,tz);
1652
1653             /**************************
1654              * CALCULATE INTERACTIONS *
1655              **************************/
1656
1657             /* REACTION-FIELD ELECTROSTATICS */
1658             felec            = _mm256_mul_pd(qq21,_mm256_sub_pd(_mm256_mul_pd(rinv21,rinvsq21),krf2));
1659
1660             fscal            = felec;
1661
1662             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1663
1664             /* Calculate temporary vectorial force */
1665             tx               = _mm256_mul_pd(fscal,dx21);
1666             ty               = _mm256_mul_pd(fscal,dy21);
1667             tz               = _mm256_mul_pd(fscal,dz21);
1668
1669             /* Update vectorial force */
1670             fix2             = _mm256_add_pd(fix2,tx);
1671             fiy2             = _mm256_add_pd(fiy2,ty);
1672             fiz2             = _mm256_add_pd(fiz2,tz);
1673
1674             fjx1             = _mm256_add_pd(fjx1,tx);
1675             fjy1             = _mm256_add_pd(fjy1,ty);
1676             fjz1             = _mm256_add_pd(fjz1,tz);
1677
1678             /**************************
1679              * CALCULATE INTERACTIONS *
1680              **************************/
1681
1682             /* REACTION-FIELD ELECTROSTATICS */
1683             felec            = _mm256_mul_pd(qq22,_mm256_sub_pd(_mm256_mul_pd(rinv22,rinvsq22),krf2));
1684
1685             fscal            = felec;
1686
1687             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
1688
1689             /* Calculate temporary vectorial force */
1690             tx               = _mm256_mul_pd(fscal,dx22);
1691             ty               = _mm256_mul_pd(fscal,dy22);
1692             tz               = _mm256_mul_pd(fscal,dz22);
1693
1694             /* Update vectorial force */
1695             fix2             = _mm256_add_pd(fix2,tx);
1696             fiy2             = _mm256_add_pd(fiy2,ty);
1697             fiz2             = _mm256_add_pd(fiz2,tz);
1698
1699             fjx2             = _mm256_add_pd(fjx2,tx);
1700             fjy2             = _mm256_add_pd(fjy2,ty);
1701             fjz2             = _mm256_add_pd(fjz2,tz);
1702
1703             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1704             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1705             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1706             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1707
1708             gmx_mm256_decrement_3rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
1709                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1710
1711             /* Inner loop uses 243 flops */
1712         }
1713
1714         /* End of innermost loop */
1715
1716         gmx_mm256_update_iforce_3atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1717                                                  f+i_coord_offset,fshift+i_shift_offset);
1718
1719         /* Increment number of inner iterations */
1720         inneriter                  += j_index_end - j_index_start;
1721
1722         /* Outer loop uses 18 flops */
1723     }
1724
1725     /* Increment number of outer iterations */
1726     outeriter        += nri;
1727
1728     /* Update outer/inner flops */
1729
1730     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3W3_F,outeriter*18 + inneriter*243);
1731 }