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