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