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