Valgrind suppression for OS X 10.9
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_c / nb_kernel_template_c.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 #ifdef HAVE_CONFIG_H
40 #include <config.h>
41 #endif
42
43 #include <math.h>
44
45 #include "../nb_kernel.h"
46 #include "types/simple.h"
47 #include "vec.h"
48 #include "nrnb.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 /* #if ('CubicSplineTable' in [KERNEL_ELEC,KERNEL_VDW]) or KERNEL_VF=='PotentialAndForce' */
85 /*     #define TABLE_POINT_SIZE 4 */
86 /* #else */
87 /*     #define TABLE_POINT_SIZE 2 */
88 /* #endif */
89
90 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
91 /*     #define TABLE_NINTERACTIONS 3 */
92 /*     #define TABLE_VDW_OFFSET TABLE_POINT_SIZE */
93 /* #elif 'Table' in KERNEL_ELEC */
94 /*     #define TABLE_NINTERACTIONS 1 */
95 /* #elif 'Table' in KERNEL_VDW */
96 /*     #define TABLE_NINTERACTIONS 2 */
97 /*     #define TABLE_VDW_OFFSET 0 */
98 /* #else */
99 /*     #define TABLE_NINTERACTIONS 0 */
100 /* #endif */
101
102 /* #if 'Buckingham' in KERNEL_VDW */
103 /*   #define NVDWPARAM 3 */
104 /* #else */
105 /*   #define NVDWPARAM 2 */
106 /* #endif */
107
108 /*
109  * Gromacs nonbonded kernel:   {KERNEL_NAME}
110  * Electrostatics interaction: {KERNEL_ELEC}
111  * VdW interaction:            {KERNEL_VDW}
112  * Geometry:                   {GEOMETRY_I}-{GEOMETRY_J}
113  * Calculate force/pot:        {KERNEL_VF}
114  */
115 void
116 {KERNEL_NAME}
117                     (t_nblist                    * gmx_restrict       nlist,
118                      rvec                        * gmx_restrict          xx,
119                      rvec                        * gmx_restrict          ff,
120                      t_forcerec                  * gmx_restrict          fr,
121                      t_mdatoms                   * gmx_restrict     mdatoms,
122                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
123                      t_nrnb                      * gmx_restrict        nrnb)
124 {
125     /* ## Not all variables are used for all kernels, but any optimizing compiler fixes that, */
126     /* ## so there is no point in going to extremes to exclude variables that are not needed. */
127     int              i_shift_offset,i_coord_offset,j_coord_offset;
128     int              j_index_start,j_index_end;
129     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
130     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
131     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
132     real             *shiftvec,*fshift,*x,*f;
133     /* #for I in PARTICLES_I */
134     int              vdwioffset{I};
135     real             ix{I},iy{I},iz{I},fix{I},fiy{I},fiz{I},iq{I},isai{I};
136     /* #endfor */
137     /* #for J in PARTICLES_J */
138     int              vdwjidx{J};
139     real             jx{J},jy{J},jz{J},fjx{J},fjy{J},fjz{J},jq{J},isaj{J};
140     /* #endfor */
141     /* #for I,J in PAIRS_IJ */
142     real             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},cexp1_{I}{J},cexp2_{I}{J};
143     /* #endfor */
144     /* #if KERNEL_ELEC != 'None' */
145     real             velec,felec,velecsum,facel,crf,krf,krf2;
146     real             *charge;
147     /* #endif */
148     /* #if 'GeneralizedBorn' in KERNEL_ELEC */
149     int              gbitab;
150     real             vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,dvdatmp;
151     real             *invsqrta,*dvda,*gbtab;
152     /* #endif */
153     /* #if KERNEL_VDW != 'None' */
154     int              nvdwtype;
155     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
156     int              *vdwtype;
157     real             *vdwparam;
158     /* #endif */
159     /* #if 'Table' in KERNEL_ELEC or 'GeneralizedBorn' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
160     int              vfitab;
161     real             rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
162     real             *vftab;
163     /* #endif */
164     /* #if 'LJEwald' in KERNEL_VDW */
165     /* #for I,J in PAIRS_IJ */
166     real             c6grid_{I}{J};
167     /* #endfor */
168     real             ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,sh_lj_ewald;
169     real             *vdwgridparam;
170     /* #endif */
171     /* #if 'Ewald' in KERNEL_ELEC */
172     int              ewitab;
173     real             ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
174     real             *ewtab;
175     /* #endif */
176     /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
177     real             rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
178     /* #endif */
179
180     x                = xx[0];
181     f                = ff[0];
182
183     nri              = nlist->nri;
184     iinr             = nlist->iinr;
185     jindex           = nlist->jindex;
186     jjnr             = nlist->jjnr;
187     shiftidx         = nlist->shift;
188     gid              = nlist->gid;
189     shiftvec         = fr->shift_vec[0];
190     fshift           = fr->fshift[0];
191     /* #if KERNEL_ELEC != 'None' */
192     facel            = fr->epsfac;
193     charge           = mdatoms->chargeA;
194     /*     #if 'ReactionField' in KERNEL_ELEC */
195     krf              = fr->ic->k_rf;
196     krf2             = krf*2.0;
197     crf              = fr->ic->c_rf;
198     /*     #endif */
199     /* #endif */
200     /* #if KERNEL_VDW != 'None' */
201     nvdwtype         = fr->ntype;
202     vdwparam         = fr->nbfp;
203     vdwtype          = mdatoms->typeA;
204     /* #endif */
205     /* #if 'LJEwald' in KERNEL_VDW */
206     vdwgridparam     = fr->ljpme_c6grid;
207     ewclj            = fr->ewaldcoeff_lj;
208     sh_lj_ewald      = fr->ic->sh_lj_ewald;
209     ewclj2           = ewclj*ewclj;
210     ewclj6           = ewclj2*ewclj2*ewclj2;
211     /* #endif */
212
213     /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
214     vftab            = kernel_data->table_elec_vdw->data;
215     vftabscale       = kernel_data->table_elec_vdw->scale;
216     /* #elif 'Table' in KERNEL_ELEC */
217     vftab            = kernel_data->table_elec->data;
218     vftabscale       = kernel_data->table_elec->scale;
219     /* #elif 'Table' in KERNEL_VDW */
220     vftab            = kernel_data->table_vdw->data;
221     vftabscale       = kernel_data->table_vdw->scale;
222     /* #endif */
223
224     /* #if 'Ewald' in KERNEL_ELEC */
225     sh_ewald         = fr->ic->sh_ewald;
226     /*     #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
227     ewtab            = fr->ic->tabq_coul_F;
228     ewtabscale       = fr->ic->tabq_scale;
229     ewtabhalfspace   = 0.5/ewtabscale;
230     /*     #else */
231     ewtab            = fr->ic->tabq_coul_FDV0;
232     ewtabscale       = fr->ic->tabq_scale;
233     ewtabhalfspace   = 0.5/ewtabscale;
234      /*     #endif */
235     /* #endif */
236
237     /* #if KERNEL_ELEC=='GeneralizedBorn' */
238     invsqrta         = fr->invsqrta;
239     dvda             = fr->dvda;
240     gbtabscale       = fr->gbtab.scale;
241     gbtab            = fr->gbtab.data;
242     gbinvepsdiff     = (1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent);
243     /* #endif */
244
245     /* #if 'Water' in GEOMETRY_I */
246     /* Setup water-specific parameters */
247     inr              = nlist->iinr[0];
248     /*     #for I in PARTICLES_ELEC_I */
249     iq{I}              = facel*charge[inr+{I}];
250     /*     #endfor */
251     /*     #for I in PARTICLES_VDW_I */
252     vdwioffset{I}      = {NVDWPARAM}*nvdwtype*vdwtype[inr+{I}];
253     /*     #endfor */
254     /* #endif */
255
256     /* #if 'Water' in GEOMETRY_J */
257     /*     #for J in PARTICLES_ELEC_J */
258     jq{J}              = charge[inr+{J}];
259     /*     #endfor */
260     /*     #for J in PARTICLES_VDW_J */
261     vdwjidx{J}         = {NVDWPARAM}*vdwtype[inr+{J}];
262     /*     #endfor */
263     /*     #for I,J in PAIRS_IJ */
264     /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
265     qq{I}{J}             = iq{I}*jq{J};
266     /*         #endif */
267     /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
268     /*             #if 'Buckingham' in KERNEL_VDW */
269     c6_{I}{J}            = vdwparam[vdwioffset{I}+vdwjidx{J}];
270     cexp1_{I}{J}         = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
271     cexp2_{I}{J}         = vdwparam[vdwioffset{I}+vdwjidx{J}+2];
272     /*             #elif 'LJEwald' in KERNEL_VDW */
273     c6_{I}{J}            = vdwparam[vdwioffset{I}+vdwjidx{J}];
274     c12_{I}{J}           = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
275     c6grid_{I}{J}        = vdwgridparam[vdwioffset{I}+vdwjidx{J}];
276     /*             #else */
277     c6_{I}{J}            = vdwparam[vdwioffset{I}+vdwjidx{J}];
278     c12_{I}{J}           = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
279     /*             #endif */
280     /*         #endif */
281     /*     #endfor */
282     /* #endif */
283
284     /* #if KERNEL_MOD_ELEC!='None' or KERNEL_MOD_VDW!='None' */
285     /*     #if KERNEL_ELEC!='None' */
286     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
287     rcutoff          = fr->rcoulomb;
288     /*     #else */
289     rcutoff          = fr->rvdw;
290     /*     #endif */
291     rcutoff2         = rcutoff*rcutoff;
292     /* #endif */
293
294     /* #if KERNEL_MOD_VDW=='PotentialShift' */
295     sh_vdw_invrcut6  = fr->ic->sh_invrc6;
296     rvdw             = fr->rvdw;
297     /* #endif */
298
299     /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
300     /*     #if KERNEL_MOD_ELEC=='PotentialSwitch'  */
301     rswitch          = fr->rcoulomb_switch;
302     /*     #else */
303     rswitch          = fr->rvdw_switch;
304     /*     #endif */
305     /* Setup switch parameters */
306     d                = rcutoff-rswitch;
307     swV3             = -10.0/(d*d*d);
308     swV4             =  15.0/(d*d*d*d);
309     swV5             =  -6.0/(d*d*d*d*d);
310     /*     #if 'Force' in KERNEL_VF */
311     swF2             = -30.0/(d*d*d);
312     swF3             =  60.0/(d*d*d*d);
313     swF4             = -30.0/(d*d*d*d*d);
314     /*     #endif */
315     /* #endif */
316
317     /* ## Keep track of the floating point operations we issue for reporting! */
318     /* #define OUTERFLOPS 0 */
319     /* #define INNERFLOPS 0 */
320     outeriter        = 0;
321     inneriter        = 0;
322
323     /* Start outer loop over neighborlists */
324     for(iidx=0; iidx<nri; iidx++)
325     {
326         /* Load shift vector for this list */
327         i_shift_offset   = DIM*shiftidx[iidx];
328         shX              = shiftvec[i_shift_offset+XX];
329         shY              = shiftvec[i_shift_offset+YY];
330         shZ              = shiftvec[i_shift_offset+ZZ];
331
332         /* Load limits for loop over neighbors */
333         j_index_start    = jindex[iidx];
334         j_index_end      = jindex[iidx+1];
335
336         /* Get outer coordinate index */
337         inr              = iinr[iidx];
338         i_coord_offset   = DIM*inr;
339
340         /* Load i particle coords and add shift vector */
341         /* ## Loop over i particles, but only include ones that we use - skip e.g. vdw-only sites for elec-only kernel */
342         /*     #for I in PARTICLES_I */
343         ix{I}              = shX + x[i_coord_offset+DIM*{I}+XX];
344         iy{I}              = shY + x[i_coord_offset+DIM*{I}+YY];
345         iz{I}              = shZ + x[i_coord_offset+DIM*{I}+ZZ];
346         /*     #define OUTERFLOPS OUTERFLOPS+3 */
347         /* #endfor */
348
349         /* #if 'Force' in KERNEL_VF */
350         /*     #for I in PARTICLES_I */
351         fix{I}             = 0.0;
352         fiy{I}             = 0.0;
353         fiz{I}             = 0.0;
354         /*     #endfor */
355         /* #endif */
356
357         /* ## For water we already preloaded parameters at the start of the kernel */
358         /* #if not 'Water' in GEOMETRY_I */
359         /* Load parameters for i particles */
360         /*     #for I in PARTICLES_ELEC_I */
361         iq{I}              = facel*charge[inr+{I}];
362         /*         #define OUTERFLOPS OUTERFLOPS+1 */
363         /*         #if KERNEL_ELEC=='GeneralizedBorn' */
364         isai{I}            = invsqrta[inr+{I}];
365         /*         #endif */
366         /*     #endfor */
367         /*     #for I in PARTICLES_VDW_I */
368         vdwioffset{I}      = {NVDWPARAM}*nvdwtype*vdwtype[inr+{I}];
369         /*     #endfor */
370         /* #endif */
371
372         /* #if 'Potential' in KERNEL_VF */
373         /* Reset potential sums */
374         /*     #if KERNEL_ELEC != 'None' */
375         velecsum         = 0.0;
376         /*     #endif */
377         /*     #if 'GeneralizedBorn' in KERNEL_ELEC */
378         vgbsum           = 0.0;
379         /*     #endif */
380         /*     #if KERNEL_VDW != 'None' */
381         vvdwsum          = 0.0;
382         /*     #endif */
383         /* #endif */
384         /*     #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
385         dvdasum          = 0.0;
386         /*     #endif */
387
388         /* Start inner kernel loop */
389         for(jidx=j_index_start; jidx<j_index_end; jidx++)
390         {
391             /* Get j neighbor index, and coordinate index */
392             jnr              = jjnr[jidx];
393             j_coord_offset   = DIM*jnr;
394
395             /* load j atom coordinates */
396             /* #for J in PARTICLES_J */
397             jx{J}              = x[j_coord_offset+DIM*{J}+XX];
398             jy{J}              = x[j_coord_offset+DIM*{J}+YY];
399             jz{J}              = x[j_coord_offset+DIM*{J}+ZZ];
400             /* #endfor */
401
402             /* Calculate displacement vector */
403             /* #for I,J in PAIRS_IJ */
404             dx{I}{J}             = ix{I} - jx{J};
405             dy{I}{J}             = iy{I} - jy{J};
406             dz{I}{J}             = iz{I} - jz{J};
407             /*     #define INNERFLOPS INNERFLOPS+3 */
408             /* #endfor */
409
410             /* Calculate squared distance and things based on it */
411             /* #for I,J in PAIRS_IJ */
412             rsq{I}{J}            = dx{I}{J}*dx{I}{J}+dy{I}{J}*dy{I}{J}+dz{I}{J}*dz{I}{J};
413             /*     #define INNERFLOPS INNERFLOPS+5 */
414             /* #endfor */
415
416             /* #for I,J in PAIRS_IJ */
417             /*     #if 'rinv' in INTERACTION_FLAGS[I][J] */
418             rinv{I}{J}           = gmx_invsqrt(rsq{I}{J});
419             /*         #define INNERFLOPS INNERFLOPS+5 */
420             /*     #endif */
421             /* #endfor */
422
423             /* #for I,J in PAIRS_IJ */
424             /*     #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
425             /*         # if 'rinv' not in INTERACTION_FLAGS[I][J] */
426             rinvsq{I}{J}         = 1.0/rsq{I}{J};
427             /*             #define INNERFLOPS INNERFLOPS+4 */
428             /*         #else */
429             rinvsq{I}{J}         = rinv{I}{J}*rinv{I}{J};
430             /*             #define INNERFLOPS INNERFLOPS+1 */
431             /*         #endif */
432             /*     #endif */
433             /* #endfor */
434
435             /* #if not 'Water' in GEOMETRY_J */
436             /* Load parameters for j particles */
437             /*     #for J in PARTICLES_ELEC_J */
438             jq{J}              = charge[jnr+{J}];
439             /*         #if KERNEL_ELEC=='GeneralizedBorn' */
440             isaj{J}           = invsqrta[jnr+{J}];
441             /*         #endif */
442             /*     #endfor */
443             /*     #for J in PARTICLES_VDW_J */
444             vdwjidx{J}         = {NVDWPARAM}*vdwtype[jnr+{J}];
445             /*     #endfor */
446             /* #endif */
447
448             /* #for I,J in PAIRS_IJ */
449
450             /**************************
451              * CALCULATE INTERACTIONS *
452              **************************/
453
454             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
455             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
456             /*         ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
457             if (rsq{I}{J}<rcutoff2)
458             {
459                 /*     #if 0    ## this and the next two lines is a hack to maintain auto-indentation in template file */
460             }
461             /*         #endif */
462             /*     #endif */
463
464             /*     #if 'r' in INTERACTION_FLAGS[I][J] */
465             r{I}{J}              = rsq{I}{J}*rinv{I}{J};
466             /*         #define INNERFLOPS INNERFLOPS+1 */
467             /*     #endif */
468
469             /*     ## For water geometries we already loaded parameters at the start of the kernel */
470             /*     #if not 'Water' in GEOMETRY_J */
471             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
472             qq{I}{J}             = iq{I}*jq{J};
473             /*             #define INNERFLOPS INNERFLOPS+1 */
474             /*         #endif */
475             /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
476             /*             #if KERNEL_VDW=='Buckingham' */
477             c6_{I}{J}            = vdwparam[vdwioffset{I}+vdwjidx{J}];
478             cexp1_{I}{J}         = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
479             cexp2_{I}{J}         = vdwparam[vdwioffset{I}+vdwjidx{J}+2];
480             /*             #elif 'LJEwald' in KERNEL_VDW */
481             c6_{I}{J}            = vdwparam[vdwioffset{I}+vdwjidx{J}];
482             c12_{I}{J}           = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
483             c6grid_{I}{J}        = vdwgridparam[vdwioffset{I}+vdwjidx{J}];
484             /*             #else */
485             c6_{I}{J}            = vdwparam[vdwioffset{I}+vdwjidx{J}];
486             c12_{I}{J}           = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
487             /*             #endif */
488             /*         #endif */
489             /*     #endif */
490
491             /*     #if 'table' in INTERACTION_FLAGS[I][J] */
492             /* Calculate table index by multiplying r with table scale and truncate to integer */
493             rt               = r{I}{J}*vftabscale;
494             vfitab           = rt;
495             vfeps            = rt-vfitab;
496             vfitab           = {TABLE_NINTERACTIONS}*{TABLE_POINT_SIZE}*vfitab;
497             /*         #define INNERFLOPS INNERFLOPS+2 */
498             /*     #endif */
499
500             /*     ## ELECTROSTATIC INTERACTIONS */
501             /*     #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
502
503             /*         #if KERNEL_ELEC=='Coulomb' */
504
505             /* COULOMB ELECTROSTATICS */
506             velec            = qq{I}{J}*rinv{I}{J};
507             /*             #define INNERFLOPS INNERFLOPS+1 */
508             /*             #if 'Force' in KERNEL_VF */
509             felec            = velec*rinvsq{I}{J};
510             /*                 #define INNERFLOPS INNERFLOPS+2 */
511             /*             #endif */
512
513             /*         #elif KERNEL_ELEC=='ReactionField' */
514
515             /* REACTION-FIELD ELECTROSTATICS */
516             /*             #if 'Potential' in KERNEL_VF */
517             velec            = qq{I}{J}*(rinv{I}{J}+krf*rsq{I}{J}-crf);
518             /*                 #define INNERFLOPS INNERFLOPS+4 */
519             /*             #endif */
520             /*             #if 'Force' in KERNEL_VF */
521             felec            = qq{I}{J}*(rinv{I}{J}*rinvsq{I}{J}-krf2);
522             /*                 #define INNERFLOPS INNERFLOPS+3 */
523             /*             #endif */
524
525             /*         #elif KERNEL_ELEC=='GeneralizedBorn' */
526
527             /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
528             isaprod          = isai{I}*isaj{J};
529             gbqqfactor       = isaprod*(-qq{I}{J})*gbinvepsdiff;
530             gbscale          = isaprod*gbtabscale;
531             dvdaj            = dvda[jnr+{J}];
532             /*             #define INNERFLOPS INNERFLOPS+5 */
533
534             /* Calculate generalized born table index - this is a separate table from the normal one,
535              * but we use the same procedure by multiplying r with scale and truncating to integer.
536              */
537             rt               = r{I}{J}*gbscale;
538             gbitab           = rt;
539             gbeps            = rt-gbitab;
540             gbitab           = 4*gbitab;
541
542             Y                = gbtab[gbitab];
543             F                = gbtab[gbitab+1];
544             Geps             = gbeps*gbtab[gbitab+2];
545             Heps2            = gbeps*gbeps*gbtab[gbitab+3];
546             Fp               = F+Geps+Heps2;
547             VV               = Y+gbeps*Fp;
548             vgb              = gbqqfactor*VV;
549             /*             #define INNERFLOPS INNERFLOPS+10 */
550
551             /*             #if 'Force' in KERNEL_VF */
552             FF               = Fp+Geps+2.0*Heps2;
553             fgb              = gbqqfactor*FF*gbscale;
554             dvdatmp          = -0.5*(vgb+fgb*r{I}{J});
555             dvdasum          = dvdasum + dvdatmp;
556             dvda[jnr]        = dvdaj+dvdatmp*isaj{J}*isaj{J};
557             /*                 #define INNERFLOPS INNERFLOPS+13 */
558             /*             #endif */
559             velec            = qq{I}{J}*rinv{I}{J};
560             /*                 #define INNERFLOPS INNERFLOPS+1 */
561             /*             #if 'Force' in KERNEL_VF */
562             felec            = (velec*rinv{I}{J}-fgb)*rinv{I}{J};
563             /*                 #define INNERFLOPS INNERFLOPS+3 */
564             /*             #endif */
565
566             /*         #elif KERNEL_ELEC=='Ewald' */
567             /* EWALD ELECTROSTATICS */
568
569             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
570             ewrt             = r{I}{J}*ewtabscale;
571             ewitab           = ewrt;
572             eweps            = ewrt-ewitab;
573             /*             #define INNERFLOPS INNERFLOPS+2 */
574             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
575             ewitab           = 4*ewitab;
576             felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
577             /*                 #define INNERFLOPS INNERFLOPS+4 */
578             /*                 #if KERNEL_MOD_ELEC=='PotentialShift' */
579             velec            = qq{I}{J}*((rinv{I}{J}-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
580             /*                     #define INNERFLOPS INNERFLOPS+7 */
581             /*                 #else */
582             velec            = qq{I}{J}*(rinv{I}{J}-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
583             /*                     #define INNERFLOPS INNERFLOPS+6 */
584             /*                 #endif */
585             /*                 #if 'Force' in KERNEL_VF */
586             felec            = qq{I}{J}*rinv{I}{J}*(rinvsq{I}{J}-felec);
587             /*                      #define INNERFLOPS INNERFLOPS+3 */
588             /*                 #endif */
589             /*             #elif KERNEL_VF=='Force' */
590             felec            = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
591             felec            = qq{I}{J}*rinv{I}{J}*(rinvsq{I}{J}-felec);
592             /*                 #define INNERFLOPS INNERFLOPS+7 */
593             /*             #endif */
594
595             /*         #elif KERNEL_ELEC=='CubicSplineTable' */
596
597             /* CUBIC SPLINE TABLE ELECTROSTATICS */
598             /*             #if 'Potential' in KERNEL_VF */
599             Y                = vftab[vfitab];
600             /*             #endif */
601             F                = vftab[vfitab+1];
602             Geps             = vfeps*vftab[vfitab+2];
603             Heps2            = vfeps*vfeps*vftab[vfitab+3];
604             Fp               = F+Geps+Heps2;
605             /*             #define INNERFLOPS INNERFLOPS+5 */
606             /*             #if 'Potential' in KERNEL_VF */
607             VV               = Y+vfeps*Fp;
608             velec            = qq{I}{J}*VV;
609             /*                 #define INNERFLOPS INNERFLOPS+3 */
610             /*             #endif */
611             /*             #if 'Force' in KERNEL_VF */
612             FF               = Fp+Geps+2.0*Heps2;
613             felec            = -qq{I}{J}*FF*vftabscale*rinv{I}{J};
614             /*                 #define INNERFLOPS INNERFLOPS+7 */
615             /*             #endif */
616             /*         #endif */
617             /*         ## End of check for electrostatics interaction forms */
618             /*     #endif */
619             /*     ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
620
621             /*     #if 'vdw' in INTERACTION_FLAGS[I][J] */
622
623             /*         #if KERNEL_VDW=='LennardJones' */
624
625             /* LENNARD-JONES DISPERSION/REPULSION */
626
627             rinvsix          = rinvsq{I}{J}*rinvsq{I}{J}*rinvsq{I}{J};
628             /*             #define INNERFLOPS INNERFLOPS+2 */
629             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
630             vvdw6            = c6_{I}{J}*rinvsix;
631             vvdw12           = c12_{I}{J}*rinvsix*rinvsix;
632             /*                 #define INNERFLOPS INNERFLOPS+3 */
633             /*                 #if KERNEL_MOD_VDW=='PotentialShift' */
634             vvdw             = (vvdw12 - c12_{I}{J}*sh_vdw_invrcut6*sh_vdw_invrcut6)*(1.0/12.0) - (vvdw6 - c6_{I}{J}*sh_vdw_invrcut6)*(1.0/6.0);
635             /*                     #define INNERFLOPS INNERFLOPS+8 */
636             /*                 #else */
637             vvdw             = vvdw12*(1.0/12.0) - vvdw6*(1.0/6.0);
638             /*                     #define INNERFLOPS INNERFLOPS+3 */
639             /*                 #endif */
640             /*                 ## Check for force inside potential check, i.e. this means we already did the potential part */
641             /*                 #if 'Force' in KERNEL_VF */
642             fvdw             = (vvdw12-vvdw6)*rinvsq{I}{J};
643             /*                     #define INNERFLOPS INNERFLOPS+2 */
644             /*                 #endif */
645             /*             #elif KERNEL_VF=='Force' */
646             /*                 ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
647             fvdw             = (c12_{I}{J}*rinvsix-c6_{I}{J})*rinvsix*rinvsq{I}{J};
648             /*                 #define INNERFLOPS INNERFLOPS+4 */
649             /*             #endif */
650
651             /*         #elif KERNEL_VDW=='Buckingham' */
652
653             /* BUCKINGHAM DISPERSION/REPULSION */
654             rinvsix          = rinvsq{I}{J}*rinvsq{I}{J}*rinvsq{I}{J};
655             vvdw6            = c6_{I}{J}*rinvsix;
656             br               = cexp2_{I}{J}*r{I}{J};
657             vvdwexp          = cexp1_{I}{J}*exp(-br);
658             /*             ## Estimate exp() to 25 flops */
659             /*             #define INNERFLOPS INNERFLOPS+31 */
660             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch'  */
661             /*                 #if KERNEL_MOD_VDW=='PotentialShift' */
662             vvdw             = (vvdwexp-cexp1_{I}{J}*exp(-cexp2_{I}{J}*rvdw)) - (vvdw6 - c6_{I}{J}*sh_vdw_invrcut6)*(1.0/6.0);
663             /*                     #define INNERFLOPS INNERFLOPS+33 */
664             /*                  #else */
665             vvdw             = vvdwexp - vvdw6*(1.0/6.0);
666             /*                     #define INNERFLOPS INNERFLOPS+2 */
667             /*                 #endif */
668             /*             #endif */
669             /*             #if 'Force' in KERNEL_VF */
670             fvdw             = (br*vvdwexp-vvdw6)*rinvsq{I}{J};
671             /*                 #define INNERFLOPS INNERFLOPS+3 */
672             /*             #endif */
673
674             /*         #elif KERNEL_VDW=='CubicSplineTable' */
675
676             /* CUBIC SPLINE TABLE DISPERSION */
677             vfitab          += {TABLE_VDW_OFFSET};
678             /*             #if 'Potential' in KERNEL_VF */
679             Y                = vftab[vfitab];
680             /*             #endif */
681             F                = vftab[vfitab+1];
682             Geps             = vfeps*vftab[vfitab+2];
683             Heps2            = vfeps*vfeps*vftab[vfitab+3];
684             Fp               = F+Geps+Heps2;
685             /*             #define INNERFLOPS INNERFLOPS+5 */
686             /*             #if 'Potential' in KERNEL_VF */
687             VV               = Y+vfeps*Fp;
688             vvdw6            = c6_{I}{J}*VV;
689             /*                 #define INNERFLOPS INNERFLOPS+3 */
690             /*             #endif */
691             /*             #if 'Force' in KERNEL_VF */
692             FF               = Fp+Geps+2.0*Heps2;
693             fvdw6            = c6_{I}{J}*FF;
694             /*                 #define INNERFLOPS INNERFLOPS+4 */
695             /*             #endif */
696
697             /* CUBIC SPLINE TABLE REPULSION */
698             /*             #if 'Potential' in KERNEL_VF */
699             Y                = vftab[vfitab+4];
700             /*             #endif */
701             F                = vftab[vfitab+5];
702             Geps             = vfeps*vftab[vfitab+6];
703             Heps2            = vfeps*vfeps*vftab[vfitab+7];
704             Fp               = F+Geps+Heps2;
705             /*             #define INNERFLOPS INNERFLOPS+5 */
706             /*             #if 'Potential' in KERNEL_VF */
707             VV               = Y+vfeps*Fp;
708             vvdw12           = c12_{I}{J}*VV;
709             /*                 #define INNERFLOPS INNERFLOPS+3 */
710             /*             #endif */
711             /*             #if 'Force' in KERNEL_VF */
712             FF               = Fp+Geps+2.0*Heps2;
713             fvdw12           = c12_{I}{J}*FF;
714             /*                 #define INNERFLOPS INNERFLOPS+4 */
715             /*             #endif */
716             /*             #if 'Potential' in KERNEL_VF */
717             vvdw             = vvdw12+vvdw6;
718             /*                 #define INNERFLOPS INNERFLOPS+1 */
719             /*             #endif */
720             /*             #if 'Force' in KERNEL_VF */
721             fvdw             = -(fvdw6+fvdw12)*vftabscale*rinv{I}{J};
722             /*                 #define INNERFLOPS INNERFLOPS+4 */
723             /*             #endif */
724
725             /*         #elif KERNEL_VDW=='LJEwald' */
726
727             rinvsix          = rinvsq{I}{J}*rinvsq{I}{J}*rinvsq{I}{J};
728             ewcljrsq         = ewclj2*rsq{I}{J};
729             exponent         = exp(-ewcljrsq);
730             poly             = exponent*(1.0 + ewcljrsq + ewcljrsq*ewcljrsq*0.5);
731             /*             #define INNERFLOPS INNERFLOPS+9 */
732             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
733             vvdw6            = (c6_{I}{J}-c6grid_{I}{J}*(1.0-poly))*rinvsix;
734             vvdw12           = c12_{I}{J}*rinvsix*rinvsix;
735             /*                 #define INNERFLOPS INNERFLOPS+6 */
736             /*                 #if KERNEL_MOD_VDW=='PotentialShift' */
737             vvdw             = (vvdw12 - c12_{I}{J}*sh_vdw_invrcut6*sh_vdw_invrcut6)*(1.0/12.0) - (vvdw6 - c6_{I}{J}*sh_vdw_invrcut6 - c6grid_{I}{J}*sh_lj_ewald)*(1.0/6.0);
738             /*                     #define INNERFLOPS INNERFLOPS+9 */
739             /*                 #else */
740             vvdw             = vvdw12*(1.0/12.0) - vvdw6*(1.0/6.0);
741             /*                     #define INNERFLOPS INNERFLOPS+3 */
742             /*                 #endif */
743             /*                 ## Check for force inside potential check, i.e. this means we already did the potential part */
744             /*                 #if 'Force' in KERNEL_VF */
745             fvdw             = (vvdw12 - vvdw6 - c6grid_{I}{J}*(1.0/6.0)*exponent*ewclj6)*rinvsq{I}{J};
746             /*                     #define INNERFLOPS INNERFLOPS+6 */
747             /*                 #endif */
748             /*             #elif KERNEL_VF=='Force' */
749             /*                 ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
750             fvdw             = (((c12_{I}{J}*rinvsix - c6_{I}{J} + c6grid_{I}{J}*(1.0-poly))*rinvsix) - c6grid_{I}{J}*(1.0/6.0)*exponent*ewclj6)*rinvsq{I}{J};
751             /*                 #define INNERFLOPS INNERFLOPS+11 */
752             /*             #endif */
753             /*         #endif */
754             /*         ## End of check for vdw interaction forms */
755             /*     #endif */
756             /*     ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
757
758             /*     #if 'switch' in INTERACTION_FLAGS[I][J] */
759             d                = r{I}{J}-rswitch;
760             d                = (d>0.0) ? d : 0.0;
761             d2               = d*d;
762             sw               = 1.0+d2*d*(swV3+d*(swV4+d*swV5));
763             /*         #define INNERFLOPS INNERFLOPS+9 */
764
765             /*         #if 'Force' in KERNEL_VF */
766             dsw              = d2*(swF2+d*(swF3+d*swF4));
767             /*             #define INNERFLOPS INNERFLOPS+5 */
768             /*         #endif */
769
770             /* Evaluate switch function */
771             /*         #if 'Force' in KERNEL_VF */
772             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
773             /*             #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
774             felec            = felec*sw - rinv{I}{J}*velec*dsw;
775             /*                 #define INNERFLOPS INNERFLOPS+3 */
776             /*             #endif */
777             /*             #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
778             fvdw             = fvdw*sw - rinv{I}{J}*vvdw*dsw;
779             /*                 #define INNERFLOPS INNERFLOPS+3 */
780             /*             #endif */
781             /*         #endif */
782             /*         #if 'Potential' in KERNEL_VF */
783             /*             #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
784             velec           *= sw;
785             /*                 #define INNERFLOPS INNERFLOPS+1 */
786             /*             #endif */
787             /*             #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
788             vvdw            *= sw;
789             /*                 #define INNERFLOPS INNERFLOPS+1 */
790             /*             #endif */
791             /*         #endif */
792             /*     #endif */
793
794             /*     #if 'Potential' in KERNEL_VF */
795             /* Update potential sums from outer loop */
796             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
797             velecsum        += velec;
798             /*             #define INNERFLOPS INNERFLOPS+1 */
799             /*             #if KERNEL_ELEC=='GeneralizedBorn' */
800             vgbsum          += vgb;
801             /*                 #define INNERFLOPS INNERFLOPS+1 */
802             /*             #endif */
803             /*         #endif */
804             /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
805             vvdwsum         += vvdw;
806             /*             #define INNERFLOPS INNERFLOPS+1 */
807             /*         #endif */
808             /*     #endif */
809
810             /*     #if 'Force' in KERNEL_VF */
811
812             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
813             fscal            = felec+fvdw;
814             /*             #define INNERFLOPS INNERFLOPS+1 */
815             /*         #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
816             fscal            = felec;
817             /*         #elif 'vdw' in INTERACTION_FLAGS[I][J] */
818             fscal            = fvdw;
819             /*        #endif */
820
821             /* Calculate temporary vectorial force */
822             tx               = fscal*dx{I}{J};
823             ty               = fscal*dy{I}{J};
824             tz               = fscal*dz{I}{J};
825
826             /* Update vectorial force */
827             fix{I}            += tx;
828             fiy{I}            += ty;
829             fiz{I}            += tz;
830             f[j_coord_offset+DIM*{J}+XX] -= tx;
831             f[j_coord_offset+DIM*{J}+YY] -= ty;
832             f[j_coord_offset+DIM*{J}+ZZ] -= tz;
833
834             /*         #define INNERFLOPS INNERFLOPS+9 */
835             /*     #endif */
836
837             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
838             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
839             /*         #if 0    ## This and next two lines is a hack to maintain indentation in template file */
840             {
841                 /*     #endif */
842             }
843             /*     #endif */
844             /*    ## End of check for the interaction being outside the cutoff */
845
846             /* #endfor */
847             /* ## End of loop over i-j interaction pairs */
848
849             /* Inner loop uses {INNERFLOPS} flops */
850         }
851         /* End of innermost loop */
852
853         /* #if 'Force' in KERNEL_VF */
854         tx = ty = tz = 0;
855         /*     #for I in PARTICLES_I */
856         f[i_coord_offset+DIM*{I}+XX] += fix{I};
857         f[i_coord_offset+DIM*{I}+YY] += fiy{I};
858         f[i_coord_offset+DIM*{I}+ZZ] += fiz{I};
859         tx                         += fix{I};
860         ty                         += fiy{I};
861         tz                         += fiz{I};
862         /*         #define OUTERFLOPS OUTERFLOPS+6 */
863         /*     #endfor */
864         fshift[i_shift_offset+XX]  += tx;
865         fshift[i_shift_offset+YY]  += ty;
866         fshift[i_shift_offset+ZZ]  += tz;
867         /*     #define OUTERFLOPS OUTERFLOPS+3 */
868         /* #endif */
869
870         /* #if 'Potential' in KERNEL_VF */
871         ggid                        = gid[iidx];
872         /* Update potential energies */
873         /*     #if KERNEL_ELEC != 'None' */
874         kernel_data->energygrp_elec[ggid] += velecsum;
875         /*         #define OUTERFLOPS OUTERFLOPS+1 */
876         /*     #endif */
877         /*     #if 'GeneralizedBorn' in KERNEL_ELEC */
878         kernel_data->energygrp_polarization[ggid] += vgbsum;
879         /*         #define OUTERFLOPS OUTERFLOPS+1 */
880         /*     #endif */
881         /*     #if KERNEL_VDW != 'None' */
882         kernel_data->energygrp_vdw[ggid] += vvdwsum;
883         /*         #define OUTERFLOPS OUTERFLOPS+1 */
884         /*     #endif */
885         /* #endif */
886         /*     #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
887         dvda[inr]                   = dvda[inr] + dvdasum*isai{I}*isai{I};
888         /*     #endif */
889
890         /* Increment number of inner iterations */
891         inneriter                  += j_index_end - j_index_start;
892
893         /* Outer loop uses {OUTERFLOPS} flops */
894     }
895
896     /* Increment number of outer iterations */
897     outeriter        += nri;
898
899     /* Update outer/inner flops */
900     /* ## NB: This is not important, it just affects the flopcount. However, since our preprocessor is */
901     /* ## primitive and replaces aggressively even in strings inside these directives, we need to      */
902     /* ## assemble the main part of the name (containing KERNEL/ELEC/VDW) directly in the source.      */
903     /* #if GEOMETRY_I == 'Water3'            */
904     /*     #define ISUFFIX '_W3'             */
905     /* #elif GEOMETRY_I == 'Water4'          */
906     /*     #define ISUFFIX '_W4'             */
907     /* #else                                 */
908     /*     #define ISUFFIX ''                */
909     /* #endif                                */
910     /* #if GEOMETRY_J == 'Water3'            */
911     /*     #define JSUFFIX 'W3'              */
912     /* #elif GEOMETRY_J == 'Water4'          */
913     /*     #define JSUFFIX 'W4'              */
914     /* #else                                 */
915     /*     #define JSUFFIX ''                */
916     /* #endif                                */
917     /* #if 'PotentialAndForce' in KERNEL_VF  */
918     /*     #define VFSUFFIX  '_VF'           */
919     /* #elif 'Potential' in KERNEL_VF        */
920     /*     #define VFSUFFIX '_V'             */
921     /* #else                                 */
922     /*     #define VFSUFFIX '_F'             */
923     /* #endif                                */
924
925     /* #if KERNEL_ELEC != 'None' and KERNEL_VDW != 'None' */
926     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
927     /* #elif KERNEL_ELEC != 'None' */
928     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
929     /* #else */
930     inc_nrnb(nrnb,eNR_NBKERNEL_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
931     /* #endif  */
932 }