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