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