Update copyright statements and change license to LGPL
[alexxy/gromacs.git] / src / gmxlib / nonbonded / nb_kernel_c / nb_kernel_template_c.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 /* #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 /* #if ('CubicSplineTable' in [KERNEL_ELEC,KERNEL_VDW]) or KERNEL_VF=='PotentialAndForce' */
85 /*     #define TABLE_POINT_SIZE 4 */
86 /* #else */
87 /*     #define TABLE_POINT_SIZE 2 */
88 /* #endif */
89
90 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
91 /*     #define TABLE_NINTERACTIONS 3 */
92 /*     #define TABLE_VDW_OFFSET TABLE_POINT_SIZE */
93 /* #elif 'Table' in KERNEL_ELEC */
94 /*     #define TABLE_NINTERACTIONS 1 */
95 /* #elif 'Table' in KERNEL_VDW */
96 /*     #define TABLE_NINTERACTIONS 2 */
97 /*     #define TABLE_VDW_OFFSET 0 */
98 /* #else */
99 /*     #define TABLE_NINTERACTIONS 0 */
100 /* #endif */
101
102 /* #if 'Buckingham' in KERNEL_VDW */
103 /*   #define NVDWPARAM 3 */
104 /* #else */
105 /*   #define NVDWPARAM 2 */
106 /* #endif */
107
108 /*
109  * GROMACS nonbonded kernel:   {KERNEL_NAME}
110  * Electrostatics interaction: {KERNEL_ELEC}
111  * VdW interaction:            {KERNEL_VDW}
112  * Geometry:                   {GEOMETRY_I}-{GEOMETRY_J}
113  * Calculate force/pot:        {KERNEL_VF}
114  */
115 void
116 {KERNEL_NAME}
117                     (t_nblist * gmx_restrict                nlist,
118                      rvec * gmx_restrict                    xx,
119                      rvec * gmx_restrict                    ff,
120                      t_forcerec * gmx_restrict              fr,
121                      t_mdatoms * gmx_restrict               mdatoms,
122                      nb_kernel_data_t * gmx_restrict        kernel_data,
123                      t_nrnb * gmx_restrict                  nrnb)
124 {
125     /* ## Not all variables are used for all kernels, but any optimizing compiler fixes that, */
126     /* ## so there is no point in going to extremes to exclude variables that are not needed. */
127     int              i_shift_offset,i_coord_offset,j_coord_offset;
128     int              j_index_start,j_index_end;
129     int              nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
130     real             shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
131     int              *iinr,*jindex,*jjnr,*shiftidx,*gid;
132     real             *shiftvec,*fshift,*x,*f;
133     /* #for I in PARTICLES_I */
134     int              vdwioffset{I};
135     real             ix{I},iy{I},iz{I},fix{I},fiy{I},fiz{I},iq{I},isai{I};
136     /* #endfor */
137     /* #for J in PARTICLES_J */
138     int              vdwjidx{J};
139     real             jx{J},jy{J},jz{J},fjx{J},fjy{J},fjz{J},jq{J},isaj{J};
140     /* #endfor */
141     /* #for I,J in PAIRS_IJ */
142     real             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},cexp1_{I}{J},cexp2_{I}{J};
143     /* #endfor */
144     /* #if KERNEL_ELEC != 'None' */
145     real             velec,felec,velecsum,facel,crf,krf,krf2;
146     real             *charge;
147     /* #endif */
148     /* #if 'GeneralizedBorn' in KERNEL_ELEC */
149     int              gbitab;
150     real             vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,dvdatmp;
151     real             *invsqrta,*dvda,*gbtab;
152     /* #endif */
153     /* #if KERNEL_VDW != 'None' */
154     int              nvdwtype;
155     real             rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
156     int              *vdwtype;
157     real             *vdwparam;
158     /* #endif */
159     /* #if 'Table' in KERNEL_ELEC or 'GeneralizedBorn' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
160     int              vfitab;
161     real             rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
162     real             *vftab;
163     /* #endif */
164     /* #if 'Ewald' in KERNEL_ELEC */
165     int              ewitab;
166     real             ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
167     real             *ewtab;
168     /* #endif */
169     /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
170     real             rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
171     /* #endif */
172
173     x                = xx[0];
174     f                = ff[0];
175
176     nri              = nlist->nri;
177     iinr             = nlist->iinr;
178     jindex           = nlist->jindex;
179     jjnr             = nlist->jjnr;
180     shiftidx         = nlist->shift;
181     gid              = nlist->gid;
182     shiftvec         = fr->shift_vec[0];
183     fshift           = fr->fshift[0];
184     /* #if KERNEL_ELEC != 'None' */
185     facel            = fr->epsfac;
186     charge           = mdatoms->chargeA;
187     /*     #if 'ReactionField' in KERNEL_ELEC */
188     krf              = fr->ic->k_rf;
189     krf2             = krf*2.0;
190     crf              = fr->ic->c_rf;
191     /*     #endif */
192     /* #endif */
193     /* #if KERNEL_VDW != 'None' */
194     nvdwtype         = fr->ntype;
195     vdwparam         = fr->nbfp;
196     vdwtype          = mdatoms->typeA;
197     /* #endif */
198
199     /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
200     vftab            = kernel_data->table_elec_vdw->data;
201     vftabscale       = kernel_data->table_elec_vdw->scale;
202     /* #elif 'Table' in KERNEL_ELEC */
203     vftab            = kernel_data->table_elec->data;
204     vftabscale       = kernel_data->table_elec->scale;
205     /* #elif 'Table' in KERNEL_VDW */
206     vftab            = kernel_data->table_vdw->data;
207     vftabscale       = kernel_data->table_vdw->scale;
208     /* #endif */
209
210     /* #if 'Ewald' in KERNEL_ELEC */
211     sh_ewald         = fr->ic->sh_ewald;
212     /*     #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
213     ewtab            = fr->ic->tabq_coul_F;
214     ewtabscale       = fr->ic->tabq_scale;
215     ewtabhalfspace   = 0.5/ewtabscale;
216     /*     #else */
217     ewtab            = fr->ic->tabq_coul_FDV0;
218     ewtabscale       = fr->ic->tabq_scale;
219     ewtabhalfspace   = 0.5/ewtabscale;
220      /*     #endif */
221     /* #endif */
222
223     /* #if KERNEL_ELEC=='GeneralizedBorn' */
224     invsqrta         = fr->invsqrta;
225     dvda             = fr->dvda;
226     gbtabscale       = fr->gbtab.scale;
227     gbtab            = fr->gbtab.data;
228     gbinvepsdiff     = (1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent);
229     /* #endif */
230
231     /* #if 'Water' in GEOMETRY_I */
232     /* Setup water-specific parameters */
233     inr              = nlist->iinr[0];
234     /*     #for I in PARTICLES_ELEC_I */
235     iq{I}              = facel*charge[inr+{I}];
236     /*     #endfor */
237     /*     #for I in PARTICLES_VDW_I */
238     vdwioffset{I}      = {NVDWPARAM}*nvdwtype*vdwtype[inr+{I}];
239     /*     #endfor */
240     /* #endif */
241
242     /* #if 'Water' in GEOMETRY_J */
243     /*     #for J in PARTICLES_ELEC_J */
244     jq{J}              = charge[inr+{J}];
245     /*     #endfor */
246     /*     #for J in PARTICLES_VDW_J */
247     vdwjidx{J}         = {NVDWPARAM}*vdwtype[inr+{J}];
248     /*     #endfor */
249     /*     #for I,J in PAIRS_IJ */
250     /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
251     qq{I}{J}             = iq{I}*jq{J};
252     /*         #endif */
253     /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
254     /*             #if 'Buckingham' in KERNEL_VDW */
255     c6_{I}{J}            = vdwparam[vdwioffset{I}+vdwjidx{J}];
256     cexp1_{I}{J}         = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
257     cexp2_{I}{J}         = vdwparam[vdwioffset{I}+vdwjidx{J}+2];
258     /*             #else */
259     c6_{I}{J}            = vdwparam[vdwioffset{I}+vdwjidx{J}];
260     c12_{I}{J}           = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
261     /*             #endif */
262     /*         #endif */
263     /*     #endfor */
264     /* #endif */
265
266     /* #if KERNEL_MOD_ELEC!='None' or KERNEL_MOD_VDW!='None' */
267     /*     #if KERNEL_ELEC!='None' */
268     /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
269     rcutoff          = fr->rcoulomb;
270     /*     #else */
271     rcutoff          = fr->rvdw;
272     /*     #endif */
273     rcutoff2         = rcutoff*rcutoff;
274     /* #endif */
275
276     /* #if KERNEL_MOD_VDW=='PotentialShift' */
277     sh_vdw_invrcut6  = fr->ic->sh_invrc6;
278     rvdw             = fr->rvdw;
279     /* #endif */
280
281     /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
282     /*     #if KERNEL_MOD_ELEC=='PotentialSwitch'  */
283     rswitch          = fr->rcoulomb_switch;
284     /*     #else */
285     rswitch          = fr->rvdw_switch;
286     /*     #endif */
287     /* Setup switch parameters */
288     d                = rcutoff-rswitch;
289     swV3             = -10.0/(d*d*d);
290     swV4             =  15.0/(d*d*d*d);
291     swV5             =  -6.0/(d*d*d*d*d);
292     /*     #if 'Force' in KERNEL_VF */
293     swF2             = -30.0/(d*d*d);
294     swF3             =  60.0/(d*d*d*d);
295     swF4             = -30.0/(d*d*d*d*d);
296     /*     #endif */
297     /* #endif */
298
299     /* ## Keep track of the floating point operations we issue for reporting! */
300     /* #define OUTERFLOPS 0 */
301     /* #define INNERFLOPS 0 */
302     outeriter        = 0;
303     inneriter        = 0;
304
305     /* Start outer loop over neighborlists */
306     for(iidx=0; iidx<nri; iidx++)
307     {
308         /* Load shift vector for this list */
309         i_shift_offset   = DIM*shiftidx[iidx];
310         shX              = shiftvec[i_shift_offset+XX];
311         shY              = shiftvec[i_shift_offset+YY];
312         shZ              = shiftvec[i_shift_offset+ZZ];
313
314         /* Load limits for loop over neighbors */
315         j_index_start    = jindex[iidx];
316         j_index_end      = jindex[iidx+1];
317
318         /* Get outer coordinate index */
319         inr              = iinr[iidx];
320         i_coord_offset   = DIM*inr;
321
322         /* Load i particle coords and add shift vector */
323         /* ## Loop over i particles, but only include ones that we use - skip e.g. vdw-only sites for elec-only kernel */
324         /*     #for I in PARTICLES_I */
325         ix{I}              = shX + x[i_coord_offset+DIM*{I}+XX];
326         iy{I}              = shY + x[i_coord_offset+DIM*{I}+YY];
327         iz{I}              = shZ + x[i_coord_offset+DIM*{I}+ZZ];
328         /*     #define OUTERFLOPS OUTERFLOPS+3 */
329         /* #endfor */
330
331         /* #if 'Force' in KERNEL_VF */
332         /*     #for I in PARTICLES_I */
333         fix{I}             = 0.0;
334         fiy{I}             = 0.0;
335         fiz{I}             = 0.0;
336         /*     #endfor */
337         /* #endif */
338
339         /* ## For water we already preloaded parameters at the start of the kernel */
340         /* #if not 'Water' in GEOMETRY_I */
341         /* Load parameters for i particles */
342         /*     #for I in PARTICLES_ELEC_I */
343         iq{I}              = facel*charge[inr+{I}];
344         /*         #define OUTERFLOPS OUTERFLOPS+1 */
345         /*         #if KERNEL_ELEC=='GeneralizedBorn' */
346         isai{I}            = invsqrta[inr+{I}];
347         /*         #endif */
348         /*     #endfor */
349         /*     #for I in PARTICLES_VDW_I */
350         vdwioffset{I}      = {NVDWPARAM}*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         = 0.0;
358         /*     #endif */
359         /*     #if 'GeneralizedBorn' in KERNEL_ELEC */
360         vgbsum           = 0.0;
361         /*     #endif */
362         /*     #if KERNEL_VDW != 'None' */
363         vvdwsum          = 0.0;
364         /*     #endif */
365         /* #endif */
366         /*     #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
367         dvdasum          = 0.0;
368         /*     #endif */
369
370         /* Start inner kernel loop */
371         for(jidx=j_index_start; jidx<j_index_end; jidx++)
372         {
373             /* Get j neighbor index, and coordinate index */
374             jnr              = jjnr[jidx];
375             j_coord_offset   = DIM*jnr;
376
377             /* load j atom coordinates */
378             /* #for J in PARTICLES_J */
379             jx{J}              = x[j_coord_offset+DIM*{J}+XX];
380             jy{J}              = x[j_coord_offset+DIM*{J}+YY];
381             jz{J}              = x[j_coord_offset+DIM*{J}+ZZ];
382             /* #endfor */
383
384             /* Calculate displacement vector */
385             /* #for I,J in PAIRS_IJ */
386             dx{I}{J}             = ix{I} - jx{J};
387             dy{I}{J}             = iy{I} - jy{J};
388             dz{I}{J}             = iz{I} - jz{J};
389             /*     #define INNERFLOPS INNERFLOPS+3 */
390             /* #endfor */
391
392             /* Calculate squared distance and things based on it */
393             /* #for I,J in PAIRS_IJ */
394             rsq{I}{J}            = dx{I}{J}*dx{I}{J}+dy{I}{J}*dy{I}{J}+dz{I}{J}*dz{I}{J};
395             /*     #define INNERFLOPS INNERFLOPS+5 */
396             /* #endfor */
397
398             /* #for I,J in PAIRS_IJ */
399             /*     #if 'rinv' in INTERACTION_FLAGS[I][J] */
400             rinv{I}{J}           = gmx_invsqrt(rsq{I}{J});
401             /*         #define INNERFLOPS INNERFLOPS+5 */
402             /*     #endif */
403             /* #endfor */
404
405             /* #for I,J in PAIRS_IJ */
406             /*     #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
407             /*         # if 'rinv' not in INTERACTION_FLAGS[I][J] */
408             rinvsq{I}{J}         = 1.0/rsq{I}{J};
409             /*             #define INNERFLOPS INNERFLOPS+4 */
410             /*         #else */
411             rinvsq{I}{J}         = rinv{I}{J}*rinv{I}{J};
412             /*             #define INNERFLOPS INNERFLOPS+1 */
413             /*         #endif */
414             /*     #endif */
415             /* #endfor */
416
417             /* #if not 'Water' in GEOMETRY_J */
418             /* Load parameters for j particles */
419             /*     #for J in PARTICLES_ELEC_J */
420             jq{J}              = charge[jnr+{J}];
421             /*         #if KERNEL_ELEC=='GeneralizedBorn' */
422             isaj{J}           = invsqrta[jnr+{J}];
423             /*         #endif */
424             /*     #endfor */
425             /*     #for J in PARTICLES_VDW_J */
426             vdwjidx{J}         = {NVDWPARAM}*vdwtype[jnr+{J}];
427             /*     #endfor */
428             /* #endif */
429
430             /* #for I,J in PAIRS_IJ */
431
432             /**************************
433              * CALCULATE INTERACTIONS *
434              **************************/
435
436             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
437             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
438             /*         ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
439             if (rsq{I}{J}<rcutoff2)
440             {
441                 /*     #if 0    ## this and the next two lines is a hack to maintain auto-indentation in template file */
442             }
443             /*         #endif */
444             /*     #endif */
445
446             /*     #if 'r' in INTERACTION_FLAGS[I][J] */
447             r{I}{J}              = rsq{I}{J}*rinv{I}{J};
448             /*         #define INNERFLOPS INNERFLOPS+1 */
449             /*     #endif */
450
451             /*     ## For water geometries we already loaded parameters at the start of the kernel */
452             /*     #if not 'Water' in GEOMETRY_J */
453             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
454             qq{I}{J}             = iq{I}*jq{J};
455             /*             #define INNERFLOPS INNERFLOPS+1 */
456             /*         #endif */
457             /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
458             /*             #if KERNEL_VDW=='Buckingham' */
459             c6_{I}{J}            = vdwparam[vdwioffset{I}+vdwjidx{J}];
460             cexp1_{I}{J}         = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
461             cexp2_{I}{J}         = vdwparam[vdwioffset{I}+vdwjidx{J}+2];
462             /*             #else */
463             c6_{I}{J}            = vdwparam[vdwioffset{I}+vdwjidx{J}];
464             c12_{I}{J}           = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
465             /*             #endif */
466             /*         #endif */
467             /*     #endif */
468
469             /*     #if 'table' in INTERACTION_FLAGS[I][J] */
470             /* Calculate table index by multiplying r with table scale and truncate to integer */
471             rt               = r{I}{J}*vftabscale;
472             vfitab           = rt;
473             vfeps            = rt-vfitab;
474             vfitab           = {TABLE_NINTERACTIONS}*{TABLE_POINT_SIZE}*vfitab;
475             /*         #define INNERFLOPS INNERFLOPS+2 */
476             /*     #endif */
477
478             /*     ## ELECTROSTATIC INTERACTIONS */
479             /*     #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
480
481             /*         #if KERNEL_ELEC=='Coulomb' */
482
483             /* COULOMB ELECTROSTATICS */
484             velec            = qq{I}{J}*rinv{I}{J};
485             /*             #define INNERFLOPS INNERFLOPS+1 */
486             /*             #if 'Force' in KERNEL_VF */
487             felec            = velec*rinvsq{I}{J};
488             /*                 #define INNERFLOPS INNERFLOPS+2 */
489             /*             #endif */
490
491             /*         #elif KERNEL_ELEC=='ReactionField' */
492
493             /* REACTION-FIELD ELECTROSTATICS */
494             /*             #if 'Potential' in KERNEL_VF */
495             velec            = qq{I}{J}*(rinv{I}{J}+krf*rsq{I}{J}-crf);
496             /*                 #define INNERFLOPS INNERFLOPS+4 */
497             /*             #endif */
498             /*             #if 'Force' in KERNEL_VF */
499             felec            = qq{I}{J}*(rinv{I}{J}*rinvsq{I}{J}-krf2);
500             /*                 #define INNERFLOPS INNERFLOPS+3 */
501             /*             #endif */
502
503             /*         #elif KERNEL_ELEC=='GeneralizedBorn' */
504
505             /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
506             isaprod          = isai{I}*isaj{J};
507             gbqqfactor       = isaprod*(-qq{I}{J})*gbinvepsdiff;
508             gbscale          = isaprod*gbtabscale;
509             dvdaj            = dvda[jnr+{J}];
510             /*             #define INNERFLOPS INNERFLOPS+5 */
511
512             /* Calculate generalized born table index - this is a separate table from the normal one,
513              * but we use the same procedure by multiplying r with scale and truncating to integer.
514              */
515             rt               = r{I}{J}*gbscale;
516             gbitab           = rt;
517             gbeps            = rt-gbitab;
518             gbitab           = 4*gbitab;
519
520             Y                = gbtab[gbitab];
521             F                = gbtab[gbitab+1];
522             Geps             = gbeps*gbtab[gbitab+2];
523             Heps2            = gbeps*gbeps*gbtab[gbitab+3];
524             Fp               = F+Geps+Heps2;
525             VV               = Y+gbeps*Fp;
526             vgb              = gbqqfactor*VV;
527             /*             #define INNERFLOPS INNERFLOPS+10 */
528
529             /*             #if 'Force' in KERNEL_VF */
530             FF               = Fp+Geps+2.0*Heps2;
531             fgb              = gbqqfactor*FF*gbscale;
532             dvdatmp          = -0.5*(vgb+fgb*r{I}{J});
533             dvdasum          = dvdasum + dvdatmp;
534             dvda[jnr]        = dvdaj+dvdatmp*isaj{J}*isaj{J};
535             /*                 #define INNERFLOPS INNERFLOPS+13 */
536             /*             #endif */
537             velec            = qq{I}{J}*rinv{I}{J};
538             /*                 #define INNERFLOPS INNERFLOPS+1 */
539             /*             #if 'Force' in KERNEL_VF */
540             felec            = (velec*rinv{I}{J}-fgb)*rinv{I}{J};
541             /*                 #define INNERFLOPS INNERFLOPS+3 */
542             /*             #endif */
543
544             /*         #elif KERNEL_ELEC=='Ewald' */
545             /* EWALD ELECTROSTATICS */
546
547             /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
548             ewrt             = r{I}{J}*ewtabscale;
549             ewitab           = ewrt;
550             eweps            = ewrt-ewitab;
551             /*             #define INNERFLOPS INNERFLOPS+2 */
552             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
553             ewitab           = 4*ewitab;
554             felec            = ewtab[ewitab]+eweps*ewtab[ewitab+1];
555             /*                 #define INNERFLOPS INNERFLOPS+4 */
556             /*                 #if KERNEL_MOD_ELEC=='PotentialShift' */
557             velec            = qq{I}{J}*((rinv{I}{J}-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
558             /*                     #define INNERFLOPS INNERFLOPS+7 */
559             /*                 #else */
560             velec            = qq{I}{J}*(rinv{I}{J}-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
561             /*                     #define INNERFLOPS INNERFLOPS+6 */
562             /*                 #endif */
563             /*                 #if 'Force' in KERNEL_VF */
564             felec            = qq{I}{J}*rinv{I}{J}*(rinvsq{I}{J}-felec);
565             /*                      #define INNERFLOPS INNERFLOPS+3 */
566             /*                 #endif */
567             /*             #elif KERNEL_VF=='Force' */
568             felec            = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
569             felec            = qq{I}{J}*rinv{I}{J}*(rinvsq{I}{J}-felec);
570             /*                 #define INNERFLOPS INNERFLOPS+7 */
571             /*             #endif */
572
573             /*         #elif KERNEL_ELEC=='CubicSplineTable' */
574
575             /* CUBIC SPLINE TABLE ELECTROSTATICS */
576             /*             #if 'Potential' in KERNEL_VF */
577             Y                = vftab[vfitab];
578             /*             #endif */
579             F                = vftab[vfitab+1];
580             Geps             = vfeps*vftab[vfitab+2];
581             Heps2            = vfeps*vfeps*vftab[vfitab+3];
582             Fp               = F+Geps+Heps2;
583             /*             #define INNERFLOPS INNERFLOPS+5 */
584             /*             #if 'Potential' in KERNEL_VF */
585             VV               = Y+vfeps*Fp;
586             velec            = qq{I}{J}*VV;
587             /*                 #define INNERFLOPS INNERFLOPS+3 */
588             /*             #endif */
589             /*             #if 'Force' in KERNEL_VF */
590             FF               = Fp+Geps+2.0*Heps2;
591             felec            = -qq{I}{J}*FF*vftabscale*rinv{I}{J};
592             /*                 #define INNERFLOPS INNERFLOPS+7 */
593             /*             #endif */
594             /*         #endif */
595             /*         ## End of check for electrostatics interaction forms */
596             /*     #endif */
597             /*     ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
598
599             /*     #if 'vdw' in INTERACTION_FLAGS[I][J] */
600
601             /*         #if KERNEL_VDW=='LennardJones' */
602
603             /* LENNARD-JONES DISPERSION/REPULSION */
604
605             rinvsix          = rinvsq{I}{J}*rinvsq{I}{J}*rinvsq{I}{J};
606             /*             #define INNERFLOPS INNERFLOPS+2 */
607             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
608             vvdw6            = c6_{I}{J}*rinvsix;
609             vvdw12           = c12_{I}{J}*rinvsix*rinvsix;
610             /*                 #define INNERFLOPS INNERFLOPS+3 */
611             /*                 #if KERNEL_MOD_VDW=='PotentialShift' */
612             vvdw             = (vvdw12 - c12_{I}{J}*sh_vdw_invrcut6*sh_vdw_invrcut6)*(1.0/12.0) - (vvdw6 - c6_{I}{J}*sh_vdw_invrcut6)*(1.0/6.0);
613             /*                     #define INNERFLOPS INNERFLOPS+8 */
614             /*                 #else */
615             vvdw             = vvdw12*(1.0/12.0) - vvdw6*(1.0/6.0);
616             /*                     #define INNERFLOPS INNERFLOPS+3 */
617             /*                 #endif */
618             /*                 ## Check for force inside potential check, i.e. this means we already did the potential part */
619             /*                 #if 'Force' in KERNEL_VF */
620             fvdw             = (vvdw12-vvdw6)*rinvsq{I}{J};
621             /*                     #define INNERFLOPS INNERFLOPS+2 */
622             /*                 #endif */
623             /*             #elif KERNEL_VF=='Force' */
624             /*                 ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
625             fvdw             = (c12_{I}{J}*rinvsix-c6_{I}{J})*rinvsix*rinvsq{I}{J};
626             /*                 #define INNERFLOPS INNERFLOPS+4 */
627             /*             #endif */
628
629             /*         #elif KERNEL_VDW=='Buckingham' */
630
631             /* BUCKINGHAM DISPERSION/REPULSION */
632             rinvsix          = rinvsq{I}{J}*rinvsq{I}{J}*rinvsq{I}{J};
633             vvdw6            = c6_{I}{J}*rinvsix;
634             br               = cexp2_{I}{J}*r{I}{J};
635             vvdwexp          = cexp1_{I}{J}*exp(-br);
636             /*             ## Estimate exp() to 25 flops */
637             /*             #define INNERFLOPS INNERFLOPS+31 */
638             /*             #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch'  */
639             /*                 #if KERNEL_MOD_VDW=='PotentialShift' */
640             vvdw             = (vvdwexp-cexp1_{I}{J}*exp(-cexp2_{I}{J}*rvdw)) - (vvdw6 - c6_{I}{J}*sh_vdw_invrcut6)*(1.0/6.0);
641             /*                     #define INNERFLOPS INNERFLOPS+33 */
642             /*                  #else */
643             vvdw             = vvdwexp - vvdw6*(1.0/6.0);
644             /*                     #define INNERFLOPS INNERFLOPS+2 */
645             /*                 #endif */
646             /*             #endif */
647             /*             #if 'Force' in KERNEL_VF */
648             fvdw             = (br*vvdwexp-vvdw6)*rinvsq{I}{J};
649             /*                 #define INNERFLOPS INNERFLOPS+3 */
650             /*             #endif */
651
652             /*         #elif KERNEL_VDW=='CubicSplineTable' */
653
654             /* CUBIC SPLINE TABLE DISPERSION */
655             vfitab          += {TABLE_VDW_OFFSET};
656             /*             #if 'Potential' in KERNEL_VF */
657             Y                = vftab[vfitab];
658             /*             #endif */
659             F                = vftab[vfitab+1];
660             Geps             = vfeps*vftab[vfitab+2];
661             Heps2            = vfeps*vfeps*vftab[vfitab+3];
662             Fp               = F+Geps+Heps2;
663             /*             #define INNERFLOPS INNERFLOPS+5 */
664             /*             #if 'Potential' in KERNEL_VF */
665             VV               = Y+vfeps*Fp;
666             vvdw6            = c6_{I}{J}*VV;
667             /*                 #define INNERFLOPS INNERFLOPS+3 */
668             /*             #endif */
669             /*             #if 'Force' in KERNEL_VF */
670             FF               = Fp+Geps+2.0*Heps2;
671             fvdw6            = c6_{I}{J}*FF;
672             /*                 #define INNERFLOPS INNERFLOPS+4 */
673             /*             #endif */
674
675             /* CUBIC SPLINE TABLE REPULSION */
676             /*             #if 'Potential' in KERNEL_VF */
677             Y                = vftab[vfitab+4];
678             /*             #endif */
679             F                = vftab[vfitab+5];
680             Geps             = vfeps*vftab[vfitab+6];
681             Heps2            = vfeps*vfeps*vftab[vfitab+7];
682             Fp               = F+Geps+Heps2;
683             /*             #define INNERFLOPS INNERFLOPS+5 */
684             /*             #if 'Potential' in KERNEL_VF */
685             VV               = Y+vfeps*Fp;
686             vvdw12           = c12_{I}{J}*VV;
687             /*                 #define INNERFLOPS INNERFLOPS+3 */
688             /*             #endif */
689             /*             #if 'Force' in KERNEL_VF */
690             FF               = Fp+Geps+2.0*Heps2;
691             fvdw12           = c12_{I}{J}*FF;
692             /*                 #define INNERFLOPS INNERFLOPS+4 */
693             /*             #endif */
694             /*             #if 'Potential' in KERNEL_VF */
695             vvdw             = vvdw12+vvdw6;
696             /*                 #define INNERFLOPS INNERFLOPS+1 */
697             /*             #endif */
698             /*             #if 'Force' in KERNEL_VF */
699             fvdw             = -(fvdw6+fvdw12)*vftabscale*rinv{I}{J};
700             /*                 #define INNERFLOPS INNERFLOPS+4 */
701             /*             #endif */
702             /*         #endif */
703             /*         ## End of check for vdw interaction forms */
704             /*     #endif */
705             /*     ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
706
707             /*     #if 'switch' in INTERACTION_FLAGS[I][J] */
708             d                = r{I}{J}-rswitch;
709             d                = (d>0.0) ? d : 0.0;
710             d2               = d*d;
711             sw               = 1.0+d2*d*(swV3+d*(swV4+d*swV5));
712             /*         #define INNERFLOPS INNERFLOPS+9 */
713
714             /*         #if 'Force' in KERNEL_VF */
715             dsw              = d2*(swF2+d*(swF3+d*swF4));
716             /*             #define INNERFLOPS INNERFLOPS+5 */
717             /*         #endif */
718
719             /* Evaluate switch function */
720             /*         #if 'Force' in KERNEL_VF */
721             /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
722             /*             #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
723             felec            = felec*sw - rinv{I}{J}*velec*dsw;
724             /*                 #define INNERFLOPS INNERFLOPS+3 */
725             /*             #endif */
726             /*             #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
727             fvdw             = fvdw*sw - rinv{I}{J}*vvdw*dsw;
728             /*                 #define INNERFLOPS INNERFLOPS+3 */
729             /*             #endif */
730             /*         #endif */
731             /*         #if 'Potential' in KERNEL_VF */
732             /*             #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
733             velec           *= sw;
734             /*                 #define INNERFLOPS INNERFLOPS+1 */
735             /*             #endif */
736             /*             #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
737             vvdw            *= sw;
738             /*                 #define INNERFLOPS INNERFLOPS+1 */
739             /*             #endif */
740             /*         #endif */
741             /*     #endif */
742
743             /*     #if 'Potential' in KERNEL_VF */
744             /* Update potential sums from outer loop */
745             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
746             velecsum        += velec;
747             /*             #define INNERFLOPS INNERFLOPS+1 */
748             /*             #if KERNEL_ELEC=='GeneralizedBorn' */
749             vgbsum          += vgb;
750             /*                 #define INNERFLOPS INNERFLOPS+1 */
751             /*             #endif */
752             /*         #endif */
753             /*         #if 'vdw' in INTERACTION_FLAGS[I][J] */
754             vvdwsum         += vvdw;
755             /*             #define INNERFLOPS INNERFLOPS+1 */
756             /*         #endif */
757             /*     #endif */
758
759             /*     #if 'Force' in KERNEL_VF */
760
761             /*         #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
762             fscal            = felec+fvdw;
763             /*             #define INNERFLOPS INNERFLOPS+1 */
764             /*         #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
765             fscal            = felec;
766             /*         #elif 'vdw' in INTERACTION_FLAGS[I][J] */
767             fscal            = fvdw;
768             /*        #endif */
769
770             /* Calculate temporary vectorial force */
771             tx               = fscal*dx{I}{J};
772             ty               = fscal*dy{I}{J};
773             tz               = fscal*dz{I}{J};
774
775             /* Update vectorial force */
776             fix{I}            += tx;
777             fiy{I}            += ty;
778             fiz{I}            += tz;
779             f[j_coord_offset+DIM*{J}+XX] -= tx;
780             f[j_coord_offset+DIM*{J}+YY] -= ty;
781             f[j_coord_offset+DIM*{J}+ZZ] -= tz;
782
783             /*         #define INNERFLOPS INNERFLOPS+9 */
784             /*     #endif */
785
786             /*     ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
787             /*     #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
788             /*         #if 0    ## This and next two lines is a hack to maintain indentation in template file */
789             {
790                 /*     #endif */
791             }
792             /*     #endif */
793             /*    ## End of check for the interaction being outside the cutoff */
794
795             /* #endfor */
796             /* ## End of loop over i-j interaction pairs */
797
798             /* Inner loop uses {INNERFLOPS} flops */
799         }
800         /* End of innermost loop */
801
802         /* #if 'Force' in KERNEL_VF */
803         tx = ty = tz = 0;
804         /*     #for I in PARTICLES_I */
805         f[i_coord_offset+DIM*{I}+XX] += fix{I};
806         f[i_coord_offset+DIM*{I}+YY] += fiy{I};
807         f[i_coord_offset+DIM*{I}+ZZ] += fiz{I};
808         tx                         += fix{I};
809         ty                         += fiy{I};
810         tz                         += fiz{I};
811         /*         #define OUTERFLOPS OUTERFLOPS+6 */
812         /*     #endfor */
813         fshift[i_shift_offset+XX]  += tx;
814         fshift[i_shift_offset+YY]  += ty;
815         fshift[i_shift_offset+ZZ]  += tz;
816         /*     #define OUTERFLOPS OUTERFLOPS+3 */
817         /* #endif */
818
819         /* #if 'Potential' in KERNEL_VF */
820         ggid                        = gid[iidx];
821         /* Update potential energies */
822         /*     #if KERNEL_ELEC != 'None' */
823         kernel_data->energygrp_elec[ggid] += velecsum;
824         /*         #define OUTERFLOPS OUTERFLOPS+1 */
825         /*     #endif */
826         /*     #if 'GeneralizedBorn' in KERNEL_ELEC */
827         kernel_data->energygrp_polarization[ggid] += vgbsum;
828         /*         #define OUTERFLOPS OUTERFLOPS+1 */
829         /*     #endif */
830         /*     #if KERNEL_VDW != 'None' */
831         kernel_data->energygrp_vdw[ggid] += vvdwsum;
832         /*         #define OUTERFLOPS OUTERFLOPS+1 */
833         /*     #endif */
834         /* #endif */
835         /*     #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
836         dvda[inr]                   = dvda[inr] + dvdasum*isai{I}*isai{I};
837         /*     #endif */
838
839         /* Increment number of inner iterations */
840         inneriter                  += j_index_end - j_index_start;
841
842         /* Outer loop uses {OUTERFLOPS} flops */
843     }
844
845     /* Increment number of outer iterations */
846     outeriter        += nri;
847
848     /* Update outer/inner flops */
849     /* ## NB: This is not important, it just affects the flopcount. However, since our preprocessor is */
850     /* ## primitive and replaces aggressively even in strings inside these directives, we need to      */
851     /* ## assemble the main part of the name (containing KERNEL/ELEC/VDW) directly in the source.      */
852     /* #if GEOMETRY_I == 'Water3'            */
853     /*     #define ISUFFIX '_W3'             */
854     /* #elif GEOMETRY_I == 'Water4'          */
855     /*     #define ISUFFIX '_W4'             */
856     /* #else                                 */
857     /*     #define ISUFFIX ''                */
858     /* #endif                                */
859     /* #if GEOMETRY_J == 'Water3'            */
860     /*     #define JSUFFIX 'W3'              */
861     /* #elif GEOMETRY_J == 'Water4'          */
862     /*     #define JSUFFIX 'W4'              */
863     /* #else                                 */
864     /*     #define JSUFFIX ''                */
865     /* #endif                                */
866     /* #if 'PotentialAndForce' in KERNEL_VF  */
867     /*     #define VFSUFFIX  '_VF'           */
868     /* #elif 'Potential' in KERNEL_VF        */
869     /*     #define VFSUFFIX '_V'             */
870     /* #else                                 */
871     /*     #define VFSUFFIX '_F'             */
872     /* #endif                                */
873
874     /* #if KERNEL_ELEC != 'None' and KERNEL_VDW != 'None' */
875     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
876     /* #elif KERNEL_ELEC != 'None' */
877     inc_nrnb(nrnb,eNR_NBKERNEL_ELEC{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
878     /* #else */
879     inc_nrnb(nrnb,eNR_NBKERNEL_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
880     /* #endif  */
881 }