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