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