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