Compile nonbonded kernels as C++
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_256_single / nb_kernel_ElecRFCut_VdwNone_GeomW3P1_avx_256_single.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_single kernel generator.
37  */
38 #include "gmxpre.h"
39
40 #include "config.h"
41
42 #include <math.h>
43
44 #include "../nb_kernel.h"
45 #include "gromacs/gmxlib/nrnb.h"
46
47 #include "kernelutil_x86_avx_256_single.h"
48
49 /*
50  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwNone_GeomW3P1_VF_avx_256_single
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_single
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,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight 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              jnrE,jnrF,jnrG,jnrH;
75     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
76     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
77     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
78     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
79     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
80     real             rcutoff_scalar;
81     real             *shiftvec,*fshift,*x,*f;
82     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
83     real             scratch[4*DIM];
84     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
85     real *           vdwioffsetptr0;
86     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
87     real *           vdwioffsetptr1;
88     __m256           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
89     real *           vdwioffsetptr2;
90     __m256           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
91     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
92     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
93     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
94     __m256           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
95     __m256           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
96     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
97     real             *charge;
98     __m256           dummy_mask,cutoff_mask;
99     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
100     __m256           one     = _mm256_set1_ps(1.0);
101     __m256           two     = _mm256_set1_ps(2.0);
102     x                = xx[0];
103     f                = ff[0];
104
105     nri              = nlist->nri;
106     iinr             = nlist->iinr;
107     jindex           = nlist->jindex;
108     jjnr             = nlist->jjnr;
109     shiftidx         = nlist->shift;
110     gid              = nlist->gid;
111     shiftvec         = fr->shift_vec[0];
112     fshift           = fr->fshift[0];
113     facel            = _mm256_set1_ps(fr->ic->epsfac);
114     charge           = mdatoms->chargeA;
115     krf              = _mm256_set1_ps(fr->ic->k_rf);
116     krf2             = _mm256_set1_ps(fr->ic->k_rf*2.0);
117     crf              = _mm256_set1_ps(fr->ic->c_rf);
118
119     /* Setup water-specific parameters */
120     inr              = nlist->iinr[0];
121     iq0              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
122     iq1              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
123     iq2              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
124
125     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
126     rcutoff_scalar   = fr->ic->rcoulomb;
127     rcutoff          = _mm256_set1_ps(rcutoff_scalar);
128     rcutoff2         = _mm256_mul_ps(rcutoff,rcutoff);
129
130     /* Avoid stupid compiler warnings */
131     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
132     j_coord_offsetA = 0;
133     j_coord_offsetB = 0;
134     j_coord_offsetC = 0;
135     j_coord_offsetD = 0;
136     j_coord_offsetE = 0;
137     j_coord_offsetF = 0;
138     j_coord_offsetG = 0;
139     j_coord_offsetH = 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_ps(shiftvec+i_shift_offset,x+i_coord_offset,
165                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
166
167         fix0             = _mm256_setzero_ps();
168         fiy0             = _mm256_setzero_ps();
169         fiz0             = _mm256_setzero_ps();
170         fix1             = _mm256_setzero_ps();
171         fiy1             = _mm256_setzero_ps();
172         fiz1             = _mm256_setzero_ps();
173         fix2             = _mm256_setzero_ps();
174         fiy2             = _mm256_setzero_ps();
175         fiz2             = _mm256_setzero_ps();
176
177         /* Reset potential sums */
178         velecsum         = _mm256_setzero_ps();
179
180         /* Start inner kernel loop */
181         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
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             jnrE             = jjnr[jidx+4];
190             jnrF             = jjnr[jidx+5];
191             jnrG             = jjnr[jidx+6];
192             jnrH             = jjnr[jidx+7];
193             j_coord_offsetA  = DIM*jnrA;
194             j_coord_offsetB  = DIM*jnrB;
195             j_coord_offsetC  = DIM*jnrC;
196             j_coord_offsetD  = DIM*jnrD;
197             j_coord_offsetE  = DIM*jnrE;
198             j_coord_offsetF  = DIM*jnrF;
199             j_coord_offsetG  = DIM*jnrG;
200             j_coord_offsetH  = DIM*jnrH;
201
202             /* load j atom coordinates */
203             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
204                                                  x+j_coord_offsetC,x+j_coord_offsetD,
205                                                  x+j_coord_offsetE,x+j_coord_offsetF,
206                                                  x+j_coord_offsetG,x+j_coord_offsetH,
207                                                  &jx0,&jy0,&jz0);
208
209             /* Calculate displacement vector */
210             dx00             = _mm256_sub_ps(ix0,jx0);
211             dy00             = _mm256_sub_ps(iy0,jy0);
212             dz00             = _mm256_sub_ps(iz0,jz0);
213             dx10             = _mm256_sub_ps(ix1,jx0);
214             dy10             = _mm256_sub_ps(iy1,jy0);
215             dz10             = _mm256_sub_ps(iz1,jz0);
216             dx20             = _mm256_sub_ps(ix2,jx0);
217             dy20             = _mm256_sub_ps(iy2,jy0);
218             dz20             = _mm256_sub_ps(iz2,jz0);
219
220             /* Calculate squared distance and things based on it */
221             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
222             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
223             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
224
225             rinv00           = avx256_invsqrt_f(rsq00);
226             rinv10           = avx256_invsqrt_f(rsq10);
227             rinv20           = avx256_invsqrt_f(rsq20);
228
229             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
230             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
231             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
232
233             /* Load parameters for j particles */
234             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
235                                                                  charge+jnrC+0,charge+jnrD+0,
236                                                                  charge+jnrE+0,charge+jnrF+0,
237                                                                  charge+jnrG+0,charge+jnrH+0);
238
239             fjx0             = _mm256_setzero_ps();
240             fjy0             = _mm256_setzero_ps();
241             fjz0             = _mm256_setzero_ps();
242
243             /**************************
244              * CALCULATE INTERACTIONS *
245              **************************/
246
247             if (gmx_mm256_any_lt(rsq00,rcutoff2))
248             {
249
250             /* Compute parameters for interactions between i and j atoms */
251             qq00             = _mm256_mul_ps(iq0,jq0);
252
253             /* REACTION-FIELD ELECTROSTATICS */
254             velec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
255             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
256
257             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
258
259             /* Update potential sum for this i atom from the interaction with this j atom. */
260             velec            = _mm256_and_ps(velec,cutoff_mask);
261             velecsum         = _mm256_add_ps(velecsum,velec);
262
263             fscal            = felec;
264
265             fscal            = _mm256_and_ps(fscal,cutoff_mask);
266
267             /* Calculate temporary vectorial force */
268             tx               = _mm256_mul_ps(fscal,dx00);
269             ty               = _mm256_mul_ps(fscal,dy00);
270             tz               = _mm256_mul_ps(fscal,dz00);
271
272             /* Update vectorial force */
273             fix0             = _mm256_add_ps(fix0,tx);
274             fiy0             = _mm256_add_ps(fiy0,ty);
275             fiz0             = _mm256_add_ps(fiz0,tz);
276
277             fjx0             = _mm256_add_ps(fjx0,tx);
278             fjy0             = _mm256_add_ps(fjy0,ty);
279             fjz0             = _mm256_add_ps(fjz0,tz);
280
281             }
282
283             /**************************
284              * CALCULATE INTERACTIONS *
285              **************************/
286
287             if (gmx_mm256_any_lt(rsq10,rcutoff2))
288             {
289
290             /* Compute parameters for interactions between i and j atoms */
291             qq10             = _mm256_mul_ps(iq1,jq0);
292
293             /* REACTION-FIELD ELECTROSTATICS */
294             velec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_add_ps(rinv10,_mm256_mul_ps(krf,rsq10)),crf));
295             felec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_mul_ps(rinv10,rinvsq10),krf2));
296
297             cutoff_mask      = _mm256_cmp_ps(rsq10,rcutoff2,_CMP_LT_OQ);
298
299             /* Update potential sum for this i atom from the interaction with this j atom. */
300             velec            = _mm256_and_ps(velec,cutoff_mask);
301             velecsum         = _mm256_add_ps(velecsum,velec);
302
303             fscal            = felec;
304
305             fscal            = _mm256_and_ps(fscal,cutoff_mask);
306
307             /* Calculate temporary vectorial force */
308             tx               = _mm256_mul_ps(fscal,dx10);
309             ty               = _mm256_mul_ps(fscal,dy10);
310             tz               = _mm256_mul_ps(fscal,dz10);
311
312             /* Update vectorial force */
313             fix1             = _mm256_add_ps(fix1,tx);
314             fiy1             = _mm256_add_ps(fiy1,ty);
315             fiz1             = _mm256_add_ps(fiz1,tz);
316
317             fjx0             = _mm256_add_ps(fjx0,tx);
318             fjy0             = _mm256_add_ps(fjy0,ty);
319             fjz0             = _mm256_add_ps(fjz0,tz);
320
321             }
322
323             /**************************
324              * CALCULATE INTERACTIONS *
325              **************************/
326
327             if (gmx_mm256_any_lt(rsq20,rcutoff2))
328             {
329
330             /* Compute parameters for interactions between i and j atoms */
331             qq20             = _mm256_mul_ps(iq2,jq0);
332
333             /* REACTION-FIELD ELECTROSTATICS */
334             velec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_add_ps(rinv20,_mm256_mul_ps(krf,rsq20)),crf));
335             felec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_mul_ps(rinv20,rinvsq20),krf2));
336
337             cutoff_mask      = _mm256_cmp_ps(rsq20,rcutoff2,_CMP_LT_OQ);
338
339             /* Update potential sum for this i atom from the interaction with this j atom. */
340             velec            = _mm256_and_ps(velec,cutoff_mask);
341             velecsum         = _mm256_add_ps(velecsum,velec);
342
343             fscal            = felec;
344
345             fscal            = _mm256_and_ps(fscal,cutoff_mask);
346
347             /* Calculate temporary vectorial force */
348             tx               = _mm256_mul_ps(fscal,dx20);
349             ty               = _mm256_mul_ps(fscal,dy20);
350             tz               = _mm256_mul_ps(fscal,dz20);
351
352             /* Update vectorial force */
353             fix2             = _mm256_add_ps(fix2,tx);
354             fiy2             = _mm256_add_ps(fiy2,ty);
355             fiz2             = _mm256_add_ps(fiz2,tz);
356
357             fjx0             = _mm256_add_ps(fjx0,tx);
358             fjy0             = _mm256_add_ps(fjy0,ty);
359             fjz0             = _mm256_add_ps(fjz0,tz);
360
361             }
362
363             fjptrA             = f+j_coord_offsetA;
364             fjptrB             = f+j_coord_offsetB;
365             fjptrC             = f+j_coord_offsetC;
366             fjptrD             = f+j_coord_offsetD;
367             fjptrE             = f+j_coord_offsetE;
368             fjptrF             = f+j_coord_offsetF;
369             fjptrG             = f+j_coord_offsetG;
370             fjptrH             = f+j_coord_offsetH;
371
372             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
373
374             /* Inner loop uses 111 flops */
375         }
376
377         if(jidx<j_index_end)
378         {
379
380             /* Get j neighbor index, and coordinate index */
381             jnrlistA         = jjnr[jidx];
382             jnrlistB         = jjnr[jidx+1];
383             jnrlistC         = jjnr[jidx+2];
384             jnrlistD         = jjnr[jidx+3];
385             jnrlistE         = jjnr[jidx+4];
386             jnrlistF         = jjnr[jidx+5];
387             jnrlistG         = jjnr[jidx+6];
388             jnrlistH         = jjnr[jidx+7];
389             /* Sign of each element will be negative for non-real atoms.
390              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
391              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
392              */
393             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
394                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
395                                             
396             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
397             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
398             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
399             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
400             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
401             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
402             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
403             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
404             j_coord_offsetA  = DIM*jnrA;
405             j_coord_offsetB  = DIM*jnrB;
406             j_coord_offsetC  = DIM*jnrC;
407             j_coord_offsetD  = DIM*jnrD;
408             j_coord_offsetE  = DIM*jnrE;
409             j_coord_offsetF  = DIM*jnrF;
410             j_coord_offsetG  = DIM*jnrG;
411             j_coord_offsetH  = DIM*jnrH;
412
413             /* load j atom coordinates */
414             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
415                                                  x+j_coord_offsetC,x+j_coord_offsetD,
416                                                  x+j_coord_offsetE,x+j_coord_offsetF,
417                                                  x+j_coord_offsetG,x+j_coord_offsetH,
418                                                  &jx0,&jy0,&jz0);
419
420             /* Calculate displacement vector */
421             dx00             = _mm256_sub_ps(ix0,jx0);
422             dy00             = _mm256_sub_ps(iy0,jy0);
423             dz00             = _mm256_sub_ps(iz0,jz0);
424             dx10             = _mm256_sub_ps(ix1,jx0);
425             dy10             = _mm256_sub_ps(iy1,jy0);
426             dz10             = _mm256_sub_ps(iz1,jz0);
427             dx20             = _mm256_sub_ps(ix2,jx0);
428             dy20             = _mm256_sub_ps(iy2,jy0);
429             dz20             = _mm256_sub_ps(iz2,jz0);
430
431             /* Calculate squared distance and things based on it */
432             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
433             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
434             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
435
436             rinv00           = avx256_invsqrt_f(rsq00);
437             rinv10           = avx256_invsqrt_f(rsq10);
438             rinv20           = avx256_invsqrt_f(rsq20);
439
440             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
441             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
442             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
443
444             /* Load parameters for j particles */
445             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
446                                                                  charge+jnrC+0,charge+jnrD+0,
447                                                                  charge+jnrE+0,charge+jnrF+0,
448                                                                  charge+jnrG+0,charge+jnrH+0);
449
450             fjx0             = _mm256_setzero_ps();
451             fjy0             = _mm256_setzero_ps();
452             fjz0             = _mm256_setzero_ps();
453
454             /**************************
455              * CALCULATE INTERACTIONS *
456              **************************/
457
458             if (gmx_mm256_any_lt(rsq00,rcutoff2))
459             {
460
461             /* Compute parameters for interactions between i and j atoms */
462             qq00             = _mm256_mul_ps(iq0,jq0);
463
464             /* REACTION-FIELD ELECTROSTATICS */
465             velec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_add_ps(rinv00,_mm256_mul_ps(krf,rsq00)),crf));
466             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
467
468             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
469
470             /* Update potential sum for this i atom from the interaction with this j atom. */
471             velec            = _mm256_and_ps(velec,cutoff_mask);
472             velec            = _mm256_andnot_ps(dummy_mask,velec);
473             velecsum         = _mm256_add_ps(velecsum,velec);
474
475             fscal            = felec;
476
477             fscal            = _mm256_and_ps(fscal,cutoff_mask);
478
479             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
480
481             /* Calculate temporary vectorial force */
482             tx               = _mm256_mul_ps(fscal,dx00);
483             ty               = _mm256_mul_ps(fscal,dy00);
484             tz               = _mm256_mul_ps(fscal,dz00);
485
486             /* Update vectorial force */
487             fix0             = _mm256_add_ps(fix0,tx);
488             fiy0             = _mm256_add_ps(fiy0,ty);
489             fiz0             = _mm256_add_ps(fiz0,tz);
490
491             fjx0             = _mm256_add_ps(fjx0,tx);
492             fjy0             = _mm256_add_ps(fjy0,ty);
493             fjz0             = _mm256_add_ps(fjz0,tz);
494
495             }
496
497             /**************************
498              * CALCULATE INTERACTIONS *
499              **************************/
500
501             if (gmx_mm256_any_lt(rsq10,rcutoff2))
502             {
503
504             /* Compute parameters for interactions between i and j atoms */
505             qq10             = _mm256_mul_ps(iq1,jq0);
506
507             /* REACTION-FIELD ELECTROSTATICS */
508             velec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_add_ps(rinv10,_mm256_mul_ps(krf,rsq10)),crf));
509             felec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_mul_ps(rinv10,rinvsq10),krf2));
510
511             cutoff_mask      = _mm256_cmp_ps(rsq10,rcutoff2,_CMP_LT_OQ);
512
513             /* Update potential sum for this i atom from the interaction with this j atom. */
514             velec            = _mm256_and_ps(velec,cutoff_mask);
515             velec            = _mm256_andnot_ps(dummy_mask,velec);
516             velecsum         = _mm256_add_ps(velecsum,velec);
517
518             fscal            = felec;
519
520             fscal            = _mm256_and_ps(fscal,cutoff_mask);
521
522             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
523
524             /* Calculate temporary vectorial force */
525             tx               = _mm256_mul_ps(fscal,dx10);
526             ty               = _mm256_mul_ps(fscal,dy10);
527             tz               = _mm256_mul_ps(fscal,dz10);
528
529             /* Update vectorial force */
530             fix1             = _mm256_add_ps(fix1,tx);
531             fiy1             = _mm256_add_ps(fiy1,ty);
532             fiz1             = _mm256_add_ps(fiz1,tz);
533
534             fjx0             = _mm256_add_ps(fjx0,tx);
535             fjy0             = _mm256_add_ps(fjy0,ty);
536             fjz0             = _mm256_add_ps(fjz0,tz);
537
538             }
539
540             /**************************
541              * CALCULATE INTERACTIONS *
542              **************************/
543
544             if (gmx_mm256_any_lt(rsq20,rcutoff2))
545             {
546
547             /* Compute parameters for interactions between i and j atoms */
548             qq20             = _mm256_mul_ps(iq2,jq0);
549
550             /* REACTION-FIELD ELECTROSTATICS */
551             velec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_add_ps(rinv20,_mm256_mul_ps(krf,rsq20)),crf));
552             felec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_mul_ps(rinv20,rinvsq20),krf2));
553
554             cutoff_mask      = _mm256_cmp_ps(rsq20,rcutoff2,_CMP_LT_OQ);
555
556             /* Update potential sum for this i atom from the interaction with this j atom. */
557             velec            = _mm256_and_ps(velec,cutoff_mask);
558             velec            = _mm256_andnot_ps(dummy_mask,velec);
559             velecsum         = _mm256_add_ps(velecsum,velec);
560
561             fscal            = felec;
562
563             fscal            = _mm256_and_ps(fscal,cutoff_mask);
564
565             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
566
567             /* Calculate temporary vectorial force */
568             tx               = _mm256_mul_ps(fscal,dx20);
569             ty               = _mm256_mul_ps(fscal,dy20);
570             tz               = _mm256_mul_ps(fscal,dz20);
571
572             /* Update vectorial force */
573             fix2             = _mm256_add_ps(fix2,tx);
574             fiy2             = _mm256_add_ps(fiy2,ty);
575             fiz2             = _mm256_add_ps(fiz2,tz);
576
577             fjx0             = _mm256_add_ps(fjx0,tx);
578             fjy0             = _mm256_add_ps(fjy0,ty);
579             fjz0             = _mm256_add_ps(fjz0,tz);
580
581             }
582
583             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
584             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
585             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
586             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
587             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
588             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
589             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
590             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
591
592             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
593
594             /* Inner loop uses 111 flops */
595         }
596
597         /* End of innermost loop */
598
599         gmx_mm256_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
600                                                  f+i_coord_offset,fshift+i_shift_offset);
601
602         ggid                        = gid[iidx];
603         /* Update potential energies */
604         gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
605
606         /* Increment number of inner iterations */
607         inneriter                  += j_index_end - j_index_start;
608
609         /* Outer loop uses 19 flops */
610     }
611
612     /* Increment number of outer iterations */
613     outeriter        += nri;
614
615     /* Update outer/inner flops */
616
617     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3_VF,outeriter*19 + inneriter*111);
618 }
619 /*
620  * Gromacs nonbonded kernel:   nb_kernel_ElecRFCut_VdwNone_GeomW3P1_F_avx_256_single
621  * Electrostatics interaction: ReactionField
622  * VdW interaction:            None
623  * Geometry:                   Water3-Particle
624  * Calculate force/pot:        Force
625  */
626 void
627 nb_kernel_ElecRFCut_VdwNone_GeomW3P1_F_avx_256_single
628                     (t_nblist                    * gmx_restrict       nlist,
629                      rvec                        * gmx_restrict          xx,
630                      rvec                        * gmx_restrict          ff,
631                      struct t_forcerec           * gmx_restrict          fr,
632                      t_mdatoms                   * gmx_restrict     mdatoms,
633                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
634                      t_nrnb                      * gmx_restrict        nrnb)
635 {
636     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
637      * just 0 for non-waters.
638      * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
639      * jnr indices corresponding to data put in the four positions in the SIMD register.
640      */
641     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
642     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
643     int              jnrA,jnrB,jnrC,jnrD;
644     int              jnrE,jnrF,jnrG,jnrH;
645     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
646     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
647     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
648     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
649     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
650     real             rcutoff_scalar;
651     real             *shiftvec,*fshift,*x,*f;
652     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
653     real             scratch[4*DIM];
654     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
655     real *           vdwioffsetptr0;
656     __m256           ix0,iy0,iz0,fix0,fiy0,fiz0,iq0,isai0;
657     real *           vdwioffsetptr1;
658     __m256           ix1,iy1,iz1,fix1,fiy1,fiz1,iq1,isai1;
659     real *           vdwioffsetptr2;
660     __m256           ix2,iy2,iz2,fix2,fiy2,fiz2,iq2,isai2;
661     int              vdwjidx0A,vdwjidx0B,vdwjidx0C,vdwjidx0D,vdwjidx0E,vdwjidx0F,vdwjidx0G,vdwjidx0H;
662     __m256           jx0,jy0,jz0,fjx0,fjy0,fjz0,jq0,isaj0;
663     __m256           dx00,dy00,dz00,rsq00,rinv00,rinvsq00,r00,qq00,c6_00,c12_00;
664     __m256           dx10,dy10,dz10,rsq10,rinv10,rinvsq10,r10,qq10,c6_10,c12_10;
665     __m256           dx20,dy20,dz20,rsq20,rinv20,rinvsq20,r20,qq20,c6_20,c12_20;
666     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
667     real             *charge;
668     __m256           dummy_mask,cutoff_mask;
669     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
670     __m256           one     = _mm256_set1_ps(1.0);
671     __m256           two     = _mm256_set1_ps(2.0);
672     x                = xx[0];
673     f                = ff[0];
674
675     nri              = nlist->nri;
676     iinr             = nlist->iinr;
677     jindex           = nlist->jindex;
678     jjnr             = nlist->jjnr;
679     shiftidx         = nlist->shift;
680     gid              = nlist->gid;
681     shiftvec         = fr->shift_vec[0];
682     fshift           = fr->fshift[0];
683     facel            = _mm256_set1_ps(fr->ic->epsfac);
684     charge           = mdatoms->chargeA;
685     krf              = _mm256_set1_ps(fr->ic->k_rf);
686     krf2             = _mm256_set1_ps(fr->ic->k_rf*2.0);
687     crf              = _mm256_set1_ps(fr->ic->c_rf);
688
689     /* Setup water-specific parameters */
690     inr              = nlist->iinr[0];
691     iq0              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+0]));
692     iq1              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+1]));
693     iq2              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+2]));
694
695     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
696     rcutoff_scalar   = fr->ic->rcoulomb;
697     rcutoff          = _mm256_set1_ps(rcutoff_scalar);
698     rcutoff2         = _mm256_mul_ps(rcutoff,rcutoff);
699
700     /* Avoid stupid compiler warnings */
701     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
702     j_coord_offsetA = 0;
703     j_coord_offsetB = 0;
704     j_coord_offsetC = 0;
705     j_coord_offsetD = 0;
706     j_coord_offsetE = 0;
707     j_coord_offsetF = 0;
708     j_coord_offsetG = 0;
709     j_coord_offsetH = 0;
710
711     outeriter        = 0;
712     inneriter        = 0;
713
714     for(iidx=0;iidx<4*DIM;iidx++)
715     {
716         scratch[iidx] = 0.0;
717     }
718
719     /* Start outer loop over neighborlists */
720     for(iidx=0; iidx<nri; iidx++)
721     {
722         /* Load shift vector for this list */
723         i_shift_offset   = DIM*shiftidx[iidx];
724
725         /* Load limits for loop over neighbors */
726         j_index_start    = jindex[iidx];
727         j_index_end      = jindex[iidx+1];
728
729         /* Get outer coordinate index */
730         inr              = iinr[iidx];
731         i_coord_offset   = DIM*inr;
732
733         /* Load i particle coords and add shift vector */
734         gmx_mm256_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
735                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
736
737         fix0             = _mm256_setzero_ps();
738         fiy0             = _mm256_setzero_ps();
739         fiz0             = _mm256_setzero_ps();
740         fix1             = _mm256_setzero_ps();
741         fiy1             = _mm256_setzero_ps();
742         fiz1             = _mm256_setzero_ps();
743         fix2             = _mm256_setzero_ps();
744         fiy2             = _mm256_setzero_ps();
745         fiz2             = _mm256_setzero_ps();
746
747         /* Start inner kernel loop */
748         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
749         {
750
751             /* Get j neighbor index, and coordinate index */
752             jnrA             = jjnr[jidx];
753             jnrB             = jjnr[jidx+1];
754             jnrC             = jjnr[jidx+2];
755             jnrD             = jjnr[jidx+3];
756             jnrE             = jjnr[jidx+4];
757             jnrF             = jjnr[jidx+5];
758             jnrG             = jjnr[jidx+6];
759             jnrH             = jjnr[jidx+7];
760             j_coord_offsetA  = DIM*jnrA;
761             j_coord_offsetB  = DIM*jnrB;
762             j_coord_offsetC  = DIM*jnrC;
763             j_coord_offsetD  = DIM*jnrD;
764             j_coord_offsetE  = DIM*jnrE;
765             j_coord_offsetF  = DIM*jnrF;
766             j_coord_offsetG  = DIM*jnrG;
767             j_coord_offsetH  = DIM*jnrH;
768
769             /* load j atom coordinates */
770             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
771                                                  x+j_coord_offsetC,x+j_coord_offsetD,
772                                                  x+j_coord_offsetE,x+j_coord_offsetF,
773                                                  x+j_coord_offsetG,x+j_coord_offsetH,
774                                                  &jx0,&jy0,&jz0);
775
776             /* Calculate displacement vector */
777             dx00             = _mm256_sub_ps(ix0,jx0);
778             dy00             = _mm256_sub_ps(iy0,jy0);
779             dz00             = _mm256_sub_ps(iz0,jz0);
780             dx10             = _mm256_sub_ps(ix1,jx0);
781             dy10             = _mm256_sub_ps(iy1,jy0);
782             dz10             = _mm256_sub_ps(iz1,jz0);
783             dx20             = _mm256_sub_ps(ix2,jx0);
784             dy20             = _mm256_sub_ps(iy2,jy0);
785             dz20             = _mm256_sub_ps(iz2,jz0);
786
787             /* Calculate squared distance and things based on it */
788             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
789             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
790             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
791
792             rinv00           = avx256_invsqrt_f(rsq00);
793             rinv10           = avx256_invsqrt_f(rsq10);
794             rinv20           = avx256_invsqrt_f(rsq20);
795
796             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
797             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
798             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
799
800             /* Load parameters for j particles */
801             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
802                                                                  charge+jnrC+0,charge+jnrD+0,
803                                                                  charge+jnrE+0,charge+jnrF+0,
804                                                                  charge+jnrG+0,charge+jnrH+0);
805
806             fjx0             = _mm256_setzero_ps();
807             fjy0             = _mm256_setzero_ps();
808             fjz0             = _mm256_setzero_ps();
809
810             /**************************
811              * CALCULATE INTERACTIONS *
812              **************************/
813
814             if (gmx_mm256_any_lt(rsq00,rcutoff2))
815             {
816
817             /* Compute parameters for interactions between i and j atoms */
818             qq00             = _mm256_mul_ps(iq0,jq0);
819
820             /* REACTION-FIELD ELECTROSTATICS */
821             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
822
823             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
824
825             fscal            = felec;
826
827             fscal            = _mm256_and_ps(fscal,cutoff_mask);
828
829             /* Calculate temporary vectorial force */
830             tx               = _mm256_mul_ps(fscal,dx00);
831             ty               = _mm256_mul_ps(fscal,dy00);
832             tz               = _mm256_mul_ps(fscal,dz00);
833
834             /* Update vectorial force */
835             fix0             = _mm256_add_ps(fix0,tx);
836             fiy0             = _mm256_add_ps(fiy0,ty);
837             fiz0             = _mm256_add_ps(fiz0,tz);
838
839             fjx0             = _mm256_add_ps(fjx0,tx);
840             fjy0             = _mm256_add_ps(fjy0,ty);
841             fjz0             = _mm256_add_ps(fjz0,tz);
842
843             }
844
845             /**************************
846              * CALCULATE INTERACTIONS *
847              **************************/
848
849             if (gmx_mm256_any_lt(rsq10,rcutoff2))
850             {
851
852             /* Compute parameters for interactions between i and j atoms */
853             qq10             = _mm256_mul_ps(iq1,jq0);
854
855             /* REACTION-FIELD ELECTROSTATICS */
856             felec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_mul_ps(rinv10,rinvsq10),krf2));
857
858             cutoff_mask      = _mm256_cmp_ps(rsq10,rcutoff2,_CMP_LT_OQ);
859
860             fscal            = felec;
861
862             fscal            = _mm256_and_ps(fscal,cutoff_mask);
863
864             /* Calculate temporary vectorial force */
865             tx               = _mm256_mul_ps(fscal,dx10);
866             ty               = _mm256_mul_ps(fscal,dy10);
867             tz               = _mm256_mul_ps(fscal,dz10);
868
869             /* Update vectorial force */
870             fix1             = _mm256_add_ps(fix1,tx);
871             fiy1             = _mm256_add_ps(fiy1,ty);
872             fiz1             = _mm256_add_ps(fiz1,tz);
873
874             fjx0             = _mm256_add_ps(fjx0,tx);
875             fjy0             = _mm256_add_ps(fjy0,ty);
876             fjz0             = _mm256_add_ps(fjz0,tz);
877
878             }
879
880             /**************************
881              * CALCULATE INTERACTIONS *
882              **************************/
883
884             if (gmx_mm256_any_lt(rsq20,rcutoff2))
885             {
886
887             /* Compute parameters for interactions between i and j atoms */
888             qq20             = _mm256_mul_ps(iq2,jq0);
889
890             /* REACTION-FIELD ELECTROSTATICS */
891             felec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_mul_ps(rinv20,rinvsq20),krf2));
892
893             cutoff_mask      = _mm256_cmp_ps(rsq20,rcutoff2,_CMP_LT_OQ);
894
895             fscal            = felec;
896
897             fscal            = _mm256_and_ps(fscal,cutoff_mask);
898
899             /* Calculate temporary vectorial force */
900             tx               = _mm256_mul_ps(fscal,dx20);
901             ty               = _mm256_mul_ps(fscal,dy20);
902             tz               = _mm256_mul_ps(fscal,dz20);
903
904             /* Update vectorial force */
905             fix2             = _mm256_add_ps(fix2,tx);
906             fiy2             = _mm256_add_ps(fiy2,ty);
907             fiz2             = _mm256_add_ps(fiz2,tz);
908
909             fjx0             = _mm256_add_ps(fjx0,tx);
910             fjy0             = _mm256_add_ps(fjy0,ty);
911             fjz0             = _mm256_add_ps(fjz0,tz);
912
913             }
914
915             fjptrA             = f+j_coord_offsetA;
916             fjptrB             = f+j_coord_offsetB;
917             fjptrC             = f+j_coord_offsetC;
918             fjptrD             = f+j_coord_offsetD;
919             fjptrE             = f+j_coord_offsetE;
920             fjptrF             = f+j_coord_offsetF;
921             fjptrG             = f+j_coord_offsetG;
922             fjptrH             = f+j_coord_offsetH;
923
924             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
925
926             /* Inner loop uses 93 flops */
927         }
928
929         if(jidx<j_index_end)
930         {
931
932             /* Get j neighbor index, and coordinate index */
933             jnrlistA         = jjnr[jidx];
934             jnrlistB         = jjnr[jidx+1];
935             jnrlistC         = jjnr[jidx+2];
936             jnrlistD         = jjnr[jidx+3];
937             jnrlistE         = jjnr[jidx+4];
938             jnrlistF         = jjnr[jidx+5];
939             jnrlistG         = jjnr[jidx+6];
940             jnrlistH         = jjnr[jidx+7];
941             /* Sign of each element will be negative for non-real atoms.
942              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
943              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
944              */
945             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
946                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
947                                             
948             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
949             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
950             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
951             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
952             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
953             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
954             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
955             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
956             j_coord_offsetA  = DIM*jnrA;
957             j_coord_offsetB  = DIM*jnrB;
958             j_coord_offsetC  = DIM*jnrC;
959             j_coord_offsetD  = DIM*jnrD;
960             j_coord_offsetE  = DIM*jnrE;
961             j_coord_offsetF  = DIM*jnrF;
962             j_coord_offsetG  = DIM*jnrG;
963             j_coord_offsetH  = DIM*jnrH;
964
965             /* load j atom coordinates */
966             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
967                                                  x+j_coord_offsetC,x+j_coord_offsetD,
968                                                  x+j_coord_offsetE,x+j_coord_offsetF,
969                                                  x+j_coord_offsetG,x+j_coord_offsetH,
970                                                  &jx0,&jy0,&jz0);
971
972             /* Calculate displacement vector */
973             dx00             = _mm256_sub_ps(ix0,jx0);
974             dy00             = _mm256_sub_ps(iy0,jy0);
975             dz00             = _mm256_sub_ps(iz0,jz0);
976             dx10             = _mm256_sub_ps(ix1,jx0);
977             dy10             = _mm256_sub_ps(iy1,jy0);
978             dz10             = _mm256_sub_ps(iz1,jz0);
979             dx20             = _mm256_sub_ps(ix2,jx0);
980             dy20             = _mm256_sub_ps(iy2,jy0);
981             dz20             = _mm256_sub_ps(iz2,jz0);
982
983             /* Calculate squared distance and things based on it */
984             rsq00            = gmx_mm256_calc_rsq_ps(dx00,dy00,dz00);
985             rsq10            = gmx_mm256_calc_rsq_ps(dx10,dy10,dz10);
986             rsq20            = gmx_mm256_calc_rsq_ps(dx20,dy20,dz20);
987
988             rinv00           = avx256_invsqrt_f(rsq00);
989             rinv10           = avx256_invsqrt_f(rsq10);
990             rinv20           = avx256_invsqrt_f(rsq20);
991
992             rinvsq00         = _mm256_mul_ps(rinv00,rinv00);
993             rinvsq10         = _mm256_mul_ps(rinv10,rinv10);
994             rinvsq20         = _mm256_mul_ps(rinv20,rinv20);
995
996             /* Load parameters for j particles */
997             jq0              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+0,charge+jnrB+0,
998                                                                  charge+jnrC+0,charge+jnrD+0,
999                                                                  charge+jnrE+0,charge+jnrF+0,
1000                                                                  charge+jnrG+0,charge+jnrH+0);
1001
1002             fjx0             = _mm256_setzero_ps();
1003             fjy0             = _mm256_setzero_ps();
1004             fjz0             = _mm256_setzero_ps();
1005
1006             /**************************
1007              * CALCULATE INTERACTIONS *
1008              **************************/
1009
1010             if (gmx_mm256_any_lt(rsq00,rcutoff2))
1011             {
1012
1013             /* Compute parameters for interactions between i and j atoms */
1014             qq00             = _mm256_mul_ps(iq0,jq0);
1015
1016             /* REACTION-FIELD ELECTROSTATICS */
1017             felec            = _mm256_mul_ps(qq00,_mm256_sub_ps(_mm256_mul_ps(rinv00,rinvsq00),krf2));
1018
1019             cutoff_mask      = _mm256_cmp_ps(rsq00,rcutoff2,_CMP_LT_OQ);
1020
1021             fscal            = felec;
1022
1023             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1024
1025             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1026
1027             /* Calculate temporary vectorial force */
1028             tx               = _mm256_mul_ps(fscal,dx00);
1029             ty               = _mm256_mul_ps(fscal,dy00);
1030             tz               = _mm256_mul_ps(fscal,dz00);
1031
1032             /* Update vectorial force */
1033             fix0             = _mm256_add_ps(fix0,tx);
1034             fiy0             = _mm256_add_ps(fiy0,ty);
1035             fiz0             = _mm256_add_ps(fiz0,tz);
1036
1037             fjx0             = _mm256_add_ps(fjx0,tx);
1038             fjy0             = _mm256_add_ps(fjy0,ty);
1039             fjz0             = _mm256_add_ps(fjz0,tz);
1040
1041             }
1042
1043             /**************************
1044              * CALCULATE INTERACTIONS *
1045              **************************/
1046
1047             if (gmx_mm256_any_lt(rsq10,rcutoff2))
1048             {
1049
1050             /* Compute parameters for interactions between i and j atoms */
1051             qq10             = _mm256_mul_ps(iq1,jq0);
1052
1053             /* REACTION-FIELD ELECTROSTATICS */
1054             felec            = _mm256_mul_ps(qq10,_mm256_sub_ps(_mm256_mul_ps(rinv10,rinvsq10),krf2));
1055
1056             cutoff_mask      = _mm256_cmp_ps(rsq10,rcutoff2,_CMP_LT_OQ);
1057
1058             fscal            = felec;
1059
1060             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1061
1062             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1063
1064             /* Calculate temporary vectorial force */
1065             tx               = _mm256_mul_ps(fscal,dx10);
1066             ty               = _mm256_mul_ps(fscal,dy10);
1067             tz               = _mm256_mul_ps(fscal,dz10);
1068
1069             /* Update vectorial force */
1070             fix1             = _mm256_add_ps(fix1,tx);
1071             fiy1             = _mm256_add_ps(fiy1,ty);
1072             fiz1             = _mm256_add_ps(fiz1,tz);
1073
1074             fjx0             = _mm256_add_ps(fjx0,tx);
1075             fjy0             = _mm256_add_ps(fjy0,ty);
1076             fjz0             = _mm256_add_ps(fjz0,tz);
1077
1078             }
1079
1080             /**************************
1081              * CALCULATE INTERACTIONS *
1082              **************************/
1083
1084             if (gmx_mm256_any_lt(rsq20,rcutoff2))
1085             {
1086
1087             /* Compute parameters for interactions between i and j atoms */
1088             qq20             = _mm256_mul_ps(iq2,jq0);
1089
1090             /* REACTION-FIELD ELECTROSTATICS */
1091             felec            = _mm256_mul_ps(qq20,_mm256_sub_ps(_mm256_mul_ps(rinv20,rinvsq20),krf2));
1092
1093             cutoff_mask      = _mm256_cmp_ps(rsq20,rcutoff2,_CMP_LT_OQ);
1094
1095             fscal            = felec;
1096
1097             fscal            = _mm256_and_ps(fscal,cutoff_mask);
1098
1099             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
1100
1101             /* Calculate temporary vectorial force */
1102             tx               = _mm256_mul_ps(fscal,dx20);
1103             ty               = _mm256_mul_ps(fscal,dy20);
1104             tz               = _mm256_mul_ps(fscal,dz20);
1105
1106             /* Update vectorial force */
1107             fix2             = _mm256_add_ps(fix2,tx);
1108             fiy2             = _mm256_add_ps(fiy2,ty);
1109             fiz2             = _mm256_add_ps(fiz2,tz);
1110
1111             fjx0             = _mm256_add_ps(fjx0,tx);
1112             fjy0             = _mm256_add_ps(fjy0,ty);
1113             fjz0             = _mm256_add_ps(fjz0,tz);
1114
1115             }
1116
1117             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1118             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1119             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1120             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1121             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1122             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1123             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1124             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1125
1126             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
1127
1128             /* Inner loop uses 93 flops */
1129         }
1130
1131         /* End of innermost loop */
1132
1133         gmx_mm256_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1134                                                  f+i_coord_offset,fshift+i_shift_offset);
1135
1136         /* Increment number of inner iterations */
1137         inneriter                  += j_index_end - j_index_start;
1138
1139         /* Outer loop uses 18 flops */
1140     }
1141
1142     /* Increment number of outer iterations */
1143     outeriter        += nri;
1144
1145     /* Update outer/inner flops */
1146
1147     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_W3_F,outeriter*18 + inneriter*93);
1148 }