3 * This file is part of the GROMACS molecular simulation package.
5 * Copyright (c) 2012, 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 */
45 #include "../nb_kernel.h"
46 #include "types/simple.h"
51 #define ALMOST_ZERO 1e-30
52 #define ALMOST_ONE 1-(1e-30)
54 /* ## List of variables set by the generating script: */
56 /* ## Setttings that apply to the entire kernel: */
57 /* ## KERNEL_ELEC: String, choice for electrostatic interactions */
58 /* ## KERNEL_VDW: String, choice for van der Waals interactions */
59 /* ## KERNEL_NAME: String, name of this kernel */
60 /* ## KERNEL_RESOLUTION: String, choice for resoultion */
61 /* ## KERNEL_VF: String telling if we calculate potential, force, or both */
62 /* ## GEOMETRY_I/GEOMETRY_J: String, name of each geometry, e.g. 'Water3' or '1Particle' */
64 /* ## Setttings that apply to particles in the outer (I) or inner (J) loops: */
65 /* ## PARTICLES_I[]/ Arrays with lists of i/j particles to use in kernel. It is */
66 /* ## PARTICLES_J[]: just [0] for particle geometry, but can be longer for water */
67 /* ## PARTICLES_ELEC_I[]/ Arrays with lists of i/j particle that have electrostatics */
68 /* ## PARTICLES_ELEC_J[]: interactions that should be calculated in this kernel. */
69 /* ## PARTICLES_VDW_I[]/ Arrays with the list of i/j particle that have VdW */
70 /* ## PARTICLES_VDW_J[]: interactions that should be calculated in this kernel. */
72 /* ## Setttings for pairs of interactions (e.g. 2nd i particle against 1st j particle) */
73 /* ## PAIRS_IJ[]: Array with (i,j) tuples of pairs for which interactions */
74 /* ## should be calculated in this kernel. Zero-charge particles */
75 /* ## do not have interactions with particles without vdw, and */
76 /* ## Vdw-only interactions are not evaluated in a no-vdw-kernel. */
77 /* ## INTERACTION_FLAGS[][]: 2D matrix, dimension e.g. 3*3 for water-water interactions. */
78 /* ## For each i-j pair, the element [I][J] is a list of strings */
79 /* ## defining properties/flags of this interaction. Examples */
80 /* ## include 'electrostatics'/'vdw' if that type of interaction */
81 /* ## should be evaluated, 'rsq'/'rinv'/'rinvsq' if those values */
82 /* ## are needed, and 'exactcutoff' or 'shift','switch' to */
83 /* ## decide if the force/potential should be modified. This way */
84 /* ## we only calculate values absolutely needed for each case. */
86 /* ## Calculate the size and offset for (merged/interleaved) table data */
88 /* #if ('CubicSplineTable' in [KERNEL_ELEC,KERNEL_VDW]) or KERNEL_VF=='PotentialAndForce' */
89 /* #define TABLE_POINT_SIZE 4 */
91 /* #define TABLE_POINT_SIZE 2 */
94 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
95 /* #define TABLE_NINTERACTIONS 3 */
96 /* #define TABLE_VDW_OFFSET TABLE_POINT_SIZE */
97 /* #elif 'Table' in KERNEL_ELEC */
98 /* #define TABLE_NINTERACTIONS 1 */
99 /* #elif 'Table' in KERNEL_VDW */
100 /* #define TABLE_NINTERACTIONS 2 */
101 /* #define TABLE_VDW_OFFSET 0 */
103 /* #define TABLE_NINTERACTIONS 0 */
106 /* #if 'Buckingham' in KERNEL_VDW */
107 /* #define NVDWPARAM 3 */
109 /* #define NVDWPARAM 2 */
113 * Gromacs nonbonded kernel: {KERNEL_NAME}
114 * Electrostatics interaction: {KERNEL_ELEC}
115 * VdW interaction: {KERNEL_VDW}
116 * Geometry: {GEOMETRY_I}-{GEOMETRY_J}
117 * Calculate force/pot: {KERNEL_VF}
121 (t_nblist * gmx_restrict nlist,
122 rvec * gmx_restrict xx,
123 rvec * gmx_restrict ff,
124 t_forcerec * gmx_restrict fr,
125 t_mdatoms * gmx_restrict mdatoms,
126 nb_kernel_data_t * gmx_restrict kernel_data,
127 t_nrnb * gmx_restrict nrnb)
129 /* ## Not all variables are used for all kernels, but any optimizing compiler fixes that, */
130 /* ## so there is no point in going to extremes to exclude variables that are not needed. */
131 int i_shift_offset,i_coord_offset,j_coord_offset;
132 int j_index_start,j_index_end;
133 int nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
134 real shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
135 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
136 real *shiftvec,*fshift,*x,*f;
137 /* #for I in PARTICLES_I */
139 real ix{I},iy{I},iz{I},fix{I},fiy{I},fiz{I},iq{I},isai{I};
141 /* #for J in PARTICLES_J */
143 real jx{J},jy{J},jz{J},fjx{J},fjy{J},fjz{J},jq{J},isaj{J};
145 /* #for I,J in PAIRS_IJ */
146 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};
148 /* #if KERNEL_ELEC != 'None' */
149 real velec,felec,velecsum,facel,crf,krf,krf2;
152 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
154 real vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,dvdatmp;
155 real *invsqrta,*dvda,*gbtab;
157 /* #if KERNEL_VDW != 'None' */
159 real rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
163 /* #if 'Table' in KERNEL_ELEC or 'GeneralizedBorn' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
165 real rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
168 /* #if 'Ewald' in KERNEL_ELEC */
170 real ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
173 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
174 real rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
181 real hybscal; /* the multiplicator to the force for hybrid interactions*/
182 gmx_bool bHybrid; /*Are we in the hybrid zone ?*/
187 force_cap = fr->adress_ex_forcecap;
194 jindex = nlist->jindex;
196 shiftidx = nlist->shift;
198 shiftvec = fr->shift_vec[0];
199 fshift = fr->fshift[0];
200 /* #if KERNEL_ELEC != 'None' */
202 charge = mdatoms->chargeA;
203 /* #if 'ReactionField' in KERNEL_ELEC */
209 /* #if KERNEL_VDW != 'None' */
210 nvdwtype = fr->ntype;
212 vdwtype = mdatoms->typeA;
215 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
216 vftab = kernel_data->table_elec_vdw->data;
217 vftabscale = kernel_data->table_elec_vdw->scale;
218 /* #elif 'Table' in KERNEL_ELEC */
219 vftab = kernel_data->table_elec->data;
220 vftabscale = kernel_data->table_elec->scale;
221 /* #elif 'Table' in KERNEL_VDW */
222 vftab = kernel_data->table_vdw->data;
223 vftabscale = kernel_data->table_vdw->scale;
226 /* #if 'Ewald' in KERNEL_ELEC */
227 sh_ewald = fr->ic->sh_ewald;
228 /* #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
229 ewtab = fr->ic->tabq_coul_F;
230 ewtabscale = fr->ic->tabq_scale;
231 ewtabhalfspace = 0.5/ewtabscale;
233 ewtab = fr->ic->tabq_coul_FDV0;
234 ewtabscale = fr->ic->tabq_scale;
235 ewtabhalfspace = 0.5/ewtabscale;
239 /* #if KERNEL_ELEC=='GeneralizedBorn' */
240 invsqrta = fr->invsqrta;
242 gbtabscale = fr->gbtab.scale;
243 gbtab = fr->gbtab.data;
244 gbinvepsdiff = (1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent);
247 /* #if 'Water' in GEOMETRY_I */
248 /* Setup water-specific parameters */
249 inr = nlist->iinr[0];
250 /* #for I in PARTICLES_ELEC_I */
251 iq{I} = facel*charge[inr+{I}];
253 /* #for I in PARTICLES_VDW_I */
254 vdwioffset{I} = {NVDWPARAM}*nvdwtype*vdwtype[inr+{I}];
258 /* #if 'Water' in GEOMETRY_J */
259 /* #for J in PARTICLES_ELEC_J */
260 jq{J} = charge[inr+{J}];
262 /* #for J in PARTICLES_VDW_J */
263 vdwjidx{J} = {NVDWPARAM}*vdwtype[inr+{J}];
265 /* #for I,J in PAIRS_IJ */
266 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
267 qq{I}{J} = iq{I}*jq{J};
269 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
270 /* #if 'Buckingham' in KERNEL_VDW */
271 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
272 cexp1_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
273 cexp2_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+2];
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 weight_cg1 = wf[inr];
357 /* ## For water we already preloaded parameters at the start of the kernel */
358 /* #if not 'Water' in GEOMETRY_I */
359 /* Load parameters for i particles */
360 /* #for I in PARTICLES_ELEC_I */
361 iq{I} = facel*charge[inr+{I}];
362 /* #define OUTERFLOPS OUTERFLOPS+1 */
363 /* #if KERNEL_ELEC=='GeneralizedBorn' */
364 isai{I} = invsqrta[inr+{I}];
367 /* #for I in PARTICLES_VDW_I */
368 vdwioffset{I} = {NVDWPARAM}*nvdwtype*vdwtype[inr+{I}];
372 /* #if 'Potential' in KERNEL_VF */
373 /* Reset potential sums */
374 /* #if KERNEL_ELEC != 'None' */
377 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
380 /* #if KERNEL_VDW != 'None' */
384 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
388 /* Start inner kernel loop */
389 for(jidx=j_index_start; jidx<j_index_end; jidx++)
391 /* Get j neighbor index, and coordinate index */
393 weight_cg2 = wf[jnr];
394 weight_product = weight_cg1*weight_cg2;
395 if (weight_product < ALMOST_ZERO) {
396 /* #if KERNEL_RESOLUTION=='CG' */
402 else if (weight_product >= ALMOST_ONE)
404 /* #if KERNEL_RESOLUTION=='CG' */
412 /* #if KERNEL_RESOLUTION=='CG' */
413 hybscal = 1.0 - weight_product;
415 hybscal = weight_product;
418 j_coord_offset = DIM*jnr;
420 /* load j atom coordinates */
421 /* #for J in PARTICLES_J */
422 jx{J} = x[j_coord_offset+DIM*{J}+XX];
423 jy{J} = x[j_coord_offset+DIM*{J}+YY];
424 jz{J} = x[j_coord_offset+DIM*{J}+ZZ];
427 /* Calculate displacement vector */
428 /* #for I,J in PAIRS_IJ */
429 dx{I}{J} = ix{I} - jx{J};
430 dy{I}{J} = iy{I} - jy{J};
431 dz{I}{J} = iz{I} - jz{J};
432 /* #define INNERFLOPS INNERFLOPS+3 */
435 /* Calculate squared distance and things based on it */
436 /* #for I,J in PAIRS_IJ */
437 rsq{I}{J} = dx{I}{J}*dx{I}{J}+dy{I}{J}*dy{I}{J}+dz{I}{J}*dz{I}{J};
438 /* #define INNERFLOPS INNERFLOPS+5 */
441 /* #for I,J in PAIRS_IJ */
442 /* #if 'rinv' in INTERACTION_FLAGS[I][J] */
443 rinv{I}{J} = gmx_invsqrt(rsq{I}{J});
444 /* #define INNERFLOPS INNERFLOPS+5 */
448 /* #for I,J in PAIRS_IJ */
449 /* #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
450 /* # if 'rinv' not in INTERACTION_FLAGS[I][J] */
451 rinvsq{I}{J} = 1.0/rsq{I}{J};
452 /* #define INNERFLOPS INNERFLOPS+4 */
454 rinvsq{I}{J} = rinv{I}{J}*rinv{I}{J};
455 /* #define INNERFLOPS INNERFLOPS+1 */
460 /* #if not 'Water' in GEOMETRY_J */
461 /* Load parameters for j particles */
462 /* #for J in PARTICLES_ELEC_J */
463 jq{J} = charge[jnr+{J}];
464 /* #if KERNEL_ELEC=='GeneralizedBorn' */
465 isaj{J} = invsqrta[jnr+{J}];
468 /* #for J in PARTICLES_VDW_J */
469 vdwjidx{J} = {NVDWPARAM}*vdwtype[jnr+{J}];
473 /* #for I,J in PAIRS_IJ */
475 /**************************
476 * CALCULATE INTERACTIONS *
477 **************************/
479 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
480 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
481 /* ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
482 if (rsq{I}{J}<rcutoff2)
484 /* #if 0 ## this and the next two lines is a hack to maintain auto-indentation in template file */
489 /* #if 'r' in INTERACTION_FLAGS[I][J] */
490 r{I}{J} = rsq{I}{J}*rinv{I}{J};
491 /* #define INNERFLOPS INNERFLOPS+1 */
494 /* ## For water geometries we already loaded parameters at the start of the kernel */
495 /* #if not 'Water' in GEOMETRY_J */
496 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
497 qq{I}{J} = iq{I}*jq{J};
498 /* #define INNERFLOPS INNERFLOPS+1 */
500 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
501 /* #if KERNEL_VDW=='Buckingham' */
502 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
503 cexp1_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
504 cexp2_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+2];
506 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
507 c12_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
512 /* #if 'table' in INTERACTION_FLAGS[I][J] */
513 /* Calculate table index by multiplying r with table scale and truncate to integer */
514 rt = r{I}{J}*vftabscale;
517 vfitab = {TABLE_NINTERACTIONS}*{TABLE_POINT_SIZE}*vfitab;
518 /* #define INNERFLOPS INNERFLOPS+2 */
521 /* ## ELECTROSTATIC INTERACTIONS */
522 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
524 /* #if KERNEL_ELEC=='Coulomb' */
526 /* COULOMB ELECTROSTATICS */
527 velec = qq{I}{J}*rinv{I}{J};
528 /* #define INNERFLOPS INNERFLOPS+1 */
529 /* #if 'Force' in KERNEL_VF */
530 felec = velec*rinvsq{I}{J};
531 /* #define INNERFLOPS INNERFLOPS+2 */
534 /* #elif KERNEL_ELEC=='ReactionField' */
536 /* REACTION-FIELD ELECTROSTATICS */
537 /* #if 'Potential' in KERNEL_VF */
538 velec = qq{I}{J}*(rinv{I}{J}+krf*rsq{I}{J}-crf);
539 /* #define INNERFLOPS INNERFLOPS+4 */
541 /* #if 'Force' in KERNEL_VF */
542 felec = qq{I}{J}*(rinv{I}{J}*rinvsq{I}{J}-krf2);
543 /* #define INNERFLOPS INNERFLOPS+3 */
546 /* #elif KERNEL_ELEC=='GeneralizedBorn' */
548 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
549 isaprod = isai{I}*isaj{J};
550 gbqqfactor = isaprod*(-qq{I}{J})*gbinvepsdiff;
551 gbscale = isaprod*gbtabscale;
552 dvdaj = dvda[jnr+{J}];
553 /* #define INNERFLOPS INNERFLOPS+5 */
555 /* Calculate generalized born table index - this is a separate table from the normal one,
556 * but we use the same procedure by multiplying r with scale and truncating to integer.
558 rt = r{I}{J}*gbscale;
565 Geps = gbeps*gbtab[gbitab+2];
566 Heps2 = gbeps*gbeps*gbtab[gbitab+3];
570 /* #define INNERFLOPS INNERFLOPS+10 */
572 /* #if 'Force' in KERNEL_VF */
573 FF = Fp+Geps+2.0*Heps2;
574 fgb = gbqqfactor*FF*gbscale;
575 dvdatmp = -0.5*(vgb+fgb*r{I}{J});
576 dvdasum = dvdasum + dvdatmp;
577 dvda[jnr] = dvdaj+dvdatmp*isaj{J}*isaj{J};
578 /* #define INNERFLOPS INNERFLOPS+13 */
580 velec = qq{I}{J}*rinv{I}{J};
581 /* #define INNERFLOPS INNERFLOPS+1 */
582 /* #if 'Force' in KERNEL_VF */
583 felec = (velec*rinv{I}{J}-fgb)*rinv{I}{J};
584 /* #define INNERFLOPS INNERFLOPS+3 */
587 /* #elif KERNEL_ELEC=='Ewald' */
588 /* EWALD ELECTROSTATICS */
590 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
591 ewrt = r{I}{J}*ewtabscale;
594 /* #define INNERFLOPS INNERFLOPS+2 */
595 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
597 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
598 /* #define INNERFLOPS INNERFLOPS+4 */
599 /* #if KERNEL_MOD_ELEC=='PotentialShift' */
600 velec = qq{I}{J}*((rinv{I}{J}-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
601 /* #define INNERFLOPS INNERFLOPS+7 */
603 velec = qq{I}{J}*(rinv{I}{J}-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
604 /* #define INNERFLOPS INNERFLOPS+6 */
606 /* #if 'Force' in KERNEL_VF */
607 felec = qq{I}{J}*rinv{I}{J}*(rinvsq{I}{J}-felec);
608 /* #define INNERFLOPS INNERFLOPS+3 */
610 /* #elif KERNEL_VF=='Force' */
611 felec = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
612 felec = qq{I}{J}*rinv{I}{J}*(rinvsq{I}{J}-felec);
613 /* #define INNERFLOPS INNERFLOPS+7 */
616 /* #elif KERNEL_ELEC=='CubicSplineTable' */
618 /* CUBIC SPLINE TABLE ELECTROSTATICS */
619 /* #if 'Potential' in KERNEL_VF */
623 Geps = vfeps*vftab[vfitab+2];
624 Heps2 = vfeps*vfeps*vftab[vfitab+3];
626 /* #define INNERFLOPS INNERFLOPS+5 */
627 /* #if 'Potential' in KERNEL_VF */
630 /* #define INNERFLOPS INNERFLOPS+3 */
632 /* #if 'Force' in KERNEL_VF */
633 FF = Fp+Geps+2.0*Heps2;
634 felec = -qq{I}{J}*FF*vftabscale*rinv{I}{J};
635 /* #define INNERFLOPS INNERFLOPS+7 */
638 /* ## End of check for electrostatics interaction forms */
640 /* ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
642 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
644 /* #if KERNEL_VDW=='LennardJones' */
646 /* LENNARD-JONES DISPERSION/REPULSION */
648 rinvsix = rinvsq{I}{J}*rinvsq{I}{J}*rinvsq{I}{J};
649 /* #define INNERFLOPS INNERFLOPS+2 */
650 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
651 vvdw6 = c6_{I}{J}*rinvsix;
652 vvdw12 = c12_{I}{J}*rinvsix*rinvsix;
653 /* #define INNERFLOPS INNERFLOPS+3 */
654 /* #if KERNEL_MOD_VDW=='PotentialShift' */
655 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);
656 /* #define INNERFLOPS INNERFLOPS+8 */
658 vvdw = vvdw12*(1.0/12.0) - vvdw6*(1.0/6.0);
659 /* #define INNERFLOPS INNERFLOPS+3 */
661 /* ## Check for force inside potential check, i.e. this means we already did the potential part */
662 /* #if 'Force' in KERNEL_VF */
663 fvdw = (vvdw12-vvdw6)*rinvsq{I}{J};
664 /* #define INNERFLOPS INNERFLOPS+2 */
666 /* #elif KERNEL_VF=='Force' */
667 /* ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
668 fvdw = (c12_{I}{J}*rinvsix-c6_{I}{J})*rinvsix*rinvsq{I}{J};
669 /* #define INNERFLOPS INNERFLOPS+4 */
672 /* #elif KERNEL_VDW=='Buckingham' */
674 /* BUCKINGHAM DISPERSION/REPULSION */
675 rinvsix = rinvsq{I}{J}*rinvsq{I}{J}*rinvsq{I}{J};
676 vvdw6 = c6_{I}{J}*rinvsix;
677 br = cexp2_{I}{J}*r{I}{J};
678 vvdwexp = cexp1_{I}{J}*exp(-br);
679 /* ## Estimate exp() to 25 flops */
680 /* #define INNERFLOPS INNERFLOPS+31 */
681 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
682 /* #if KERNEL_MOD_VDW=='PotentialShift' */
683 vvdw = (vvdwexp-cexp1_{I}{J}*exp(-cexp2_{I}{J}*rvdw)) - (vvdw6 - c6_{I}{J}*sh_vdw_invrcut6)*(1.0/6.0);
684 /* #define INNERFLOPS INNERFLOPS+33 */
686 vvdw = vvdwexp - vvdw6*(1.0/6.0);
687 /* #define INNERFLOPS INNERFLOPS+2 */
690 /* #if 'Force' in KERNEL_VF */
691 fvdw = (br*vvdwexp-vvdw6)*rinvsq{I}{J};
692 /* #define INNERFLOPS INNERFLOPS+3 */
695 /* #elif KERNEL_VDW=='CubicSplineTable' */
697 /* CUBIC SPLINE TABLE DISPERSION */
698 vfitab += {TABLE_VDW_OFFSET};
699 /* #if 'Potential' in KERNEL_VF */
703 Geps = vfeps*vftab[vfitab+2];
704 Heps2 = vfeps*vfeps*vftab[vfitab+3];
706 /* #define INNERFLOPS INNERFLOPS+5 */
707 /* #if 'Potential' in KERNEL_VF */
709 vvdw6 = c6_{I}{J}*VV;
710 /* #define INNERFLOPS INNERFLOPS+3 */
712 /* #if 'Force' in KERNEL_VF */
713 FF = Fp+Geps+2.0*Heps2;
714 fvdw6 = c6_{I}{J}*FF;
715 /* #define INNERFLOPS INNERFLOPS+4 */
718 /* CUBIC SPLINE TABLE REPULSION */
719 /* #if 'Potential' in KERNEL_VF */
723 Geps = vfeps*vftab[vfitab+6];
724 Heps2 = vfeps*vfeps*vftab[vfitab+7];
726 /* #define INNERFLOPS INNERFLOPS+5 */
727 /* #if 'Potential' in KERNEL_VF */
729 vvdw12 = c12_{I}{J}*VV;
730 /* #define INNERFLOPS INNERFLOPS+3 */
732 /* #if 'Force' in KERNEL_VF */
733 FF = Fp+Geps+2.0*Heps2;
734 fvdw12 = c12_{I}{J}*FF;
735 /* #define INNERFLOPS INNERFLOPS+4 */
737 /* #if 'Potential' in KERNEL_VF */
739 /* #define INNERFLOPS INNERFLOPS+1 */
741 /* #if 'Force' in KERNEL_VF */
742 fvdw = -(fvdw6+fvdw12)*vftabscale*rinv{I}{J};
743 /* #define INNERFLOPS INNERFLOPS+4 */
746 /* ## End of check for vdw interaction forms */
748 /* ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
750 /* #if 'switch' in INTERACTION_FLAGS[I][J] */
752 d = (d>0.0) ? d : 0.0;
754 sw = 1.0+d2*d*(swV3+d*(swV4+d*swV5));
755 /* #define INNERFLOPS INNERFLOPS+9 */
757 /* #if 'Force' in KERNEL_VF */
758 dsw = d2*(swF2+d*(swF3+d*swF4));
759 /* #define INNERFLOPS INNERFLOPS+5 */
762 /* Evaluate switch function */
763 /* #if 'Force' in KERNEL_VF */
764 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
765 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
766 felec = felec*sw - rinv{I}{J}*velec*dsw;
767 /* #define INNERFLOPS INNERFLOPS+3 */
769 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
770 fvdw = fvdw*sw - rinv{I}{J}*vvdw*dsw;
771 /* #define INNERFLOPS INNERFLOPS+3 */
774 /* #if 'Potential' in KERNEL_VF */
775 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
777 /* #define INNERFLOPS INNERFLOPS+1 */
779 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
781 /* #define INNERFLOPS INNERFLOPS+1 */
786 /* #if 'Potential' in KERNEL_VF */
787 /* Update potential sums from outer loop */
788 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
790 /* #define INNERFLOPS INNERFLOPS+1 */
791 /* #if KERNEL_ELEC=='GeneralizedBorn' */
793 /* #define INNERFLOPS INNERFLOPS+1 */
796 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
798 /* #define INNERFLOPS INNERFLOPS+1 */
802 /* #if 'Force' in KERNEL_VF */
804 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
806 /* #define INNERFLOPS INNERFLOPS+1 */
807 /* #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
809 /* #elif 'vdw' in INTERACTION_FLAGS[I][J] */
813 /* #if KERNEL_RESOLUTION=='CG' */
814 if(force_cap>0 && (fabs(fscal)> force_cap))
816 fscal=force_cap*fscal/fabs(fscal);
821 /* Calculate temporary vectorial force */
826 /* Update vectorial force */
830 f[j_coord_offset+DIM*{J}+XX] -= tx;
831 f[j_coord_offset+DIM*{J}+YY] -= ty;
832 f[j_coord_offset+DIM*{J}+ZZ] -= tz;
834 /* #define INNERFLOPS INNERFLOPS+9 */
837 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
838 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
839 /* #if 0 ## This and next two lines is a hack to maintain indentation in template file */
844 /* ## End of check for the interaction being outside the cutoff */
847 /* ## End of loop over i-j interaction pairs */
849 /* Inner loop uses {INNERFLOPS} flops */
851 /* End of innermost loop */
853 /* #if 'Force' in KERNEL_VF */
855 /* #for I in PARTICLES_I */
856 f[i_coord_offset+DIM*{I}+XX] += fix{I};
857 f[i_coord_offset+DIM*{I}+YY] += fiy{I};
858 f[i_coord_offset+DIM*{I}+ZZ] += fiz{I};
862 /* #define OUTERFLOPS OUTERFLOPS+6 */
864 fshift[i_shift_offset+XX] += tx;
865 fshift[i_shift_offset+YY] += ty;
866 fshift[i_shift_offset+ZZ] += tz;
867 /* #define OUTERFLOPS OUTERFLOPS+3 */
870 /* #if 'Potential' in KERNEL_VF */
872 /* Update potential energies */
873 /* #if KERNEL_ELEC != 'None' */
874 kernel_data->energygrp_elec[ggid] += velecsum;
875 /* #define OUTERFLOPS OUTERFLOPS+1 */
877 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
878 kernel_data->energygrp_polarization[ggid] += vgbsum;
879 /* #define OUTERFLOPS OUTERFLOPS+1 */
881 /* #if KERNEL_VDW != 'None' */
882 kernel_data->energygrp_vdw[ggid] += vvdwsum;
883 /* #define OUTERFLOPS OUTERFLOPS+1 */
886 /* #if 'GeneralizedBorn' in KERNEL_ELEC and 'Force' in KERNEL_VF */
887 dvda[inr] = dvda[inr] + dvdasum*isai{I}*isai{I};
890 /* Increment number of inner iterations */
891 inneriter += j_index_end - j_index_start;
893 /* Outer loop uses {OUTERFLOPS} flops */
896 /* Increment number of outer iterations */
899 /* Update outer/inner flops */
900 /* ## NB: This is not important, it just affects the flopcount. However, since our preprocessor is */
901 /* ## primitive and replaces aggressively even in strings inside these directives, we need to */
902 /* ## assemble the main part of the name (containing KERNEL/ELEC/VDW) directly in the source. */
903 /* #if GEOMETRY_I == 'Water3' */
904 /* #define ISUFFIX '_W3' */
905 /* #elif GEOMETRY_I == 'Water4' */
906 /* #define ISUFFIX '_W4' */
908 /* #define ISUFFIX '' */
910 /* #if GEOMETRY_J == 'Water3' */
911 /* #define JSUFFIX 'W3' */
912 /* #elif GEOMETRY_J == 'Water4' */
913 /* #define JSUFFIX 'W4' */
915 /* #define JSUFFIX '' */
917 /* #if 'PotentialAndForce' in KERNEL_VF */
918 /* #define VFSUFFIX '_VF' */
919 /* #elif 'Potential' in KERNEL_VF */
920 /* #define VFSUFFIX '_V' */
922 /* #define VFSUFFIX '_F' */
925 /* #if KERNEL_ELEC != 'None' and KERNEL_VDW != 'None' */
926 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
927 /* #elif KERNEL_ELEC != 'None' */
928 inc_nrnb(nrnb,eNR_NBKERNEL_ELEC{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});
930 inc_nrnb(nrnb,eNR_NBKERNEL_VDW{ISUFFIX}{JSUFFIX}{VFSUFFIX},outeriter*{OUTERFLOPS} + inneriter*{INNERFLOPS});