Remove support for implicit solvation
[alexxy/gromacs.git] / src / gromacs / gmxlib / nonbonded / nb_kernel_avx_128_fma_single / nb_kernel_template_avx_128_fma_single.pre
1 /* #if 0 */
2 /*
3  * This file is part of the GROMACS molecular simulation package.
4  *
5  * Copyright (c) 2012,2013,2014,2015,2017,2018, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 #error This file must be processed with the Gromacs pre-preprocessor
37 /* #endif */
38 /* #if INCLUDE_HEADER */
39 #include "gmxpre.h"
40
41 #include "config.h"
42
43 #include <math.h>
44
45 #include "../nb_kernel.h"
46 #include "gromacs/gmxlib/nrnb.h"
47
48 #include "kernelutil_x86_avx_128_fma_single.h"
49 /* #endif */
50
51 /* ## List of variables set by the generating script:                                    */
52 /* ##                                                                                    */
53 /* ## Setttings that apply to the entire kernel:                                         */
54 /* ## KERNEL_ELEC:           String, choice for electrostatic interactions               */
55 /* ## KERNEL_VDW:            String, choice for van der Waals interactions               */
56 /* ## KERNEL_NAME:           String, name of this kernel                                 */
57 /* ## KERNEL_VF:             String telling if we calculate potential, force, or both    */
58 /* ## GEOMETRY_I/GEOMETRY_J: String, name of each geometry, e.g. 'Water3' or '1Particle' */
59 /* ##                                                                                    */
60 /* ## Setttings that apply to particles in the outer (I) or inner (J) loops:             */
61 /* ## PARTICLES_I[]/         Arrays with lists of i/j particles to use in kernel. It is  */
62 /* ## PARTICLES_J[]:         just [0] for particle geometry, but can be longer for water */
63 /* ## PARTICLES_ELEC_I[]/    Arrays with lists of i/j particle that have electrostatics  */
64 /* ## PARTICLES_ELEC_J[]:    interactions that should be calculated in this kernel.      */
65 /* ## PARTICLES_VDW_I[]/     Arrays with the list of i/j particle that have VdW          */
66 /* ## PARTICLES_VDW_J[]:     interactions that should be calculated in this kernel.      */
67 /* ##                                                                                    */
68 /* ## Setttings for pairs of interactions (e.g. 2nd i particle against 1st j particle)   */
69 /* ## PAIRS_IJ[]:            Array with (i,j) tuples of pairs for which interactions     */
70 /* ##                        should be calculated in this kernel. Zero-charge particles  */
71 /* ##                        do not have interactions with particles without vdw, and    */
72 /* ##                        Vdw-only interactions are not evaluated in a no-vdw-kernel. */
73 /* ## INTERACTION_FLAGS[][]: 2D matrix, dimension e.g. 3*3 for water-water interactions. */
74 /* ##                        For each i-j pair, the element [I][J] is a list of strings  */
75 /* ##                        defining properties/flags of this interaction. Examples     */
76 /* ##                        include 'electrostatics'/'vdw' if that type of interaction  */
77 /* ##                        should be evaluated, 'rsq'/'rinv'/'rinvsq' if those values  */
78 /* ##                        are needed, and 'exactcutoff' or 'shift','switch' to        */
79 /* ##                        decide if the force/potential should be modified. This way  */
80 /* ##                        we only calculate values absolutely needed for each case.   */
81
82 /* ## Calculate the size and offset for (merged/interleaved) table data */
83
84 /*
85  * Gromacs nonbonded kernel:   {KERNEL_NAME}
86  * Electrostatics interaction: {KERNEL_ELEC}
87  * VdW interaction:            {KERNEL_VDW}
88  * Geometry:                   {GEOMETRY_I}-{GEOMETRY_J}
89  * Calculate force/pot:        {KERNEL_VF}
90  */
91 void
92 {KERNEL_NAME}
93                     (t_nblist                    * gmx_restrict       nlist,
94                      rvec                        * gmx_restrict          xx,
95                      rvec                        * gmx_restrict          ff,
96                      struct t_forcerec           * gmx_restrict          fr,
97                      t_mdatoms                   * gmx_restrict     mdatoms,
98                      nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
99                      t_nrnb                      * gmx_restrict        nrnb)
100 {
101     /* ## Not all variables are used for all kernels, but any optimizing compiler fixes that, */
102     /* ## so there is no point in going to extremes to exclude variables that are not needed. */
103     /* Suffixes 0,1,2,3 refer to particle indices for waters in the inner or outer loop, or
104      * just 0 for non-waters.
105      * Suffixes A,B,C,D refer to j loop unrolling done with AVX_128, e.g. for the four different
106      * jnr indices corresponding to data put in the four positions in the SIMD register.
107      */
108     int              i_shift_offset,i_coord_offset,outeriter,inneriter;
109     int              j_index_start,j_index_end,jidx,nri,inr,ggid,iidx;
110     int              jnrA,jnrB,jnrC,jnrD;
111     int              jnrlistA,jnrlistB,jnrlistC,jnrlistD;
112     int              j_coord_offsetA,j_coord_offsetB,j_coord_offsetC,j_coord_offsetD;
113     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
114     real             rcutoff_scalar;
115     real             *shiftvec,*fshift,*x,*f;
116     real             *fjptrA,*fjptrB,*fjptrC,*fjptrD;
117     real             scratch[4*DIM];
118     __m128           fscal,rcutoff,rcutoff2,jidxall;
119     /* #for I in PARTICLES_I */
120     int              vdwioffset{I};
121     __m128           ix{I},iy{I},iz{I},fix{I},fiy{I},fiz{I},iq{I},isai{I};
122     /* #endfor */
123     /* #for J in PARTICLES_J */
124     int              vdwjidx{J}A,vdwjidx{J}B,vdwjidx{J}C,vdwjidx{J}D;
125     __m128           jx{J},jy{J},jz{J},fjx{J},fjy{J},fjz{J},jq{J},isaj{J};
126     /* #endfor */
127     /* #for I,J in PAIRS_IJ */
128     __m128           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};
129     /* #endfor */
130     /* #if KERNEL_ELEC != 'None' */
131     __m128           velec,felec,velecsum,facel,crf,krf,krf2;
132     real             *charge;
133     /* #endif */
134     /* #if KERNEL_VDW != 'None' */
135     int              nvdwtype;
136     __m128           rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,sh_vdw_invrcut6;
137     int              *vdwtype;
138     real             *vdwparam;
139     __m128           one_sixth   = _mm_set1_ps(1.0/6.0);
140     __m128           one_twelfth = _mm_set1_ps(1.0/12.0);
141     /* #endif */
142     /* #if 'Table' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
143     __m128i          vfitab;
144     __m128i          ifour       = _mm_set1_epi32(4);
145     __m128           rt,vfeps,twovfeps,vftabscale,Y,F,G,H,Fp,VV,FF;
146     real             *vftab;
147     /* #endif */
148     /* #if 'LJEwald' in KERNEL_VDW */
149     /* #for I,J in PAIRS_IJ */
150     __m128           c6grid_{I}{J};
151     /* #endfor */
152     real             *vdwgridparam;
153     __m128           ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,f6A,f6B,sh_lj_ewald;
154     __m128           one_half = _mm_set1_ps(0.5);
155     __m128           minus_one = _mm_set1_ps(-1.0);
156     /* #endif */
157     /* #if 'Ewald' in KERNEL_ELEC */
158     __m128i          ewitab;
159     __m128           ewtabscale,eweps,twoeweps,sh_ewald,ewrt,ewtabhalfspace,ewtabF,ewtabFn,ewtabD,ewtabV;
160     __m128           beta,beta2,beta3,zeta2,pmecorrF,pmecorrV,rinv3;
161     real             *ewtab;
162     /* #endif */
163     /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
164     __m128           rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
165     real             rswitch_scalar,d_scalar;
166     /* #endif */
167     __m128           dummy_mask,cutoff_mask;
168     __m128           signbit = _mm_castsi128_ps( _mm_set1_epi32(0x80000000) );
169     __m128           one     = _mm_set1_ps(1.0);
170     __m128           two     = _mm_set1_ps(2.0);
171     x                = xx[0];
172     f                = ff[0];
173
174     nri              = nlist->nri;
175     iinr             = nlist->iinr;
176     jindex           = nlist->jindex;
177     jjnr             = nlist->jjnr;
178     shiftidx         = nlist->shift;
179     gid              = nlist->gid;
180     shiftvec         = fr->shift_vec[0];
181     fshift           = fr->fshift[0];
182     /* #if KERNEL_ELEC != 'None' */
183     facel            = _mm_set1_ps(fr->ic->epsfac);
184     charge           = mdatoms->chargeA;
185     /*     #if 'ReactionField' in KERNEL_ELEC */
186     krf              = _mm_set1_ps(fr->ic->k_rf);
187     krf2             = _mm_set1_ps(fr->ic->k_rf*2.0);
188     crf              = _mm_set1_ps(fr->ic->c_rf);
189     /*     #endif */
190     /* #endif */
191     /* #if KERNEL_VDW != 'None' */
192     nvdwtype         = fr->ntype;
193     vdwparam         = fr->nbfp;
194     vdwtype          = mdatoms->typeA;
195     /* #endif */
196     /* #if 'LJEwald' in KERNEL_VDW */
197     vdwgridparam     = fr->ljpme_c6grid;
198     sh_lj_ewald      = _mm_set1_ps(fr->ic->sh_lj_ewald);
199     ewclj            = _mm_set1_ps(fr->ic->ewaldcoeff_lj);
200     ewclj2           = _mm_mul_ps(minus_one,_mm_mul_ps(ewclj,ewclj));
201     /* #endif */
202
203     /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
204     vftab            = kernel_data->table_elec_vdw->data;
205     vftabscale       = _mm_set1_ps(kernel_data->table_elec_vdw->scale);
206     /* #elif 'Table' in KERNEL_ELEC */
207     vftab            = kernel_data->table_elec->data;
208     vftabscale       = _mm_set1_ps(kernel_data->table_elec->scale);
209     /* #elif 'Table' in KERNEL_VDW */
210     vftab            = kernel_data->table_vdw->data;
211     vftabscale       = _mm_set1_ps(kernel_data->table_vdw->scale);
212     /* #endif */
213
214     /* #if 'Ewald' in KERNEL_ELEC */
215     sh_ewald         = _mm_set1_ps(fr->ic->sh_ewald);
216     beta             = _mm_set1_ps(fr->ic->ewaldcoeff_q);
217     beta2            = _mm_mul_ps(beta,beta);
218     beta3            = _mm_mul_ps(beta,beta2);
219     /*     #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
220     ewtab            = fr->ic->tabq_coul_F;
221     ewtabscale       = _mm_set1_ps(fr->ic->tabq_scale);
222     ewtabhalfspace   = _mm_set1_ps(0.5/fr->ic->tabq_scale);
223     /*     #else */
224     ewtab            = fr->ic->tabq_coul_FDV0;
225     ewtabscale       = _mm_set1_ps(fr->ic->tabq_scale);
226     ewtabhalfspace   = _mm_set1_ps(0.5/fr->ic->tabq_scale);
227      /*     #endif */
228     /* #endif */
229
230     /* #if 'Water' in GEOMETRY_I */
231     /* Setup water-specific parameters */
232     inr              = nlist->iinr[0];
233     /*     #for I in PARTICLES_ELEC_I */
234     iq{I}              = _mm_mul_ps(facel,_mm_set1_ps(charge[inr+{I}]));
235     /*     #endfor */
236     /*     #for I in PARTICLES_VDW_I */
237     vdwioffset{I}      = 2*nvdwtype*vdwtype[inr+{I}];
238     /*     #endfor */
239     /* #endif */
240
241     /* #if 'Water' in GEOMETRY_J */
242     /*     #for J in PARTICLES_ELEC_J */
243     jq{J}              = _mm_set1_ps(charge[inr+{J}]);
244     /*     #endfor */
245     /*     #for J in PARTICLES_VDW_J */
246     vdwjidx{J}A        = 2*vdwtype[inr+{J}];
247     /*     #endfor */
248     /*     #for I,J in PAIRS_IJ */
249     /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
250     qq{I}{J}             = _mm_mul_ps(iq{I},jq{J});
251     /*         #endif */
252     /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
253     /*             #if 'LJEwald' in KERNEL_VDW */
254     c6_{I}{J}            = _mm_set1_ps(vdwparam[vdwioffset{I}+vdwjidx{J}A]);
255     c12_{I}{J}           = _mm_set1_ps(vdwparam[vdwioffset{I}+vdwjidx{J}A+1]);
256     c6grid_{I}{J}        = _mm_set1_ps(vdwgridparam[vdwioffset{I}+vdwjidx{J}A]);
257     /*             #else */
258     c6_{I}{J}            = _mm_set1_ps(vdwparam[vdwioffset{I}+vdwjidx{J}A]);
259     c12_{I}{J}           = _mm_set1_ps(vdwparam[vdwioffset{I}+vdwjidx{J}A+1]);
260     /*             #endif */
261     /*         #endif */
262     /*     #endfor */
263     /* #endif */
264
265     /* #if KERNEL_MOD_ELEC!='None' or KERNEL_MOD_VDW!='None' */
266     /*     #if KERNEL_ELEC!='None' */
267     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
268     rcutoff_scalar   = fr->ic->rcoulomb;
269     /*     #else */
270     rcutoff_scalar   = fr->ic->rvdw;
271     /*     #endif */
272     rcutoff          = _mm_set1_ps(rcutoff_scalar);
273     rcutoff2         = _mm_mul_ps(rcutoff,rcutoff);
274     /* #endif */
275
276     /* #if KERNEL_MOD_VDW=='PotentialShift' */
277     sh_vdw_invrcut6  = _mm_set1_ps(fr->ic->sh_invrc6);
278     rvdw             = _mm_set1_ps(fr->ic->rvdw);
279     /* #endif */
280
281     /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
282     /*     #if KERNEL_MOD_ELEC=='PotentialSwitch'  */
283     rswitch_scalar   = fr->ic->rcoulomb_switch;
284     rswitch          = _mm_set1_ps(rswitch_scalar);
285     /*     #else */
286     rswitch_scalar   = fr->ic->rvdw_switch;
287     rswitch          = _mm_set1_ps(rswitch_scalar);
288     /*     #endif */
289     /* Setup switch parameters */
290     d_scalar         = rcutoff_scalar-rswitch_scalar;
291     d                = _mm_set1_ps(d_scalar);
292     swV3             = _mm_set1_ps(-10.0/(d_scalar*d_scalar*d_scalar));
293     swV4             = _mm_set1_ps( 15.0/(d_scalar*d_scalar*d_scalar*d_scalar));
294     swV5             = _mm_set1_ps( -6.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
295     /*     #if 'Force' in KERNEL_VF */
296     swF2             = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar));
297     swF3             = _mm_set1_ps( 60.0/(d_scalar*d_scalar*d_scalar*d_scalar));
298     swF4             = _mm_set1_ps(-30.0/(d_scalar*d_scalar*d_scalar*d_scalar*d_scalar));
299     /*     #endif */
300     /* #endif */
301
302     /* Avoid stupid compiler warnings */
303     jnrA = jnrB = jnrC = jnrD = 0;
304     j_coord_offsetA = 0;
305     j_coord_offsetB = 0;
306     j_coord_offsetC = 0;
307     j_coord_offsetD = 0;
308
309     /* ## Keep track of the floating point operations we issue for reporting! */
310     /* #define OUTERFLOPS 0 */
311     outeriter        = 0;
312     inneriter        = 0;
313
314     for(iidx=0;iidx<4*DIM;iidx++)
315     {
316         scratch[iidx] = 0.0;
317     }
318
319     /* Start outer loop over neighborlists */
320     for(iidx=0; iidx<nri; iidx++)
321     {
322         /* Load shift vector for this list */
323         i_shift_offset   = DIM*shiftidx[iidx];
324
325         /* Load limits for loop over neighbors */
326         j_index_start    = jindex[iidx];
327         j_index_end      = jindex[iidx+1];
328
329         /* Get outer coordinate index */
330         inr              = iinr[iidx];
331         i_coord_offset   = DIM*inr;
332
333         /* Load i particle coords and add shift vector */
334         /* #if GEOMETRY_I == 'Particle' */
335         gmx_mm_load_shift_and_1rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,&ix0,&iy0,&iz0);
336         /* #elif GEOMETRY_I == 'Water3' */
337         gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
338                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2);
339         /* #elif GEOMETRY_I == 'Water4' */
340         /*     #if 0 in PARTICLES_I                 */
341         gmx_mm_load_shift_and_4rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset,
342                                                  &ix0,&iy0,&iz0,&ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
343         /*     #else                                */
344         gmx_mm_load_shift_and_3rvec_broadcast_ps(shiftvec+i_shift_offset,x+i_coord_offset+DIM,
345                                                  &ix1,&iy1,&iz1,&ix2,&iy2,&iz2,&ix3,&iy3,&iz3);
346         /*     #endif                               */
347         /* #endif                                   */
348
349         /* #if 'Force' in KERNEL_VF */
350         /*     #for I in PARTICLES_I */
351         fix{I}             = _mm_setzero_ps();
352         fiy{I}             = _mm_setzero_ps();
353         fiz{I}             = _mm_setzero_ps();
354         /*     #endfor */
355         /* #endif */
356
357         /* ## For water we already preloaded parameters at the start of the kernel */
358         /* #if not 'Water' in GEOMETRY_I */
359         /* Load parameters for i particles */
360         /*     #for I in PARTICLES_ELEC_I */
361         iq{I}              = _mm_mul_ps(facel,_mm_load1_ps(charge+inr+{I}));
362         /*         #define OUTERFLOPS OUTERFLOPS+1 */
363         /*     #endfor */
364         /*     #for I in PARTICLES_VDW_I */
365         vdwioffset{I}      = 2*nvdwtype*vdwtype[inr+{I}];
366         /*     #endfor */
367         /* #endif */
368
369         /* #if 'Potential' in KERNEL_VF */
370         /* Reset potential sums */
371         /*     #if KERNEL_ELEC != 'None' */
372         velecsum         = _mm_setzero_ps();
373         /*     #endif */
374         /*     #if KERNEL_VDW != 'None' */
375         vvdwsum          = _mm_setzero_ps();
376         /*     #endif */
377         /* #endif */
378
379         /* #for ROUND in ['Loop','Epilogue'] */
380
381         /* #if ROUND =='Loop' */
382         /* Start inner kernel loop */
383         for(jidx=j_index_start; jidx<j_index_end && jjnr[jidx+3]>=0; jidx+=4)
384         {
385         /* ## First round is normal loop (next statement resets indentation) */
386         /*     #if 0 */
387         }
388         /*     #endif */
389         /* #else */
390         if(jidx<j_index_end)
391         {
392         /* ## Second round is epilogue */
393         /* #endif */
394         /* #define INNERFLOPS 0 */
395
396             /* Get j neighbor index, and coordinate index */
397             /* #if ROUND =='Loop' */
398             jnrA             = jjnr[jidx];
399             jnrB             = jjnr[jidx+1];
400             jnrC             = jjnr[jidx+2];
401             jnrD             = jjnr[jidx+3];
402             /* #else */
403             jnrlistA         = jjnr[jidx];
404             jnrlistB         = jjnr[jidx+1];
405             jnrlistC         = jjnr[jidx+2];
406             jnrlistD         = jjnr[jidx+3];
407             /* Sign of each element will be negative for non-real atoms.
408              * This mask will be 0xFFFFFFFF for dummy entries and 0x0 for real ones,
409              * so use it as val = _mm_andnot_ps(mask,val) to clear dummy entries.
410              */
411             dummy_mask = gmx_mm_castsi128_ps(_mm_cmplt_epi32(_mm_loadu_si128((const __m128i *)(jjnr+jidx)),_mm_setzero_si128()));
412             jnrA       = (jnrlistA>=0) ? jnrlistA : 0;
413             jnrB       = (jnrlistB>=0) ? jnrlistB : 0;
414             jnrC       = (jnrlistC>=0) ? jnrlistC : 0;
415             jnrD       = (jnrlistD>=0) ? jnrlistD : 0;
416             /* #endif */
417             j_coord_offsetA  = DIM*jnrA;
418             j_coord_offsetB  = DIM*jnrB;
419             j_coord_offsetC  = DIM*jnrC;
420             j_coord_offsetD  = DIM*jnrD;
421
422             /* load j atom coordinates */
423             /* #if GEOMETRY_J == 'Particle'             */
424             gmx_mm_load_1rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
425                                               x+j_coord_offsetC,x+j_coord_offsetD,
426                                               &jx0,&jy0,&jz0);
427             /* #elif GEOMETRY_J == 'Water3'             */
428             gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
429                                               x+j_coord_offsetC,x+j_coord_offsetD,
430                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,&jy2,&jz2);
431             /* #elif GEOMETRY_J == 'Water4'             */
432             /*     #if 0 in PARTICLES_J                 */
433             gmx_mm_load_4rvec_4ptr_swizzle_ps(x+j_coord_offsetA,x+j_coord_offsetB,
434                                               x+j_coord_offsetC,x+j_coord_offsetD,
435                                               &jx0,&jy0,&jz0,&jx1,&jy1,&jz1,&jx2,
436                                               &jy2,&jz2,&jx3,&jy3,&jz3);
437             /*     #else                                */
438             gmx_mm_load_3rvec_4ptr_swizzle_ps(x+j_coord_offsetA+DIM,x+j_coord_offsetB+DIM,
439                                               x+j_coord_offsetC+DIM,x+j_coord_offsetD+DIM,
440                                               &jx1,&jy1,&jz1,&jx2,&jy2,&jz2,&jx3,&jy3,&jz3);
441             /*     #endif                               */
442             /* #endif                                   */
443
444             /* Calculate displacement vector */
445             /* #for I,J in PAIRS_IJ */
446             dx{I}{J}             = _mm_sub_ps(ix{I},jx{J});
447             dy{I}{J}             = _mm_sub_ps(iy{I},jy{J});
448             dz{I}{J}             = _mm_sub_ps(iz{I},jz{J});
449             /*     #define INNERFLOPS INNERFLOPS+3 */
450             /* #endfor */
451
452             /* Calculate squared distance and things based on it */
453             /* #for I,J in PAIRS_IJ */
454             rsq{I}{J}            = gmx_mm_calc_rsq_ps(dx{I}{J},dy{I}{J},dz{I}{J});
455             /*     #define INNERFLOPS INNERFLOPS+5 */
456             /* #endfor */
457
458             /* #for I,J in PAIRS_IJ */
459             /*     #if 'rinv' in INTERACTION_FLAGS[I][J] */
460             rinv{I}{J}           = avx128fma_invsqrt_f(rsq{I}{J});
461             /*         #define INNERFLOPS INNERFLOPS+5 */
462             /*     #endif */
463             /* #endfor */
464
465             /* #for I,J in PAIRS_IJ */
466             /*     #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
467             /*         # if 'rinv' not in INTERACTION_FLAGS[I][J] */
468             rinvsq{I}{J}         = avx128fma_inv_f(rsq{I}{J});
469             /*             #define INNERFLOPS INNERFLOPS+4 */
470             /*         #else */
471             rinvsq{I}{J}         = _mm_mul_ps(rinv{I}{J},rinv{I}{J});
472             /*             #define INNERFLOPS INNERFLOPS+1 */
473             /*         #endif */
474             /*     #endif */
475             /* #endfor */
476
477             /* #if not 'Water' in GEOMETRY_J */
478             /* Load parameters for j particles */
479             /*     #for J in PARTICLES_ELEC_J */
480             jq{J}              = gmx_mm_load_4real_swizzle_ps(charge+jnrA+{J},charge+jnrB+{J},
481                                                               charge+jnrC+{J},charge+jnrD+{J});
482             /*     #endfor */
483             /*     #for J in PARTICLES_VDW_J */
484             vdwjidx{J}A        = 2*vdwtype[jnrA+{J}];
485             vdwjidx{J}B        = 2*vdwtype[jnrB+{J}];
486             vdwjidx{J}C        = 2*vdwtype[jnrC+{J}];
487             vdwjidx{J}D        = 2*vdwtype[jnrD+{J}];
488             /*     #endfor */
489             /* #endif */
490
491             /* #if 'Force' in KERNEL_VF and not 'Particle' in GEOMETRY_I */
492             /*     #for J in PARTICLES_J */
493             fjx{J}             = _mm_setzero_ps();
494             fjy{J}             = _mm_setzero_ps();
495             fjz{J}             = _mm_setzero_ps();
496             /*     #endfor */
497             /* #endif */
498
499             /* #for I,J in PAIRS_IJ */
500
501             /**************************
502              * CALCULATE INTERACTIONS *
503              **************************/
504
505             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
506             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
507             /*         ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
508             if (gmx_mm_any_lt(rsq{I}{J},rcutoff2))
509             {
510                 /*     #if 0    ## this and the next two lines is a hack to maintain auto-indentation in template file */
511             }
512             /*         #endif */
513             /*         #define INNERFLOPS INNERFLOPS+1 */
514             /*     #endif */
515
516             /*     #if 'r' in INTERACTION_FLAGS[I][J] */
517             r{I}{J}              = _mm_mul_ps(rsq{I}{J},rinv{I}{J});
518             /*         #if ROUND == 'Epilogue' */
519             r{I}{J}              = _mm_andnot_ps(dummy_mask,r{I}{J});
520             /*             #define INNERFLOPS INNERFLOPS+1 */
521             /*         #endif */
522             /*         #define INNERFLOPS INNERFLOPS+1 */
523             /*     #endif */
524
525             /*     ## For water geometries we already loaded parameters at the start of the kernel */
526             /*     #if not 'Water' in GEOMETRY_J */
527             /* Compute parameters for interactions between i and j atoms */
528             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
529             qq{I}{J}             = _mm_mul_ps(iq{I},jq{J});
530             /*             #define INNERFLOPS INNERFLOPS+1 */
531             /*         #endif */
532             /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
533             gmx_mm_load_4pair_swizzle_ps(vdwparam+vdwioffset{I}+vdwjidx{J}A,
534                                          vdwparam+vdwioffset{I}+vdwjidx{J}B,
535                                          vdwparam+vdwioffset{I}+vdwjidx{J}C,
536                                          vdwparam+vdwioffset{I}+vdwjidx{J}D,
537                                          &c6_{I}{J},&c12_{I}{J});
538
539             /*           #if 'LJEwald' in KERNEL_VDW */
540             c6grid_{I}{J}       = gmx_mm_load_4real_swizzle_ps(vdwgridparam+vdwioffset{I}+vdwjidx{J}A,
541                                                                vdwgridparam+vdwioffset{I}+vdwjidx{J}B,
542                                                                vdwgridparam+vdwioffset{I}+vdwjidx{J}C,
543                                                                vdwgridparam+vdwioffset{I}+vdwjidx{J}D);
544             /*           #endif */
545
546             /*         #endif */
547             /*     #endif */
548
549             /*     #if 'table' in INTERACTION_FLAGS[I][J] */
550             /* Calculate table index by multiplying r with table scale and truncate to integer */
551             rt               = _mm_mul_ps(r{I}{J},vftabscale);
552             vfitab           = _mm_cvttps_epi32(rt);
553 #ifdef __XOP__
554             vfeps            = _mm_frcz_ps(rt);
555 #else
556             vfeps            = _mm_sub_ps(rt,_mm_round_ps(rt, _MM_FROUND_FLOOR));
557 #endif
558             twovfeps         = _mm_add_ps(vfeps,vfeps);
559             /*         #define INNERFLOPS INNERFLOPS+4                          */
560             /*         #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW     */
561             /*             ## 3 tables, 4 bytes per point: multiply index by 12 */
562             vfitab           = _mm_slli_epi32(_mm_add_epi32(vfitab,_mm_slli_epi32(vfitab,1)),2);
563             /*         #elif 'Table' in KERNEL_ELEC                             */
564             /*             ## 1 table, 4 bytes per point: multiply index by 4   */
565             vfitab           = _mm_slli_epi32(vfitab,2);
566             /*         #elif 'Table' in KERNEL_VDW                              */
567             /*             ## 2 tables, 4 bytes per point: multiply index by 8  */
568             vfitab           = _mm_slli_epi32(vfitab,3);
569             /*         #endif                                                   */
570             /*     #endif */
571
572             /*     ## ELECTROSTATIC INTERACTIONS */
573             /*     #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
574
575             /*         #if KERNEL_ELEC=='Coulomb' */
576
577             /* COULOMB ELECTROSTATICS */
578             velec            = _mm_mul_ps(qq{I}{J},rinv{I}{J});
579             /*             #define INNERFLOPS INNERFLOPS+1 */
580             /*             #if 'Force' in KERNEL_VF */
581             felec            = _mm_mul_ps(velec,rinvsq{I}{J});
582             /*                 #define INNERFLOPS INNERFLOPS+2 */
583             /*             #endif */
584
585             /*         #elif KERNEL_ELEC=='ReactionField' */
586
587             /* REACTION-FIELD ELECTROSTATICS */
588             /*             #if 'Potential' in KERNEL_VF */
589             velec            = _mm_mul_ps(qq{I}{J},_mm_sub_ps(_mm_macc_ps(krf,rsq{I}{J},rinv{I}{J}),crf));
590             /*                 #define INNERFLOPS INNERFLOPS+4 */
591             /*             #endif */
592             /*             #if 'Force' in KERNEL_VF */
593             felec            = _mm_mul_ps(qq{I}{J},_mm_msub_ps(rinv{I}{J},rinvsq{I}{J},krf2));
594             /*                 #define INNERFLOPS INNERFLOPS+3 */
595             /*             #endif */
596
597             /*         #elif KERNEL_ELEC=='Ewald' */
598             /* EWALD ELECTROSTATICS */
599
600             /* Analytical PME correction */
601             zeta2            = _mm_mul_ps(beta2,rsq{I}{J});
602             /*             #if 'Force' in KERNEL_VF */
603             rinv3            = _mm_mul_ps(rinvsq{I}{J},rinv{I}{J});
604             pmecorrF         = avx128fma_pmecorrF_f(zeta2);
605             felec            = _mm_macc_ps(pmecorrF,beta3,rinv3);
606             felec            = _mm_mul_ps(qq{I}{J},felec);
607             /*             #endif */
608             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
609             pmecorrV         = avx128fma_pmecorrV_f(zeta2);
610             /*                 #if KERNEL_MOD_ELEC=='PotentialShift' */
611             velec            = _mm_nmacc_ps(pmecorrV,beta,_mm_sub_ps(rinv{I}{J},sh_ewald));
612             /*                 #else */
613             velec            = _mm_nmacc_ps(pmecorrV,beta,rinv{I}{J});
614             /*                 #endif */
615             velec            = _mm_mul_ps(qq{I}{J},velec);
616             /*             #endif */
617
618             /*         #elif KERNEL_ELEC=='CubicSplineTable' */
619
620             /* CUBIC SPLINE TABLE ELECTROSTATICS */
621             Y                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
622             F                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
623             G                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
624             H                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
625             _MM_TRANSPOSE4_PS(Y,F,G,H);
626             Fp               = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
627             /*             #define INNERFLOPS INNERFLOPS+4 */
628             /*             #if 'Potential' in KERNEL_VF */
629             VV               = _mm_macc_ps(vfeps,Fp,Y);
630             velec            = _mm_mul_ps(qq{I}{J},VV);
631             /*                 #define INNERFLOPS INNERFLOPS+3 */
632             /*             #endif */
633             /*             #if 'Force' in KERNEL_VF */
634             FF               = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
635             felec            = _mm_xor_ps(signbit,_mm_mul_ps(_mm_mul_ps(qq{I}{J},FF),_mm_mul_ps(vftabscale,rinv{I}{J})));
636             /*                 #define INNERFLOPS INNERFLOPS+7 */
637             /*             #endif */
638             /*         #endif */
639             /*         ## End of check for electrostatics interaction forms */
640             /*     #endif */
641             /*     ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
642
643             /*     #if 'vdw' in INTERACTION_FLAGS[I][J] */
644
645             /*         #if KERNEL_VDW=='LennardJones' */
646
647             /* LENNARD-JONES DISPERSION/REPULSION */
648
649             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq{I}{J},rinvsq{I}{J}),rinvsq{I}{J});
650             /*             #define INNERFLOPS INNERFLOPS+2 */
651             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
652             vvdw6            = _mm_mul_ps(c6_{I}{J},rinvsix);
653             vvdw12           = _mm_mul_ps(c12_{I}{J},_mm_mul_ps(rinvsix,rinvsix));
654             /*                 #define INNERFLOPS INNERFLOPS+3 */
655             /*                 #if KERNEL_MOD_VDW=='PotentialShift' */
656             vvdw             = _mm_msub_ps(_mm_nmacc_ps(c12_{I}{J},_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
657                                           _mm_mul_ps( _mm_nmacc_ps(c6_{I}{J},sh_vdw_invrcut6,vvdw6),one_sixth));
658             /*                     #define INNERFLOPS INNERFLOPS+8 */
659             /*                 #else */
660             vvdw             = _mm_msub_ps(vvdw12,one_twelfth,_mm_mul_ps(vvdw6,one_sixth));
661             /*                     #define INNERFLOPS INNERFLOPS+3 */
662             /*                 #endif */
663             /*                 ## Check for force inside potential check, i.e. this means we already did the potential part */
664             /*                 #if 'Force' in KERNEL_VF */
665             fvdw             = _mm_mul_ps(_mm_sub_ps(vvdw12,vvdw6),rinvsq{I}{J});
666             /*                     #define INNERFLOPS INNERFLOPS+2 */
667             /*                 #endif */
668             /*             #elif KERNEL_VF=='Force' */
669             /*                 ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
670             fvdw             = _mm_mul_ps(_mm_msub_ps(c12_{I}{J},rinvsix,c6_{I}{J}),_mm_mul_ps(rinvsix,rinvsq{I}{J}));
671             /*                 #define INNERFLOPS INNERFLOPS+4 */
672             /*             #endif */
673
674             /*         #elif KERNEL_VDW=='CubicSplineTable' */
675
676             /* CUBIC SPLINE TABLE DISPERSION */
677             /*             #if 'Table' in KERNEL_ELEC */
678             vfitab           = _mm_add_epi32(vfitab,ifour);
679             /*             #endif                     */
680             Y                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
681             F                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
682             G                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
683             H                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
684             _MM_TRANSPOSE4_PS(Y,F,G,H);
685             Fp               = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
686             /*             #define INNERFLOPS INNERFLOPS+4 */
687             /*             #if 'Potential' in KERNEL_VF */
688             VV               = _mm_macc_ps(vfeps,Fp,Y);
689             vvdw6            = _mm_mul_ps(c6_{I}{J},VV);
690             /*                 #define INNERFLOPS INNERFLOPS+3 */
691             /*             #endif */
692             /*             #if 'Force' in KERNEL_VF */
693             FF               = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
694             fvdw6            = _mm_mul_ps(c6_{I}{J},FF);
695             /*                 #define INNERFLOPS INNERFLOPS+4 */
696             /*             #endif */
697
698             /* CUBIC SPLINE TABLE REPULSION */
699             vfitab           = _mm_add_epi32(vfitab,ifour);
700             Y                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,0) );
701             F                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,1) );
702             G                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,2) );
703             H                = _mm_load_ps( vftab + _mm_extract_epi32(vfitab,3) );
704             _MM_TRANSPOSE4_PS(Y,F,G,H);
705             Fp               = _mm_macc_ps(vfeps,_mm_macc_ps(H,vfeps,G),F);
706             /*             #define INNERFLOPS INNERFLOPS+4 */
707             /*             #if 'Potential' in KERNEL_VF */
708             VV               = _mm_macc_ps(vfeps,Fp,Y);
709             vvdw12           = _mm_mul_ps(c12_{I}{J},VV);
710             /*                 #define INNERFLOPS INNERFLOPS+3 */
711             /*             #endif */
712             /*             #if 'Force' in KERNEL_VF */
713             FF               = _mm_macc_ps(vfeps,_mm_macc_ps(twovfeps,H,G),Fp);
714             fvdw12           = _mm_mul_ps(c12_{I}{J},FF);
715             /*                 #define INNERFLOPS INNERFLOPS+5 */
716             /*             #endif */
717             /*             #if 'Potential' in KERNEL_VF */
718             vvdw             = _mm_add_ps(vvdw12,vvdw6);
719             /*                 #define INNERFLOPS INNERFLOPS+1 */
720             /*             #endif */
721             /*             #if 'Force' in KERNEL_VF */
722             fvdw             = _mm_xor_ps(signbit,_mm_mul_ps(_mm_add_ps(fvdw6,fvdw12),_mm_mul_ps(vftabscale,rinv{I}{J})));
723             /*                 #define INNERFLOPS INNERFLOPS+4 */
724             /*             #endif */
725
726             /*         #elif KERNEL_VDW=='LJEwald' */
727
728             /* Analytical LJ-PME */
729             rinvsix          = _mm_mul_ps(_mm_mul_ps(rinvsq{I}{J},rinvsq{I}{J}),rinvsq{I}{J});
730             ewcljrsq         = _mm_mul_ps(ewclj2,rsq{I}{J});
731             ewclj6           = _mm_mul_ps(ewclj2,_mm_mul_ps(ewclj2,ewclj2));
732             exponent         = avx128fma_exp_f(ewcljrsq);
733             /* poly = exp(-(beta*r)^2) * (1 + (beta*r)^2 + (beta*r)^4 /2) */
734             poly             = _mm_mul_ps(exponent,_mm_macc_ps(_mm_mul_ps(ewcljrsq,ewcljrsq),one_half,_mm_sub_ps(one,ewcljrsq)));
735             /*                 #define INNERFLOPS INNERFLOPS+10 */
736             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
737             /* vvdw6 = [C6 - C6grid * (1-poly)]/r6 */
738             vvdw6            = _mm_mul_ps(_mm_macc_ps(-c6grid_{I}{J},_mm_sub_ps(one,poly),c6_{I}{J}),rinvsix);
739             vvdw12           = _mm_mul_ps(c12_{I}{J},_mm_mul_ps(rinvsix,rinvsix));
740             /*                     #define INNERFLOPS INNERFLOPS+5 */
741             /*                 #if KERNEL_MOD_VDW=='PotentialShift' */
742             vvdw             = _mm_msub_ps(_mm_nmacc_ps(c12_{I}{J},_mm_mul_ps(sh_vdw_invrcut6,sh_vdw_invrcut6),vvdw12),one_twelfth,
743                                           _mm_mul_ps(_mm_sub_ps(vvdw6,_mm_macc_ps(c6grid_{I}{J},sh_lj_ewald,_mm_mul_ps(c6_{I}{J},sh_vdw_invrcut6))),one_sixth));
744             /*                     #define INNERFLOPS INNERFLOPS+7 */
745             /*                 #else */
746             vvdw             = _mm_msub_ps(vvdw12,one_twelfth,_mm_mul_ps(vvdw6,one_sixth));
747             /*                     #define INNERFLOPS INNERFLOPS+2 */
748             /*                 #endif */
749             /*                  ## Check for force inside potential check, i.e. this means we already did the potential part */
750             /*                  #if 'Force' in KERNEL_VF */
751             /* fvdw = vvdw12/r - (vvdw6/r + (C6grid * exponent * beta^6)/r) */
752             fvdw             = _mm_mul_ps(_mm_add_ps(vvdw12,_mm_msub_ps(_mm_mul_ps(c6grid_{I}{J},one_sixth),_mm_mul_ps(exponent,ewclj6),vvdw6)),rinvsq{I}{J});
753             /*                 #define INNERFLOPS INNERFLOPS+5 */
754             /*                  #endif */
755             /*              #elif KERNEL_VF=='Force' */
756             /* f6A = 6 * C6grid * (1 - poly) */
757             f6A              = _mm_mul_ps(c6grid_{I}{J},_mm_sub_ps(one,poly));
758             /* f6B = C6grid * exponent * beta^6 */
759             f6B              = _mm_mul_ps(_mm_mul_ps(c6grid_{I}{J},one_sixth),_mm_mul_ps(exponent,ewclj6));
760             /* fvdw = 12*C12/r13 - ((6*C6 - f6A)/r6 + f6B)/r */
761             fvdw              = _mm_mul_ps(_mm_macc_ps(_mm_msub_ps(c12_{I}{J},rinvsix,_mm_sub_ps(c6_{I}{J},f6A)),rinvsix,f6B),rinvsq{I}{J});
762             /*                 #define INNERFLOPS INNERFLOPS+10 */
763             /*              #endif */
764             /*         #endif */
765             /*         ## End of check for vdw interaction forms */
766             /*     #endif */
767             /*     ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
768
769             /*     #if 'switch' in INTERACTION_FLAGS[I][J] */
770             d                = _mm_sub_ps(r{I}{J},rswitch);
771             d                = _mm_max_ps(d,_mm_setzero_ps());
772             d2               = _mm_mul_ps(d,d);
773             sw               = _mm_add_ps(one,_mm_mul_ps(d2,_mm_mul_ps(d,_mm_macc_ps(d,_mm_macc_ps(d,swV5,swV4),swV3))));
774             /*         #define INNERFLOPS INNERFLOPS+10 */
775
776             /*         #if 'Force' in KERNEL_VF */
777             dsw              = _mm_mul_ps(d2,_mm_macc_ps(d,_mm_macc_ps(d,swF4,swF3),swF2));
778             /*             #define INNERFLOPS INNERFLOPS+5 */
779             /*         #endif */
780
781             /* Evaluate switch function */
782             /*         #if 'Force' in KERNEL_VF */
783             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
784             /*             #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
785             felec            = _mm_msub_ps( felec,sw , _mm_mul_ps(rinv{I}{J},_mm_mul_ps(velec,dsw)) );
786             /*                 #define INNERFLOPS INNERFLOPS+4 */
787             /*             #endif */
788             /*             #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
789             fvdw             = _mm_msub_ps( fvdw,sw , _mm_mul_ps(rinv{I}{J},_mm_mul_ps(vvdw,dsw)) );
790             /*                 #define INNERFLOPS INNERFLOPS+4 */
791             /*             #endif */
792             /*         #endif */
793             /*         #if 'Potential' in KERNEL_VF */
794             /*             #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
795             velec            = _mm_mul_ps(velec,sw);
796             /*                 #define INNERFLOPS INNERFLOPS+1 */
797             /*             #endif */
798             /*             #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
799             vvdw             = _mm_mul_ps(vvdw,sw);
800             /*                 #define INNERFLOPS INNERFLOPS+1 */
801             /*             #endif */
802             /*         #endif */
803             /*     #endif */
804             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
805             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
806             cutoff_mask      = _mm_cmplt_ps(rsq{I}{J},rcutoff2);
807             /*         #define INNERFLOPS INNERFLOPS+1 */
808             /*     #endif */
809
810             /*     #if 'Potential' in KERNEL_VF */
811             /* Update potential sum for this i atom from the interaction with this j atom. */
812             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
813             /*             #if 'exactcutoff' in INTERACTION_FLAGS[I][J] */
814             velec            = _mm_and_ps(velec,cutoff_mask);
815             /*                 #define INNERFLOPS INNERFLOPS+1 */
816             /*             #endif                                       */
817             /*             #if ROUND == 'Epilogue' */
818             velec            = _mm_andnot_ps(dummy_mask,velec);
819             /*             #endif */
820             velecsum         = _mm_add_ps(velecsum,velec);
821             /*             #define INNERFLOPS INNERFLOPS+1 */
822             /*         #endif */
823             /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
824             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
825             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
826             vvdw             = _mm_and_ps(vvdw,cutoff_mask);
827             /*                 #define INNERFLOPS INNERFLOPS+1 */
828             /*             #endif                                       */
829             /*             #if ROUND == 'Epilogue' */
830             vvdw             = _mm_andnot_ps(dummy_mask,vvdw);
831             /*             #endif */
832             vvdwsum          = _mm_add_ps(vvdwsum,vvdw);
833             /*             #define INNERFLOPS INNERFLOPS+1 */
834             /*         #endif */
835             /*     #endif */
836
837             /*     #if 'Force' in KERNEL_VF */
838
839             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
840             fscal            = _mm_add_ps(felec,fvdw);
841             /*             #define INNERFLOPS INNERFLOPS+1 */
842             /*         #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
843             fscal            = felec;
844             /*         #elif 'vdw' in INTERACTION_FLAGS[I][J] */
845             fscal            = fvdw;
846             /*        #endif */
847
848             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
849             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
850             fscal            = _mm_and_ps(fscal,cutoff_mask);
851             /*                 #define INNERFLOPS INNERFLOPS+1 */
852             /*             #endif                                       */
853
854             /*             #if ROUND == 'Epilogue' */
855             fscal            = _mm_andnot_ps(dummy_mask,fscal);
856             /*             #endif */
857
858             /* ## Construction of vectorial force built into FMA instructions now */
859             /* #define INNERFLOPS INNERFLOPS+3      */
860
861              /* Update vectorial force */
862             fix{I}             = _mm_macc_ps(dx{I}{J},fscal,fix{I});
863             fiy{I}             = _mm_macc_ps(dy{I}{J},fscal,fiy{I});
864             fiz{I}             = _mm_macc_ps(dz{I}{J},fscal,fiz{I});
865             /*             #define INNERFLOPS INNERFLOPS+6 */
866
867             /* #if GEOMETRY_I == 'Particle'             */
868             /*     #if ROUND == 'Loop' */
869             fjptrA             = f+j_coord_offsetA;
870             fjptrB             = f+j_coord_offsetB;
871             fjptrC             = f+j_coord_offsetC;
872             fjptrD             = f+j_coord_offsetD;
873             /*     #else */
874             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
875             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
876             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
877             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
878             /*     #endif */
879             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
880                                                    _mm_mul_ps(dx{I}{J},fscal),
881                                                    _mm_mul_ps(dy{I}{J},fscal),
882                                                    _mm_mul_ps(dz{I}{J},fscal));
883             /*     #define INNERFLOPS INNERFLOPS+3      */
884             /* #else                                    */
885             fjx{J}             = _mm_macc_ps(dx{I}{J},fscal,fjx{J});
886             fjy{J}             = _mm_macc_ps(dy{I}{J},fscal,fjy{J});
887             fjz{J}             = _mm_macc_ps(dz{I}{J},fscal,fjz{J});
888             /*     #define INNERFLOPS INNERFLOPS+3      */
889             /* #endif                                   */
890
891             /*     #endif */
892
893             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
894             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
895             /*         #if 0    ## This and next two lines is a hack to maintain indentation in template file */
896             {
897                 /*     #endif */
898             }
899             /*     #endif */
900             /*    ## End of check for the interaction being outside the cutoff */
901
902             /* #endfor */
903             /* ## End of loop over i-j interaction pairs */
904
905             /* #if GEOMETRY_I != 'Particle' */
906             /*     #if ROUND == 'Loop' */
907             fjptrA             = f+j_coord_offsetA;
908             fjptrB             = f+j_coord_offsetB;
909             fjptrC             = f+j_coord_offsetC;
910             fjptrD             = f+j_coord_offsetD;
911             /*     #else */
912             fjptrA             = (jnrlistA>=0) ? f+j_coord_offsetA : scratch;
913             fjptrB             = (jnrlistB>=0) ? f+j_coord_offsetB : scratch;
914             fjptrC             = (jnrlistC>=0) ? f+j_coord_offsetC : scratch;
915             fjptrD             = (jnrlistD>=0) ? f+j_coord_offsetD : scratch;
916             /*     #endif */
917             /* #endif */
918
919             /* #if 'Water' in GEOMETRY_I and GEOMETRY_J == 'Particle' */
920             gmx_mm_decrement_1rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,fjx0,fjy0,fjz0);
921             /* #elif GEOMETRY_J == 'Water3'               */
922             gmx_mm_decrement_3rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
923                                                    fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,fjx2,fjy2,fjz2);
924             /*     #define INNERFLOPS INNERFLOPS+9      */
925             /* #elif GEOMETRY_J == 'Water4'             */
926             /*     #if 0 in PARTICLES_J                 */
927             gmx_mm_decrement_4rvec_4ptr_swizzle_ps(fjptrA,fjptrB,fjptrC,fjptrD,
928                                                    fjx0,fjy0,fjz0,fjx1,fjy1,fjz1,
929                                                    fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
930             /*     #define INNERFLOPS INNERFLOPS+12     */
931             /*     #else                                */
932             gmx_mm_decrement_3rvec_4ptr_swizzle_ps(fjptrA+DIM,fjptrB+DIM,fjptrC+DIM,fjptrD+DIM,
933                                                    fjx1,fjy1,fjz1,fjx2,fjy2,fjz2,fjx3,fjy3,fjz3);
934             /*     #define INNERFLOPS INNERFLOPS+9      */
935             /*     #endif                               */
936             /* #endif                                   */
937
938             /* Inner loop uses {INNERFLOPS} flops */
939         }
940
941         /* #endfor */
942
943         /* End of innermost loop */
944
945         /* #if 'Force' in KERNEL_VF */
946         /*     #if GEOMETRY_I == 'Particle'            */
947         gmx_mm_update_iforce_1atom_swizzle_ps(fix0,fiy0,fiz0,
948                                               f+i_coord_offset,fshift+i_shift_offset);
949         /*         #define OUTERFLOPS OUTERFLOPS+6     */
950         /*     #elif GEOMETRY_I == 'Water3'            */
951         gmx_mm_update_iforce_3atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,
952                                               f+i_coord_offset,fshift+i_shift_offset);
953         /*         #define OUTERFLOPS OUTERFLOPS+18    */
954         /*     #elif GEOMETRY_I == 'Water4'            */
955         /*         #if 0 in PARTICLES_I                */
956         gmx_mm_update_iforce_4atom_swizzle_ps(fix0,fiy0,fiz0,fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
957                                               f+i_coord_offset,fshift+i_shift_offset);
958         /*             #define OUTERFLOPS OUTERFLOPS+24    */
959         /*         #else                               */
960         gmx_mm_update_iforce_3atom_swizzle_ps(fix1,fiy1,fiz1,fix2,fiy2,fiz2,fix3,fiy3,fiz3,
961                                               f+i_coord_offset+DIM,fshift+i_shift_offset);
962         /*             #define OUTERFLOPS OUTERFLOPS+18    */
963         /*         #endif                              */
964         /*     #endif                                  */
965         /* #endif                                      */
966
967         /* #if 'Potential' in KERNEL_VF */
968         ggid                        = gid[iidx];
969         /* Update potential energies */
970         /*     #if KERNEL_ELEC != 'None' */
971         gmx_mm_update_1pot_ps(velecsum,kernel_data->energygrp_elec+ggid);
972         /*         #define OUTERFLOPS OUTERFLOPS+1 */
973         /*     #endif */
974         /*     #if KERNEL_VDW != 'None' */
975         gmx_mm_update_1pot_ps(vvdwsum,kernel_data->energygrp_vdw+ggid);
976         /*         #define OUTERFLOPS OUTERFLOPS+1 */
977         /*     #endif */
978         /* #endif */
979
980         /* Increment number of inner iterations */
981         inneriter                  += j_index_end - j_index_start;
982
983         /* Outer loop uses {OUTERFLOPS} flops */
984     }
985
986     /* Increment number of outer iterations */
987     outeriter        += nri;
988
989     /* Update outer/inner flops */
990     /* ## NB: This is not important, it just affects the flopcount. However, since our preprocessor is */
991     /* ## primitive and replaces aggressively even in strings inside these directives, we need to      */
992     /* ## assemble the main part of the name (containing KERNEL/ELEC/VDW) directly in the source.      */
993     /* #if GEOMETRY_I == 'Water3'            */
994     /*     #define ISUFFIX '_W3'             */
995     /* #elif GEOMETRY_I == 'Water4'          */
996     /*     #define ISUFFIX '_W4'             */
997     /* #else                                 */
998     /*     #define ISUFFIX ''                */
999     /* #endif                                */
1000     /* #if GEOMETRY_J == 'Water3'            */
1001     /*     #define JSUFFIX 'W3'              */
1002     /* #elif GEOMETRY_J == 'Water4'          */
1003     /*     #define JSUFFIX 'W4'              */
1004     /* #else                                 */
1005     /*     #define JSUFFIX ''                */
1006     /* #endif                                */
1007     /* #if 'PotentialAndForce' in KERNEL_VF  */
1008     /*     #define VFSUFFIX  '_VF'           */
1009     /* #elif 'Potential' in KERNEL_VF        */
1010     /*     #define VFSUFFIX '_V'             */
1011     /* #else                                 */
1012     /*     #define VFSUFFIX '_F'             */
1013     /* #endif                                */
1014
1015     /* #if KERNEL_ELEC != 'None' and KERNEL_VDW != 'None' */
1016     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1017     /* #elif KERNEL_ELEC != 'None' */
1018     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1019     /* #else */
1020     inc_nrnb(nrnb,eNR_NBKERNEL_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
1021     /* #endif  */
1022 }