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