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