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