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