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