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