3 * This file is part of the GROMACS molecular simulation package.
5 * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
36 #error This file must be processed with the Gromacs pre-preprocessor
38 /* #if INCLUDE_HEADER */
43 #include "../nb_kernel.h"
44 #include "types/simple.h"
45 #include "gromacs/math/vec.h"
49 /* ## List of variables set by the generating script: */
51 /* ## Setttings that apply to the entire kernel: */
52 /* ## KERNEL_ELEC: String, choice for electrostatic interactions */
53 /* ## KERNEL_VDW: String, choice for van der Waals interactions */
54 /* ## KERNEL_NAME: String, name of this kernel */
55 /* ## KERNEL_VF: String telling if we calculate potential, force, or both */
56 /* ## GEOMETRY_I/GEOMETRY_J: String, name of each geometry, e.g. 'Water3' or '1Particle' */
58 /* ## Setttings that apply to particles in the outer (I) or inner (J) loops: */
59 /* ## PARTICLES_I[]/ Arrays with lists of i/j particles to use in kernel. It is */
60 /* ## PARTICLES_J[]: just [0] for particle geometry, but can be longer for water */
61 /* ## PARTICLES_ELEC_I[]/ Arrays with lists of i/j particle that have electrostatics */
62 /* ## PARTICLES_ELEC_J[]: interactions that should be calculated in this kernel. */
63 /* ## PARTICLES_VDW_I[]/ Arrays with the list of i/j particle that have VdW */
64 /* ## PARTICLES_VDW_J[]: interactions that should be calculated in this kernel. */
66 /* ## Setttings for pairs of interactions (e.g. 2nd i particle against 1st j particle) */
67 /* ## PAIRS_IJ[]: Array with (i,j) tuples of pairs for which interactions */
68 /* ## should be calculated in this kernel. Zero-charge particles */
69 /* ## do not have interactions with particles without vdw, and */
70 /* ## Vdw-only interactions are not evaluated in a no-vdw-kernel. */
71 /* ## INTERACTION_FLAGS[][]: 2D matrix, dimension e.g. 3*3 for water-water interactions. */
72 /* ## For each i-j pair, the element [I][J] is a list of strings */
73 /* ## defining properties/flags of this interaction. Examples */
74 /* ## include 'electrostatics'/'vdw' if that type of interaction */
75 /* ## should be evaluated, 'rsq'/'rinv'/'rinvsq' if those values */
76 /* ## are needed, and 'exactcutoff' or 'shift','switch' to */
77 /* ## decide if the force/potential should be modified. This way */
78 /* ## we only calculate values absolutely needed for each case. */
80 /* ## Calculate the size and offset for (merged/interleaved) table data */
82 /* #if ('CubicSplineTable' in [KERNEL_ELEC,KERNEL_VDW]) or KERNEL_VF=='PotentialAndForce' */
83 /* #define TABLE_POINT_SIZE 4 */
85 /* #define TABLE_POINT_SIZE 2 */
88 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
89 /* #define TABLE_NINTERACTIONS 3 */
90 /* #define TABLE_VDW_OFFSET TABLE_POINT_SIZE */
91 /* #elif 'Table' in KERNEL_ELEC */
92 /* #define TABLE_NINTERACTIONS 1 */
93 /* #elif 'Table' in KERNEL_VDW */
94 /* #define TABLE_NINTERACTIONS 2 */
95 /* #define TABLE_VDW_OFFSET 0 */
97 /* #define TABLE_NINTERACTIONS 0 */
100 /* #if 'Buckingham' in KERNEL_VDW */
101 /* #define NVDWPARAM 3 */
103 /* #define NVDWPARAM 2 */
107 * Gromacs nonbonded kernel: {KERNEL_NAME}
108 * Electrostatics interaction: {KERNEL_ELEC}
109 * VdW interaction: {KERNEL_VDW}
110 * Geometry: {GEOMETRY_I}-{GEOMETRY_J}
111 * Calculate force/pot: {KERNEL_VF}
115 (t_nblist * gmx_restrict nlist,
116 rvec * gmx_restrict xx,
117 rvec * gmx_restrict ff,
118 t_forcerec * gmx_restrict fr,
119 t_mdatoms * gmx_restrict mdatoms,
120 nb_kernel_data_t gmx_unused * gmx_restrict kernel_data,
121 t_nrnb * gmx_restrict nrnb)
123 /* ## Not all variables are used for all kernels, but any optimizing compiler fixes that, */
124 /* ## so there is no point in going to extremes to exclude variables that are not needed. */
125 int i_shift_offset,i_coord_offset,j_coord_offset;
126 int j_index_start,j_index_end;
127 int nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
128 real shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
129 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
130 real *shiftvec,*fshift,*x,*f;
131 /* #for I in PARTICLES_I */
133 real ix{I},iy{I},iz{I},fix{I},fiy{I},fiz{I},iq{I},isai{I};
135 /* #for J in PARTICLES_J */
137 real jx{J},jy{J},jz{J},fjx{J},fjy{J},fjz{J},jq{J},isaj{J};
139 /* #for I,J in PAIRS_IJ */
140 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};
142 /* #if KERNEL_ELEC != 'None' */
143 real velec,felec,velecsum,facel,crf,krf,krf2;
146 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
148 real vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,dvdatmp;
149 real *invsqrta,*dvda,*gbtab;
151 /* #if KERNEL_VDW != 'None' */
153 real rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
157 /* #if 'Table' in KERNEL_ELEC or 'GeneralizedBorn' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
159 real rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
162 /* #if 'LJEwald' in KERNEL_VDW */
163 /* #for I,J in PAIRS_IJ */
166 real ewclj,ewclj2,ewclj6,ewcljrsq,poly,exponent,sh_lj_ewald;
169 /* #if 'Ewald' in KERNEL_ELEC */
171 real ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
174 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
175 real rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
183 jindex = nlist->jindex;
185 shiftidx = nlist->shift;
187 shiftvec = fr->shift_vec[0];
188 fshift = fr->fshift[0];
189 /* #if KERNEL_ELEC != 'None' */
191 charge = mdatoms->chargeA;
192 /* #if 'ReactionField' in KERNEL_ELEC */
198 /* #if KERNEL_VDW != 'None' */
199 nvdwtype = fr->ntype;
201 vdwtype = mdatoms->typeA;
203 /* #if 'LJEwald' in KERNEL_VDW */
204 vdwgridparam = fr->ljpme_c6grid;
205 ewclj = fr->ewaldcoeff_lj;
206 sh_lj_ewald = fr->ic->sh_lj_ewald;
207 ewclj2 = ewclj*ewclj;
208 ewclj6 = ewclj2*ewclj2*ewclj2;
211 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
212 vftab = kernel_data->table_elec_vdw->data;
213 vftabscale = kernel_data->table_elec_vdw->scale;
214 /* #elif 'Table' in KERNEL_ELEC */
215 vftab = kernel_data->table_elec->data;
216 vftabscale = kernel_data->table_elec->scale;
217 /* #elif 'Table' in KERNEL_VDW */
218 vftab = kernel_data->table_vdw->data;
219 vftabscale = kernel_data->table_vdw->scale;
222 /* #if 'Ewald' in KERNEL_ELEC */
223 sh_ewald = fr->ic->sh_ewald;
224 /* #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
225 ewtab = fr->ic->tabq_coul_F;
226 ewtabscale = fr->ic->tabq_scale;
227 ewtabhalfspace = 0.5/ewtabscale;
229 ewtab = fr->ic->tabq_coul_FDV0;
230 ewtabscale = fr->ic->tabq_scale;
231 ewtabhalfspace = 0.5/ewtabscale;
235 /* #if KERNEL_ELEC=='GeneralizedBorn' */
236 invsqrta = fr->invsqrta;
238 gbtabscale = fr->gbtab.scale;
239 gbtab = fr->gbtab.data;
240 gbinvepsdiff = (1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent);
243 /* #if 'Water' in GEOMETRY_I */
244 /* Setup water-specific parameters */
245 inr = nlist->iinr[0];
246 /* #for I in PARTICLES_ELEC_I */
247 iq{I} = facel*charge[inr+{I}];
249 /* #for I in PARTICLES_VDW_I */
250 vdwioffset{I} = {NVDWPARAM}*nvdwtype*vdwtype[inr+{I}];
254 /* #if 'Water' in GEOMETRY_J */
255 /* #for J in PARTICLES_ELEC_J */
256 jq{J} = charge[inr+{J}];
258 /* #for J in PARTICLES_VDW_J */
259 vdwjidx{J} = {NVDWPARAM}*vdwtype[inr+{J}];
261 /* #for I,J in PAIRS_IJ */
262 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
263 qq{I}{J} = iq{I}*jq{J};
265 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
266 /* #if 'Buckingham' in KERNEL_VDW */
267 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
268 cexp1_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
269 cexp2_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+2];
270 /* #elif 'LJEwald' in KERNEL_VDW */
271 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
272 c12_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
273 c6grid_{I}{J} = vdwgridparam[vdwioffset{I}+vdwjidx{J}];
275 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
276 c12_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
282 /* #if KERNEL_MOD_ELEC!='None' or KERNEL_MOD_VDW!='None' */
283 /* #if KERNEL_ELEC!='None' */
284 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
285 rcutoff = fr->rcoulomb;
289 rcutoff2 = rcutoff*rcutoff;
292 /* #if KERNEL_MOD_VDW=='PotentialShift' */
293 sh_vdw_invrcut6 = fr->ic->sh_invrc6;
297 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
298 /* #if KERNEL_MOD_ELEC=='PotentialSwitch' */
299 rswitch = fr->rcoulomb_switch;
301 rswitch = fr->rvdw_switch;
303 /* Setup switch parameters */
305 swV3 = -10.0/(d*d*d);
306 swV4 = 15.0/(d*d*d*d);
307 swV5 = -6.0/(d*d*d*d*d);
308 /* #if 'Force' in KERNEL_VF */
309 swF2 = -30.0/(d*d*d);
310 swF3 = 60.0/(d*d*d*d);
311 swF4 = -30.0/(d*d*d*d*d);
315 /* ## Keep track of the floating point operations we issue for reporting! */
316 /* #define OUTERFLOPS 0 */
317 /* #define INNERFLOPS 0 */
321 /* Start outer loop over neighborlists */
322 for(iidx=0; iidx<nri; iidx++)
324 /* Load shift vector for this list */
325 i_shift_offset = DIM*shiftidx[iidx];
326 shX = shiftvec[i_shift_offset+XX];
327 shY = shiftvec[i_shift_offset+YY];
328 shZ = shiftvec[i_shift_offset+ZZ];
330 /* Load limits for loop over neighbors */
331 j_index_start = jindex[iidx];
332 j_index_end = jindex[iidx+1];
334 /* Get outer coordinate index */
336 i_coord_offset = DIM*inr;
338 /* Load i particle coords and add shift vector */
339 /* ## Loop over i particles, but only include ones that we use - skip e.g. vdw-only sites for elec-only kernel */
340 /* #for I in PARTICLES_I */
341 ix{I} = shX + x[i_coord_offset+DIM*{I}+XX];
342 iy{I} = shY + x[i_coord_offset+DIM*{I}+YY];
343 iz{I} = shZ + x[i_coord_offset+DIM*{I}+ZZ];
344 /* #define OUTERFLOPS OUTERFLOPS+3 */
347 /* #if 'Force' in KERNEL_VF */
348 /* #for I in PARTICLES_I */
355 /* ## For water we already preloaded parameters at the start of the kernel */
356 /* #if not 'Water' in GEOMETRY_I */
357 /* Load parameters for i particles */
358 /* #for I in PARTICLES_ELEC_I */
359 iq{I} = facel*charge[inr+{I}];
360 /* #define OUTERFLOPS OUTERFLOPS+1 */
361 /* #if KERNEL_ELEC=='GeneralizedBorn' */
362 isai{I} = invsqrta[inr+{I}];
365 /* #for I in PARTICLES_VDW_I */
366 vdwioffset{I} = {NVDWPARAM}*nvdwtype*vdwtype[inr+{I}];
370 /* #if 'Potential' in KERNEL_VF */
371 /* Reset potential sums */
372 /* #if KERNEL_ELEC != 'None' */
375 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
378 /* #if KERNEL_VDW != 'None' */
382 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
386 /* Start inner kernel loop */
387 for(jidx=j_index_start; jidx<j_index_end; jidx++)
389 /* Get j neighbor index, and coordinate index */
391 j_coord_offset = DIM*jnr;
393 /* load j atom coordinates */
394 /* #for J in PARTICLES_J */
395 jx{J} = x[j_coord_offset+DIM*{J}+XX];
396 jy{J} = x[j_coord_offset+DIM*{J}+YY];
397 jz{J} = x[j_coord_offset+DIM*{J}+ZZ];
400 /* Calculate displacement vector */
401 /* #for I,J in PAIRS_IJ */
402 dx{I}{J} = ix{I} - jx{J};
403 dy{I}{J} = iy{I} - jy{J};
404 dz{I}{J} = iz{I} - jz{J};
405 /* #define INNERFLOPS INNERFLOPS+3 */
408 /* Calculate squared distance and things based on it */
409 /* #for I,J in PAIRS_IJ */
410 rsq{I}{J} = dx{I}{J}*dx{I}{J}+dy{I}{J}*dy{I}{J}+dz{I}{J}*dz{I}{J};
411 /* #define INNERFLOPS INNERFLOPS+5 */
414 /* #for I,J in PAIRS_IJ */
415 /* #if 'rinv' in INTERACTION_FLAGS[I][J] */
416 rinv{I}{J} = gmx_invsqrt(rsq{I}{J});
417 /* #define INNERFLOPS INNERFLOPS+5 */
421 /* #for I,J in PAIRS_IJ */
422 /* #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
423 /* # if 'rinv' not in INTERACTION_FLAGS[I][J] */
424 rinvsq{I}{J} = 1.0/rsq{I}{J};
425 /* #define INNERFLOPS INNERFLOPS+4 */
427 rinvsq{I}{J} = rinv{I}{J}*rinv{I}{J};
428 /* #define INNERFLOPS INNERFLOPS+1 */
433 /* #if not 'Water' in GEOMETRY_J */
434 /* Load parameters for j particles */
435 /* #for J in PARTICLES_ELEC_J */
436 jq{J} = charge[jnr+{J}];
437 /* #if KERNEL_ELEC=='GeneralizedBorn' */
438 isaj{J} = invsqrta[jnr+{J}];
441 /* #for J in PARTICLES_VDW_J */
442 vdwjidx{J} = {NVDWPARAM}*vdwtype[jnr+{J}];
446 /* #for I,J in PAIRS_IJ */
448 /**************************
449 * CALCULATE INTERACTIONS *
450 **************************/
452 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
453 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
454 /* ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
455 if (rsq{I}{J}<rcutoff2)
457 /* #if 0 ## this and the next two lines is a hack to maintain auto-indentation in template file */
462 /* #if 'r' in INTERACTION_FLAGS[I][J] */
463 r{I}{J} = rsq{I}{J}*rinv{I}{J};
464 /* #define INNERFLOPS INNERFLOPS+1 */
467 /* ## For water geometries we already loaded parameters at the start of the kernel */
468 /* #if not 'Water' in GEOMETRY_J */
469 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
470 qq{I}{J} = iq{I}*jq{J};
471 /* #define INNERFLOPS INNERFLOPS+1 */
473 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
474 /* #if KERNEL_VDW=='Buckingham' */
475 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
476 cexp1_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
477 cexp2_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+2];
478 /* #elif 'LJEwald' in KERNEL_VDW */
479 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
480 c12_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
481 c6grid_{I}{J} = vdwgridparam[vdwioffset{I}+vdwjidx{J}];
483 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
484 c12_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
489 /* #if 'table' in INTERACTION_FLAGS[I][J] */
490 /* Calculate table index by multiplying r with table scale and truncate to integer */
491 rt = r{I}{J}*vftabscale;
494 vfitab = {TABLE_NINTERACTIONS}*{TABLE_POINT_SIZE}*vfitab;
495 /* #define INNERFLOPS INNERFLOPS+2 */
498 /* ## ELECTROSTATIC INTERACTIONS */
499 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
501 /* #if KERNEL_ELEC=='Coulomb' */
503 /* COULOMB ELECTROSTATICS */
504 velec = qq{I}{J}*rinv{I}{J};
505 /* #define INNERFLOPS INNERFLOPS+1 */
506 /* #if 'Force' in KERNEL_VF */
507 felec = velec*rinvsq{I}{J};
508 /* #define INNERFLOPS INNERFLOPS+2 */
511 /* #elif KERNEL_ELEC=='ReactionField' */
513 /* REACTION-FIELD ELECTROSTATICS */
514 /* #if 'Potential' in KERNEL_VF */
515 velec = qq{I}{J}*(rinv{I}{J}+krf*rsq{I}{J}-crf);
516 /* #define INNERFLOPS INNERFLOPS+4 */
518 /* #if 'Force' in KERNEL_VF */
519 felec = qq{I}{J}*(rinv{I}{J}*rinvsq{I}{J}-krf2);
520 /* #define INNERFLOPS INNERFLOPS+3 */
523 /* #elif KERNEL_ELEC=='GeneralizedBorn' */
525 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
526 isaprod = isai{I}*isaj{J};
527 gbqqfactor = isaprod*(-qq{I}{J})*gbinvepsdiff;
528 gbscale = isaprod*gbtabscale;
529 dvdaj = dvda[jnr+{J}];
530 /* #define INNERFLOPS INNERFLOPS+5 */
532 /* Calculate generalized born table index - this is a separate table from the normal one,
533 * but we use the same procedure by multiplying r with scale and truncating to integer.
535 rt = r{I}{J}*gbscale;
542 Geps = gbeps*gbtab[gbitab+2];
543 Heps2 = gbeps*gbeps*gbtab[gbitab+3];
547 /* #define INNERFLOPS INNERFLOPS+10 */
549 /* #if 'Force' in KERNEL_VF */
550 FF = Fp+Geps+2.0*Heps2;
551 fgb = gbqqfactor*FF*gbscale;
552 dvdatmp = -0.5*(vgb+fgb*r{I}{J});
553 dvdasum = dvdasum + dvdatmp;
554 dvda[jnr] = dvdaj+dvdatmp*isaj{J}*isaj{J};
555 /* #define INNERFLOPS INNERFLOPS+13 */
557 velec = qq{I}{J}*rinv{I}{J};
558 /* #define INNERFLOPS INNERFLOPS+1 */
559 /* #if 'Force' in KERNEL_VF */
560 felec = (velec*rinv{I}{J}-fgb)*rinv{I}{J};
561 /* #define INNERFLOPS INNERFLOPS+3 */
564 /* #elif KERNEL_ELEC=='Ewald' */
565 /* EWALD ELECTROSTATICS */
567 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
568 ewrt = r{I}{J}*ewtabscale;
571 /* #define INNERFLOPS INNERFLOPS+2 */
572 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
574 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
575 /* #define INNERFLOPS INNERFLOPS+4 */
576 /* #if KERNEL_MOD_ELEC=='PotentialShift' */
577 velec = qq{I}{J}*((rinv{I}{J}-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
578 /* #define INNERFLOPS INNERFLOPS+7 */
580 velec = qq{I}{J}*(rinv{I}{J}-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
581 /* #define INNERFLOPS INNERFLOPS+6 */
583 /* #if 'Force' in KERNEL_VF */
584 felec = qq{I}{J}*rinv{I}{J}*(rinvsq{I}{J}-felec);
585 /* #define INNERFLOPS INNERFLOPS+3 */
587 /* #elif KERNEL_VF=='Force' */
588 felec = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
589 felec = qq{I}{J}*rinv{I}{J}*(rinvsq{I}{J}-felec);
590 /* #define INNERFLOPS INNERFLOPS+7 */
593 /* #elif KERNEL_ELEC=='CubicSplineTable' */
595 /* CUBIC SPLINE TABLE ELECTROSTATICS */
596 /* #if 'Potential' in KERNEL_VF */
600 Geps = vfeps*vftab[vfitab+2];
601 Heps2 = vfeps*vfeps*vftab[vfitab+3];
603 /* #define INNERFLOPS INNERFLOPS+5 */
604 /* #if 'Potential' in KERNEL_VF */
607 /* #define INNERFLOPS INNERFLOPS+3 */
609 /* #if 'Force' in KERNEL_VF */
610 FF = Fp+Geps+2.0*Heps2;
611 felec = -qq{I}{J}*FF*vftabscale*rinv{I}{J};
612 /* #define INNERFLOPS INNERFLOPS+7 */
615 /* ## End of check for electrostatics interaction forms */
617 /* ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
619 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
621 /* #if KERNEL_VDW=='LennardJones' */
623 /* LENNARD-JONES DISPERSION/REPULSION */
625 rinvsix = rinvsq{I}{J}*rinvsq{I}{J}*rinvsq{I}{J};
626 /* #define INNERFLOPS INNERFLOPS+2 */
627 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
628 vvdw6 = c6_{I}{J}*rinvsix;
629 vvdw12 = c12_{I}{J}*rinvsix*rinvsix;
630 /* #define INNERFLOPS INNERFLOPS+3 */
631 /* #if KERNEL_MOD_VDW=='PotentialShift' */
632 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);
633 /* #define INNERFLOPS INNERFLOPS+8 */
635 vvdw = vvdw12*(1.0/12.0) - vvdw6*(1.0/6.0);
636 /* #define INNERFLOPS INNERFLOPS+3 */
638 /* ## Check for force inside potential check, i.e. this means we already did the potential part */
639 /* #if 'Force' in KERNEL_VF */
640 fvdw = (vvdw12-vvdw6)*rinvsq{I}{J};
641 /* #define INNERFLOPS INNERFLOPS+2 */
643 /* #elif KERNEL_VF=='Force' */
644 /* ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
645 fvdw = (c12_{I}{J}*rinvsix-c6_{I}{J})*rinvsix*rinvsq{I}{J};
646 /* #define INNERFLOPS INNERFLOPS+4 */
649 /* #elif KERNEL_VDW=='Buckingham' */
651 /* BUCKINGHAM DISPERSION/REPULSION */
652 rinvsix = rinvsq{I}{J}*rinvsq{I}{J}*rinvsq{I}{J};
653 vvdw6 = c6_{I}{J}*rinvsix;
654 br = cexp2_{I}{J}*r{I}{J};
655 vvdwexp = cexp1_{I}{J}*exp(-br);
656 /* ## Estimate exp() to 25 flops */
657 /* #define INNERFLOPS INNERFLOPS+31 */
658 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
659 /* #if KERNEL_MOD_VDW=='PotentialShift' */
660 vvdw = (vvdwexp-cexp1_{I}{J}*exp(-cexp2_{I}{J}*rvdw)) - (vvdw6 - c6_{I}{J}*sh_vdw_invrcut6)*(1.0/6.0);
661 /* #define INNERFLOPS INNERFLOPS+33 */
663 vvdw = vvdwexp - vvdw6*(1.0/6.0);
664 /* #define INNERFLOPS INNERFLOPS+2 */
667 /* #if 'Force' in KERNEL_VF */
668 fvdw = (br*vvdwexp-vvdw6)*rinvsq{I}{J};
669 /* #define INNERFLOPS INNERFLOPS+3 */
672 /* #elif KERNEL_VDW=='CubicSplineTable' */
674 /* CUBIC SPLINE TABLE DISPERSION */
675 vfitab += {TABLE_VDW_OFFSET};
676 /* #if 'Potential' in KERNEL_VF */
680 Geps = vfeps*vftab[vfitab+2];
681 Heps2 = vfeps*vfeps*vftab[vfitab+3];
683 /* #define INNERFLOPS INNERFLOPS+5 */
684 /* #if 'Potential' in KERNEL_VF */
686 vvdw6 = c6_{I}{J}*VV;
687 /* #define INNERFLOPS INNERFLOPS+3 */
689 /* #if 'Force' in KERNEL_VF */
690 FF = Fp+Geps+2.0*Heps2;
691 fvdw6 = c6_{I}{J}*FF;
692 /* #define INNERFLOPS INNERFLOPS+4 */
695 /* CUBIC SPLINE TABLE REPULSION */
696 /* #if 'Potential' in KERNEL_VF */
700 Geps = vfeps*vftab[vfitab+6];
701 Heps2 = vfeps*vfeps*vftab[vfitab+7];
703 /* #define INNERFLOPS INNERFLOPS+5 */
704 /* #if 'Potential' in KERNEL_VF */
706 vvdw12 = c12_{I}{J}*VV;
707 /* #define INNERFLOPS INNERFLOPS+3 */
709 /* #if 'Force' in KERNEL_VF */
710 FF = Fp+Geps+2.0*Heps2;
711 fvdw12 = c12_{I}{J}*FF;
712 /* #define INNERFLOPS INNERFLOPS+4 */
714 /* #if 'Potential' in KERNEL_VF */
716 /* #define INNERFLOPS INNERFLOPS+1 */
718 /* #if 'Force' in KERNEL_VF */
719 fvdw = -(fvdw6+fvdw12)*vftabscale*rinv{I}{J};
720 /* #define INNERFLOPS INNERFLOPS+4 */
723 /* #elif KERNEL_VDW=='LJEwald' */
725 rinvsix = rinvsq{I}{J}*rinvsq{I}{J}*rinvsq{I}{J};
726 ewcljrsq = ewclj2*rsq{I}{J};
727 exponent = exp(-ewcljrsq);
728 poly = exponent*(1.0 + ewcljrsq + ewcljrsq*ewcljrsq*0.5);
729 /* #define INNERFLOPS INNERFLOPS+9 */
730 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
731 vvdw6 = (c6_{I}{J}-c6grid_{I}{J}*(1.0-poly))*rinvsix;
732 vvdw12 = c12_{I}{J}*rinvsix*rinvsix;
733 /* #define INNERFLOPS INNERFLOPS+6 */
734 /* #if KERNEL_MOD_VDW=='PotentialShift' */
735 vvdw = (vvdw12 - c12_{I}{J}*sh_vdw_invrcut6*sh_vdw_invrcut6)*(1.0/12.0) - (vvdw6 - c6_{I}{J}*sh_vdw_invrcut6 - c6grid_{I}{J}*sh_lj_ewald)*(1.0/6.0);
736 /* #define INNERFLOPS INNERFLOPS+9 */
738 vvdw = vvdw12*(1.0/12.0) - vvdw6*(1.0/6.0);
739 /* #define INNERFLOPS INNERFLOPS+3 */
741 /* ## Check for force inside potential check, i.e. this means we already did the potential part */
742 /* #if 'Force' in KERNEL_VF */
743 fvdw = (vvdw12 - vvdw6 - c6grid_{I}{J}*(1.0/6.0)*exponent*ewclj6)*rinvsq{I}{J};
744 /* #define INNERFLOPS INNERFLOPS+6 */
746 /* #elif KERNEL_VF=='Force' */
747 /* ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
748 fvdw = (((c12_{I}{J}*rinvsix - c6_{I}{J} + c6grid_{I}{J}*(1.0-poly))*rinvsix) - c6grid_{I}{J}*(1.0/6.0)*exponent*ewclj6)*rinvsq{I}{J};
749 /* #define INNERFLOPS INNERFLOPS+11 */
752 /* ## End of check for vdw interaction forms */
754 /* ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
756 /* #if 'switch' in INTERACTION_FLAGS[I][J] */
758 d = (d>0.0) ? d : 0.0;
760 sw = 1.0+d2*d*(swV3+d*(swV4+d*swV5));
761 /* #define INNERFLOPS INNERFLOPS+9 */
763 /* #if 'Force' in KERNEL_VF */
764 dsw = d2*(swF2+d*(swF3+d*swF4));
765 /* #define INNERFLOPS INNERFLOPS+5 */
768 /* Evaluate switch function */
769 /* #if 'Force' in KERNEL_VF */
770 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
771 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
772 felec = felec*sw - rinv{I}{J}*velec*dsw;
773 /* #define INNERFLOPS INNERFLOPS+3 */
775 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
776 fvdw = fvdw*sw - rinv{I}{J}*vvdw*dsw;
777 /* #define INNERFLOPS INNERFLOPS+3 */
780 /* #if 'Potential' in KERNEL_VF */
781 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
783 /* #define INNERFLOPS INNERFLOPS+1 */
785 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
787 /* #define INNERFLOPS INNERFLOPS+1 */
792 /* #if 'Potential' in KERNEL_VF */
793 /* Update potential sums from outer loop */
794 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
796 /* #define INNERFLOPS INNERFLOPS+1 */
797 /* #if KERNEL_ELEC=='GeneralizedBorn' */
799 /* #define INNERFLOPS INNERFLOPS+1 */
802 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
804 /* #define INNERFLOPS INNERFLOPS+1 */
808 /* #if 'Force' in KERNEL_VF */
810 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
812 /* #define INNERFLOPS INNERFLOPS+1 */
813 /* #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
815 /* #elif 'vdw' in INTERACTION_FLAGS[I][J] */
819 /* Calculate temporary vectorial force */
824 /* Update vectorial force */
828 f[j_coord_offset+DIM*{J}+XX] -= tx;
829 f[j_coord_offset+DIM*{J}+YY] -= ty;
830 f[j_coord_offset+DIM*{J}+ZZ] -= tz;
832 /* #define INNERFLOPS INNERFLOPS+9 */
835 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
836 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
837 /* #if 0 ## This and next two lines is a hack to maintain indentation in template file */
842 /* ## End of check for the interaction being outside the cutoff */
845 /* ## End of loop over i-j interaction pairs */
847 /* Inner loop uses {INNERFLOPS} flops */
849 /* End of innermost loop */
851 /* #if 'Force' in KERNEL_VF */
853 /* #for I in PARTICLES_I */
854 f[i_coord_offset+DIM*{I}+XX] += fix{I};
855 f[i_coord_offset+DIM*{I}+YY] += fiy{I};
856 f[i_coord_offset+DIM*{I}+ZZ] += fiz{I};
860 /* #define OUTERFLOPS OUTERFLOPS+6 */
862 fshift[i_shift_offset+XX] += tx;
863 fshift[i_shift_offset+YY] += ty;
864 fshift[i_shift_offset+ZZ] += tz;
865 /* #define OUTERFLOPS OUTERFLOPS+3 */
868 /* #if 'Potential' in KERNEL_VF */
870 /* Update potential energies */
871 /* #if KERNEL_ELEC != 'None' */
872 kernel_data->energygrp_elec[ggid] += velecsum;
873 /* #define OUTERFLOPS OUTERFLOPS+1 */
875 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
876 kernel_data->energygrp_polarization[ggid] += vgbsum;
877 /* #define OUTERFLOPS OUTERFLOPS+1 */
879 /* #if KERNEL_VDW != 'None' */
880 kernel_data->energygrp_vdw[ggid] += vvdwsum;
881 /* #define OUTERFLOPS OUTERFLOPS+1 */
884 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
885 dvda[inr] = dvda[inr] + dvdasum*isai{I}*isai{I};
888 /* Increment number of inner iterations */
889 inneriter += j_index_end - j_index_start;
891 /* Outer loop uses {OUTERFLOPS} flops */
894 /* Increment number of outer iterations */
897 /* Update outer/inner flops */
898 /* ## NB: This is not important, it just affects the flopcount. However, since our preprocessor is */
899 /* ## primitive and replaces aggressively even in strings inside these directives, we need to */
900 /* ## assemble the main part of the name (containing KERNEL/ELEC/VDW) directly in the source. */
901 /* #if GEOMETRY_I == 'Water3' */
902 /* #define ISUFFIX '_W3' */
903 /* #elif GEOMETRY_I == 'Water4' */
904 /* #define ISUFFIX '_W4' */
906 /* #define ISUFFIX '' */
908 /* #if GEOMETRY_J == 'Water3' */
909 /* #define JSUFFIX 'W3' */
910 /* #elif GEOMETRY_J == 'Water4' */
911 /* #define JSUFFIX 'W4' */
913 /* #define JSUFFIX '' */
915 /* #if 'PotentialAndForce' in KERNEL_VF */
916 /* #define VFSUFFIX '_VF' */
917 /* #elif 'Potential' in KERNEL_VF */
918 /* #define VFSUFFIX '_V' */
920 /* #define VFSUFFIX '_F' */
923 /* #if KERNEL_ELEC != 'None' and KERNEL_VDW != 'None' */
924 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
925 /* #elif KERNEL_ELEC != 'None' */
926 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
928 inc_nrnb(nrnb,eNR_NBKERNEL_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});