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