89bfdd3514646bab46276d15bdaa493ffd873f85
[alexxy/gromacs.git] / src / gmxlib / nonbonded / nb_kernel_avx_256_single / nb_kernel_template_avx_256_single.pre
1 /* #if 0 */
2 #error This file must be processed with the Gromacs pre-preprocessor
3 /* #endif */
4 /* #if INCLUDE_HEADER */
5 #ifdef HAVE_CONFIG_H
6 #include <config.h>
7 #endif
8
9 #include <math.h>
10
11 #include "../nb_kernel.h"
12 #include "types/simple.h"
13 #include "vec.h"
14 #include "nrnb.h"
15
16 #include "gmx_math_x86_avx_256_single.h"
17 #include "kernelutil_x86_avx_256_single.h"
18 /* #endif */
19
20 /* ## List of variables set by the generating script:                                    */
21 /* ##                                                                                    */
22 /* ## Setttings that apply to the entire kernel:                                         */
23 /* ## KERNEL_ELEC:           String, choice for electrostatic interactions               */
24 /* ## KERNEL_VDW:            String, choice for van der Waals interactions               */
25 /* ## KERNEL_NAME:           String, name of this kernel                                 */
26 /* ## KERNEL_VF:             String telling if we calculate potential, force, or both    */
27 /* ## GEOMETRY_I/GEOMETRY_J: String, name of each geometry, e.g. 'Water3' or '1Particle' */
28 /* ##                                                                                    */
29 /* ## Setttings that apply to particles in the outer (I) or inner (J) loops:             */
30 /* ## PARTICLES_I[]/         Arrays with lists of i/j particles to use in kernel. It is  */
31 /* ## PARTICLES_J[]:         just [0] for particle geometry, but can be longer for water */
32 /* ## PARTICLES_ELEC_I[]/    Arrays with lists of i/j particle that have electrostatics  */
33 /* ## PARTICLES_ELEC_J[]:    interactions that should be calculated in this kernel.      */
34 /* ## PARTICLES_VDW_I[]/     Arrays with the list of i/j particle that have VdW          */
35 /* ## PARTICLES_VDW_J[]:     interactions that should be calculated in this kernel.      */
36 /* ##                                                                                    */
37 /* ## Setttings for pairs of interactions (e.g. 2nd i particle against 1st j particle)   */
38 /* ## PAIRS_IJ[]:            Array with (i,j) tuples of pairs for which interactions     */
39 /* ##                        should be calculated in this kernel. Zero-charge particles  */
40 /* ##                        do not have interactions with particles without vdw, and    */
41 /* ##                        Vdw-only interactions are not evaluated in a no-vdw-kernel. */
42 /* ## INTERACTION_FLAGS[][]: 2D matrix, dimension e.g. 3*3 for water-water interactions. */
43 /* ##                        For each i-j pair, the element [I][J] is a list of strings  */
44 /* ##                        defining properties/flags of this interaction. Examples     */
45 /* ##                        include 'electrostatics'/'vdw' if that type of interaction  */
46 /* ##                        should be evaluated, 'rsq'/'rinv'/'rinvsq' if those values  */
47 /* ##                        are needed, and 'exactcutoff' or 'shift','switch' to        */
48 /* ##                        decide if the force/potential should be modified. This way  */
49 /* ##                        we only calculate values absolutely needed for each case.   */
50
51 /* ## Calculate the size and offset for (merged/interleaved) table data */
52
53 /*
54  * Gromacs nonbonded kernel:   {KERNEL_NAME}
55  * Electrostatics interaction: {KERNEL_ELEC}
56  * VdW interaction:            {KERNEL_VDW}
57  * Geometry:                   {GEOMETRY_I}-{GEOMETRY_J}
58  * Calculate force/pot:        {KERNEL_VF}
59  */
60 void
61 {KERNEL_NAME}
62                     (t_nblist * gmx_restrict                nlist,
63                      rvec * gmx_restrict                    xx,
64                      rvec * gmx_restrict                    ff,
65                      t_forcerec * gmx_restrict              fr,
66                      t_mdatoms * gmx_restrict               mdatoms,
67                      nb_kernel_data_t * gmx_restrict        kernel_data,
68                      t_nrnb * gmx_restrict                  nrnb)
69 {
70     /* ## Not all variables are used for all kernels, but any optimizing compiler fixes that, */
71     /* ## so there is no point in going to extremes to exclude variables that are not needed. */
72     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or 
73      * just 0 for non-waters.
74      * Suffixes A,B,C,D,E,F,G,H refer to j loop unrolling done with AVX, e.g. for the eight different
75      * jnr indices corresponding to data put in the four positions in the SIMD register.
76      */
77     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
78     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
79     int              jnrA,jnrB,jnrC,jnrD;
80     int              jnrE,jnrF,jnrG,jnrH;
81     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
82     int              jnrlistE,jnrlistF,jnrlistG,jnrlistH;
83     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
84     int              j_coord_offsetE,j_coord_offsetF,j_coord_offsetG,j_coord_offsetH;
85     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
86     real             rcutoff_scalar;
87     real             *shiftvec,*fshift,*x,*f;
88     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD,*fjptrE,*fjptrF,*fjptrG,*fjptrH;
89     real             scratch[4*DIM];
90     __m256           tx,ty,tz,fscal,rcutoff,rcutoff2,jidxall;
91     /* #for I in PARTICLES_I */
92     real *           vdwioffsetptr{I};
93     __m256           ix{I},iy{I},iz{I},fix{I},fiy{I},fiz{I},iq{I},isai{I};
94     /* #endfor */
95     /* #for J in PARTICLES_J */
96     int              vdwjidx{J}A,vdwjidx{J}B,vdwjidx{J}C,vdwjidx{J}D,vdwjidx{J}E,vdwjidx{J}F,vdwjidx{J}G,vdwjidx{J}H;
97     __m256           jx{J},jy{J},jz{J},fjx{J},fjy{J},fjz{J},jq{J},isaj{J};
98     /* #endfor */
99     /* #for I,J in PAIRS_IJ */
100     __m256           dx{I}{J},dy{I}{J},dz{I}{J},rsq{I}{J},rinv{I}{J},rinvsq{I}{J},r{I}{J},qq{I}{J},c6_{I}{J},c12_{I}{J};
101     /* #endfor */
102     /* #if KERNEL_ELEC != 'None' */
103     __m256           velec,felec,velecsum,facel,crf,krf,krf2;
104     real             *charge;
105     /* #endif */
106     /* #if 'GeneralizedBorn' in KERNEL_ELEC */
107     __m256i          gbitab;
108     __m128i          gbitab_lo,gbitab_hi;
109     __m256           vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
110     __m256           minushalf = _mm256_set1_ps(-0.5);
111     real             *invsqrta,*dvda,*gbtab;
112     /* #endif */
113     /* #if KERNEL_VDW != 'None' */
114     int              nvdwtype;
115     __m256           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
116     int              *vdwtype;
117     real             *vdwparam;
118     __m256           one_sixth   = _mm256_set1_ps(1.0/6.0);
119     __m256           one_twelfth = _mm256_set1_ps(1.0/12.0);
120     /* #endif */
121     /* #if 'Table' in KERNEL_ELEC or 'GeneralizedBorn' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
122     __m256i          vfitab;
123     __m128i          vfitab_lo,vfitab_hi;
124     __m128i          ifour       = _mm_set1_epi32(4);
125     __m256           rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
126     real             *vftab;
127     /* #endif */
128     /* #if 'Ewald' in KERNEL_ELEC */
129     __m256i          ewitab;
130     __m128i          ewitab_lo,ewitab_hi;
131     __m256           ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
132     __m256           beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
133     real             *ewtab;
134     /* #endif */
135     /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
136     __m256           rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
137     real             rswitch_scalar,d_scalar;
138     /* #endif */
139     __m256           dummy_mask,cutoff_mask;
140     __m256           signbit = _mm256_castsi256_ps( _mm256_set1_epi32(0x80000000) );
141     __m256           one     = _mm256_set1_ps(1.0);
142     __m256           two     = _mm256_set1_ps(2.0);
143     x                = xx[0];
144     f                = ff[0];
145
146     nri              = nlist->nri;
147     iinr             = nlist->iinr;
148     jindex           = nlist->jindex;
149     jjnr             = nlist->jjnr;
150     shiftidx         = nlist->shift;
151     gid              = nlist->gid;
152     shiftvec         = fr->shift_vec[0];
153     fshift           = fr->fshift[0];
154     /* #if KERNEL_ELEC != 'None' */
155     facel            = _mm256_set1_ps(fr->epsfac);
156     charge           = mdatoms->chargeA;
157     /*     #if 'ReactionField' in KERNEL_ELEC */
158     krf              = _mm256_set1_ps(fr->ic->k_rf);
159     krf2             = _mm256_set1_ps(fr->ic->k_rf*2.0);
160     crf              = _mm256_set1_ps(fr->ic->c_rf);
161     /*     #endif */
162     /* #endif */
163     /* #if KERNEL_VDW != 'None' */
164     nvdwtype         = fr->ntype;
165     vdwparam         = fr->nbfp;
166     vdwtype          = mdatoms->typeA;
167     /* #endif */
168
169     /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
170     vftab            = kernel_data->table_elec_vdw->data;
171     vftabscale       = _mm256_set1_ps(kernel_data->table_elec_vdw->scale);
172     /* #elif 'Table' in KERNEL_ELEC */
173     vftab            = kernel_data->table_elec->data;
174     vftabscale       = _mm256_set1_ps(kernel_data->table_elec->scale);
175     /* #elif 'Table' in KERNEL_VDW */
176     vftab            = kernel_data->table_vdw->data;
177     vftabscale       = _mm256_set1_ps(kernel_data->table_vdw->scale);
178     /* #endif */
179
180     /* #if 'Ewald' in KERNEL_ELEC */
181     sh_ewald         = _mm256_set1_ps(fr->ic->sh_ewald);
182     beta             = _mm256_set1_ps(fr->ic->ewaldcoeff);
183     beta2            = _mm256_mul_ps(beta,beta);
184     beta3            = _mm256_mul_ps(beta,beta2);
185
186     /*     #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
187     ewtab            = fr->ic->tabq_coul_F;
188     ewtabscale       = _mm256_set1_ps(fr->ic->tabq_scale);
189     ewtabhalfspace   = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
190     /*     #else */
191     ewtab            = fr->ic->tabq_coul_FDV0;
192     ewtabscale       = _mm256_set1_ps(fr->ic->tabq_scale);
193     ewtabhalfspace   = _mm256_set1_ps(0.5/fr->ic->tabq_scale);
194      /*     #endif */
195     /* #endif */
196
197     /* #if KERNEL_ELEC=='GeneralizedBorn' */
198     invsqrta         = fr->invsqrta;
199     dvda             = fr->dvda;
200     gbtabscale       = _mm256_set1_ps(fr->gbtab.scale);
201     gbtab            = fr->gbtab.data;
202     gbinvepsdiff     = _mm256_set1_ps((1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
203     /* #endif */
204
205     /* #if 'Water' in GEOMETRY_I */
206     /* Setup water-specific parameters */
207     inr              = nlist->iinr[0];
208     /*     #for I in PARTICLES_ELEC_I */
209     iq{I}              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+{I}]));
210     /*     #endfor */
211     /*     #for I in PARTICLES_VDW_I */
212     vdwioffsetptr{I}   = vdwparam+2*nvdwtype*vdwtype[inr+{I}];
213     /*     #endfor */
214     /* #endif */
215
216     /* #if 'Water' in GEOMETRY_J */
217     /*     #for J in PARTICLES_ELEC_J */
218     jq{J}              = _mm256_set1_ps(charge[inr+{J}]);
219     /*     #endfor */
220     /*     #for J in PARTICLES_VDW_J */
221     vdwjidx{J}A        = 2*vdwtype[inr+{J}];
222     /*     #endfor */
223     /*     #for I,J in PAIRS_IJ */
224     /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
225     qq{I}{J}             = _mm256_mul_ps(iq{I},jq{J});
226     /*         #endif */
227     /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
228     c6_{I}{J}            = _mm256_set1_ps(vdwioffsetptr{I}[vdwjidx{J}A]);
229     c12_{I}{J}           = _mm256_set1_ps(vdwioffsetptr{I}[vdwjidx{J}A+1]);
230     /*         #endif */
231     /*     #endfor */
232     /* #endif */
233
234     /* #if KERNEL_MOD_ELEC!='None' or KERNEL_MOD_VDW!='None' */
235     /*     #if KERNEL_ELEC!='None' */
236     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
237     rcutoff_scalar   = fr->rcoulomb;
238     /*     #else */
239     rcutoff_scalar   = fr->rvdw;
240     /*     #endif */
241     rcutoff          = _mm256_set1_ps(rcutoff_scalar);
242     rcutoff2         = _mm256_mul_ps(rcutoff,rcutoff);
243     /* #endif */
244
245     /* #if KERNEL_MOD_VDW=='PotentialShift' */
246     sh_vdw_invrcut6  = _mm256_set1_ps(fr->ic->sh_invrc6);
247     rvdw             = _mm256_set1_ps(fr->rvdw);
248     /* #endif */
249
250     /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
251     /*     #if KERNEL_MOD_ELEC=='PotentialSwitch'  */
252     rswitch_scalar   = fr->rcoulomb_switch;
253     rswitch          = _mm256_set1_ps(rswitch_scalar);
254     /*     #else */
255     rswitch_scalar   = fr->rvdw_switch;
256     rswitch          = _mm256_set1_ps(rswitch_scalar);
257     /*     #endif */
258     /* Setup switch parameters */
259     d_scalar         = rcutoff_scalar-rswitch_scalar;
260     d                = _mm256_set1_ps(d_scalar);
261     swV3             = _mm256_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
262     swV4             = _mm256_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
263     swV5             = _mm256_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
264     /*     #if 'Force' in KERNEL_VF */
265     swF2             = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
266     swF3             = _mm256_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
267     swF4             = _mm256_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
268     /*     #endif */
269     /* #endif */
270
271     /* Avoid stupid compiler warnings */
272     jnrA = jnrB = jnrC = jnrD = jnrE = jnrF = jnrG = jnrH = 0;
273     j_coord_offsetA = 0;
274     j_coord_offsetB = 0;
275     j_coord_offsetC = 0;
276     j_coord_offsetD = 0;
277     j_coord_offsetE = 0;
278     j_coord_offsetF = 0;
279     j_coord_offsetG = 0;
280     j_coord_offsetH = 0;
281
282     /* ## Keep track of the floating point operations we issue for reporting! */
283     /* #define OUTERFLOPS 0 */
284     outeriter        = 0;
285     inneriter        = 0;
286
287     for(iidx=0;iidx<4*DIM;iidx++)
288     {
289         scratch[iidx] = 0.0;
290     }
291
292     /* Start outer loop over neighborlists */
293     for(iidx=0; iidx<nri; iidx++)
294     {
295         /* Load shift vector for this list */
296         i_shift_offset   = DIM*shiftidx[iidx];
297
298         /* Load limits for loop over neighbors */
299         j_index_start    = jindex[iidx];
300         j_index_end      = jindex[iidx+1];
301
302         /* Get outer coordinate index */
303         inr              = iinr[iidx];
304         i_coord_offset   = DIM*inr;
305
306         /* Load i particle coords and add shift vector */
307         /* #if GEOMETRY_I == 'Particle' */
308         gmx_mm256_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
309         /* #elif GEOMETRY_I == 'Water3' */
310         gmx_mm256_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
311                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
312         /* #elif GEOMETRY_I == 'Water4' */
313         /*     #if 0 in PARTICLES_I                 */
314         gmx_mm256_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
315                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
316         /*     #else                                */
317         gmx_mm256_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
318                                                     &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
319         /*     #endif                               */
320         /* #endif                                   */
321
322         /* #if 'Force' in KERNEL_VF */
323         /*     #for I in PARTICLES_I */
324         fix{I}             = _mm256_setzero_ps();
325         fiy{I}             = _mm256_setzero_ps();
326         fiz{I}             = _mm256_setzero_ps();
327         /*     #endfor */
328         /* #endif */
329
330         /* ## For water we already preloaded parameters at the start of the kernel */
331         /* #if not 'Water' in GEOMETRY_I */
332         /* Load parameters for i particles */
333         /*     #for I in PARTICLES_ELEC_I */
334         iq{I}              = _mm256_mul_ps(facel,_mm256_set1_ps(charge[inr+{I}]));
335         /*         #define OUTERFLOPS OUTERFLOPS+1 */
336         /*         #if KERNEL_ELEC=='GeneralizedBorn' */
337         isai{I}            = _mm256_set1_ps(invsqrta[inr+{I}]);
338         /*         #endif */
339         /*     #endfor */
340         /*     #for I in PARTICLES_VDW_I */
341         vdwioffsetptr{I}   = vdwparam+2*nvdwtype*vdwtype[inr+{I}];
342         /*     #endfor */
343         /* #endif */
344
345         /* #if 'Potential' in KERNEL_VF */
346         /* Reset potential sums */
347         /*     #if KERNEL_ELEC != 'None' */
348         velecsum         = _mm256_setzero_ps();
349         /*     #endif */
350         /*     #if 'GeneralizedBorn' in KERNEL_ELEC */
351         vgbsum           = _mm256_setzero_ps();
352         /*     #endif */
353         /*     #if KERNEL_VDW != 'None' */
354         vvdwsum          = _mm256_setzero_ps();
355         /*     #endif */
356         /* #endif */
357         /*     #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
358         dvdasum          = _mm256_setzero_ps();
359         /*     #endif */
360
361         /* #for ROUND in ['Loop','Epilogue'] */
362
363         /* #if ROUND =='Loop' */
364         /* Start inner kernel loop */
365         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+7]>=0; jidx+=8)
366         {
367         /* ## First round is normal loop (next statement resets indentation) */
368         /*     #if 0 */
369         }
370         /*     #endif */
371         /* #else */
372         if(jidx<j_index_end)
373         {
374         /* ## Second round is epilogue */
375         /* #endif */
376         /* #define INNERFLOPS 0 */
377
378             /* Get j neighbor index, and coordinate index */
379             /* #if ROUND =='Loop' */
380             jnrA             = jjnr[jidx];
381             jnrB             = jjnr[jidx+1];
382             jnrC             = jjnr[jidx+2];
383             jnrD             = jjnr[jidx+3];
384             jnrE             = jjnr[jidx+4];
385             jnrF             = jjnr[jidx+5];
386             jnrG             = jjnr[jidx+6];
387             jnrH             = jjnr[jidx+7];
388             /* #else */
389             jnrlistA         = jjnr[jidx];
390             jnrlistB         = jjnr[jidx+1];
391             jnrlistC         = jjnr[jidx+2];
392             jnrlistD         = jjnr[jidx+3];
393             jnrlistE         = jjnr[jidx+4];
394             jnrlistF         = jjnr[jidx+5];
395             jnrlistG         = jjnr[jidx+6];
396             jnrlistH         = jjnr[jidx+7];
397             /* Sign of each element will be negative for non-real atoms.
398              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
399              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
400              */
401             dummy_mask = gmx_mm256_set_m128(gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx+4)),_mm_setzero_si128())),
402                                             gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128())));
403                                             
404             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
405             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
406             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
407             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
408             jnrE       = (jnrlistE>=0) ? jnrlistE : 0;
409             jnrF       = (jnrlistF>=0) ? jnrlistF : 0;
410             jnrG       = (jnrlistG>=0) ? jnrlistG : 0;
411             jnrH       = (jnrlistH>=0) ? jnrlistH : 0;
412             /* #endif */
413             j_coord_offsetA  = DIM*jnrA;
414             j_coord_offsetB  = DIM*jnrB;
415             j_coord_offsetC  = DIM*jnrC;
416             j_coord_offsetD  = DIM*jnrD;
417             j_coord_offsetE  = DIM*jnrE;
418             j_coord_offsetF  = DIM*jnrF;
419             j_coord_offsetG  = DIM*jnrG;
420             j_coord_offsetH  = DIM*jnrH;
421
422             /* load j atom coordinates */
423             /* #if GEOMETRY_J == 'Particle'             */
424             gmx_mm256_load_1rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
425                                                  x+j_coord_offsetC,x+j_coord_offsetD,
426                                                  x+j_coord_offsetE,x+j_coord_offsetF,
427                                                  x+j_coord_offsetG,x+j_coord_offsetH,
428                                                  &jx0,&jy0,&jz0);
429             /* #elif GEOMETRY_J == 'Water3'             */
430             gmx_mm256_load_3rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
431                                                  x+j_coord_offsetC,x+j_coord_offsetD,
432                                                  x+j_coord_offsetE,x+j_coord_offsetF,
433                                                  x+j_coord_offsetG,x+j_coord_offsetH,
434                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
435             /* #elif GEOMETRY_J == 'Water4'             */
436             /*     #if 0 in PARTICLES_J                 */
437             gmx_mm256_load_4rvec_8ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
438                                                  x+j_coord_offsetC,x+j_coord_offsetD,
439                                                  x+j_coord_offsetE,x+j_coord_offsetF,
440                                                  x+j_coord_offsetG,x+j_coord_offsetH,
441                                                  &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
442                                                  &jy2,&jz2,&jx3,&jy3,&jz3);
443             /*     #else                                */
444             gmx_mm256_load_3rvec_8ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
445                                                  x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
446                                                  x+j_coord_offsetE+DIM,x+j_coord_offsetF+DIM,
447                                                  x+j_coord_offsetG+DIM,x+j_coord_offsetH+DIM,
448                                                  &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
449             /*     #endif                               */
450             /* #endif                                   */
451
452             /* Calculate displacement vector */
453             /* #for I,J in PAIRS_IJ */
454             dx{I}{J}             = _mm256_sub_ps(ix{I},jx{J});
455             dy{I}{J}             = _mm256_sub_ps(iy{I},jy{J});
456             dz{I}{J}             = _mm256_sub_ps(iz{I},jz{J});
457             /*     #define INNERFLOPS INNERFLOPS+3 */
458             /* #endfor */
459
460             /* Calculate squared distance and things based on it */
461             /* #for I,J in PAIRS_IJ */
462             rsq{I}{J}            = gmx_mm256_calc_rsq_ps(dx{I}{J},dy{I}{J},dz{I}{J});
463             /*     #define INNERFLOPS INNERFLOPS+5 */
464             /* #endfor */
465
466             /* #for I,J in PAIRS_IJ */
467             /*     #if 'rinv' in INTERACTION_FLAGS[I][J] */
468             rinv{I}{J}           = gmx_mm256_invsqrt_ps(rsq{I}{J});
469             /*         #define INNERFLOPS INNERFLOPS+5 */
470             /*     #endif */
471             /* #endfor */
472
473             /* #for I,J in PAIRS_IJ */
474             /*     #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
475             /*         # if 'rinv' not in INTERACTION_FLAGS[I][J] */
476             rinvsq{I}{J}         = gmx_mm256_inv_ps(rsq{I}{J});
477             /*             #define INNERFLOPS INNERFLOPS+4 */
478             /*         #else */
479             rinvsq{I}{J}         = _mm256_mul_ps(rinv{I}{J},rinv{I}{J});
480             /*             #define INNERFLOPS INNERFLOPS+1 */
481             /*         #endif */
482             /*     #endif */
483             /* #endfor */
484
485             /* #if not 'Water' in GEOMETRY_J */
486             /* Load parameters for j particles */
487             /*     #for J in PARTICLES_ELEC_J */
488             jq{J}              = gmx_mm256_load_8real_swizzle_ps(charge+jnrA+{J},charge+jnrB+{J},
489                                                                  charge+jnrC+{J},charge+jnrD+{J},
490                                                                  charge+jnrE+{J},charge+jnrF+{J},
491                                                                  charge+jnrG+{J},charge+jnrH+{J});
492             /*         #if KERNEL_ELEC=='GeneralizedBorn' */
493             isaj{J}            = gmx_mm256_load_8real_swizzle_ps(invsqrta+jnrA+{J},invsqrta+jnrB+{J},
494                                                                  invsqrta+jnrC+{J},invsqrta+jnrD+{J},
495                                                                  invsqrta+jnrE+{J},invsqrta+jnrF+{J},
496                                                                  invsqrta+jnrG+{J},invsqrta+jnrH+{J});
497             /*         #endif */
498             /*     #endfor */
499             /*     #for J in PARTICLES_VDW_J */
500             vdwjidx{J}A        = 2*vdwtype[jnrA+{J}];
501             vdwjidx{J}B        = 2*vdwtype[jnrB+{J}];
502             vdwjidx{J}C        = 2*vdwtype[jnrC+{J}];
503             vdwjidx{J}D        = 2*vdwtype[jnrD+{J}];
504             vdwjidx{J}E        = 2*vdwtype[jnrE+{J}];
505             vdwjidx{J}F        = 2*vdwtype[jnrF+{J}];
506             vdwjidx{J}G        = 2*vdwtype[jnrG+{J}];
507             vdwjidx{J}H        = 2*vdwtype[jnrH+{J}];
508             /*     #endfor */
509             /* #endif */
510
511             /* #if 'Force' in KERNEL_VF and not 'Particle' in GEOMETRY_I */
512             /*     #for J in PARTICLES_J */
513             fjx{J}             = _mm256_setzero_ps();
514             fjy{J}             = _mm256_setzero_ps();
515             fjz{J}             = _mm256_setzero_ps();
516             /*     #endfor */
517             /* #endif */
518
519             /* #for I,J in PAIRS_IJ */
520
521             /**************************
522              * CALCULATE INTERACTIONS *
523              **************************/
524
525             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
526             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
527             /*         ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
528             if (gmx_mm256_any_lt(rsq{I}{J},rcutoff2))
529             {
530                 /*     #if 0    ## this and the next two lines is a hack to maintain auto-indentation in template file */
531             }
532             /*         #endif */
533             /*         #define INNERFLOPS INNERFLOPS+1 */
534             /*     #endif */
535
536             /*     #if 'r' in INTERACTION_FLAGS[I][J] */
537             r{I}{J}              = _mm256_mul_ps(rsq{I}{J},rinv{I}{J});
538             /*         #if ROUND == 'Epilogue' */
539             r{I}{J}              = _mm256_andnot_ps(dummy_mask,r{I}{J});
540             /*             #define INNERFLOPS INNERFLOPS+1 */
541             /*         #endif */
542             /*         #define INNERFLOPS INNERFLOPS+1 */
543             /*     #endif */
544
545             /*     ## For water geometries we already loaded parameters at the start of the kernel */
546             /*     #if not 'Water' in GEOMETRY_J */
547             /* Compute parameters for interactions between i and j atoms */
548             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
549             qq{I}{J}             = _mm256_mul_ps(iq{I},jq{J});
550             /*             #define INNERFLOPS INNERFLOPS+1 */
551             /*         #endif */
552             /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
553             gmx_mm256_load_8pair_swizzle_ps(vdwioffsetptr{I}+vdwjidx{J}A,
554                                             vdwioffsetptr{I}+vdwjidx{J}B,
555                                             vdwioffsetptr{I}+vdwjidx{J}C,
556                                             vdwioffsetptr{I}+vdwjidx{J}D,
557                                             vdwioffsetptr{I}+vdwjidx{J}E,
558                                             vdwioffsetptr{I}+vdwjidx{J}F,
559                                             vdwioffsetptr{I}+vdwjidx{J}G,
560                                             vdwioffsetptr{I}+vdwjidx{J}H,
561                                             &c6_{I}{J},&c12_{I}{J});
562             /*         #endif */
563             /*     #endif */
564
565             /*     #if 'table' in INTERACTION_FLAGS[I][J] */
566             /* Calculate table index by multiplying r with table scale and truncate to integer */
567             rt               = _mm256_mul_ps(r{I}{J},vftabscale);
568             vfitab           = _mm256_cvttps_epi32(rt);
569             vfeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
570             /*         #define INNERFLOPS INNERFLOPS+4                          */
571             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
572             vfitab_lo        = _mm256_extractf128_si256(vfitab,0x0);
573             vfitab_hi        = _mm256_extractf128_si256(vfitab,0x1);
574             /*         #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW     */
575             /*             ## 3 tables, 4 bytes per point: multiply index by 12 */
576             vfitab_lo        = _mm_slli_epi32(_mm_add_epi32(vfitab_lo,_mm_slli_epi32(vfitab_lo,1)),2);
577             vfitab_hi        = _mm_slli_epi32(_mm_add_epi32(vfitab_hi,_mm_slli_epi32(vfitab_hi,1)),2);
578             /*         #elif 'Table' in KERNEL_ELEC                             */
579             /*             ## 1 table, 4 bytes per point: multiply index by 4   */
580             vfitab_lo        = _mm_slli_epi32(vfitab_lo,2);
581             vfitab_hi        = _mm_slli_epi32(vfitab_hi,2);
582             /*         #elif 'Table' in KERNEL_VDW                              */
583             /*             ## 2 tables, 4 bytes per point: multiply index by 8  */
584             vfitab_lo        = _mm_slli_epi32(vfitab_lo,3);
585             vfitab_hi        = _mm_slli_epi32(vfitab_hi,3);
586             /*         #endif                                                   */
587             /*     #endif */
588
589             /*     ## ELECTROSTATIC INTERACTIONS */
590             /*     #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
591
592             /*         #if KERNEL_ELEC=='Coulomb' */
593
594             /* COULOMB ELECTROSTATICS */
595             velec            = _mm256_mul_ps(qq{I}{J},rinv{I}{J});
596             /*             #define INNERFLOPS INNERFLOPS+1 */
597             /*             #if 'Force' in KERNEL_VF */
598             felec            = _mm256_mul_ps(velec,rinvsq{I}{J});
599             /*                 #define INNERFLOPS INNERFLOPS+1 */
600             /*             #endif */
601
602             /*         #elif KERNEL_ELEC=='ReactionField' */
603
604             /* REACTION-FIELD ELECTROSTATICS */
605             /*             #if 'Potential' in KERNEL_VF */
606             velec            = _mm256_mul_ps(qq{I}{J},_mm256_sub_ps(_mm256_add_ps(rinv{I}{J},_mm256_mul_ps(krf,rsq{I}{J})),crf));
607             /*                 #define INNERFLOPS INNERFLOPS+4 */
608             /*             #endif */
609             /*             #if 'Force' in KERNEL_VF */
610             felec            = _mm256_mul_ps(qq{I}{J},_mm256_sub_ps(_mm256_mul_ps(rinv{I}{J},rinvsq{I}{J}),krf2));
611             /*                 #define INNERFLOPS INNERFLOPS+3 */
612             /*             #endif */
613
614             /*         #elif KERNEL_ELEC=='GeneralizedBorn' */
615
616             /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
617             isaprod          = _mm256_mul_ps(isai{I},isaj{J});
618             gbqqfactor       = _mm256_xor_ps(signbit,_mm256_mul_ps(qq{I}{J},_mm256_mul_ps(isaprod,gbinvepsdiff)));
619             gbscale          = _mm256_mul_ps(isaprod,gbtabscale);
620             /*             #define INNERFLOPS INNERFLOPS+5 */
621
622             /* Calculate generalized born table index - this is a separate table from the normal one,
623              * but we use the same procedure by multiplying r with scale and truncating to integer.
624              */
625             rt               = _mm256_mul_ps(r{I}{J},gbscale);
626             gbitab           = _mm256_cvttps_epi32(rt);
627             gbeps            = _mm256_sub_ps(rt,_mm256_round_ps(rt, _MM_FROUND_FLOOR));
628             /*         AVX1 does not support 256-bit integer operations, so now we go to 128-bit mode... */
629             gbitab_lo        = _mm256_extractf128_si256(gbitab,0x0);
630             gbitab_hi        = _mm256_extractf128_si256(gbitab,0x1);
631             gbitab_lo        = _mm_slli_epi32(gbitab_lo,2);
632             gbitab_hi        = _mm_slli_epi32(gbitab_hi,2);
633             Y                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,0)),
634                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,0)));
635             F                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,1)),
636                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,1)));
637             G                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,2)),
638                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,2)));
639             H                = gmx_mm256_set_m128(_mm_load_ps(gbtab + _mm_extract_epi32(gbitab_hi,3)),
640                                                   _mm_load_ps(gbtab + _mm_extract_epi32(gbitab_lo,3)));
641             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
642             Heps             = _mm256_mul_ps(gbeps,H);
643             Fp               = _mm256_add_ps(F,_mm256_mul_ps(gbeps,_mm256_add_ps(G,Heps)));
644             VV               = _mm256_add_ps(Y,_mm256_mul_ps(gbeps,Fp));
645             vgb              = _mm256_mul_ps(gbqqfactor,VV);
646             /*             #define INNERFLOPS INNERFLOPS+10 */
647
648             /*             #if 'Force' in KERNEL_VF */
649             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(gbeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
650             fgb              = _mm256_mul_ps(gbqqfactor,_mm256_mul_ps(FF,gbscale));
651             dvdatmp          = _mm256_mul_ps(minushalf,_mm256_add_ps(vgb,_mm256_mul_ps(fgb,r{I}{J})));
652             dvdasum          = _mm256_add_ps(dvdasum,dvdatmp);
653             /*                 #if ROUND == 'Loop' */
654             fjptrA           = dvda+jnrA;
655             fjptrB           = dvda+jnrB;
656             fjptrC           = dvda+jnrC;
657             fjptrD           = dvda+jnrD;
658             fjptrE           = dvda+jnrE;
659             fjptrF           = dvda+jnrF;
660             fjptrG           = dvda+jnrG;
661             fjptrH           = dvda+jnrH;
662             /*                 #else */
663             /* The pointers to scratch make sure that this code with compilers that take gmx_restrict seriously (e.g. icc 13) really can't screw things up. */
664             fjptrA             = (jnrlistA>=0) ? dvda+jnrA : scratch;
665             fjptrB             = (jnrlistB>=0) ? dvda+jnrB : scratch;
666             fjptrC             = (jnrlistC>=0) ? dvda+jnrC : scratch;
667             fjptrD             = (jnrlistD>=0) ? dvda+jnrD : scratch;
668             fjptrE             = (jnrlistE>=0) ? dvda+jnrE : scratch;
669             fjptrF             = (jnrlistF>=0) ? dvda+jnrF : scratch;
670             fjptrG             = (jnrlistG>=0) ? dvda+jnrG : scratch;
671             fjptrH             = (jnrlistH>=0) ? dvda+jnrH : scratch;
672             /*                 #endif */
673             gmx_mm256_increment_8real_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
674                                                  _mm256_mul_ps(dvdatmp,_mm256_mul_ps(isaj{J},isaj{J})));
675             /*                 #define INNERFLOPS INNERFLOPS+12 */
676             /*             #endif */
677             velec            = _mm256_mul_ps(qq{I}{J},rinv{I}{J});
678             /*                 #define INNERFLOPS INNERFLOPS+1 */
679             /*             #if 'Force' in KERNEL_VF */
680             felec            = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(velec,rinv{I}{J}),fgb),rinv{I}{J});
681             /*                 #define INNERFLOPS INNERFLOPS+3 */
682             /*             #endif */
683
684             /*         #elif KERNEL_ELEC=='Ewald' */
685             /* EWALD ELECTROSTATICS */
686             
687             /* Analytical PME correction */
688             zeta2            = _mm256_mul_ps(beta2,rsq{I}{J});
689             /*             #if 'Force' in KERNEL_VF */
690             rinv3            = _mm256_mul_ps(rinvsq{I}{J},rinv{I}{J});
691             pmecorrF         = gmx_mm256_pmecorrF_ps(zeta2);
692             felec            = _mm256_add_ps( _mm256_mul_ps(pmecorrF,beta3), rinv3);
693             felec            = _mm256_mul_ps(qq{I}{J},felec);
694             /*                 #define INNERFLOPS INNERFLOPS+31 */
695             /*             #endif */
696             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
697             pmecorrV         = gmx_mm256_pmecorrV_ps(zeta2);
698             pmecorrV         = _mm256_mul_ps(pmecorrV,beta);
699             /*                 #define INNERFLOPS INNERFLOPS+27 */
700             /*                 #if KERNEL_MOD_ELEC=='PotentialShift' */
701             velec            = _mm256_sub_ps(_mm256_sub_ps(rinv{I}{J},sh_ewald),pmecorrV);
702             /*                     #define INNERFLOPS INNERFLOPS+21 */
703             /*                 #else */
704             velec            = _mm256_sub_ps(rinv{I}{J},pmecorrV);
705             /*                 #endif */
706             velec            = _mm256_mul_ps(qq{I}{J},velec);
707             /*             #endif */
708             
709             /*         #elif KERNEL_ELEC=='CubicSplineTable' */
710
711             /* CUBIC SPLINE TABLE ELECTROSTATICS */
712             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
713                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
714             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
715                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
716             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
717                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
718             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
719                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
720             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
721             Heps             = _mm256_mul_ps(vfeps,H);
722             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
723             /*             #define INNERFLOPS INNERFLOPS+4 */
724             /*             #if 'Potential' in KERNEL_VF */
725             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
726             velec            = _mm256_mul_ps(qq{I}{J},VV);
727             /*                 #define INNERFLOPS INNERFLOPS+3 */
728             /*             #endif */
729             /*             #if 'Force' in KERNEL_VF */
730             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
731             felec            = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_mul_ps(qq{I}{J},FF),_mm256_mul_ps(vftabscale,rinv{I}{J})));
732             /*                 #define INNERFLOPS INNERFLOPS+7 */
733             /*             #endif */
734             /*         #endif */
735             /*         ## End of check for electrostatics interaction forms */
736             /*     #endif */
737             /*     ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
738
739             /*     #if 'vdw' in INTERACTION_FLAGS[I][J] */
740
741             /*         #if KERNEL_VDW=='LennardJones' */
742
743             /* LENNARD-JONES DISPERSION/REPULSION */
744
745             rinvsix          = _mm256_mul_ps(_mm256_mul_ps(rinvsq{I}{J},rinvsq{I}{J}),rinvsq{I}{J});
746             /*             #define INNERFLOPS INNERFLOPS+2 */
747             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
748             vvdw6            = _mm256_mul_ps(c6_{I}{J},rinvsix);
749             vvdw12           = _mm256_mul_ps(c12_{I}{J},_mm256_mul_ps(rinvsix,rinvsix));
750             /*                 #define INNERFLOPS INNERFLOPS+3 */
751             /*                 #if KERNEL_MOD_VDW=='PotentialShift' */
752             vvdw             = _mm256_sub_ps(_mm256_mul_ps( _mm256_sub_ps(vvdw12 , _mm256_mul_ps(c12_{I}{J},_mm256_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6))), one_twelfth) ,
753                                           _mm256_mul_ps( _mm256_sub_ps(vvdw6,_mm256_mul_ps(c6_{I}{J},sh_vdw_invrcut6)),one_sixth));
754             /*                     #define INNERFLOPS INNERFLOPS+8 */
755             /*                 #else */
756             vvdw             = _mm256_sub_ps( _mm256_mul_ps(vvdw12,one_twelfth) , _mm256_mul_ps(vvdw6,one_sixth) );
757             /*                     #define INNERFLOPS INNERFLOPS+3 */
758             /*                 #endif */
759             /*                 ## Check for force inside potential check, i.e. this means we already did the potential part */
760             /*                 #if 'Force' in KERNEL_VF */
761             fvdw             = _mm256_mul_ps(_mm256_sub_ps(vvdw12,vvdw6),rinvsq{I}{J});
762             /*                     #define INNERFLOPS INNERFLOPS+2 */
763             /*                 #endif */
764             /*             #elif KERNEL_VF=='Force' */
765             /*                 ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
766             fvdw             = _mm256_mul_ps(_mm256_sub_ps(_mm256_mul_ps(c12_{I}{J},rinvsix),c6_{I}{J}),_mm256_mul_ps(rinvsix,rinvsq{I}{J}));
767             /*                 #define INNERFLOPS INNERFLOPS+4 */
768             /*             #endif */
769
770             /*         #elif KERNEL_VDW=='CubicSplineTable' */
771
772             /* CUBIC SPLINE TABLE DISPERSION */
773             /*             #if 'Table' in KERNEL_ELEC */
774             vfitab_lo        = _mm_add_epi32(vfitab_lo,ifour);
775             vfitab_hi        = _mm_add_epi32(vfitab_hi,ifour);
776             /*             #endif                     */
777             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
778                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
779             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
780                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
781             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
782                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
783             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
784                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
785             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
786             Heps             = _mm256_mul_ps(vfeps,H);
787             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
788             /*             #define INNERFLOPS INNERFLOPS+4 */
789             /*             #if 'Potential' in KERNEL_VF */
790             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
791             vvdw6            = _mm256_mul_ps(c6_{I}{J},VV);
792             /*                 #define INNERFLOPS INNERFLOPS+3 */
793             /*             #endif */
794             /*             #if 'Force' in KERNEL_VF */
795             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
796             fvdw6            = _mm256_mul_ps(c6_{I}{J},FF);
797             /*                 #define INNERFLOPS INNERFLOPS+4 */
798             /*             #endif */
799
800             /* CUBIC SPLINE TABLE REPULSION */
801             vfitab_lo        = _mm_add_epi32(vfitab_lo,ifour);
802             vfitab_hi        = _mm_add_epi32(vfitab_hi,ifour);
803             Y                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,0)),
804                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,0)));
805             F                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,1)),
806                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,1)));
807             G                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,2)),
808                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,2)));
809             H                = gmx_mm256_set_m128(_mm_load_ps(vftab + _mm_extract_epi32(vfitab_hi,3)),
810                                                   _mm_load_ps(vftab + _mm_extract_epi32(vfitab_lo,3)));
811             GMX_MM256_HALFTRANSPOSE4_PS(Y,F,G,H);
812             Heps             = _mm256_mul_ps(vfeps,H);
813             Fp               = _mm256_add_ps(F,_mm256_mul_ps(vfeps,_mm256_add_ps(G,Heps)));
814             /*             #define INNERFLOPS INNERFLOPS+4 */
815             /*             #if 'Potential' in KERNEL_VF */
816             VV               = _mm256_add_ps(Y,_mm256_mul_ps(vfeps,Fp));
817             vvdw12           = _mm256_mul_ps(c12_{I}{J},VV);
818             /*                 #define INNERFLOPS INNERFLOPS+3 */
819             /*             #endif */
820             /*             #if 'Force' in KERNEL_VF */
821             FF               = _mm256_add_ps(Fp,_mm256_mul_ps(vfeps,_mm256_add_ps(G,_mm256_add_ps(Heps,Heps))));
822             fvdw12           = _mm256_mul_ps(c12_{I}{J},FF);
823             /*                 #define INNERFLOPS INNERFLOPS+5 */
824             /*             #endif */
825             /*             #if 'Potential' in KERNEL_VF */
826             vvdw             = _mm256_add_ps(vvdw12,vvdw6);
827             /*                 #define INNERFLOPS INNERFLOPS+1 */
828             /*             #endif */
829             /*             #if 'Force' in KERNEL_VF */
830             fvdw             = _mm256_xor_ps(signbit,_mm256_mul_ps(_mm256_add_ps(fvdw6,fvdw12),_mm256_mul_ps(vftabscale,rinv{I}{J})));
831             /*                 #define INNERFLOPS INNERFLOPS+4 */
832             /*             #endif */
833             /*         #endif */
834             /*         ## End of check for vdw interaction forms */
835             /*     #endif */
836             /*     ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
837
838             /*     #if 'switch' in INTERACTION_FLAGS[I][J] */
839             d                = _mm256_sub_ps(r{I}{J},rswitch);
840             d                = _mm256_max_ps(d,_mm256_setzero_ps());
841             d2               = _mm256_mul_ps(d,d);
842             sw               = _mm256_add_ps(one,_mm256_mul_ps(d2,_mm256_mul_ps(d,_mm256_add_ps(swV3,_mm256_mul_ps(d,_mm256_add_ps(swV4,_mm256_mul_ps(d,swV5)))))));
843             /*         #define INNERFLOPS INNERFLOPS+10 */
844
845             /*         #if 'Force' in KERNEL_VF */
846             dsw              = _mm256_mul_ps(d2,_mm256_add_ps(swF2,_mm256_mul_ps(d,_mm256_add_ps(swF3,_mm256_mul_ps(d,swF4)))));
847             /*             #define INNERFLOPS INNERFLOPS+5 */
848             /*         #endif */
849
850             /* Evaluate switch function */
851             /*         #if 'Force' in KERNEL_VF */
852             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
853             /*             #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
854             felec            = _mm256_sub_ps( _mm256_mul_ps(felec,sw) , _mm256_mul_ps(rinv{I}{J},_mm256_mul_ps(velec,dsw)) );
855             /*                 #define INNERFLOPS INNERFLOPS+4 */
856             /*             #endif */
857             /*             #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
858             fvdw             = _mm256_sub_ps( _mm256_mul_ps(fvdw,sw) , _mm256_mul_ps(rinv{I}{J},_mm256_mul_ps(vvdw,dsw)) );
859             /*                 #define INNERFLOPS INNERFLOPS+4 */
860             /*             #endif */
861             /*         #endif */
862             /*         #if 'Potential' in KERNEL_VF */
863             /*             #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
864             velec            = _mm256_mul_ps(velec,sw);
865             /*                 #define INNERFLOPS INNERFLOPS+1 */
866             /*             #endif */
867             /*             #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
868             vvdw             = _mm256_mul_ps(vvdw,sw);
869             /*                 #define INNERFLOPS INNERFLOPS+1 */
870             /*             #endif */
871             /*         #endif */
872             /*     #endif */
873             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
874             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
875             cutoff_mask      = _mm256_cmp_ps(rsq{I}{J},rcutoff2,_CMP_LT_OQ);
876             /*         #define INNERFLOPS INNERFLOPS+1 */
877             /*     #endif */
878
879             /*     #if 'Potential' in KERNEL_VF */
880             /* Update potential sum for this i atom from the interaction with this j atom. */
881             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
882             /*             #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
883             velec            = _mm256_and_ps(velec,cutoff_mask);
884             /*                 #define INNERFLOPS INNERFLOPS+1 */
885             /*             #endif                                       */
886             /*             #if ROUND == 'Epilogue' */
887             velec            = _mm256_andnot_ps(dummy_mask,velec);
888             /*             #endif */
889             velecsum         = _mm256_add_ps(velecsum,velec);
890             /*             #define INNERFLOPS INNERFLOPS+1 */
891             /*             #if KERNEL_ELEC=='GeneralizedBorn' */
892             /*             #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
893             vgb              = _mm256_and_ps(vgb,cutoff_mask);
894             /*                 #define INNERFLOPS INNERFLOPS+1 */
895             /*             #endif                                       */
896             /*             #if ROUND == 'Epilogue' */
897             vgb              = _mm256_andnot_ps(dummy_mask,vgb);
898             /*             #endif */
899             vgbsum           = _mm256_add_ps(vgbsum,vgb);
900             /*                 #define INNERFLOPS INNERFLOPS+1 */
901             /*             #endif */
902             /*         #endif */
903             /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
904             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
905             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
906             vvdw             = _mm256_and_ps(vvdw,cutoff_mask);
907             /*                 #define INNERFLOPS INNERFLOPS+1 */
908             /*             #endif                                       */
909             /*             #if ROUND == 'Epilogue' */
910             vvdw             = _mm256_andnot_ps(dummy_mask,vvdw);
911             /*             #endif */
912             vvdwsum          = _mm256_add_ps(vvdwsum,vvdw);
913             /*             #define INNERFLOPS INNERFLOPS+1 */
914             /*         #endif */
915             /*     #endif */
916
917             /*     #if 'Force' in KERNEL_VF */
918
919             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
920             fscal            = _mm256_add_ps(felec,fvdw);
921             /*             #define INNERFLOPS INNERFLOPS+1 */
922             /*         #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
923             fscal            = felec;
924             /*         #elif 'vdw' in INTERACTION_FLAGS[I][J] */
925             fscal            = fvdw;
926             /*        #endif */
927
928             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
929             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
930             fscal            = _mm256_and_ps(fscal,cutoff_mask);
931             /*                 #define INNERFLOPS INNERFLOPS+1 */
932             /*             #endif                                       */
933
934             /*             #if ROUND == 'Epilogue' */
935             fscal            = _mm256_andnot_ps(dummy_mask,fscal);
936             /*             #endif */
937
938             /* Calculate temporary vectorial force */
939             tx               = _mm256_mul_ps(fscal,dx{I}{J});
940             ty               = _mm256_mul_ps(fscal,dy{I}{J});
941             tz               = _mm256_mul_ps(fscal,dz{I}{J});
942
943             /* Update vectorial force */
944             fix{I}             = _mm256_add_ps(fix{I},tx);
945             fiy{I}             = _mm256_add_ps(fiy{I},ty);
946             fiz{I}             = _mm256_add_ps(fiz{I},tz);
947             /*             #define INNERFLOPS INNERFLOPS+6 */
948
949             /* #if GEOMETRY_I == 'Particle'             */
950             /*     #if ROUND == 'Loop' */
951             fjptrA             = f+j_coord_offsetA;
952             fjptrB             = f+j_coord_offsetB;
953             fjptrC             = f+j_coord_offsetC;
954             fjptrD             = f+j_coord_offsetD;
955             fjptrE             = f+j_coord_offsetE;
956             fjptrF             = f+j_coord_offsetF;
957             fjptrG             = f+j_coord_offsetG;
958             fjptrH             = f+j_coord_offsetH;
959             /*     #else */
960             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
961             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
962             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
963             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
964             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
965             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
966             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
967             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
968             /*     #endif */
969             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,tx,ty,tz);
970             /*     #define INNERFLOPS INNERFLOPS+3      */
971             /* #else                                    */
972             fjx{J}             = _mm256_add_ps(fjx{J},tx);
973             fjy{J}             = _mm256_add_ps(fjy{J},ty);
974             fjz{J}             = _mm256_add_ps(fjz{J},tz);
975             /*     #define INNERFLOPS INNERFLOPS+3      */
976             /* #endif                                   */
977
978             /*     #endif */
979
980             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
981             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
982             /*         #if 0    ## This and next two lines is a hack to maintain indentation in template file */
983             {
984                 /*     #endif */
985             }
986             /*     #endif */
987             /*    ## End of check for the interaction being outside the cutoff */
988
989             /* #endfor */
990             /* ## End of loop over i-j interaction pairs */
991
992             /* #if GEOMETRY_I != 'Particle' */
993             /*     #if ROUND == 'Loop' */
994             fjptrA             = f+j_coord_offsetA;
995             fjptrB             = f+j_coord_offsetB;
996             fjptrC             = f+j_coord_offsetC;
997             fjptrD             = f+j_coord_offsetD;
998             fjptrE             = f+j_coord_offsetE;
999             fjptrF             = f+j_coord_offsetF;
1000             fjptrG             = f+j_coord_offsetG;
1001             fjptrH             = f+j_coord_offsetH;
1002             /*     #else */
1003             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1004             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1005             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1006             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1007             fjptrE             = (jnrlistE>=0) ? f+j_coord_offsetE : scratch;
1008             fjptrF             = (jnrlistF>=0) ? f+j_coord_offsetF : scratch;
1009             fjptrG             = (jnrlistG>=0) ? f+j_coord_offsetG : scratch;
1010             fjptrH             = (jnrlistH>=0) ? f+j_coord_offsetH : scratch;
1011             /*     #endif */
1012             /* #endif */
1013
1014             /* #if 'Water' in GEOMETRY_I and GEOMETRY_J == 'Particle' */
1015             gmx_mm256_decrement_1rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,fjx0,fjy0,fjz0);
1016             /*     #define INNERFLOPS INNERFLOPS+3      */
1017             /* #elif GEOMETRY_J == 'Water3'             */
1018             gmx_mm256_decrement_3rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
1019                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1020             /*     #define INNERFLOPS INNERFLOPS+9      */
1021             /* #elif GEOMETRY_J == 'Water4'             */
1022             /*     #if 0 in PARTICLES_J                 */
1023             gmx_mm256_decrement_4rvec_8ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjptrE,fjptrF,fjptrG,fjptrH,
1024                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
1025                                                       fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1026             /*     #define INNERFLOPS INNERFLOPS+12     */
1027             /*     #else                                */
1028             gmx_mm256_decrement_3rvec_8ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
1029                                                       fjptrE+DIM,fjptrF+DIM,fjptrG+DIM,fjptrH+DIM,
1030                                                       fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1031             /*     #define INNERFLOPS INNERFLOPS+9      */
1032             /*     #endif                               */
1033             /* #endif                                   */
1034
1035             /* Inner loop uses {INNERFLOPS} flops */
1036         }
1037
1038         /* #endfor */
1039
1040         /* End of innermost loop */
1041
1042         /* #if 'Force' in KERNEL_VF */
1043         /*     #if GEOMETRY_I == 'Particle'            */
1044         gmx_mm256_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
1045                                                  f+i_coord_offset,fshift+i_shift_offset);
1046         /*         #define OUTERFLOPS OUTERFLOPS+6     */
1047         /*     #elif GEOMETRY_I == 'Water3'            */
1048         gmx_mm256_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1049                                                  f+i_coord_offset,fshift+i_shift_offset);
1050         /*         #define OUTERFLOPS OUTERFLOPS+18    */
1051         /*     #elif GEOMETRY_I == 'Water4'            */
1052         /*         #if 0 in PARTICLES_I                */
1053         gmx_mm256_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1054                                                  f+i_coord_offset,fshift+i_shift_offset);
1055         /*             #define OUTERFLOPS OUTERFLOPS+24    */
1056         /*         #else                               */
1057         gmx_mm256_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1058                                                  f+i_coord_offset+DIM,fshift+i_shift_offset);
1059         /*             #define OUTERFLOPS OUTERFLOPS+18    */
1060         /*         #endif                              */
1061         /*     #endif                                  */
1062         /* #endif                                      */
1063
1064         /* #if 'Potential' in KERNEL_VF */
1065         ggid                        = gid[iidx];
1066         /* Update potential energies */
1067         /*     #if KERNEL_ELEC != 'None' */
1068         gmx_mm256_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
1069         /*         #define OUTERFLOPS OUTERFLOPS+1 */
1070         /*     #endif */
1071         /*     #if 'GeneralizedBorn' in KERNEL_ELEC */
1072         gmx_mm256_update_1pot_ps(vgbsum,kernel_data->energygrp_polarization+ggid);
1073         /*         #define OUTERFLOPS OUTERFLOPS+1 */
1074         /*     #endif */
1075         /*     #if KERNEL_VDW != 'None' */
1076         gmx_mm256_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
1077         /*         #define OUTERFLOPS OUTERFLOPS+1 */
1078         /*     #endif */
1079         /* #endif */
1080         /*     #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
1081         dvdasum = _mm256_mul_ps(dvdasum, _mm256_mul_ps(isai{I},isai{I}));
1082         gmx_mm256_update_1pot_ps(dvdasum,dvda+inr);
1083         /*     #endif */
1084
1085         /* Increment number of inner iterations */
1086         inneriter                  += j_index_end - j_index_start;
1087
1088         /* Outer loop uses {OUTERFLOPS} flops */
1089     }
1090
1091     /* Increment number of outer iterations */
1092     outeriter        += nri;
1093
1094     /* Update outer/inner flops */
1095     /* ## NB: This is not important, it just affects the flopcount. However, since our preprocessor is */
1096     /* ## primitive and replaces aggressively even in strings inside these directives, we need to      */
1097     /* ## assemble the main part of the name (containing KERNEL/ELEC/VDW) directly in the source.      */
1098     /* #if GEOMETRY_I == 'Water3'            */
1099     /*     #define ISUFFIX '_W3'             */
1100     /* #elif GEOMETRY_I == 'Water4'          */
1101     /*     #define ISUFFIX '_W4'             */
1102     /* #else                                 */
1103     /*     #define ISUFFIX ''                */
1104     /* #endif                                */
1105     /* #if GEOMETRY_J == 'Water3'            */
1106     /*     #define JSUFFIX 'W3'              */
1107     /* #elif GEOMETRY_J == 'Water4'          */
1108     /*     #define JSUFFIX 'W4'              */
1109     /* #else                                 */
1110     /*     #define JSUFFIX ''                */
1111     /* #endif                                */
1112     /* #if 'PotentialAndForce' in KERNEL_VF  */
1113     /*     #define VFSUFFIX  '_VF'           */
1114     /* #elif 'Potential' in KERNEL_VF        */
1115     /*     #define VFSUFFIX '_V'             */
1116     /* #else                                 */
1117     /*     #define VFSUFFIX '_F'             */
1118     /* #endif                                */
1119
1120     /* #if KERNEL_ELEC != 'None' and KERNEL_VDW != 'None' */
1121     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1122     /* #elif KERNEL_ELEC != 'None' */
1123     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1124     /* #else */
1125     inc_nrnb(nrnb,eNR_NBKERNEL_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1126     /* #endif  */
1127 }