2 #error This file must be processed with the Gromacs pre-preprocessor
4 /* #if INCLUDE_HEADER */
11 #include "../nb_kernel.h"
12 #include "types/simple.h"
17 /* ## List of variables set by the generating script: */
19 /* ## Setttings that apply to the entire kernel: */
20 /* ## KERNEL_ELEC: String, choice for electrostatic interactions */
21 /* ## KERNEL_VDW: String, choice for van der Waals interactions */
22 /* ## KERNEL_NAME: String, name of this kernel */
23 /* ## KERNEL_VF: String telling if we calculate potential, force, or both */
24 /* ## GEOMETRY_I/GEOMETRY_J: String, name of each geometry, e.g. 'Water3' or '1Particle' */
26 /* ## Setttings that apply to particles in the outer (I) or inner (J) loops: */
27 /* ## PARTICLES_I[]/ Arrays with lists of i/j particles to use in kernel. It is */
28 /* ## PARTICLES_J[]: just [0] for particle geometry, but can be longer for water */
29 /* ## PARTICLES_ELEC_I[]/ Arrays with lists of i/j particle that have electrostatics */
30 /* ## PARTICLES_ELEC_J[]: interactions that should be calculated in this kernel. */
31 /* ## PARTICLES_VDW_I[]/ Arrays with the list of i/j particle that have VdW */
32 /* ## PARTICLES_VDW_J[]: interactions that should be calculated in this kernel. */
34 /* ## Setttings for pairs of interactions (e.g. 2nd i particle against 1st j particle) */
35 /* ## PAIRS_IJ[]: Array with (i,j) tuples of pairs for which interactions */
36 /* ## should be calculated in this kernel. Zero-charge particles */
37 /* ## do not have interactions with particles without vdw, and */
38 /* ## Vdw-only interactions are not evaluated in a no-vdw-kernel. */
39 /* ## INTERACTION_FLAGS[][]: 2D matrix, dimension e.g. 3*3 for water-water interactions. */
40 /* ## For each i-j pair, the element [I][J] is a list of strings */
41 /* ## defining properties/flags of this interaction. Examples */
42 /* ## include 'electrostatics'/'vdw' if that type of interaction */
43 /* ## should be evaluated, 'rsq'/'rinv'/'rinvsq' if those values */
44 /* ## are needed, and 'exactcutoff' or 'shift','switch' to */
45 /* ## decide if the force/potential should be modified. This way */
46 /* ## we only calculate values absolutely needed for each case. */
48 /* ## Calculate the size and offset for (merged/interleaved) table data */
50 /* #if ('CubicSplineTable' in [KERNEL_ELEC,KERNEL_VDW]) or KERNEL_VF=='PotentialAndForce' */
51 /* #define TABLE_POINT_SIZE 4 */
53 /* #define TABLE_POINT_SIZE 2 */
56 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
57 /* #define TABLE_NINTERACTIONS 3 */
58 /* #define TABLE_VDW_OFFSET TABLE_POINT_SIZE */
59 /* #elif 'Table' in KERNEL_ELEC */
60 /* #define TABLE_NINTERACTIONS 1 */
61 /* #elif 'Table' in KERNEL_VDW */
62 /* #define TABLE_NINTERACTIONS 2 */
63 /* #define TABLE_VDW_OFFSET 0 */
65 /* #define TABLE_NINTERACTIONS 0 */
68 /* #if 'Buckingham' in KERNEL_VDW */
69 /* #define NVDWPARAM 3 */
71 /* #define NVDWPARAM 2 */
75 * Gromacs nonbonded kernel: {KERNEL_NAME}
76 * Electrostatics interaction: {KERNEL_ELEC}
77 * VdW interaction: {KERNEL_VDW}
78 * Geometry: {GEOMETRY_I}-{GEOMETRY_J}
79 * Calculate force/pot: {KERNEL_VF}
83 (t_nblist * gmx_restrict nlist,
84 rvec * gmx_restrict xx,
85 rvec * gmx_restrict ff,
86 t_forcerec * gmx_restrict fr,
87 t_mdatoms * gmx_restrict mdatoms,
88 nb_kernel_data_t * gmx_restrict kernel_data,
89 t_nrnb * gmx_restrict nrnb)
91 /* ## Not all variables are used for all kernels, but any optimizing compiler fixes that, */
92 /* ## so there is no point in going to extremes to exclude variables that are not needed. */
93 int i_shift_offset,i_coord_offset,j_coord_offset;
94 int j_index_start,j_index_end;
95 int nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
96 real shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
97 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
98 real *shiftvec,*fshift,*x,*f;
99 /* #for I in PARTICLES_I */
101 real ix{I},iy{I},iz{I},fix{I},fiy{I},fiz{I},iq{I},isai{I};
103 /* #for J in PARTICLES_J */
105 real jx{J},jy{J},jz{J},fjx{J},fjy{J},fjz{J},jq{J},isaj{J};
107 /* #for I,J in PAIRS_IJ */
108 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};
110 /* #if KERNEL_ELEC != 'None' */
111 real velec,felec,velecsum,facel,crf,krf,krf2;
114 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
116 real vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,dvdatmp;
117 real *invsqrta,*dvda,*gbtab;
119 /* #if KERNEL_VDW != 'None' */
121 real rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
125 /* #if 'Table' in KERNEL_ELEC or 'GeneralizedBorn' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
127 real rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
130 /* #if 'Ewald' in KERNEL_ELEC */
132 real ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
135 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
136 real rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
144 jindex = nlist->jindex;
146 shiftidx = nlist->shift;
148 shiftvec = fr->shift_vec[0];
149 fshift = fr->fshift[0];
150 /* #if KERNEL_ELEC != 'None' */
152 charge = mdatoms->chargeA;
153 /* #if 'ReactionField' in KERNEL_ELEC */
159 /* #if KERNEL_VDW != 'None' */
160 nvdwtype = fr->ntype;
162 vdwtype = mdatoms->typeA;
165 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
166 vftab = kernel_data->table_elec_vdw->data;
167 vftabscale = kernel_data->table_elec_vdw->scale;
168 /* #elif 'Table' in KERNEL_ELEC */
169 vftab = kernel_data->table_elec->data;
170 vftabscale = kernel_data->table_elec->scale;
171 /* #elif 'Table' in KERNEL_VDW */
172 vftab = kernel_data->table_vdw->data;
173 vftabscale = kernel_data->table_vdw->scale;
176 /* #if 'Ewald' in KERNEL_ELEC */
177 sh_ewald = fr->ic->sh_ewald;
178 /* #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
179 ewtab = fr->ic->tabq_coul_F;
180 ewtabscale = fr->ic->tabq_scale;
181 ewtabhalfspace = 0.5/ewtabscale;
183 ewtab = fr->ic->tabq_coul_FDV0;
184 ewtabscale = fr->ic->tabq_scale;
185 ewtabhalfspace = 0.5/ewtabscale;
189 /* #if KERNEL_ELEC=='GeneralizedBorn' */
190 invsqrta = fr->invsqrta;
192 gbtabscale = fr->gbtab.scale;
193 gbtab = fr->gbtab.data;
194 gbinvepsdiff = (1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent);
197 /* #if 'Water' in GEOMETRY_I */
198 /* Setup water-specific parameters */
199 inr = nlist->iinr[0];
200 /* #for I in PARTICLES_ELEC_I */
201 iq{I} = facel*charge[inr+{I}];
203 /* #for I in PARTICLES_VDW_I */
204 vdwioffset{I} = {NVDWPARAM}*nvdwtype*vdwtype[inr+{I}];
208 /* #if 'Water' in GEOMETRY_J */
209 /* #for J in PARTICLES_ELEC_J */
210 jq{J} = charge[inr+{J}];
212 /* #for J in PARTICLES_VDW_J */
213 vdwjidx{J} = {NVDWPARAM}*vdwtype[inr+{J}];
215 /* #for I,J in PAIRS_IJ */
216 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
217 qq{I}{J} = iq{I}*jq{J};
219 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
220 /* #if 'Buckingham' in KERNEL_VDW */
221 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
222 cexp1_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
223 cexp2_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+2];
225 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
226 c12_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
232 /* #if KERNEL_MOD_ELEC!='None' or KERNEL_MOD_VDW!='None' */
233 /* #if KERNEL_ELEC!='None' */
234 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
235 rcutoff = fr->rcoulomb;
239 rcutoff2 = rcutoff*rcutoff;
242 /* #if KERNEL_MOD_VDW=='PotentialShift' */
243 sh_vdw_invrcut6 = fr->ic->sh_invrc6;
247 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
248 /* #if KERNEL_MOD_ELEC=='PotentialSwitch' */
249 rswitch = fr->rcoulomb_switch;
251 rswitch = fr->rvdw_switch;
253 /* Setup switch parameters */
255 swV3 = -10.0/(d*d*d);
256 swV4 = 15.0/(d*d*d*d);
257 swV5 = -6.0/(d*d*d*d*d);
258 /* #if 'Force' in KERNEL_VF */
259 swF2 = -30.0/(d*d*d);
260 swF3 = 60.0/(d*d*d*d);
261 swF4 = -30.0/(d*d*d*d*d);
265 /* ## Keep track of the floating point operations we issue for reporting! */
266 /* #define OUTERFLOPS 0 */
267 /* #define INNERFLOPS 0 */
271 /* Start outer loop over neighborlists */
272 for(iidx=0; iidx<nri; iidx++)
274 /* Load shift vector for this list */
275 i_shift_offset = DIM*shiftidx[iidx];
276 shX = shiftvec[i_shift_offset+XX];
277 shY = shiftvec[i_shift_offset+YY];
278 shZ = shiftvec[i_shift_offset+ZZ];
280 /* Load limits for loop over neighbors */
281 j_index_start = jindex[iidx];
282 j_index_end = jindex[iidx+1];
284 /* Get outer coordinate index */
286 i_coord_offset = DIM*inr;
288 /* Load i particle coords and add shift vector */
289 /* ## Loop over i particles, but only include ones that we use - skip e.g. vdw-only sites for elec-only kernel */
290 /* #for I in PARTICLES_I */
291 ix{I} = shX + x[i_coord_offset+DIM*{I}+XX];
292 iy{I} = shY + x[i_coord_offset+DIM*{I}+YY];
293 iz{I} = shZ + x[i_coord_offset+DIM*{I}+ZZ];
294 /* #define OUTERFLOPS OUTERFLOPS+3 */
297 /* #if 'Force' in KERNEL_VF */
298 /* #for I in PARTICLES_I */
305 /* ## For water we already preloaded parameters at the start of the kernel */
306 /* #if not 'Water' in GEOMETRY_I */
307 /* Load parameters for i particles */
308 /* #for I in PARTICLES_ELEC_I */
309 iq{I} = facel*charge[inr+{I}];
310 /* #define OUTERFLOPS OUTERFLOPS+1 */
311 /* #if KERNEL_ELEC=='GeneralizedBorn' */
312 isai{I} = invsqrta[inr+{I}];
315 /* #for I in PARTICLES_VDW_I */
316 vdwioffset{I} = {NVDWPARAM}*nvdwtype*vdwtype[inr+{I}];
320 /* #if 'Potential' in KERNEL_VF */
321 /* Reset potential sums */
322 /* #if KERNEL_ELEC != 'None' */
325 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
328 /* #if KERNEL_VDW != 'None' */
332 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
336 /* Start inner kernel loop */
337 for(jidx=j_index_start; jidx<j_index_end; jidx++)
339 /* Get j neighbor index, and coordinate index */
341 j_coord_offset = DIM*jnr;
343 /* load j atom coordinates */
344 /* #for J in PARTICLES_J */
345 jx{J} = x[j_coord_offset+DIM*{J}+XX];
346 jy{J} = x[j_coord_offset+DIM*{J}+YY];
347 jz{J} = x[j_coord_offset+DIM*{J}+ZZ];
350 /* Calculate displacement vector */
351 /* #for I,J in PAIRS_IJ */
352 dx{I}{J} = ix{I} - jx{J};
353 dy{I}{J} = iy{I} - jy{J};
354 dz{I}{J} = iz{I} - jz{J};
355 /* #define INNERFLOPS INNERFLOPS+3 */
358 /* Calculate squared distance and things based on it */
359 /* #for I,J in PAIRS_IJ */
360 rsq{I}{J} = dx{I}{J}*dx{I}{J}+dy{I}{J}*dy{I}{J}+dz{I}{J}*dz{I}{J};
361 /* #define INNERFLOPS INNERFLOPS+5 */
364 /* #for I,J in PAIRS_IJ */
365 /* #if 'rinv' in INTERACTION_FLAGS[I][J] */
366 rinv{I}{J} = gmx_invsqrt(rsq{I}{J});
367 /* #define INNERFLOPS INNERFLOPS+5 */
371 /* #for I,J in PAIRS_IJ */
372 /* #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
373 /* # if 'rinv' not in INTERACTION_FLAGS[I][J] */
374 rinvsq{I}{J} = 1.0/rsq{I}{J};
375 /* #define INNERFLOPS INNERFLOPS+4 */
377 rinvsq{I}{J} = rinv{I}{J}*rinv{I}{J};
378 /* #define INNERFLOPS INNERFLOPS+1 */
383 /* #if not 'Water' in GEOMETRY_J */
384 /* Load parameters for j particles */
385 /* #for J in PARTICLES_ELEC_J */
386 jq{J} = charge[jnr+{J}];
387 /* #if KERNEL_ELEC=='GeneralizedBorn' */
388 isaj{J} = invsqrta[jnr+{J}];
391 /* #for J in PARTICLES_VDW_J */
392 vdwjidx{J} = {NVDWPARAM}*vdwtype[jnr+{J}];
396 /* #for I,J in PAIRS_IJ */
398 /**************************
399 * CALCULATE INTERACTIONS *
400 **************************/
402 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
403 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
404 /* ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
405 if (rsq{I}{J}<rcutoff2)
407 /* #if 0 ## this and the next two lines is a hack to maintain auto-indentation in template file */
412 /* #if 'r' in INTERACTION_FLAGS[I][J] */
413 r{I}{J} = rsq{I}{J}*rinv{I}{J};
414 /* #define INNERFLOPS INNERFLOPS+1 */
417 /* ## For water geometries we already loaded parameters at the start of the kernel */
418 /* #if not 'Water' in GEOMETRY_J */
419 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
420 qq{I}{J} = iq{I}*jq{J};
421 /* #define INNERFLOPS INNERFLOPS+1 */
423 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
424 /* #if KERNEL_VDW=='Buckingham' */
425 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
426 cexp1_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
427 cexp2_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+2];
429 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
430 c12_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
435 /* #if 'table' in INTERACTION_FLAGS[I][J] */
436 /* Calculate table index by multiplying r with table scale and truncate to integer */
437 rt = r{I}{J}*vftabscale;
440 vfitab = {TABLE_NINTERACTIONS}*{TABLE_POINT_SIZE}*vfitab;
441 /* #define INNERFLOPS INNERFLOPS+2 */
444 /* ## ELECTROSTATIC INTERACTIONS */
445 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
447 /* #if KERNEL_ELEC=='Coulomb' */
449 /* COULOMB ELECTROSTATICS */
450 velec = qq{I}{J}*rinv{I}{J};
451 /* #define INNERFLOPS INNERFLOPS+1 */
452 /* #if 'Force' in KERNEL_VF */
453 felec = velec*rinvsq{I}{J};
454 /* #define INNERFLOPS INNERFLOPS+2 */
457 /* #elif KERNEL_ELEC=='ReactionField' */
459 /* REACTION-FIELD ELECTROSTATICS */
460 /* #if 'Potential' in KERNEL_VF */
461 velec = qq{I}{J}*(rinv{I}{J}+krf*rsq{I}{J}-crf);
462 /* #define INNERFLOPS INNERFLOPS+4 */
464 /* #if 'Force' in KERNEL_VF */
465 felec = qq{I}{J}*(rinv{I}{J}*rinvsq{I}{J}-krf2);
466 /* #define INNERFLOPS INNERFLOPS+3 */
469 /* #elif KERNEL_ELEC=='GeneralizedBorn' */
471 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
472 isaprod = isai{I}*isaj{J};
473 gbqqfactor = isaprod*(-qq{I}{J})*gbinvepsdiff;
474 gbscale = isaprod*gbtabscale;
475 dvdaj = dvda[jnr+{J}];
476 /* #define INNERFLOPS INNERFLOPS+5 */
478 /* Calculate generalized born table index - this is a separate table from the normal one,
479 * but we use the same procedure by multiplying r with scale and truncating to integer.
481 rt = r{I}{J}*gbscale;
488 Geps = gbeps*gbtab[gbitab+2];
489 Heps2 = gbeps*gbeps*gbtab[gbitab+3];
493 /* #define INNERFLOPS INNERFLOPS+10 */
495 /* #if 'Force' in KERNEL_VF */
496 FF = Fp+Geps+2.0*Heps2;
497 fgb = gbqqfactor*FF*gbscale;
498 dvdatmp = -0.5*(vgb+fgb*r{I}{J});
499 dvdasum = dvdasum + dvdatmp;
500 dvda[jnr] = dvdaj+dvdatmp*isaj{J}*isaj{J};
501 /* #define INNERFLOPS INNERFLOPS+13 */
503 velec = qq{I}{J}*rinv{I}{J};
504 /* #define INNERFLOPS INNERFLOPS+1 */
505 /* #if 'Force' in KERNEL_VF */
506 felec = (velec*rinv{I}{J}-fgb)*rinv{I}{J};
507 /* #define INNERFLOPS INNERFLOPS+3 */
510 /* #elif KERNEL_ELEC=='Ewald' */
511 /* EWALD ELECTROSTATICS */
513 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
514 ewrt = r{I}{J}*ewtabscale;
517 /* #define INNERFLOPS INNERFLOPS+2 */
518 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
520 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
521 /* #define INNERFLOPS INNERFLOPS+4 */
522 /* #if KERNEL_MOD_ELEC=='PotentialShift' */
523 velec = qq{I}{J}*((rinv{I}{J}-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
524 /* #define INNERFLOPS INNERFLOPS+7 */
526 velec = qq{I}{J}*(rinv{I}{J}-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
527 /* #define INNERFLOPS INNERFLOPS+6 */
529 /* #if 'Force' in KERNEL_VF */
530 felec = qq{I}{J}*rinv{I}{J}*(rinvsq{I}{J}-felec);
531 /* #define INNERFLOPS INNERFLOPS+3 */
533 /* #elif KERNEL_VF=='Force' */
534 felec = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
535 felec = qq{I}{J}*rinv{I}{J}*(rinvsq{I}{J}-felec);
536 /* #define INNERFLOPS INNERFLOPS+7 */
539 /* #elif KERNEL_ELEC=='CubicSplineTable' */
541 /* CUBIC SPLINE TABLE ELECTROSTATICS */
542 /* #if 'Potential' in KERNEL_VF */
546 Geps = vfeps*vftab[vfitab+2];
547 Heps2 = vfeps*vfeps*vftab[vfitab+3];
549 /* #define INNERFLOPS INNERFLOPS+5 */
550 /* #if 'Potential' in KERNEL_VF */
553 /* #define INNERFLOPS INNERFLOPS+3 */
555 /* #if 'Force' in KERNEL_VF */
556 FF = Fp+Geps+2.0*Heps2;
557 felec = -qq{I}{J}*FF*vftabscale*rinv{I}{J};
558 /* #define INNERFLOPS INNERFLOPS+7 */
561 /* ## End of check for electrostatics interaction forms */
563 /* ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
565 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
567 /* #if KERNEL_VDW=='LennardJones' */
569 /* LENNARD-JONES DISPERSION/REPULSION */
571 rinvsix = rinvsq{I}{J}*rinvsq{I}{J}*rinvsq{I}{J};
572 /* #define INNERFLOPS INNERFLOPS+2 */
573 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
574 vvdw6 = c6_{I}{J}*rinvsix;
575 vvdw12 = c12_{I}{J}*rinvsix*rinvsix;
576 /* #define INNERFLOPS INNERFLOPS+3 */
577 /* #if KERNEL_MOD_VDW=='PotentialShift' */
578 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);
579 /* #define INNERFLOPS INNERFLOPS+8 */
581 vvdw = vvdw12*(1.0/12.0) - vvdw6*(1.0/6.0);
582 /* #define INNERFLOPS INNERFLOPS+3 */
584 /* ## Check for force inside potential check, i.e. this means we already did the potential part */
585 /* #if 'Force' in KERNEL_VF */
586 fvdw = (vvdw12-vvdw6)*rinvsq{I}{J};
587 /* #define INNERFLOPS INNERFLOPS+2 */
589 /* #elif KERNEL_VF=='Force' */
590 /* ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
591 fvdw = (c12_{I}{J}*rinvsix-c6_{I}{J})*rinvsix*rinvsq{I}{J};
592 /* #define INNERFLOPS INNERFLOPS+4 */
595 /* #elif KERNEL_VDW=='Buckingham' */
597 /* BUCKINGHAM DISPERSION/REPULSION */
598 rinvsix = rinvsq{I}{J}*rinvsq{I}{J}*rinvsq{I}{J};
599 vvdw6 = c6_{I}{J}*rinvsix;
600 br = cexp2_{I}{J}*r{I}{J};
601 vvdwexp = cexp1_{I}{J}*exp(-br);
602 /* ## Estimate exp() to 25 flops */
603 /* #define INNERFLOPS INNERFLOPS+31 */
604 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
605 /* #if KERNEL_MOD_VDW=='PotentialShift' */
606 vvdw = (vvdwexp-cexp1_{I}{J}*exp(-cexp2_{I}{J}*rvdw)) - (vvdw6 - c6_{I}{J}*sh_vdw_invrcut6)*(1.0/6.0);
607 /* #define INNERFLOPS INNERFLOPS+33 */
609 vvdw = vvdwexp - vvdw6*(1.0/6.0);
610 /* #define INNERFLOPS INNERFLOPS+2 */
613 /* #if 'Force' in KERNEL_VF */
614 fvdw = (br*vvdwexp-vvdw6)*rinvsq{I}{J};
615 /* #define INNERFLOPS INNERFLOPS+3 */
618 /* #elif KERNEL_VDW=='CubicSplineTable' */
620 /* CUBIC SPLINE TABLE DISPERSION */
621 vfitab += {TABLE_VDW_OFFSET};
622 /* #if 'Potential' in KERNEL_VF */
626 Geps = vfeps*vftab[vfitab+2];
627 Heps2 = vfeps*vfeps*vftab[vfitab+3];
629 /* #define INNERFLOPS INNERFLOPS+5 */
630 /* #if 'Potential' in KERNEL_VF */
632 vvdw6 = c6_{I}{J}*VV;
633 /* #define INNERFLOPS INNERFLOPS+3 */
635 /* #if 'Force' in KERNEL_VF */
636 FF = Fp+Geps+2.0*Heps2;
637 fvdw6 = c6_{I}{J}*FF;
638 /* #define INNERFLOPS INNERFLOPS+4 */
641 /* CUBIC SPLINE TABLE REPULSION */
642 /* #if 'Potential' in KERNEL_VF */
646 Geps = vfeps*vftab[vfitab+6];
647 Heps2 = vfeps*vfeps*vftab[vfitab+7];
649 /* #define INNERFLOPS INNERFLOPS+5 */
650 /* #if 'Potential' in KERNEL_VF */
652 vvdw12 = c12_{I}{J}*VV;
653 /* #define INNERFLOPS INNERFLOPS+3 */
655 /* #if 'Force' in KERNEL_VF */
656 FF = Fp+Geps+2.0*Heps2;
657 fvdw12 = c12_{I}{J}*FF;
658 /* #define INNERFLOPS INNERFLOPS+4 */
660 /* #if 'Potential' in KERNEL_VF */
662 /* #define INNERFLOPS INNERFLOPS+1 */
664 /* #if 'Force' in KERNEL_VF */
665 fvdw = -(fvdw6+fvdw12)*vftabscale*rinv{I}{J};
666 /* #define INNERFLOPS INNERFLOPS+4 */
669 /* ## End of check for vdw interaction forms */
671 /* ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
673 /* #if 'switch' in INTERACTION_FLAGS[I][J] */
675 d = (d>0.0) ? d : 0.0;
677 sw = 1.0+d2*d*(swV3+d*(swV4+d*swV5));
678 /* #define INNERFLOPS INNERFLOPS+9 */
680 /* #if 'Force' in KERNEL_VF */
681 dsw = d2*(swF2+d*(swF3+d*swF4));
682 /* #define INNERFLOPS INNERFLOPS+5 */
685 /* Evaluate switch function */
686 /* #if 'Force' in KERNEL_VF */
687 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
688 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
689 felec = felec*sw - rinv{I}{J}*velec*dsw;
690 /* #define INNERFLOPS INNERFLOPS+3 */
692 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
693 fvdw = fvdw*sw - rinv{I}{J}*vvdw*dsw;
694 /* #define INNERFLOPS INNERFLOPS+3 */
697 /* #if 'Potential' in KERNEL_VF */
698 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
700 /* #define INNERFLOPS INNERFLOPS+1 */
702 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
704 /* #define INNERFLOPS INNERFLOPS+1 */
709 /* #if 'Potential' in KERNEL_VF */
710 /* Update potential sums from outer loop */
711 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
713 /* #define INNERFLOPS INNERFLOPS+1 */
714 /* #if KERNEL_ELEC=='GeneralizedBorn' */
716 /* #define INNERFLOPS INNERFLOPS+1 */
719 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
721 /* #define INNERFLOPS INNERFLOPS+1 */
725 /* #if 'Force' in KERNEL_VF */
727 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
729 /* #define INNERFLOPS INNERFLOPS+1 */
730 /* #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
732 /* #elif 'vdw' in INTERACTION_FLAGS[I][J] */
736 /* Calculate temporary vectorial force */
741 /* Update vectorial force */
745 f[j_coord_offset+DIM*{J}+XX] -= tx;
746 f[j_coord_offset+DIM*{J}+YY] -= ty;
747 f[j_coord_offset+DIM*{J}+ZZ] -= tz;
749 /* #define INNERFLOPS INNERFLOPS+9 */
752 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
753 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
754 /* #if 0 ## This and next two lines is a hack to maintain indentation in template file */
759 /* ## End of check for the interaction being outside the cutoff */
762 /* ## End of loop over i-j interaction pairs */
764 /* Inner loop uses {INNERFLOPS} flops */
766 /* End of innermost loop */
768 /* #if 'Force' in KERNEL_VF */
770 /* #for I in PARTICLES_I */
771 f[i_coord_offset+DIM*{I}+XX] += fix{I};
772 f[i_coord_offset+DIM*{I}+YY] += fiy{I};
773 f[i_coord_offset+DIM*{I}+ZZ] += fiz{I};
777 /* #define OUTERFLOPS OUTERFLOPS+6 */
779 fshift[i_shift_offset+XX] += tx;
780 fshift[i_shift_offset+YY] += ty;
781 fshift[i_shift_offset+ZZ] += tz;
782 /* #define OUTERFLOPS OUTERFLOPS+3 */
785 /* #if 'Potential' in KERNEL_VF */
787 /* Update potential energies */
788 /* #if KERNEL_ELEC != 'None' */
789 kernel_data->energygrp_elec[ggid] += velecsum;
790 /* #define OUTERFLOPS OUTERFLOPS+1 */
792 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
793 kernel_data->energygrp_polarization[ggid] += vgbsum;
794 /* #define OUTERFLOPS OUTERFLOPS+1 */
796 /* #if KERNEL_VDW != 'None' */
797 kernel_data->energygrp_vdw[ggid] += vvdwsum;
798 /* #define OUTERFLOPS OUTERFLOPS+1 */
801 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
802 dvda[inr] = dvda[inr] + dvdasum*isai{I}*isai{I};
805 /* Increment number of inner iterations */
806 inneriter += j_index_end - j_index_start;
808 /* Outer loop uses {OUTERFLOPS} flops */
811 /* Increment number of outer iterations */
814 /* Update outer/inner flops */
815 /* ## NB: This is not important, it just affects the flopcount. However, since our preprocessor is */
816 /* ## primitive and replaces aggressively even in strings inside these directives, we need to */
817 /* ## assemble the main part of the name (containing KERNEL/ELEC/VDW) directly in the source. */
818 /* #if GEOMETRY_I == 'Water3' */
819 /* #define ISUFFIX '_W3' */
820 /* #elif GEOMETRY_I == 'Water4' */
821 /* #define ISUFFIX '_W4' */
823 /* #define ISUFFIX '' */
825 /* #if GEOMETRY_J == 'Water3' */
826 /* #define JSUFFIX 'W3' */
827 /* #elif GEOMETRY_J == 'Water4' */
828 /* #define JSUFFIX 'W4' */
830 /* #define JSUFFIX '' */
832 /* #if 'PotentialAndForce' in KERNEL_VF */
833 /* #define VFSUFFIX '_VF' */
834 /* #elif 'Potential' in KERNEL_VF */
835 /* #define VFSUFFIX '_V' */
837 /* #define VFSUFFIX '_F' */
840 /* #if KERNEL_ELEC != 'None' and KERNEL_VDW != 'None' */
841 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
842 /* #elif KERNEL_ELEC != 'None' */
843 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
845 inc_nrnb(nrnb,eNR_NBKERNEL_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});