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