a25ec0e3908bfd8de78d3bffcb773a09b289c95d
[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, 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 'GeneralizedBorn' in KERNEL_ELEC */
139     __m128i          gbitab;
140     __m256d          vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,gbeps,dvdatmp;
141     __m256d          minushalf = _mm256_set1_pd(-0.5);
142     real             *invsqrta,*dvda,*gbtab;
143     /* #endif */
144     /* #if KERNEL_VDW != 'None' */
145     int              nvdwtype;
146     __m256d          rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
147     int              *vdwtype;
148     real             *vdwparam;
149     __m256d          one_sixth   = _mm256_set1_pd(1.0/6.0);
150     __m256d          one_twelfth = _mm256_set1_pd(1.0/12.0);
151     /* #endif */
152     /* #if 'Table' in KERNEL_ELEC or 'GeneralizedBorn' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
153     __m128i          vfitab;
154     __m128i          ifour       = _mm_set1_epi32(4);
155     __m256d          rt,vfeps,vftabscale,Y,F,G,H,Heps,Fp,VV,FF;
156     real             *vftab;
157     /* #endif */
158     /* #if 'LJEwald' in KERNEL_VDW */
159     /* #for I,J in PAIRS_IJ */
160     __m256d           c6grid_{I}{J};
161     /* #endfor */
162     real             *vdwgridparam;
163     __m256d           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
164     __m256d           one_half  = _mm256_set1_pd(0.5);
165     __m256d           minus_one = _mm256_set1_pd(-1.0);
166     /* #endif */
167     /* #if 'Ewald' in KERNEL_ELEC */
168     __m128i          ewitab;
169     __m256d          ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
170     __m256d          beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
171     real             *ewtab;
172     /* #endif */
173     /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
174     __m256d          rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
175     real             rswitch_scalar,d_scalar;
176     /* #endif */
177     __m256d          dummy_mask,cutoff_mask;
178     __m128           tmpmask0,tmpmask1;
179     __m256d          signbit = _mm256_castsi256_pd( _mm256_set1_epi32(0x80000000) );
180     __m256d          one     = _mm256_set1_pd(1.0);
181     __m256d          two     = _mm256_set1_pd(2.0);
182     x                = xx[0];
183     f                = ff[0];
184
185     nri              = nlist->nri;
186     iinr             = nlist->iinr;
187     jindex           = nlist->jindex;
188     jjnr             = nlist->jjnr;
189     shiftidx         = nlist->shift;
190     gid              = nlist->gid;
191     shiftvec         = fr->shift_vec[0];
192     fshift           = fr->fshift[0];
193     /* #if KERNEL_ELEC != 'None' */
194     facel            = _mm256_set1_pd(fr->ic->epsfac);
195     charge           = mdatoms->chargeA;
196     /*     #if 'ReactionField' in KERNEL_ELEC */
197     krf              = _mm256_set1_pd(fr->ic->k_rf);
198     krf2             = _mm256_set1_pd(fr->ic->k_rf*2.0);
199     crf              = _mm256_set1_pd(fr->ic->c_rf);
200     /*     #endif */
201     /* #endif */
202     /* #if KERNEL_VDW != 'None' */
203     nvdwtype         = fr->ntype;
204     vdwparam         = fr->nbfp;
205     vdwtype          = mdatoms->typeA;
206     /* #endif */
207     /* #if 'LJEwald' in KERNEL_VDW */
208     vdwgridparam     = fr->ljpme_c6grid;
209     sh_lj_ewald      = _mm256_set1_pd(fr->ic->sh_lj_ewald);
210     ewclj            = _mm256_set1_pd(fr->ic->ewaldcoeff_lj);
211     ewclj2           = _mm256_mul_pd(minus_one,_mm256_mul_pd(ewclj,ewclj));
212     /* #endif */
213
214     /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
215     vftab            = kernel_data->table_elec_vdw->data;
216     vftabscale       = _mm256_set1_pd(kernel_data->table_elec_vdw->scale);
217     /* #elif 'Table' in KERNEL_ELEC */
218     vftab            = kernel_data->table_elec->data;
219     vftabscale       = _mm256_set1_pd(kernel_data->table_elec->scale);
220     /* #elif 'Table' in KERNEL_VDW */
221     vftab            = kernel_data->table_vdw->data;
222     vftabscale       = _mm256_set1_pd(kernel_data->table_vdw->scale);
223     /* #endif */
224
225     /* #if 'Ewald' in KERNEL_ELEC */
226     sh_ewald         = _mm256_set1_pd(fr->ic->sh_ewald);
227     beta             = _mm256_set1_pd(fr->ic->ewaldcoeff_q);
228     beta2            = _mm256_mul_pd(beta,beta);
229     beta3            = _mm256_mul_pd(beta,beta2);
230
231     /*     #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
232     ewtab            = fr->ic->tabq_coul_F;
233     ewtabscale       = _mm256_set1_pd(fr->ic->tabq_scale);
234     ewtabhalfspace   = _mm256_set1_pd(0.5/fr->ic->tabq_scale);
235     /*     #else */
236     ewtab            = fr->ic->tabq_coul_FDV0;
237     ewtabscale       = _mm256_set1_pd(fr->ic->tabq_scale);
238     ewtabhalfspace   = _mm256_set1_pd(0.5/fr->ic->tabq_scale);
239      /*     #endif */
240     /* #endif */
241
242     /* #if KERNEL_ELEC=='GeneralizedBorn' */
243     invsqrta         = fr->invsqrta;
244     dvda             = fr->dvda;
245     gbtabscale       = _mm256_set1_pd(fr->gbtab->scale);
246     gbtab            = fr->gbtab->data;
247     gbinvepsdiff     = _mm256_set1_pd((1.0/fr->ic->epsilon_r) - (1.0/fr->gb_epsilon_solvent));
248     /* #endif */
249
250     /* #if 'Water' in GEOMETRY_I */
251     /* Setup water-specific parameters */
252     inr              = nlist->iinr[0];
253     /*     #for I in PARTICLES_ELEC_I */
254     iq{I}              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+{I}]));
255     /*     #endfor */
256     /*     #for I in PARTICLES_VDW_I */
257     vdwioffsetptr{I}   = vdwparam+2*nvdwtype*vdwtype[inr+{I}];
258     /*         #if 'LJEwald' in KERNEL_VDW */
259     vdwgridioffsetptr{I} = vdwgridparam+2*nvdwtype*vdwtype[inr+{I}];
260     /*         #endif */
261     /*     #endfor */
262     /* #endif */
263
264     /* #if 'Water' in GEOMETRY_J */
265     /*     #for J in PARTICLES_ELEC_J */
266     jq{J}              = _mm256_set1_pd(charge[inr+{J}]);
267     /*     #endfor */
268     /*     #for J in PARTICLES_VDW_J */
269     vdwjidx{J}A        = 2*vdwtype[inr+{J}];
270     /*     #endfor */
271     /*     #for I,J in PAIRS_IJ */
272     /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
273     qq{I}{J}             = _mm256_mul_pd(iq{I},jq{J});
274     /*         #endif */
275     /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
276     /*             #if 'LJEwald' in KERNEL_VDW */
277     c6_{I}{J}            = _mm256_set1_pd(vdwioffsetptr{I}[vdwjidx{J}A]);
278     c12_{I}{J}           = _mm256_set1_pd(vdwioffsetptr{I}[vdwjidx{J}A+1]);
279     c6grid_{I}{J}        = _mm256_set1_pd(vdwgridioffsetptr{I}[vdwjidx{J}A]);
280     /*             #else */
281     c6_{I}{J}            = _mm256_set1_pd(vdwioffsetptr{I}[vdwjidx{J}A]);
282     c12_{I}{J}           = _mm256_set1_pd(vdwioffsetptr{I}[vdwjidx{J}A+1]);
283     /*             #endif */
284     /*         #endif */
285     /*     #endfor */
286     /* #endif */
287
288     /* #if KERNEL_MOD_ELEC!='None' or KERNEL_MOD_VDW!='None' */
289     /*     #if KERNEL_ELEC!='None' */
290     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
291     rcutoff_scalar   = fr->ic->rcoulomb;
292     /*     #else */
293     rcutoff_scalar   = fr->ic->rvdw;
294     /*     #endif */
295     rcutoff          = _mm256_set1_pd(rcutoff_scalar);
296     rcutoff2         = _mm256_mul_pd(rcutoff,rcutoff);
297     /* #endif */
298
299     /* #if KERNEL_MOD_VDW=='PotentialShift' */
300     sh_vdw_invrcut6  = _mm256_set1_pd(fr->ic->sh_invrc6);
301     rvdw             = _mm256_set1_pd(fr->ic->rvdw);
302     /* #endif */
303
304     /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
305     /*     #if KERNEL_MOD_ELEC=='PotentialSwitch'  */
306     rswitch_scalar   = fr->ic->rcoulomb_switch;
307     rswitch          = _mm256_set1_pd(rswitch_scalar);
308     /*     #else */
309     rswitch_scalar   = fr->ic->rvdw_switch;
310     rswitch          = _mm256_set1_pd(rswitch_scalar);
311     /*     #endif */
312     /* Setup switch parameters */
313     d_scalar         = rcutoff_scalar-rswitch_scalar;
314     d                = _mm256_set1_pd(d_scalar);
315     swV3             = _mm256_set1_pd(-10.0/(d_scalar*d_scalar*d_scalar));
316     swV4             = _mm256_set1_pd( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
317     swV5             = _mm256_set1_pd( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
318     /*     #if 'Force' in KERNEL_VF */
319     swF2             = _mm256_set1_pd(-30.0/(d_scalar*d_scalar*d_scalar));
320     swF3             = _mm256_set1_pd( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
321     swF4             = _mm256_set1_pd(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
322     /*     #endif */
323     /* #endif */
324
325     /* Avoid stupid compiler warnings */
326     jnrA = jnrB = jnrC = jnrD = 0;
327     j_coord_offsetA = 0;
328     j_coord_offsetB = 0;
329     j_coord_offsetC = 0;
330     j_coord_offsetD = 0;
331
332     /* ## Keep track of the floating point operations we issue for reporting! */
333     /* #define OUTERFLOPS 0 */
334     outeriter        = 0;
335     inneriter        = 0;
336
337     for(iidx=0;iidx<4*DIM;iidx++)
338     {
339         scratch[iidx] = 0.0;
340     }
341
342     /* Start outer loop over neighborlists */
343     for(iidx=0; iidx<nri; iidx++)
344     {
345         /* Load shift vector for this list */
346         i_shift_offset   = DIM*shiftidx[iidx];
347
348         /* Load limits for loop over neighbors */
349         j_index_start    = jindex[iidx];
350         j_index_end      = jindex[iidx+1];
351
352         /* Get outer coordinate index */
353         inr              = iinr[iidx];
354         i_coord_offset   = DIM*inr;
355
356         /* Load i particle coords and add shift vector */
357         /* #if GEOMETRY_I == 'Particle' */
358         gmx_mm256_load_shift_and_1rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
359         /* #elif GEOMETRY_I == 'Water3' */
360         gmx_mm256_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
361                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
362         /* #elif GEOMETRY_I == 'Water4' */
363         /*     #if 0 in PARTICLES_I                 */
364         gmx_mm256_load_shift_and_4rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset,
365                                                     &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
366         /*     #else                                */
367         gmx_mm256_load_shift_and_3rvec_broadcast_pd(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
368                                                     &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
369         /*     #endif                               */
370         /* #endif                                   */
371
372         /* #if 'Force' in KERNEL_VF */
373         /*     #for I in PARTICLES_I */
374         fix{I}             = _mm256_setzero_pd();
375         fiy{I}             = _mm256_setzero_pd();
376         fiz{I}             = _mm256_setzero_pd();
377         /*     #endfor */
378         /* #endif */
379
380         /* ## For water we already preloaded parameters at the start of the kernel */
381         /* #if not 'Water' in GEOMETRY_I */
382         /* Load parameters for i particles */
383         /*     #for I in PARTICLES_ELEC_I */
384         iq{I}              = _mm256_mul_pd(facel,_mm256_set1_pd(charge[inr+{I}]));
385         /*         #define OUTERFLOPS OUTERFLOPS+1 */
386         /*         #if KERNEL_ELEC=='GeneralizedBorn' */
387         isai{I}            = _mm256_set1_pd(invsqrta[inr+{I}]);
388         /*         #endif */
389         /*     #endfor */
390         /*     #for I in PARTICLES_VDW_I */
391         vdwioffsetptr{I}   = vdwparam+2*nvdwtype*vdwtype[inr+{I}];
392         /*         #if 'LJEwald' in KERNEL_VDW */
393         vdwgridioffsetptr{I} = vdwgridparam+2*nvdwtype*vdwtype[inr+{I}];
394         /*         #endif */
395         /*     #endfor */
396         /* #endif */
397
398         /* #if 'Potential' in KERNEL_VF */
399         /* Reset potential sums */
400         /*     #if KERNEL_ELEC != 'None' */
401         velecsum         = _mm256_setzero_pd();
402         /*     #endif */
403         /*     #if 'GeneralizedBorn' in KERNEL_ELEC */
404         vgbsum           = _mm256_setzero_pd();
405         /*     #endif */
406         /*     #if KERNEL_VDW != 'None' */
407         vvdwsum          = _mm256_setzero_pd();
408         /*     #endif */
409         /* #endif */
410         /*     #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
411         dvdasum          = _mm256_setzero_pd();
412         /*     #endif */
413
414         /* #for ROUND in ['Loop','Epilogue'] */
415
416         /* #if ROUND =='Loop' */
417         /* Start inner kernel loop */
418         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
419         {
420         /* ## First round is normal loop (next statement resets indentation) */
421         /*     #if 0 */
422         }
423         /*     #endif */
424         /* #else */
425         if(jidx<j_index_end)
426         {
427         /* ## Second round is epilogue */
428         /* #endif */
429         /* #define INNERFLOPS 0 */
430
431             /* Get j neighbor index, and coordinate index */
432             /* #if ROUND =='Loop' */
433             jnrA             = jjnr[jidx];
434             jnrB             = jjnr[jidx+1];
435             jnrC             = jjnr[jidx+2];
436             jnrD             = jjnr[jidx+3];
437             /* #else */
438             jnrlistA         = jjnr[jidx];
439             jnrlistB         = jjnr[jidx+1];
440             jnrlistC         = jjnr[jidx+2];
441             jnrlistD         = jjnr[jidx+3];
442             /* Sign of each element will be negative for non-real atoms.
443              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
444              * so use it as val = _mm_andnot_pd(mask,val) to clear dummy entries.
445              */
446             tmpmask0 = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
447
448             tmpmask1 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(3,3,2,2));
449             tmpmask0 = _mm_permute_ps(tmpmask0,_GMX_MM_PERMUTE(1,1,0,0));
450             dummy_mask = _mm256_castps_pd(gmx_mm256_set_m128(tmpmask1,tmpmask0));
451
452             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
453             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
454             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
455             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
456             /* #endif */
457             j_coord_offsetA  = DIM*jnrA;
458             j_coord_offsetB  = DIM*jnrB;
459             j_coord_offsetC  = DIM*jnrC;
460             j_coord_offsetD  = DIM*jnrD;
461
462             /* load j atom coordinates */
463             /* #if GEOMETRY_J == 'Particle'             */
464             gmx_mm256_load_1rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
465                                                  x+j_coord_offsetC,x+j_coord_offsetD,
466                                                  &jx0,&jy0,&jz0);
467             /* #elif GEOMETRY_J == 'Water3'             */
468             gmx_mm256_load_3rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
469                                                  x+j_coord_offsetC,x+j_coord_offsetD,
470                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
471             /* #elif GEOMETRY_J == 'Water4'             */
472             /*     #if 0 in PARTICLES_J                 */
473             gmx_mm256_load_4rvec_4ptr_swizzle_pd(x+j_coord_offsetA,x+j_coord_offsetB,
474                                                  x+j_coord_offsetC,x+j_coord_offsetD,
475                                                  &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
476                                                  &jy2,&jz2,&jx3,&jy3,&jz3);
477             /*     #else                                */
478             gmx_mm256_load_3rvec_4ptr_swizzle_pd(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
479                                                  x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
480                                                  &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
481             /*     #endif                               */
482             /* #endif                                   */
483
484             /* Calculate displacement vector */
485             /* #for I,J in PAIRS_IJ */
486             dx{I}{J}             = _mm256_sub_pd(ix{I},jx{J});
487             dy{I}{J}             = _mm256_sub_pd(iy{I},jy{J});
488             dz{I}{J}             = _mm256_sub_pd(iz{I},jz{J});
489             /*     #define INNERFLOPS INNERFLOPS+3 */
490             /* #endfor */
491
492             /* Calculate squared distance and things based on it */
493             /* #for I,J in PAIRS_IJ */
494             rsq{I}{J}            = gmx_mm256_calc_rsq_pd(dx{I}{J},dy{I}{J},dz{I}{J});
495             /*     #define INNERFLOPS INNERFLOPS+5 */
496             /* #endfor */
497
498             /* #for I,J in PAIRS_IJ */
499             /*     #if 'rinv' in INTERACTION_FLAGS[I][J] */
500             rinv{I}{J}           = avx256_invsqrt_d(rsq{I}{J});
501             /*         #define INNERFLOPS INNERFLOPS+5 */
502             /*     #endif */
503             /* #endfor */
504
505             /* #for I,J in PAIRS_IJ */
506             /*     #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
507             /*         # if 'rinv' not in INTERACTION_FLAGS[I][J] */
508             rinvsq{I}{J}         = avx256_inv_d(rsq{I}{J});
509             /*             #define INNERFLOPS INNERFLOPS+4 */
510             /*         #else */
511             rinvsq{I}{J}         = _mm256_mul_pd(rinv{I}{J},rinv{I}{J});
512             /*             #define INNERFLOPS INNERFLOPS+1 */
513             /*         #endif */
514             /*     #endif */
515             /* #endfor */
516
517             /* #if not 'Water' in GEOMETRY_J */
518             /* Load parameters for j particles */
519             /*     #for J in PARTICLES_ELEC_J */
520             jq{J}              = gmx_mm256_load_4real_swizzle_pd(charge+jnrA+{J},charge+jnrB+{J},
521                                                                  charge+jnrC+{J},charge+jnrD+{J});
522             /*         #if KERNEL_ELEC=='GeneralizedBorn' */
523             isaj{J}            = gmx_mm256_load_4real_swizzle_pd(invsqrta+jnrA+{J},invsqrta+jnrB+{J},
524                                                                  invsqrta+jnrC+{J},invsqrta+jnrD+{J});
525             /*         #endif */
526             /*     #endfor */
527             /*     #for J in PARTICLES_VDW_J */
528             vdwjidx{J}A        = 2*vdwtype[jnrA+{J}];
529             vdwjidx{J}B        = 2*vdwtype[jnrB+{J}];
530             vdwjidx{J}C        = 2*vdwtype[jnrC+{J}];
531             vdwjidx{J}D        = 2*vdwtype[jnrD+{J}];
532             /*     #endfor */
533             /* #endif */
534
535             /* #if 'Force' in KERNEL_VF and not 'Particle' in GEOMETRY_I */
536             /*     #for J in PARTICLES_J */
537             fjx{J}             = _mm256_setzero_pd();
538             fjy{J}             = _mm256_setzero_pd();
539             fjz{J}             = _mm256_setzero_pd();
540             /*     #endfor */
541             /* #endif */
542
543             /* #for I,J in PAIRS_IJ */
544
545             /**************************
546              * CALCULATE INTERACTIONS *
547              **************************/
548
549             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
550             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
551             /*         ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
552             if (gmx_mm256_any_lt(rsq{I}{J},rcutoff2))
553             {
554                 /*     #if 0    ## this and the next two lines is a hack to maintain auto-indentation in template file */
555             }
556             /*         #endif */
557             /*         #define INNERFLOPS INNERFLOPS+1 */
558             /*     #endif */
559
560             /*     #if 'r' in INTERACTION_FLAGS[I][J] */
561             r{I}{J}              = _mm256_mul_pd(rsq{I}{J},rinv{I}{J});
562             /*         #if ROUND == 'Epilogue' */
563             r{I}{J}              = _mm256_andnot_pd(dummy_mask,r{I}{J});
564             /*             #define INNERFLOPS INNERFLOPS+1 */
565             /*         #endif */
566             /*         #define INNERFLOPS INNERFLOPS+1 */
567             /*     #endif */
568
569             /*     ## For water geometries we already loaded parameters at the start of the kernel */
570             /*     #if not 'Water' in GEOMETRY_J */
571             /* Compute parameters for interactions between i and j atoms */
572             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
573             qq{I}{J}             = _mm256_mul_pd(iq{I},jq{J});
574             /*             #define INNERFLOPS INNERFLOPS+1 */
575             /*         #endif */
576             /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
577             gmx_mm256_load_4pair_swizzle_pd(vdwioffsetptr{I}+vdwjidx{J}A,
578                                             vdwioffsetptr{I}+vdwjidx{J}B,
579                                             vdwioffsetptr{I}+vdwjidx{J}C,
580                                             vdwioffsetptr{I}+vdwjidx{J}D,
581                                             &c6_{I}{J},&c12_{I}{J});
582
583             /*           #if 'LJEwald' in KERNEL_VDW */
584             c6grid_{I}{J}       = gmx_mm256_load_4real_swizzle_pd(vdwgridioffsetptr{I}+vdwjidx{J}A,
585                                                                   vdwgridioffsetptr{I}+vdwjidx{J}B,
586                                                                   vdwgridioffsetptr{I}+vdwjidx{J}C,
587                                                                   vdwgridioffsetptr{I}+vdwjidx{J}D);
588             /*           #endif */
589             /*         #endif */
590             /*     #endif */
591
592             /*     #if 'table' in INTERACTION_FLAGS[I][J] */
593             /* Calculate table index by multiplying r with table scale and truncate to integer */
594             rt               = _mm256_mul_pd(r{I}{J},vftabscale);
595             vfitab           = _mm256_cvttpd_epi32(rt);
596             vfeps            = _mm256_sub_pd(rt,_mm256_round_pd(rt, _MM_FROUND_FLOOR));
597             /*         #define INNERFLOPS INNERFLOPS+4                          */
598             /*         #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW     */
599             /*             ## 3 tables, 4 bytes per point: multiply index by 12 */
600             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
601             /*         #elif 'Table' in KERNEL_ELEC                             */
602             /*             ## 1 table, 4 bytes per point: multiply index by 4   */
603             vfitab           = _mm_slli_epi32(vfitab,2);
604             /*         #elif 'Table' in KERNEL_VDW                              */
605             /*             ## 2 tables, 4 bytes per point: multiply index by 8  */
606             vfitab           = _mm_slli_epi32(vfitab,3);
607             /*         #endif                                                   */
608             /*     #endif */
609
610             /*     ## ELECTROSTATIC INTERACTIONS */
611             /*     #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
612
613             /*         #if KERNEL_ELEC=='Coulomb' */
614
615             /* COULOMB ELECTROSTATICS */
616             velec            = _mm256_mul_pd(qq{I}{J},rinv{I}{J});
617             /*             #define INNERFLOPS INNERFLOPS+1 */
618             /*             #if 'Force' in KERNEL_VF */
619             felec            = _mm256_mul_pd(velec,rinvsq{I}{J});
620             /*                 #define INNERFLOPS INNERFLOPS+1 */
621             /*             #endif */
622
623             /*         #elif KERNEL_ELEC=='ReactionField' */
624
625             /* REACTION-FIELD ELECTROSTATICS */
626             /*             #if 'Potential' in KERNEL_VF */
627             velec            = _mm256_mul_pd(qq{I}{J},_mm256_sub_pd(_mm256_add_pd(rinv{I}{J},_mm256_mul_pd(krf,rsq{I}{J})),crf));
628             /*                 #define INNERFLOPS INNERFLOPS+4 */
629             /*             #endif */
630             /*             #if 'Force' in KERNEL_VF */
631             felec            = _mm256_mul_pd(qq{I}{J},_mm256_sub_pd(_mm256_mul_pd(rinv{I}{J},rinvsq{I}{J}),krf2));
632             /*                 #define INNERFLOPS INNERFLOPS+3 */
633             /*             #endif */
634
635             /*         #elif KERNEL_ELEC=='GeneralizedBorn' */
636
637             /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
638             isaprod          = _mm256_mul_pd(isai{I},isaj{J});
639             gbqqfactor       = _mm256_xor_pd(signbit,_mm256_mul_pd(qq{I}{J},_mm256_mul_pd(isaprod,gbinvepsdiff)));
640             gbscale          = _mm256_mul_pd(isaprod,gbtabscale);
641             /*             #define INNERFLOPS INNERFLOPS+5 */
642
643             /* Calculate generalized born table index - this is a separate table from the normal one,
644              * but we use the same procedure by multiplying r with scale and truncating to integer.
645              */
646             rt               = _mm256_mul_pd(r{I}{J},gbscale);
647             gbitab           = _mm256_cvttpd_epi32(rt);
648             gbeps            = _mm256_sub_pd(rt,_mm256_round_pd(rt, _MM_FROUND_FLOOR));
649             gbitab           = _mm_slli_epi32(gbitab,2);
650             Y                = _mm256_load_pd( gbtab + _mm_extract_epi32(gbitab,0) );
651             F                = _mm256_load_pd( gbtab + _mm_extract_epi32(gbitab,1) );
652             G                = _mm256_load_pd( gbtab + _mm_extract_epi32(gbitab,2) );
653             H                = _mm256_load_pd( gbtab + _mm_extract_epi32(gbitab,3) );
654             GMX_MM256_FULLTRANSPOSE4_PD(Y,F,G,H);
655             Heps             = _mm256_mul_pd(gbeps,H);
656             Fp               = _mm256_add_pd(F,_mm256_mul_pd(gbeps,_mm256_add_pd(G,Heps)));
657             VV               = _mm256_add_pd(Y,_mm256_mul_pd(gbeps,Fp));
658             vgb              = _mm256_mul_pd(gbqqfactor,VV);
659             /*             #define INNERFLOPS INNERFLOPS+10 */
660
661             /*             #if 'Force' in KERNEL_VF */
662             FF               = _mm256_add_pd(Fp,_mm256_mul_pd(gbeps,_mm256_add_pd(G,_mm256_add_pd(Heps,Heps))));
663             fgb              = _mm256_mul_pd(gbqqfactor,_mm256_mul_pd(FF,gbscale));
664             dvdatmp          = _mm256_mul_pd(minushalf,_mm256_add_pd(vgb,_mm256_mul_pd(fgb,r{I}{J})));
665             /*                 #if ROUND == 'Epilogue' */
666             dvdatmp          = _mm256_andnot_pd(dummy_mask,dvdatmp);
667             /*                 #endif */
668             dvdasum          = _mm256_add_pd(dvdasum,dvdatmp);
669             /*                 #if ROUND == 'Loop' */
670             fjptrA           = dvda+jnrA;
671             fjptrB           = dvda+jnrB;
672             fjptrC           = dvda+jnrC;
673             fjptrD           = dvda+jnrD;
674             /*                 #else */
675             /* 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. */
676             fjptrA             = (jnrlistA>=0) ? dvda+jnrA : scratch;
677             fjptrB             = (jnrlistB>=0) ? dvda+jnrB : scratch;
678             fjptrC             = (jnrlistC>=0) ? dvda+jnrC : scratch;
679             fjptrD             = (jnrlistD>=0) ? dvda+jnrD : scratch;
680             /*                 #endif */
681             gmx_mm256_increment_4real_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
682                                                  _mm256_mul_pd(dvdatmp,_mm256_mul_pd(isaj{J},isaj{J})));
683             /*                 #define INNERFLOPS INNERFLOPS+12 */
684             /*             #endif */
685             velec            = _mm256_mul_pd(qq{I}{J},rinv{I}{J});
686             /*                 #define INNERFLOPS INNERFLOPS+1 */
687             /*             #if 'Force' in KERNEL_VF */
688             felec            = _mm256_mul_pd(_mm256_sub_pd(_mm256_mul_pd(velec,rinv{I}{J}),fgb),rinv{I}{J});
689             /*                 #define INNERFLOPS INNERFLOPS+3 */
690             /*             #endif */
691
692             /*         #elif KERNEL_ELEC=='Ewald' */
693             /* EWALD ELECTROSTATICS */
694
695             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
696             ewrt             = _mm256_mul_pd(r{I}{J},ewtabscale);
697             ewitab           = _mm256_cvttpd_epi32(ewrt);
698             eweps            = _mm256_sub_pd(ewrt,_mm256_round_pd(ewrt, _MM_FROUND_FLOOR));
699             /*             #define INNERFLOPS INNERFLOPS+4 */
700             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
701             ewitab           = _mm_slli_epi32(ewitab,2);
702             ewtabF           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,0) );
703             ewtabD           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,1) );
704             ewtabV           = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,2) );
705             ewtabFn          = _mm256_load_pd( ewtab + _mm_extract_epi32(ewitab,3) );
706             GMX_MM256_FULLTRANSPOSE4_PD(ewtabF,ewtabD,ewtabV,ewtabFn);
707             felec            = _mm256_add_pd(ewtabF,_mm256_mul_pd(eweps,ewtabD));
708             /*                 #define INNERFLOPS INNERFLOPS+2 */
709             /*                 #if KERNEL_MOD_ELEC=='PotentialShift' */
710             velec            = _mm256_sub_pd(ewtabV,_mm256_mul_pd(_mm256_mul_pd(ewtabhalfspace,eweps),_mm256_add_pd(ewtabF,felec)));
711             velec            = _mm256_mul_pd(qq{I}{J},_mm256_sub_pd(_mm256_sub_pd(rinv{I}{J},sh_ewald),velec));
712             /*                     #define INNERFLOPS INNERFLOPS+7 */
713             /*                 #else */
714             velec            = _mm256_sub_pd(ewtabV,_mm256_mul_pd(_mm256_mul_pd(ewtabhalfspace,eweps),_mm256_add_pd(ewtabF,felec)));
715             velec            = _mm256_mul_pd(qq{I}{J},_mm256_sub_pd(rinv{I}{J},velec));
716             /*                     #define INNERFLOPS INNERFLOPS+6 */
717             /*                 #endif */
718             /*                 #if 'Force' in KERNEL_VF */
719             felec            = _mm256_mul_pd(_mm256_mul_pd(qq{I}{J},rinv{I}{J}),_mm256_sub_pd(rinvsq{I}{J},felec));
720             /*                      #define INNERFLOPS INNERFLOPS+3 */
721             /*                 #endif */
722             /*             #elif KERNEL_VF=='Force' */
723             gmx_mm256_load_4pair_swizzle_pd(ewtab + _mm_extract_epi32(ewitab,0),ewtab + _mm_extract_epi32(ewitab,1),
724                                             ewtab + _mm_extract_epi32(ewitab,2),ewtab + _mm_extract_epi32(ewitab,3),
725                                             &ewtabF,&ewtabFn);
726             felec            = _mm256_add_pd(_mm256_mul_pd( _mm256_sub_pd(one,eweps),ewtabF),_mm256_mul_pd(eweps,ewtabFn));
727             felec            = _mm256_mul_pd(_mm256_mul_pd(qq{I}{J},rinv{I}{J}),_mm256_sub_pd(rinvsq{I}{J},felec));
728             /*                 #define INNERFLOPS INNERFLOPS+7 */
729             /*             #endif */
730
731             /*         #elif KERNEL_ELEC=='CubicSplineTable' */
732
733             /* CUBIC SPLINE TABLE ELECTROSTATICS */
734             Y                = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
735             F                = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
736             G                = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,2) );
737             H                = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,3) );
738             GMX_MM256_FULLTRANSPOSE4_PD(Y,F,G,H);
739             Heps             = _mm256_mul_pd(vfeps,H);
740             Fp               = _mm256_add_pd(F,_mm256_mul_pd(vfeps,_mm256_add_pd(G,Heps)));
741             /*             #define INNERFLOPS INNERFLOPS+4 */
742             /*             #if 'Potential' in KERNEL_VF */
743             VV               = _mm256_add_pd(Y,_mm256_mul_pd(vfeps,Fp));
744             velec            = _mm256_mul_pd(qq{I}{J},VV);
745             /*                 #define INNERFLOPS INNERFLOPS+3 */
746             /*             #endif */
747             /*             #if 'Force' in KERNEL_VF */
748             FF               = _mm256_add_pd(Fp,_mm256_mul_pd(vfeps,_mm256_add_pd(G,_mm256_add_pd(Heps,Heps))));
749             felec            = _mm256_xor_pd(signbit,_mm256_mul_pd(_mm256_mul_pd(qq{I}{J},FF),_mm256_mul_pd(vftabscale,rinv{I}{J})));
750             /*                 #define INNERFLOPS INNERFLOPS+7 */
751             /*             #endif */
752             /*         #endif */
753             /*         ## End of check for electrostatics interaction forms */
754             /*     #endif */
755             /*     ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
756
757             /*     #if 'vdw' in INTERACTION_FLAGS[I][J] */
758
759             /*         #if KERNEL_VDW=='LennardJones' */
760
761             /* LENNARD-JONES DISPERSION/REPULSION */
762
763             rinvsix          = _mm256_mul_pd(_mm256_mul_pd(rinvsq{I}{J},rinvsq{I}{J}),rinvsq{I}{J});
764             /*             #define INNERFLOPS INNERFLOPS+2 */
765             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
766             vvdw6            = _mm256_mul_pd(c6_{I}{J},rinvsix);
767             vvdw12           = _mm256_mul_pd(c12_{I}{J},_mm256_mul_pd(rinvsix,rinvsix));
768             /*                 #define INNERFLOPS INNERFLOPS+3 */
769             /*                 #if KERNEL_MOD_VDW=='PotentialShift' */
770             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) ,
771                                           _mm256_mul_pd( _mm256_sub_pd(vvdw6,_mm256_mul_pd(c6_{I}{J},sh_vdw_invrcut6)),one_sixth));
772             /*                     #define INNERFLOPS INNERFLOPS+8 */
773             /*                 #else */
774             vvdw             = _mm256_sub_pd( _mm256_mul_pd(vvdw12,one_twelfth) , _mm256_mul_pd(vvdw6,one_sixth) );
775             /*                     #define INNERFLOPS INNERFLOPS+3 */
776             /*                 #endif */
777             /*                 ## Check for force inside potential check, i.e. this means we already did the potential part */
778             /*                 #if 'Force' in KERNEL_VF */
779             fvdw             = _mm256_mul_pd(_mm256_sub_pd(vvdw12,vvdw6),rinvsq{I}{J});
780             /*                     #define INNERFLOPS INNERFLOPS+2 */
781             /*                 #endif */
782             /*             #elif KERNEL_VF=='Force' */
783             /*                 ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
784             fvdw             = _mm256_mul_pd(_mm256_sub_pd(_mm256_mul_pd(c12_{I}{J},rinvsix),c6_{I}{J}),_mm256_mul_pd(rinvsix,rinvsq{I}{J}));
785             /*                 #define INNERFLOPS INNERFLOPS+4 */
786             /*             #endif */
787
788             /*         #elif KERNEL_VDW=='CubicSplineTable' */
789
790             /* CUBIC SPLINE TABLE DISPERSION */
791             /*             #if 'Table' in KERNEL_ELEC */
792             vfitab           = _mm_add_epi32(vfitab,ifour);
793             /*             #endif                     */
794             Y                = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
795             F                = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
796             G                = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,2) );
797             H                = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,3) );
798             GMX_MM256_FULLTRANSPOSE4_PD(Y,F,G,H);
799             Heps             = _mm256_mul_pd(vfeps,H);
800             Fp               = _mm256_add_pd(F,_mm256_mul_pd(vfeps,_mm256_add_pd(G,Heps)));
801             /*             #define INNERFLOPS INNERFLOPS+4 */
802             /*             #if 'Potential' in KERNEL_VF */
803             VV               = _mm256_add_pd(Y,_mm256_mul_pd(vfeps,Fp));
804             vvdw6            = _mm256_mul_pd(c6_{I}{J},VV);
805             /*                 #define INNERFLOPS INNERFLOPS+3 */
806             /*             #endif */
807             /*             #if 'Force' in KERNEL_VF */
808             FF               = _mm256_add_pd(Fp,_mm256_mul_pd(vfeps,_mm256_add_pd(G,_mm256_add_pd(Heps,Heps))));
809             fvdw6            = _mm256_mul_pd(c6_{I}{J},FF);
810             /*                 #define INNERFLOPS INNERFLOPS+4 */
811             /*             #endif */
812
813             /* CUBIC SPLINE TABLE REPULSION */
814             vfitab           = _mm_add_epi32(vfitab,ifour);
815             Y                = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,0) );
816             F                = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,1) );
817             G                = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,2) );
818             H                = _mm256_load_pd( vftab + _mm_extract_epi32(vfitab,3) );
819             GMX_MM256_FULLTRANSPOSE4_PD(Y,F,G,H);
820             Heps             = _mm256_mul_pd(vfeps,H);
821             Fp               = _mm256_add_pd(F,_mm256_mul_pd(vfeps,_mm256_add_pd(G,Heps)));
822             /*             #define INNERFLOPS INNERFLOPS+4 */
823             /*             #if 'Potential' in KERNEL_VF */
824             VV               = _mm256_add_pd(Y,_mm256_mul_pd(vfeps,Fp));
825             vvdw12           = _mm256_mul_pd(c12_{I}{J},VV);
826             /*                 #define INNERFLOPS INNERFLOPS+3 */
827             /*             #endif */
828             /*             #if 'Force' in KERNEL_VF */
829             FF               = _mm256_add_pd(Fp,_mm256_mul_pd(vfeps,_mm256_add_pd(G,_mm256_add_pd(Heps,Heps))));
830             fvdw12           = _mm256_mul_pd(c12_{I}{J},FF);
831             /*                 #define INNERFLOPS INNERFLOPS+5 */
832             /*             #endif */
833             /*             #if 'Potential' in KERNEL_VF */
834             vvdw             = _mm256_add_pd(vvdw12,vvdw6);
835             /*                 #define INNERFLOPS INNERFLOPS+1 */
836             /*             #endif */
837             /*             #if 'Force' in KERNEL_VF */
838             fvdw             = _mm256_xor_pd(signbit,_mm256_mul_pd(_mm256_add_pd(fvdw6,fvdw12),_mm256_mul_pd(vftabscale,rinv{I}{J})));
839             /*                 #define INNERFLOPS INNERFLOPS+4 */
840             /*             #endif */
841
842             /*         #elif KERNEL_VDW=='LJEwald' */
843
844             /* Analytical LJ-PME */
845             rinvsix          = _mm256_mul_pd(_mm256_mul_pd(rinvsq{I}{J},rinvsq{I}{J}),rinvsq{I}{J});
846             ewcljrsq         = _mm256_mul_pd(ewclj2,rsq{I}{J});
847             ewclj6           = _mm256_mul_pd(ewclj2,_mm256_mul_pd(ewclj2,ewclj2));
848             exponent         = avx256_exp_d(ewcljrsq);
849             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
850             poly             = _mm256_mul_pd(exponent,_mm256_add_pd(_mm256_sub_pd(one,ewcljrsq),_mm256_mul_pd(_mm256_mul_pd(ewcljrsq,ewcljrsq),one_half)));
851             /*                 #define INNERFLOPS INNERFLOPS+11 */
852             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
853             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
854             vvdw6            = _mm256_mul_pd(_mm256_sub_pd(c6_{I}{J},_mm256_mul_pd(c6grid_{I}{J},_mm256_sub_pd(one,poly))),rinvsix);
855             vvdw12           = _mm256_mul_pd(c12_{I}{J},_mm256_mul_pd(rinvsix,rinvsix));
856             /*                     #define INNERFLOPS INNERFLOPS+6 */
857             /*                 #if KERNEL_MOD_VDW=='PotentialShift' */
858             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) ,
859                                           _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));
860             /*                     #define INNERFLOPS INNERFLOPS+10 */
861             /*                 #else */
862             vvdw             = _mm256_sub_pd(_mm256_mul_pd(vvdw12,one_twelfth),_mm256_mul_pd(vvdw6,one_sixth));
863             /*                 #define INNERFLOPS INNERFLOPS+6 */
864             /*                 #endif */
865             /*                  ## Check for force inside potential check, i.e. this means we already did the potential part */
866             /*                  #if 'Force' in KERNEL_VF */
867             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
868             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});
869             /*                 #define INNERFLOPS INNERFLOPS+6 */
870             /*                  #endif */
871             /*              #elif KERNEL_VF=='Force' */
872             /* f6A = 6 * C6grid * (1 - poly) */
873             f6A              = _mm256_mul_pd(c6grid_{I}{J},_mm256_sub_pd(one,poly));
874             /* f6B = C6grid * exponent * beta^6 */
875             f6B              = _mm256_mul_pd(_mm256_mul_pd(c6grid_{I}{J},one_sixth),_mm256_mul_pd(exponent,ewclj6));
876             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
877             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});
878             /*                 #define INNERFLOPS INNERFLOPS+11 */
879             /*              #endif */
880             /*         #endif */
881             /*         ## End of check for vdw interaction forms */
882             /*     #endif */
883             /*     ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
884
885             /*     #if 'switch' in INTERACTION_FLAGS[I][J] */
886             d                = _mm256_sub_pd(r{I}{J},rswitch);
887             d                = _mm256_max_pd(d,_mm256_setzero_pd());
888             d2               = _mm256_mul_pd(d,d);
889             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)))))));
890             /*         #define INNERFLOPS INNERFLOPS+10 */
891
892             /*         #if 'Force' in KERNEL_VF */
893             dsw              = _mm256_mul_pd(d2,_mm256_add_pd(swF2,_mm256_mul_pd(d,_mm256_add_pd(swF3,_mm256_mul_pd(d,swF4)))));
894             /*             #define INNERFLOPS INNERFLOPS+5 */
895             /*         #endif */
896
897             /* Evaluate switch function */
898             /*         #if 'Force' in KERNEL_VF */
899             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
900             /*             #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
901             felec            = _mm256_sub_pd( _mm256_mul_pd(felec,sw) , _mm256_mul_pd(rinv{I}{J},_mm256_mul_pd(velec,dsw)) );
902             /*                 #define INNERFLOPS INNERFLOPS+4 */
903             /*             #endif */
904             /*             #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
905             fvdw             = _mm256_sub_pd( _mm256_mul_pd(fvdw,sw) , _mm256_mul_pd(rinv{I}{J},_mm256_mul_pd(vvdw,dsw)) );
906             /*                 #define INNERFLOPS INNERFLOPS+4 */
907             /*             #endif */
908             /*         #endif */
909             /*         #if 'Potential' in KERNEL_VF */
910             /*             #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
911             velec            = _mm256_mul_pd(velec,sw);
912             /*                 #define INNERFLOPS INNERFLOPS+1 */
913             /*             #endif */
914             /*             #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
915             vvdw             = _mm256_mul_pd(vvdw,sw);
916             /*                 #define INNERFLOPS INNERFLOPS+1 */
917             /*             #endif */
918             /*         #endif */
919             /*     #endif */
920             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
921             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
922             cutoff_mask      = _mm256_cmp_pd(rsq{I}{J},rcutoff2,_CMP_LT_OQ);
923             /*         #define INNERFLOPS INNERFLOPS+1 */
924             /*     #endif */
925
926             /*     #if 'Potential' in KERNEL_VF */
927             /* Update potential sum for this i atom from the interaction with this j atom. */
928             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
929             /*             #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
930             velec            = _mm256_and_pd(velec,cutoff_mask);
931             /*                 #define INNERFLOPS INNERFLOPS+1 */
932             /*             #endif                                       */
933             /*             #if ROUND == 'Epilogue' */
934             velec            = _mm256_andnot_pd(dummy_mask,velec);
935             /*             #endif */
936             velecsum         = _mm256_add_pd(velecsum,velec);
937             /*             #define INNERFLOPS INNERFLOPS+1 */
938             /*             #if KERNEL_ELEC=='GeneralizedBorn' */
939             /*             #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
940             vgb              = _mm256_and_pd(vgb,cutoff_mask);
941             /*                 #define INNERFLOPS INNERFLOPS+1 */
942             /*             #endif                                       */
943             /*             #if ROUND == 'Epilogue' */
944             vgb              = _mm256_andnot_pd(dummy_mask,vgb);
945             /*             #endif */
946             vgbsum           = _mm256_add_pd(vgbsum,vgb);
947             /*                 #define INNERFLOPS INNERFLOPS+1 */
948             /*             #endif */
949             /*         #endif */
950             /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
951             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
952             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
953             vvdw             = _mm256_and_pd(vvdw,cutoff_mask);
954             /*                 #define INNERFLOPS INNERFLOPS+1 */
955             /*             #endif                                       */
956             /*             #if ROUND == 'Epilogue' */
957             vvdw             = _mm256_andnot_pd(dummy_mask,vvdw);
958             /*             #endif */
959             vvdwsum          = _mm256_add_pd(vvdwsum,vvdw);
960             /*             #define INNERFLOPS INNERFLOPS+1 */
961             /*         #endif */
962             /*     #endif */
963
964             /*     #if 'Force' in KERNEL_VF */
965
966             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
967             fscal            = _mm256_add_pd(felec,fvdw);
968             /*             #define INNERFLOPS INNERFLOPS+1 */
969             /*         #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
970             fscal            = felec;
971             /*         #elif 'vdw' in INTERACTION_FLAGS[I][J] */
972             fscal            = fvdw;
973             /*        #endif */
974
975             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
976             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
977             fscal            = _mm256_and_pd(fscal,cutoff_mask);
978             /*                 #define INNERFLOPS INNERFLOPS+1 */
979             /*             #endif                                       */
980
981             /*             #if ROUND == 'Epilogue' */
982             fscal            = _mm256_andnot_pd(dummy_mask,fscal);
983             /*             #endif */
984
985             /* Calculate temporary vectorial force */
986             tx               = _mm256_mul_pd(fscal,dx{I}{J});
987             ty               = _mm256_mul_pd(fscal,dy{I}{J});
988             tz               = _mm256_mul_pd(fscal,dz{I}{J});
989
990             /* Update vectorial force */
991             fix{I}             = _mm256_add_pd(fix{I},tx);
992             fiy{I}             = _mm256_add_pd(fiy{I},ty);
993             fiz{I}             = _mm256_add_pd(fiz{I},tz);
994             /*             #define INNERFLOPS INNERFLOPS+6 */
995
996             /* #if GEOMETRY_I == 'Particle'             */
997             /*     #if ROUND == 'Loop' */
998             fjptrA             = f+j_coord_offsetA;
999             fjptrB             = f+j_coord_offsetB;
1000             fjptrC             = f+j_coord_offsetC;
1001             fjptrD             = f+j_coord_offsetD;
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             /*     #endif */
1008             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,tx,ty,tz);
1009             /*     #define INNERFLOPS INNERFLOPS+3      */
1010             /* #else                                    */
1011             fjx{J}             = _mm256_add_pd(fjx{J},tx);
1012             fjy{J}             = _mm256_add_pd(fjy{J},ty);
1013             fjz{J}             = _mm256_add_pd(fjz{J},tz);
1014             /*     #define INNERFLOPS INNERFLOPS+3      */
1015             /* #endif                                   */
1016
1017             /*     #endif */
1018
1019             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
1020             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
1021             /*         #if 0    ## This and next two lines is a hack to maintain indentation in template file */
1022             {
1023                 /*     #endif */
1024             }
1025             /*     #endif */
1026             /*    ## End of check for the interaction being outside the cutoff */
1027
1028             /* #endfor */
1029             /* ## End of loop over i-j interaction pairs */
1030
1031             /* #if GEOMETRY_I != 'Particle' */
1032             /*     #if ROUND == 'Loop' */
1033             fjptrA             = f+j_coord_offsetA;
1034             fjptrB             = f+j_coord_offsetB;
1035             fjptrC             = f+j_coord_offsetC;
1036             fjptrD             = f+j_coord_offsetD;
1037             /*     #else */
1038             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
1039             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
1040             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
1041             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
1042             /*     #endif */
1043             /* #endif */
1044
1045             /* #if 'Water' in GEOMETRY_I and GEOMETRY_J == 'Particle' */
1046             gmx_mm256_decrement_1rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
1047             /*     #define INNERFLOPS INNERFLOPS+3      */
1048             /* #elif GEOMETRY_J == 'Water3'             */
1049             gmx_mm256_decrement_3rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
1050                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
1051             /*     #define INNERFLOPS INNERFLOPS+9      */
1052             /* #elif GEOMETRY_J == 'Water4'             */
1053             /*     #if 0 in PARTICLES_J                 */
1054             gmx_mm256_decrement_4rvec_4ptr_swizzle_pd(fjptrA,fjptrB,fjptrC,fjptrD,
1055                                                       fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
1056                                                       fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1057             /*     #define INNERFLOPS INNERFLOPS+12     */
1058             /*     #else                                */
1059             gmx_mm256_decrement_3rvec_4ptr_swizzle_pd(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
1060                                                       fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
1061             /*     #define INNERFLOPS INNERFLOPS+9      */
1062             /*     #endif                               */
1063             /* #endif                                   */
1064
1065             /* Inner loop uses {INNERFLOPS} flops */
1066         }
1067
1068         /* #endfor */
1069
1070         /* End of innermost loop */
1071
1072         /* #if 'Force' in KERNEL_VF */
1073         /*     #if GEOMETRY_I == 'Particle'            */
1074         gmx_mm256_update_iforce_1atom_swizzle_pd(fix0,fiy0,fiz0,
1075                                                  f+i_coord_offset,fshift+i_shift_offset);
1076         /*         #define OUTERFLOPS OUTERFLOPS+6     */
1077         /*     #elif GEOMETRY_I == 'Water3'            */
1078         gmx_mm256_update_iforce_3atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
1079                                                  f+i_coord_offset,fshift+i_shift_offset);
1080         /*         #define OUTERFLOPS OUTERFLOPS+18    */
1081         /*     #elif GEOMETRY_I == 'Water4'            */
1082         /*         #if 0 in PARTICLES_I                */
1083         gmx_mm256_update_iforce_4atom_swizzle_pd(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1084                                                  f+i_coord_offset,fshift+i_shift_offset);
1085         /*             #define OUTERFLOPS OUTERFLOPS+24    */
1086         /*         #else                               */
1087         gmx_mm256_update_iforce_3atom_swizzle_pd(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
1088                                                  f+i_coord_offset+DIM,fshift+i_shift_offset);
1089         /*             #define OUTERFLOPS OUTERFLOPS+18    */
1090         /*         #endif                              */
1091         /*     #endif                                  */
1092         /* #endif                                      */
1093
1094         /* #if 'Potential' in KERNEL_VF */
1095         ggid                        = gid[iidx];
1096         /* Update potential energies */
1097         /*     #if KERNEL_ELEC != 'None' */
1098         gmx_mm256_update_1pot_pd(velecsum,kernel_data->energygrp_elec+ggid);
1099         /*         #define OUTERFLOPS OUTERFLOPS+1 */
1100         /*     #endif */
1101         /*     #if 'GeneralizedBorn' in KERNEL_ELEC */
1102         gmx_mm256_update_1pot_pd(vgbsum,kernel_data->energygrp_polarization+ggid);
1103         /*         #define OUTERFLOPS OUTERFLOPS+1 */
1104         /*     #endif */
1105         /*     #if KERNEL_VDW != 'None' */
1106         gmx_mm256_update_1pot_pd(vvdwsum,kernel_data->energygrp_vdw+ggid);
1107         /*         #define OUTERFLOPS OUTERFLOPS+1 */
1108         /*     #endif */
1109         /* #endif */
1110         /*     #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
1111         dvdasum = _mm256_mul_pd(dvdasum, _mm256_mul_pd(isai{I},isai{I}));
1112         gmx_mm256_update_1pot_pd(dvdasum,dvda+inr);
1113         /*     #endif */
1114
1115         /* Increment number of inner iterations */
1116         inneriter                  += j_index_end - j_index_start;
1117
1118         /* Outer loop uses {OUTERFLOPS} flops */
1119     }
1120
1121     /* Increment number of outer iterations */
1122     outeriter        += nri;
1123
1124     /* Update outer/inner flops */
1125     /* ## NB: This is not important, it just affects the flopcount. However, since our preprocessor is */
1126     /* ## primitive and replaces aggressively even in strings inside these directives, we need to      */
1127     /* ## assemble the main part of the name (containing KERNEL/ELEC/VDW) directly in the source.      */
1128     /* #if GEOMETRY_I == 'Water3'            */
1129     /*     #define ISUFFIX '_W3'             */
1130     /* #elif GEOMETRY_I == 'Water4'          */
1131     /*     #define ISUFFIX '_W4'             */
1132     /* #else                                 */
1133     /*     #define ISUFFIX ''                */
1134     /* #endif                                */
1135     /* #if GEOMETRY_J == 'Water3'            */
1136     /*     #define JSUFFIX 'W3'              */
1137     /* #elif GEOMETRY_J == 'Water4'          */
1138     /*     #define JSUFFIX 'W4'              */
1139     /* #else                                 */
1140     /*     #define JSUFFIX ''                */
1141     /* #endif                                */
1142     /* #if 'PotentialAndForce' in KERNEL_VF  */
1143     /*     #define VFSUFFIX  '_VF'           */
1144     /* #elif 'Potential' in KERNEL_VF        */
1145     /*     #define VFSUFFIX '_V'             */
1146     /* #else                                 */
1147     /*     #define VFSUFFIX '_F'             */
1148     /* #endif                                */
1149
1150     /* #if KERNEL_ELEC != 'None' and KERNEL_VDW != 'None' */
1151     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1152     /* #elif KERNEL_ELEC != 'None' */
1153     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1154     /* #else */
1155     inc_nrnb(nrnb,eNR_NBKERNEL_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1156     /* #endif  */
1157 }