3 * This file is part of the GROMACS molecular simulation package.
5 * Copyright (c) 2012,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 #define ALMOST_ZERO 1e-30
50 #define ALMOST_ONE 1-(1e-30)
52 /* ## List of variables set by the generating script: */
54 /* ## Setttings that apply to the entire kernel: */
55 /* ## KERNEL_ELEC: String, choice for electrostatic interactions */
56 /* ## KERNEL_VDW: String, choice for van der Waals interactions */
57 /* ## KERNEL_NAME: String, name of this kernel */
58 /* ## KERNEL_RESOLUTION: String, choice for resoultion */
59 /* ## KERNEL_VF: String telling if we calculate potential, force, or both */
60 /* ## GEOMETRY_I/GEOMETRY_J: String, name of each geometry, e.g. 'Water3' or '1Particle' */
62 /* ## Setttings that apply to particles in the outer (I) or inner (J) loops: */
63 /* ## PARTICLES_I[]/ Arrays with lists of i/j particles to use in kernel. It is */
64 /* ## PARTICLES_J[]: just [0] for particle geometry, but can be longer for water */
65 /* ## PARTICLES_ELEC_I[]/ Arrays with lists of i/j particle that have electrostatics */
66 /* ## PARTICLES_ELEC_J[]: interactions that should be calculated in this kernel. */
67 /* ## PARTICLES_VDW_I[]/ Arrays with the list of i/j particle that have VdW */
68 /* ## PARTICLES_VDW_J[]: interactions that should be calculated in this kernel. */
70 /* ## Setttings for pairs of interactions (e.g. 2nd i particle against 1st j particle) */
71 /* ## PAIRS_IJ[]: Array with (i,j) tuples of pairs for which interactions */
72 /* ## should be calculated in this kernel. Zero-charge particles */
73 /* ## do not have interactions with particles without vdw, and */
74 /* ## Vdw-only interactions are not evaluated in a no-vdw-kernel. */
75 /* ## INTERACTION_FLAGS[][]: 2D matrix, dimension e.g. 3*3 for water-water interactions. */
76 /* ## For each i-j pair, the element [I][J] is a list of strings */
77 /* ## defining properties/flags of this interaction. Examples */
78 /* ## include 'electrostatics'/'vdw' if that type of interaction */
79 /* ## should be evaluated, 'rsq'/'rinv'/'rinvsq' if those values */
80 /* ## are needed, and 'exactcutoff' or 'shift','switch' to */
81 /* ## decide if the force/potential should be modified. This way */
82 /* ## we only calculate values absolutely needed for each case. */
84 /* ## Calculate the size and offset for (merged/interleaved) table data */
86 /* #if ('CubicSplineTable' in [KERNEL_ELEC,KERNEL_VDW]) or KERNEL_VF=='PotentialAndForce' */
87 /* #define TABLE_POINT_SIZE 4 */
89 /* #define TABLE_POINT_SIZE 2 */
92 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
93 /* #define TABLE_NINTERACTIONS 3 */
94 /* #define TABLE_VDW_OFFSET TABLE_POINT_SIZE */
95 /* #elif 'Table' in KERNEL_ELEC */
96 /* #define TABLE_NINTERACTIONS 1 */
97 /* #elif 'Table' in KERNEL_VDW */
98 /* #define TABLE_NINTERACTIONS 2 */
99 /* #define TABLE_VDW_OFFSET 0 */
101 /* #define TABLE_NINTERACTIONS 0 */
104 /* #if 'Buckingham' in KERNEL_VDW */
105 /* #define NVDWPARAM 3 */
107 /* #define NVDWPARAM 2 */
111 * Gromacs nonbonded kernel: {KERNEL_NAME}
112 * Electrostatics interaction: {KERNEL_ELEC}
113 * VdW interaction: {KERNEL_VDW}
114 * Geometry: {GEOMETRY_I}-{GEOMETRY_J}
115 * Calculate force/pot: {KERNEL_VF}
119 (t_nblist * gmx_restrict nlist,
120 rvec * gmx_restrict xx,
121 rvec * gmx_restrict ff,
122 t_forcerec * gmx_restrict fr,
123 t_mdatoms * gmx_restrict mdatoms,
124 nb_kernel_data_t * gmx_restrict kernel_data,
125 t_nrnb * gmx_restrict nrnb)
127 /* ## Not all variables are used for all kernels, but any optimizing compiler fixes that, */
128 /* ## so there is no point in going to extremes to exclude variables that are not needed. */
129 int i_shift_offset,i_coord_offset,j_coord_offset;
130 int j_index_start,j_index_end;
131 int nri,inr,ggid,iidx,jidx,jnr,outeriter,inneriter;
132 real shX,shY,shZ,tx,ty,tz,fscal,rcutoff,rcutoff2;
133 int *iinr,*jindex,*jjnr,*shiftidx,*gid;
134 real *shiftvec,*fshift,*x,*f;
135 /* #for I in PARTICLES_I */
137 real ix{I},iy{I},iz{I},fix{I},fiy{I},fiz{I},iq{I},isai{I};
139 /* #for J in PARTICLES_J */
141 real jx{J},jy{J},jz{J},fjx{J},fjy{J},fjz{J},jq{J},isaj{J};
143 /* #for I,J in PAIRS_IJ */
144 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};
146 /* #if KERNEL_ELEC != 'None' */
147 real velec,felec,velecsum,facel,crf,krf,krf2;
150 /* #if 'GeneralizedBorn' in KERNEL_ELEC */
152 real vgb,fgb,vgbsum,dvdasum,gbscale,gbtabscale,isaprod,gbqqfactor,gbinvepsdiff,dvdaj,gbeps,dvdatmp;
153 real *invsqrta,*dvda,*gbtab;
155 /* #if KERNEL_VDW != 'None' */
157 real rinvsix,rvdw,vvdw,vvdw6,vvdw12,fvdw,fvdw6,fvdw12,vvdwsum,br,vvdwexp,sh_vdw_invrcut6;
161 /* #if 'Table' in KERNEL_ELEC or 'GeneralizedBorn' in KERNEL_ELEC or 'Table' in KERNEL_VDW */
163 real rt,vfeps,vftabscale,Y,F,Geps,Heps2,Fp,VV,FF;
166 /* #if 'Ewald' in KERNEL_ELEC */
168 real ewtabscale,eweps,sh_ewald,ewrt,ewtabhalfspace;
171 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
172 real rswitch,swV3,swV4,swV5,swF2,swF3,swF4,d,d2,sw,dsw;
179 real hybscal; /* the multiplicator to the force for hybrid interactions*/
180 gmx_bool bHybrid; /*Are we in the hybrid zone ?*/
185 force_cap = fr->adress_ex_forcecap;
192 jindex = nlist->jindex;
194 shiftidx = nlist->shift;
196 shiftvec = fr->shift_vec[0];
197 fshift = fr->fshift[0];
198 /* #if KERNEL_ELEC != 'None' */
200 charge = mdatoms->chargeA;
201 /* #if 'ReactionField' in KERNEL_ELEC */
207 /* #if KERNEL_VDW != 'None' */
208 nvdwtype = fr->ntype;
210 vdwtype = mdatoms->typeA;
213 /* #if 'Table' in KERNEL_ELEC and 'Table' in KERNEL_VDW */
214 vftab = kernel_data->table_elec_vdw->data;
215 vftabscale = kernel_data->table_elec_vdw->scale;
216 /* #elif 'Table' in KERNEL_ELEC */
217 vftab = kernel_data->table_elec->data;
218 vftabscale = kernel_data->table_elec->scale;
219 /* #elif 'Table' in KERNEL_VDW */
220 vftab = kernel_data->table_vdw->data;
221 vftabscale = kernel_data->table_vdw->scale;
224 /* #if 'Ewald' in KERNEL_ELEC */
225 sh_ewald = fr->ic->sh_ewald;
226 /* #if KERNEL_VF=='Force' and KERNEL_MOD_ELEC!='PotentialSwitch' */
227 ewtab = fr->ic->tabq_coul_F;
228 ewtabscale = fr->ic->tabq_scale;
229 ewtabhalfspace = 0.5/ewtabscale;
231 ewtab = fr->ic->tabq_coul_FDV0;
232 ewtabscale = fr->ic->tabq_scale;
233 ewtabhalfspace = 0.5/ewtabscale;
237 /* #if KERNEL_ELEC=='GeneralizedBorn' */
238 invsqrta = fr->invsqrta;
240 gbtabscale = fr->gbtab.scale;
241 gbtab = fr->gbtab.data;
242 gbinvepsdiff = (1.0/fr->epsilon_r) - (1.0/fr->gb_epsilon_solvent);
245 /* #if 'Water' in GEOMETRY_I */
246 /* Setup water-specific parameters */
247 inr = nlist->iinr[0];
248 /* #for I in PARTICLES_ELEC_I */
249 iq{I} = facel*charge[inr+{I}];
251 /* #for I in PARTICLES_VDW_I */
252 vdwioffset{I} = {NVDWPARAM}*nvdwtype*vdwtype[inr+{I}];
256 /* #if 'Water' in GEOMETRY_J */
257 /* #for J in PARTICLES_ELEC_J */
258 jq{J} = charge[inr+{J}];
260 /* #for J in PARTICLES_VDW_J */
261 vdwjidx{J} = {NVDWPARAM}*vdwtype[inr+{J}];
263 /* #for I,J in PAIRS_IJ */
264 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
265 qq{I}{J} = iq{I}*jq{J};
267 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
268 /* #if 'Buckingham' in KERNEL_VDW */
269 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
270 cexp1_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
271 cexp2_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+2];
273 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
274 c12_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
280 /* #if KERNEL_MOD_ELEC!='None' or KERNEL_MOD_VDW!='None' */
281 /* #if KERNEL_ELEC!='None' */
282 /* When we use explicit cutoffs the value must be identical for elec and VdW, so use elec as an arbitrary choice */
283 rcutoff = fr->rcoulomb;
287 rcutoff2 = rcutoff*rcutoff;
290 /* #if KERNEL_MOD_VDW=='PotentialShift' */
291 sh_vdw_invrcut6 = fr->ic->sh_invrc6;
295 /* #if 'PotentialSwitch' in [KERNEL_MOD_ELEC,KERNEL_MOD_VDW] */
296 /* #if KERNEL_MOD_ELEC=='PotentialSwitch' */
297 rswitch = fr->rcoulomb_switch;
299 rswitch = fr->rvdw_switch;
301 /* Setup switch parameters */
303 swV3 = -10.0/(d*d*d);
304 swV4 = 15.0/(d*d*d*d);
305 swV5 = -6.0/(d*d*d*d*d);
306 /* #if 'Force' in KERNEL_VF */
307 swF2 = -30.0/(d*d*d);
308 swF3 = 60.0/(d*d*d*d);
309 swF4 = -30.0/(d*d*d*d*d);
313 /* ## Keep track of the floating point operations we issue for reporting! */
314 /* #define OUTERFLOPS 0 */
315 /* #define INNERFLOPS 0 */
319 /* Start outer loop over neighborlists */
320 for(iidx=0; iidx<nri; iidx++)
322 /* Load shift vector for this list */
323 i_shift_offset = DIM*shiftidx[iidx];
324 shX = shiftvec[i_shift_offset+XX];
325 shY = shiftvec[i_shift_offset+YY];
326 shZ = shiftvec[i_shift_offset+ZZ];
328 /* Load limits for loop over neighbors */
329 j_index_start = jindex[iidx];
330 j_index_end = jindex[iidx+1];
332 /* Get outer coordinate index */
334 i_coord_offset = DIM*inr;
336 /* Load i particle coords and add shift vector */
337 /* ## Loop over i particles, but only include ones that we use - skip e.g. vdw-only sites for elec-only kernel */
338 /* #for I in PARTICLES_I */
339 ix{I} = shX + x[i_coord_offset+DIM*{I}+XX];
340 iy{I} = shY + x[i_coord_offset+DIM*{I}+YY];
341 iz{I} = shZ + x[i_coord_offset+DIM*{I}+ZZ];
342 /* #define OUTERFLOPS OUTERFLOPS+3 */
345 /* #if 'Force' in KERNEL_VF */
346 /* #for I in PARTICLES_I */
353 weight_cg1 = wf[inr];
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 weight_cg2 = wf[jnr];
392 weight_product = weight_cg1*weight_cg2;
393 if (weight_product < ALMOST_ZERO) {
394 /* #if KERNEL_RESOLUTION=='CG' */
400 else if (weight_product >= ALMOST_ONE)
402 /* #if KERNEL_RESOLUTION=='CG' */
410 /* #if KERNEL_RESOLUTION=='CG' */
411 hybscal = 1.0 - weight_product;
413 hybscal = weight_product;
416 j_coord_offset = DIM*jnr;
418 /* load j atom coordinates */
419 /* #for J in PARTICLES_J */
420 jx{J} = x[j_coord_offset+DIM*{J}+XX];
421 jy{J} = x[j_coord_offset+DIM*{J}+YY];
422 jz{J} = x[j_coord_offset+DIM*{J}+ZZ];
425 /* Calculate displacement vector */
426 /* #for I,J in PAIRS_IJ */
427 dx{I}{J} = ix{I} - jx{J};
428 dy{I}{J} = iy{I} - jy{J};
429 dz{I}{J} = iz{I} - jz{J};
430 /* #define INNERFLOPS INNERFLOPS+3 */
433 /* Calculate squared distance and things based on it */
434 /* #for I,J in PAIRS_IJ */
435 rsq{I}{J} = dx{I}{J}*dx{I}{J}+dy{I}{J}*dy{I}{J}+dz{I}{J}*dz{I}{J};
436 /* #define INNERFLOPS INNERFLOPS+5 */
439 /* #for I,J in PAIRS_IJ */
440 /* #if 'rinv' in INTERACTION_FLAGS[I][J] */
441 rinv{I}{J} = gmx_invsqrt(rsq{I}{J});
442 /* #define INNERFLOPS INNERFLOPS+5 */
446 /* #for I,J in PAIRS_IJ */
447 /* #if 'rinvsq' in INTERACTION_FLAGS[I][J] */
448 /* # if 'rinv' not in INTERACTION_FLAGS[I][J] */
449 rinvsq{I}{J} = 1.0/rsq{I}{J};
450 /* #define INNERFLOPS INNERFLOPS+4 */
452 rinvsq{I}{J} = rinv{I}{J}*rinv{I}{J};
453 /* #define INNERFLOPS INNERFLOPS+1 */
458 /* #if not 'Water' in GEOMETRY_J */
459 /* Load parameters for j particles */
460 /* #for J in PARTICLES_ELEC_J */
461 jq{J} = charge[jnr+{J}];
462 /* #if KERNEL_ELEC=='GeneralizedBorn' */
463 isaj{J} = invsqrta[jnr+{J}];
466 /* #for J in PARTICLES_VDW_J */
467 vdwjidx{J} = {NVDWPARAM}*vdwtype[jnr+{J}];
471 /* #for I,J in PAIRS_IJ */
473 /**************************
474 * CALCULATE INTERACTIONS *
475 **************************/
477 /* ## Note special check for TIP4P-TIP4P. Since we are cutting of all hydrogen interactions we also cut the LJ-only O-O interaction */
478 /* #if 'exactcutoff' in INTERACTION_FLAGS[I][J] or (GEOMETRY_I=='Water4' and GEOMETRY_J=='Water4' and 'exactcutoff' in INTERACTION_FLAGS[1][1]) */
479 /* ## We always calculate rinv/rinvsq above to enable pipelineing in compilers (performance tested on x86) */
480 if (rsq{I}{J}<rcutoff2)
482 /* #if 0 ## this and the next two lines is a hack to maintain auto-indentation in template file */
487 /* #if 'r' in INTERACTION_FLAGS[I][J] */
488 r{I}{J} = rsq{I}{J}*rinv{I}{J};
489 /* #define INNERFLOPS INNERFLOPS+1 */
492 /* ## For water geometries we already loaded parameters at the start of the kernel */
493 /* #if not 'Water' in GEOMETRY_J */
494 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
495 qq{I}{J} = iq{I}*jq{J};
496 /* #define INNERFLOPS INNERFLOPS+1 */
498 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
499 /* #if KERNEL_VDW=='Buckingham' */
500 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
501 cexp1_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
502 cexp2_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+2];
504 c6_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}];
505 c12_{I}{J} = vdwparam[vdwioffset{I}+vdwjidx{J}+1];
510 /* #if 'table' in INTERACTION_FLAGS[I][J] */
511 /* Calculate table index by multiplying r with table scale and truncate to integer */
512 rt = r{I}{J}*vftabscale;
515 vfitab = {TABLE_NINTERACTIONS}*{TABLE_POINT_SIZE}*vfitab;
516 /* #define INNERFLOPS INNERFLOPS+2 */
519 /* ## ELECTROSTATIC INTERACTIONS */
520 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
522 /* #if KERNEL_ELEC=='Coulomb' */
524 /* COULOMB ELECTROSTATICS */
525 velec = qq{I}{J}*rinv{I}{J};
526 /* #define INNERFLOPS INNERFLOPS+1 */
527 /* #if 'Force' in KERNEL_VF */
528 felec = velec*rinvsq{I}{J};
529 /* #define INNERFLOPS INNERFLOPS+2 */
532 /* #elif KERNEL_ELEC=='ReactionField' */
534 /* REACTION-FIELD ELECTROSTATICS */
535 /* #if 'Potential' in KERNEL_VF */
536 velec = qq{I}{J}*(rinv{I}{J}+krf*rsq{I}{J}-crf);
537 /* #define INNERFLOPS INNERFLOPS+4 */
539 /* #if 'Force' in KERNEL_VF */
540 felec = qq{I}{J}*(rinv{I}{J}*rinvsq{I}{J}-krf2);
541 /* #define INNERFLOPS INNERFLOPS+3 */
544 /* #elif KERNEL_ELEC=='GeneralizedBorn' */
546 /* GENERALIZED BORN AND COULOMB ELECTROSTATICS */
547 isaprod = isai{I}*isaj{J};
548 gbqqfactor = isaprod*(-qq{I}{J})*gbinvepsdiff;
549 gbscale = isaprod*gbtabscale;
550 dvdaj = dvda[jnr+{J}];
551 /* #define INNERFLOPS INNERFLOPS+5 */
553 /* Calculate generalized born table index - this is a separate table from the normal one,
554 * but we use the same procedure by multiplying r with scale and truncating to integer.
556 rt = r{I}{J}*gbscale;
563 Geps = gbeps*gbtab[gbitab+2];
564 Heps2 = gbeps*gbeps*gbtab[gbitab+3];
568 /* #define INNERFLOPS INNERFLOPS+10 */
570 /* #if 'Force' in KERNEL_VF */
571 FF = Fp+Geps+2.0*Heps2;
572 fgb = gbqqfactor*FF*gbscale;
573 dvdatmp = -0.5*(vgb+fgb*r{I}{J});
574 dvdasum = dvdasum + dvdatmp;
575 dvda[jnr] = dvdaj+dvdatmp*isaj{J}*isaj{J};
576 /* #define INNERFLOPS INNERFLOPS+13 */
578 velec = qq{I}{J}*rinv{I}{J};
579 /* #define INNERFLOPS INNERFLOPS+1 */
580 /* #if 'Force' in KERNEL_VF */
581 felec = (velec*rinv{I}{J}-fgb)*rinv{I}{J};
582 /* #define INNERFLOPS INNERFLOPS+3 */
585 /* #elif KERNEL_ELEC=='Ewald' */
586 /* EWALD ELECTROSTATICS */
588 /* Calculate Ewald table index by multiplying r with scale and truncate to integer */
589 ewrt = r{I}{J}*ewtabscale;
592 /* #define INNERFLOPS INNERFLOPS+2 */
593 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_ELEC=='PotentialSwitch' */
595 felec = ewtab[ewitab]+eweps*ewtab[ewitab+1];
596 /* #define INNERFLOPS INNERFLOPS+4 */
597 /* #if KERNEL_MOD_ELEC=='PotentialShift' */
598 velec = qq{I}{J}*((rinv{I}{J}-sh_ewald)-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
599 /* #define INNERFLOPS INNERFLOPS+7 */
601 velec = qq{I}{J}*(rinv{I}{J}-(ewtab[ewitab+2]-ewtabhalfspace*eweps*(ewtab[ewitab]+felec)));
602 /* #define INNERFLOPS INNERFLOPS+6 */
604 /* #if 'Force' in KERNEL_VF */
605 felec = qq{I}{J}*rinv{I}{J}*(rinvsq{I}{J}-felec);
606 /* #define INNERFLOPS INNERFLOPS+3 */
608 /* #elif KERNEL_VF=='Force' */
609 felec = (1.0-eweps)*ewtab[ewitab]+eweps*ewtab[ewitab+1];
610 felec = qq{I}{J}*rinv{I}{J}*(rinvsq{I}{J}-felec);
611 /* #define INNERFLOPS INNERFLOPS+7 */
614 /* #elif KERNEL_ELEC=='CubicSplineTable' */
616 /* CUBIC SPLINE TABLE ELECTROSTATICS */
617 /* #if 'Potential' in KERNEL_VF */
621 Geps = vfeps*vftab[vfitab+2];
622 Heps2 = vfeps*vfeps*vftab[vfitab+3];
624 /* #define INNERFLOPS INNERFLOPS+5 */
625 /* #if 'Potential' in KERNEL_VF */
628 /* #define INNERFLOPS INNERFLOPS+3 */
630 /* #if 'Force' in KERNEL_VF */
631 FF = Fp+Geps+2.0*Heps2;
632 felec = -qq{I}{J}*FF*vftabscale*rinv{I}{J};
633 /* #define INNERFLOPS INNERFLOPS+7 */
636 /* ## End of check for electrostatics interaction forms */
638 /* ## END OF ELECTROSTATIC INTERACTION CHECK FOR PAIR I-J */
640 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
642 /* #if KERNEL_VDW=='LennardJones' */
644 /* LENNARD-JONES DISPERSION/REPULSION */
646 rinvsix = rinvsq{I}{J}*rinvsq{I}{J}*rinvsq{I}{J};
647 /* #define INNERFLOPS INNERFLOPS+2 */
648 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
649 vvdw6 = c6_{I}{J}*rinvsix;
650 vvdw12 = c12_{I}{J}*rinvsix*rinvsix;
651 /* #define INNERFLOPS INNERFLOPS+3 */
652 /* #if KERNEL_MOD_VDW=='PotentialShift' */
653 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);
654 /* #define INNERFLOPS INNERFLOPS+8 */
656 vvdw = vvdw12*(1.0/12.0) - vvdw6*(1.0/6.0);
657 /* #define INNERFLOPS INNERFLOPS+3 */
659 /* ## Check for force inside potential check, i.e. this means we already did the potential part */
660 /* #if 'Force' in KERNEL_VF */
661 fvdw = (vvdw12-vvdw6)*rinvsq{I}{J};
662 /* #define INNERFLOPS INNERFLOPS+2 */
664 /* #elif KERNEL_VF=='Force' */
665 /* ## Force-only LennardJones makes it possible to save 1 flop (they do add up...) */
666 fvdw = (c12_{I}{J}*rinvsix-c6_{I}{J})*rinvsix*rinvsq{I}{J};
667 /* #define INNERFLOPS INNERFLOPS+4 */
670 /* #elif KERNEL_VDW=='Buckingham' */
672 /* BUCKINGHAM DISPERSION/REPULSION */
673 rinvsix = rinvsq{I}{J}*rinvsq{I}{J}*rinvsq{I}{J};
674 vvdw6 = c6_{I}{J}*rinvsix;
675 br = cexp2_{I}{J}*r{I}{J};
676 vvdwexp = cexp1_{I}{J}*exp(-br);
677 /* ## Estimate exp() to 25 flops */
678 /* #define INNERFLOPS INNERFLOPS+31 */
679 /* #if 'Potential' in KERNEL_VF or KERNEL_MOD_VDW=='PotentialSwitch' */
680 /* #if KERNEL_MOD_VDW=='PotentialShift' */
681 vvdw = (vvdwexp-cexp1_{I}{J}*exp(-cexp2_{I}{J}*rvdw)) - (vvdw6 - c6_{I}{J}*sh_vdw_invrcut6)*(1.0/6.0);
682 /* #define INNERFLOPS INNERFLOPS+33 */
684 vvdw = vvdwexp - vvdw6*(1.0/6.0);
685 /* #define INNERFLOPS INNERFLOPS+2 */
688 /* #if 'Force' in KERNEL_VF */
689 fvdw = (br*vvdwexp-vvdw6)*rinvsq{I}{J};
690 /* #define INNERFLOPS INNERFLOPS+3 */
693 /* #elif KERNEL_VDW=='CubicSplineTable' */
695 /* CUBIC SPLINE TABLE DISPERSION */
696 vfitab += {TABLE_VDW_OFFSET};
697 /* #if 'Potential' in KERNEL_VF */
701 Geps = vfeps*vftab[vfitab+2];
702 Heps2 = vfeps*vfeps*vftab[vfitab+3];
704 /* #define INNERFLOPS INNERFLOPS+5 */
705 /* #if 'Potential' in KERNEL_VF */
707 vvdw6 = c6_{I}{J}*VV;
708 /* #define INNERFLOPS INNERFLOPS+3 */
710 /* #if 'Force' in KERNEL_VF */
711 FF = Fp+Geps+2.0*Heps2;
712 fvdw6 = c6_{I}{J}*FF;
713 /* #define INNERFLOPS INNERFLOPS+4 */
716 /* CUBIC SPLINE TABLE REPULSION */
717 /* #if 'Potential' in KERNEL_VF */
721 Geps = vfeps*vftab[vfitab+6];
722 Heps2 = vfeps*vfeps*vftab[vfitab+7];
724 /* #define INNERFLOPS INNERFLOPS+5 */
725 /* #if 'Potential' in KERNEL_VF */
727 vvdw12 = c12_{I}{J}*VV;
728 /* #define INNERFLOPS INNERFLOPS+3 */
730 /* #if 'Force' in KERNEL_VF */
731 FF = Fp+Geps+2.0*Heps2;
732 fvdw12 = c12_{I}{J}*FF;
733 /* #define INNERFLOPS INNERFLOPS+4 */
735 /* #if 'Potential' in KERNEL_VF */
737 /* #define INNERFLOPS INNERFLOPS+1 */
739 /* #if 'Force' in KERNEL_VF */
740 fvdw = -(fvdw6+fvdw12)*vftabscale*rinv{I}{J};
741 /* #define INNERFLOPS INNERFLOPS+4 */
744 /* ## End of check for vdw interaction forms */
746 /* ## END OF VDW INTERACTION CHECK FOR PAIR I-J */
748 /* #if 'switch' in INTERACTION_FLAGS[I][J] */
750 d = (d>0.0) ? d : 0.0;
752 sw = 1.0+d2*d*(swV3+d*(swV4+d*swV5));
753 /* #define INNERFLOPS INNERFLOPS+9 */
755 /* #if 'Force' in KERNEL_VF */
756 dsw = d2*(swF2+d*(swF3+d*swF4));
757 /* #define INNERFLOPS INNERFLOPS+5 */
760 /* Evaluate switch function */
761 /* #if 'Force' in KERNEL_VF */
762 /* fscal'=f'/r=-(v*sw)'/r=-(v'*sw+v*dsw)/r=-v'*sw/r-v*dsw/r=fscal*sw-v*dsw/r */
763 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
764 felec = felec*sw - rinv{I}{J}*velec*dsw;
765 /* #define INNERFLOPS INNERFLOPS+3 */
767 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
768 fvdw = fvdw*sw - rinv{I}{J}*vvdw*dsw;
769 /* #define INNERFLOPS INNERFLOPS+3 */
772 /* #if 'Potential' in KERNEL_VF */
773 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_ELEC=='PotentialSwitch' */
775 /* #define INNERFLOPS INNERFLOPS+1 */
777 /* #if 'vdw' in INTERACTION_FLAGS[I][J] and KERNEL_MOD_VDW=='PotentialSwitch' */
779 /* #define INNERFLOPS INNERFLOPS+1 */
784 /* #if 'Potential' in KERNEL_VF */
785 /* Update potential sums from outer loop */
786 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] */
788 /* #define INNERFLOPS INNERFLOPS+1 */
789 /* #if KERNEL_ELEC=='GeneralizedBorn' */
791 /* #define INNERFLOPS INNERFLOPS+1 */
794 /* #if 'vdw' in INTERACTION_FLAGS[I][J] */
796 /* #define INNERFLOPS INNERFLOPS+1 */
800 /* #if 'Force' in KERNEL_VF */
802 /* #if 'electrostatics' in INTERACTION_FLAGS[I][J] and 'vdw' in INTERACTION_FLAGS[I][J] */
804 /* #define INNERFLOPS INNERFLOPS+1 */
805 /* #elif 'electrostatics' in INTERACTION_FLAGS[I][J] */
807 /* #elif 'vdw' in INTERACTION_FLAGS[I][J] */
811 /* #if KERNEL_RESOLUTION=='CG' */
812 if(force_cap>0 && (fabs(fscal)> force_cap))
814 fscal=force_cap*fscal/fabs(fscal);
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});