2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
37 #include "calc_verletbuf.h"
43 #include <sys/types.h>
45 #include "gromacs/legacyheaders/macros.h"
46 #include "gromacs/legacyheaders/typedefs.h"
47 #include "gromacs/math/calculate-ewald-splitting-coefficient.h"
48 #include "gromacs/math/units.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/mdlib/nb_verlet.h"
51 #include "gromacs/mdlib/nbnxn_search.h"
52 #include "gromacs/mdlib/nbnxn_simd.h"
53 #include "gromacs/utility/fatalerror.h"
54 #include "gromacs/utility/smalloc.h"
56 /* The code in this file estimates a pairlist buffer length
57 * given a target energy drift per atom per picosecond.
58 * This is done by estimating the drift given a buffer length.
59 * Ideally we would like to have a tight overestimate of the drift,
60 * but that can be difficult to achieve.
62 * Significant approximations used:
64 * Uniform particle density. UNDERESTIMATES the drift by rho_global/rho_local.
66 * Interactions don't affect particle motion. OVERESTIMATES the drift on longer
67 * time scales. This approximation probably introduces the largest errors.
69 * Only take one constraint per particle into account: OVERESTIMATES the drift.
71 * For rotating constraints assume the same functional shape for time scales
72 * where the constraints rotate significantly as the exact expression for
73 * short time scales. OVERESTIMATES the drift on long time scales.
75 * For non-linear virtual sites use the mass of the lightest constructing atom
76 * to determine the displacement. OVER/UNDERESTIMATES the drift, depending on
77 * the geometry and masses of constructing atoms.
79 * Note that the formulas for normal atoms and linear virtual sites are exact,
80 * apart from the first two approximations.
82 * Note that apart from the effect of the above approximations, the actual
83 * drift of the total energy of a system can be order of magnitude smaller
84 * due to cancellation of positive and negative drift for different pairs.
88 /* Struct for unique atom type for calculating the energy drift.
89 * The atom displacement depends on mass and constraints.
90 * The energy jump for given distance depend on LJ type and q.
95 int type; /* type (used for LJ parameters) */
97 gmx_bool bConstr; /* constrained, if TRUE, use #DOF=2 iso 3 */
98 real con_mass; /* mass of heaviest atom connected by constraints */
99 real con_len; /* constraint length to the heaviest atom */
100 } atom_nonbonded_kinetic_prop_t;
102 /* Struct for unique atom type for calculating the energy drift.
103 * The atom displacement depends on mass and constraints.
104 * The energy jump for given distance depend on LJ type and q.
108 atom_nonbonded_kinetic_prop_t prop; /* non-bonded and kinetic atom prop. */
109 int n; /* #atoms of this type in the system */
110 } verletbuf_atomtype_t;
112 void verletbuf_get_list_setup(gmx_bool gmx_unused bSIMD,
114 verletbuf_list_setup_t *list_setup)
116 /* When calling this function we often don't know which kernel type we
117 * are going to use. W choose the kernel type with the smallest possible
118 * i- and j-cluster sizes, so we potentially overestimate, but never
119 * underestimate, the buffer drift.
120 * Note that the current buffer estimation code only handles clusters
121 * of size 1, 2 or 4, so for 4x8 or 8x8 we use the estimate for 4x4.
126 /* The CUDA kernels split the j-clusters in two halves */
127 list_setup->cluster_size_i = nbnxn_kernel_to_ci_size(nbnxnk8x8x8_GPU);
128 list_setup->cluster_size_j = nbnxn_kernel_to_cj_size(nbnxnk8x8x8_GPU)/2;
134 kernel_type = nbnxnk4x4_PlainC;
136 #ifdef GMX_NBNXN_SIMD
139 #ifdef GMX_NBNXN_SIMD_2XNN
140 /* We use the smallest cluster size to be on the safe side */
141 kernel_type = nbnxnk4xN_SIMD_2xNN;
143 kernel_type = nbnxnk4xN_SIMD_4xN;
148 list_setup->cluster_size_i = nbnxn_kernel_to_ci_size(kernel_type);
149 list_setup->cluster_size_j = nbnxn_kernel_to_cj_size(kernel_type);
154 atom_nonbonded_kinetic_prop_equal(const atom_nonbonded_kinetic_prop_t *prop1,
155 const atom_nonbonded_kinetic_prop_t *prop2)
157 return (prop1->mass == prop2->mass &&
158 prop1->type == prop2->type &&
159 prop1->q == prop2->q &&
160 prop1->bConstr == prop2->bConstr &&
161 prop1->con_mass == prop2->con_mass &&
162 prop1->con_len == prop2->con_len);
165 static void add_at(verletbuf_atomtype_t **att_p, int *natt_p,
166 const atom_nonbonded_kinetic_prop_t *prop,
169 verletbuf_atomtype_t *att;
174 /* Ignore massless particles */
182 while (i < natt && !atom_nonbonded_kinetic_prop_equal(prop, &att[i].prop))
194 srenew(*att_p, *natt_p);
195 (*att_p)[i].prop = *prop;
196 (*att_p)[i].n = nmol;
200 static void get_vsite_masses(const gmx_moltype_t *moltype,
201 const gmx_ffparams_t *ffparams,
210 /* Check for virtual sites, determine mass from constructing atoms */
211 for (ft = 0; ft < F_NRE; ft++)
215 il = &moltype->ilist[ft];
217 for (i = 0; i < il->nr; i += 1+NRAL(ft))
220 real inv_mass, coeff, m_aj;
223 ip = &ffparams->iparams[il->iatoms[i]];
225 a1 = il->iatoms[i+1];
229 /* Only vsiten can have more than four
230 constructing atoms, so NRAL(ft) <= 5 */
233 const int maxj = NRAL(ft);
237 for (j = 1; j < maxj; j++)
239 cam[j] = moltype->atoms.atom[il->iatoms[i+1+j]].m;
242 cam[j] = vsite_m[il->iatoms[i+1+j]];
246 gmx_fatal(FARGS, "In molecule type '%s' %s construction involves atom %d, which is a virtual site of equal or high complexity. This is not supported.",
248 interaction_function[ft].longname,
249 il->iatoms[i+1+j]+1);
257 vsite_m[a1] = (cam[1]*cam[2])/(cam[2]*sqr(1-ip->vsite.a) + cam[1]*sqr(ip->vsite.a));
261 vsite_m[a1] = (cam[1]*cam[2]*cam[3])/(cam[2]*cam[3]*sqr(1-ip->vsite.a-ip->vsite.b) + cam[1]*cam[3]*sqr(ip->vsite.a) + cam[1]*cam[2]*sqr(ip->vsite.b));
264 gmx_incons("Invalid vsite type");
267 /* Use the mass of the lightest constructing atom.
268 * This is an approximation.
269 * If the distance of the virtual site to the
270 * constructing atom is less than all distances
271 * between constructing atoms, this is a safe
272 * over-estimate of the displacement of the vsite.
273 * This condition holds for all H mass replacement
274 * vsite constructions, except for SP2/3 groups.
275 * In SP3 groups one H will have a F_VSITE3
276 * construction, so even there the total drift
277 * estimate shouldn't be far off.
279 vsite_m[a1] = cam[1];
280 for (j = 2; j < maxj; j++)
282 vsite_m[a1] = min(vsite_m[a1], cam[j]);
295 for (j = 0; j < 3*ffparams->iparams[il->iatoms[i]].vsiten.n; j += 3)
297 aj = il->iatoms[i+j+2];
298 coeff = ffparams->iparams[il->iatoms[i+j]].vsiten.a;
299 if (moltype->atoms.atom[aj].ptype == eptVSite)
305 m_aj = moltype->atoms.atom[aj].m;
309 gmx_incons("The mass of a vsiten constructing atom is <= 0");
311 inv_mass += coeff*coeff/m_aj;
313 vsite_m[a1] = 1/inv_mass;
314 /* Correct for loop increment of i */
315 i += j - 1 - NRAL(ft);
319 fprintf(debug, "atom %4d %-20s mass %6.3f\n",
320 a1, interaction_function[ft].longname, vsite_m[a1]);
327 static void get_verlet_buffer_atomtypes(const gmx_mtop_t *mtop,
328 verletbuf_atomtype_t **att_p,
332 verletbuf_atomtype_t *att;
334 int mb, nmol, ft, i, a1, a2, a3, a;
335 const t_atoms *atoms;
338 atom_nonbonded_kinetic_prop_t *prop;
340 int n_nonlin_vsite_mol;
345 if (n_nonlin_vsite != NULL)
350 for (mb = 0; mb < mtop->nmolblock; mb++)
352 nmol = mtop->molblock[mb].nmol;
354 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
356 /* Check for constraints, as they affect the kinetic energy.
357 * For virtual sites we need the masses and geometry of
358 * the constructing atoms to determine their velocity distribution.
360 snew(prop, atoms->nr);
361 snew(vsite_m, atoms->nr);
363 for (ft = F_CONSTR; ft <= F_CONSTRNC; ft++)
365 il = &mtop->moltype[mtop->molblock[mb].type].ilist[ft];
367 for (i = 0; i < il->nr; i += 1+NRAL(ft))
369 ip = &mtop->ffparams.iparams[il->iatoms[i]];
370 a1 = il->iatoms[i+1];
371 a2 = il->iatoms[i+2];
372 if (atoms->atom[a2].m > prop[a1].con_mass)
374 prop[a1].con_mass = atoms->atom[a2].m;
375 prop[a1].con_len = ip->constr.dA;
377 if (atoms->atom[a1].m > prop[a2].con_mass)
379 prop[a2].con_mass = atoms->atom[a1].m;
380 prop[a2].con_len = ip->constr.dA;
385 il = &mtop->moltype[mtop->molblock[mb].type].ilist[F_SETTLE];
387 for (i = 0; i < il->nr; i += 1+NRAL(F_SETTLE))
389 ip = &mtop->ffparams.iparams[il->iatoms[i]];
390 a1 = il->iatoms[i+1];
391 a2 = il->iatoms[i+2];
392 a3 = il->iatoms[i+3];
393 /* Usually the mass of a1 (usually oxygen) is larger than a2/a3.
394 * If this is not the case, we overestimate the displacement,
395 * which leads to a larger buffer (ok since this is an exotic case).
397 prop[a1].con_mass = atoms->atom[a2].m;
398 prop[a1].con_len = ip->settle.doh;
400 prop[a2].con_mass = atoms->atom[a1].m;
401 prop[a2].con_len = ip->settle.doh;
403 prop[a3].con_mass = atoms->atom[a1].m;
404 prop[a3].con_len = ip->settle.doh;
407 get_vsite_masses(&mtop->moltype[mtop->molblock[mb].type],
410 &n_nonlin_vsite_mol);
411 if (n_nonlin_vsite != NULL)
413 *n_nonlin_vsite += nmol*n_nonlin_vsite_mol;
416 for (a = 0; a < atoms->nr; a++)
418 if (atoms->atom[a].ptype == eptVSite)
420 prop[a].mass = vsite_m[a];
424 prop[a].mass = atoms->atom[a].m;
426 prop[a].type = atoms->atom[a].type;
427 prop[a].q = atoms->atom[a].q;
428 /* We consider an atom constrained, #DOF=2, when it is
429 * connected with constraints to (at least one) atom with
430 * a mass of more than 0.4x its own mass. This is not a critical
431 * parameter, since with roughly equal masses the unconstrained
432 * and constrained displacement will not differ much (and both
433 * overestimate the displacement).
435 prop[a].bConstr = (prop[a].con_mass > 0.4*prop[a].mass);
437 add_at(&att, &natt, &prop[a], nmol);
446 for (a = 0; a < natt; a++)
448 fprintf(debug, "type %d: m %5.2f t %d q %6.3f con %d con_m %5.3f con_l %5.3f n %d\n",
449 a, att[a].prop.mass, att[a].prop.type, att[a].prop.q,
450 att[a].prop.bConstr, att[a].prop.con_mass, att[a].prop.con_len,
459 /* This function computes two components of the estimate of the variance
460 * in the displacement of one atom in a system of two constrained atoms.
461 * Returns in sigma2_2d the variance due to rotation of the constrained
462 * atom around the atom to which it constrained.
463 * Returns in sigma2_3d the variance due to displacement of the COM
464 * of the whole system of the two constrained atoms.
466 * Note that we only take a single constraint (the one to the heaviest atom)
467 * into account. If an atom has multiple constraints, this will result in
468 * an overestimate of the displacement, which gives a larger drift and buffer.
470 static void constrained_atom_sigma2(real kT_fac,
471 const atom_nonbonded_kinetic_prop_t *prop,
480 /* Here we decompose the motion of a constrained atom into two
481 * components: rotation around the COM and translation of the COM.
484 /* Determine the variance for the displacement of the rotational mode */
485 sigma2_rot = kT_fac/(prop->mass*(prop->mass + prop->con_mass)/prop->con_mass);
487 /* The distance from the atom to the COM, i.e. the rotational arm */
488 com_dist = prop->con_len*prop->con_mass/(prop->mass + prop->con_mass);
490 /* The variance relative to the arm */
491 sigma2_rel = sigma2_rot/(com_dist*com_dist);
492 /* At 6 the scaling formula has slope 0,
493 * so we keep sigma2_2d constant after that.
497 /* A constrained atom rotates around the atom it is constrained to.
498 * This results in a smaller linear displacement than for a free atom.
499 * For a perfectly circular displacement, this lowers the displacement
500 * by: 1/arcsin(arc_length)
501 * and arcsin(x) = 1 + x^2/6 + ...
502 * For sigma2_rel<<1 the displacement distribution is erfc
503 * (exact formula is provided below). For larger sigma, it is clear
504 * that the displacement can't be larger than 2*com_dist.
505 * It turns out that the distribution becomes nearly uniform.
506 * For intermediate sigma2_rel, scaling down sigma with the third
507 * order expansion of arcsin with argument sigma_rel turns out
508 * to give a very good approximation of the distribution and variance.
509 * Even for larger values, the variance is only slightly overestimated.
510 * Note that the most relevant displacements are in the long tail.
511 * This rotation approximation always overestimates the tail (which
512 * runs to infinity, whereas it should be <= 2*com_dist).
513 * Thus we always overestimate the drift and the buffer size.
515 scale = 1/(1 + sigma2_rel/6);
516 *sigma2_2d = sigma2_rot*scale*scale;
520 /* sigma_2d is set to the maximum given by the scaling above.
521 * For large sigma2 the real displacement distribution is close
522 * to uniform over -2*con_len to 2*com_dist.
523 * Our erfc with sigma_2d=sqrt(1.5)*com_dist (which means the sigma
524 * of the erfc output distribution is con_dist) overestimates
525 * the variance and additionally has a long tail. This means
526 * we have a (safe) overestimation of the drift.
528 *sigma2_2d = 1.5*com_dist*com_dist;
531 /* The constrained atom also moves (in 3D) with the COM of both atoms */
532 *sigma2_3d = kT_fac/(prop->mass + prop->con_mass);
535 static void get_atom_sigma2(real kT_fac,
536 const atom_nonbonded_kinetic_prop_t *prop,
542 /* Complicated constraint calculation in a separate function */
543 constrained_atom_sigma2(kT_fac, prop, sigma2_2d, sigma2_3d);
547 /* Unconstrained atom: trivial */
549 *sigma2_3d = kT_fac/prop->mass;
553 static void approx_2dof(real s2, real x, real *shift, real *scale)
555 /* A particle with 1 DOF constrained has 2 DOFs instead of 3.
556 * This code is also used for particles with multiple constraints,
557 * in which case we overestimate the displacement.
558 * The 2DOF distribution is sqrt(pi/2)*erfc(r/(sqrt(2)*s))/(2*s).
559 * We approximate this with scale*Gaussian(s,r+shift),
560 * by matching the distribution value and derivative at x.
561 * This is a tight overestimate for all r>=0 at any s and x.
565 ex = exp(-x*x/(2*s2));
566 er = gmx_erfc(x/sqrt(2*s2));
568 *shift = -x + sqrt(2*s2/M_PI)*ex/er;
569 *scale = 0.5*M_PI*exp(ex*ex/(M_PI*er*er))*er;
572 static real ener_drift(const verletbuf_atomtype_t *att, int natt,
573 const gmx_ffparams_t *ffp,
575 real md1_ljd, real d2_ljd, real md3_ljd,
576 real md1_ljr, real d2_ljr, real md3_ljr,
577 real md1_el, real d2_el,
579 real rlist, real boxvol)
581 /* Erfc(8)=1e-29, use this limit so we have some space for arithmetic
582 * on the result when using float precision.
584 const real erfc_arg_max = 8.0;
586 double drift_tot, pot1, pot2, pot3, pot;
588 real s2i_2d, s2i_3d, s2j_2d, s2j_3d, s2, s;
591 real sc_fac, rsh, rsh2;
592 double c_exp, c_erfc;
596 /* Loop over the different atom type pairs */
597 for (i = 0; i < natt; i++)
599 get_atom_sigma2(kT_fac, &att[i].prop, &s2i_2d, &s2i_3d);
600 ti = att[i].prop.type;
602 for (j = i; j < natt; j++)
604 get_atom_sigma2(kT_fac, &att[j].prop, &s2j_2d, &s2j_3d);
605 tj = att[j].prop.type;
607 /* Add up the up to four independent variances */
608 s2 = s2i_2d + s2i_3d + s2j_2d + s2j_3d;
610 /* Note that attractive and repulsive potentials for individual
611 * pairs will partially cancel.
613 /* -dV/dr at the cut-off for LJ + Coulomb */
615 md1_ljd*ffp->iparams[ti*ffp->atnr+tj].lj.c6 +
616 md1_ljr*ffp->iparams[ti*ffp->atnr+tj].lj.c12 +
617 md1_el*att[i].prop.q*att[j].prop.q;
619 /* d2V/dr2 at the cut-off for LJ + Coulomb */
621 d2_ljd*ffp->iparams[ti*ffp->atnr+tj].lj.c6 +
622 d2_ljr*ffp->iparams[ti*ffp->atnr+tj].lj.c12 +
623 d2_el*att[i].prop.q*att[j].prop.q;
625 /* -d3V/dr3 at the cut-off for LJ, we neglect Coulomb */
627 md3_ljd*ffp->iparams[ti*ffp->atnr+tj].lj.c6 +
628 md3_ljr*ffp->iparams[ti*ffp->atnr+tj].lj.c12;
633 if (rsh*rsh > 2*s2*erfc_arg_max*erfc_arg_max)
635 /* Erfc might run out of float and become 0, somewhat before
636 * c_exp becomes 0. To avoid this and to avoid NaN in
637 * approx_2dof, we set both c_expc and c_erfc to zero.
638 * In any relevant case this has no effect on the results,
639 * since c_exp < 6e-29, so the displacement is completely
640 * negligible for such atom pairs (and an overestimate).
641 * In nearly all use cases, there will be other atom
642 * pairs that contribute much more to the total, so zeroing
643 * this particular contribution has no effect at all.
650 /* For constraints: adapt r and scaling for the Gaussian */
651 if (att[i].prop.bConstr)
655 approx_2dof(s2i_2d, r_buffer*s2i_2d/s2, &sh, &sc);
659 if (att[j].prop.bConstr)
663 approx_2dof(s2j_2d, r_buffer*s2j_2d/s2, &sh, &sc);
668 /* Exact contribution of an atom pair with Gaussian displacement
669 * with sigma s to the energy drift for a potential with
670 * derivative -md and second derivative dd at the cut-off.
671 * The only catch is that for potentials that change sign
672 * near the cut-off there could be an unlucky compensation
673 * of positive and negative energy drift.
674 * Such potentials are extremely rare though.
676 * Note that pot has unit energy*length, as the linear
677 * atom density still needs to be put in.
679 c_exp = exp(-rsh*rsh/(2*s2))/sqrt(2*M_PI);
680 c_erfc = 0.5*gmx_erfc(rsh/(sqrt(2*s2)));
686 md1/2*((rsh2 + s2)*c_erfc - rsh*s*c_exp);
688 d2/6*(s*(rsh2 + 2*s2)*c_exp - rsh*(rsh2 + 3*s2)*c_erfc);
690 md3/24*((rsh2*rsh2 + 6*rsh2*s2 + 3*s2*s2)*c_erfc - rsh*s*(rsh2 + 5*s2)*c_exp);
691 pot = pot1 + pot2 + pot3;
695 fprintf(debug, "n %d %d d s %.3f %.3f %.3f %.3f con %d -d1 %8.1e d2 %8.1e -d3 %8.1e pot1 %8.1e pot2 %8.1e pot3 %8.1e pot %8.1e\n",
697 sqrt(s2i_2d), sqrt(s2i_3d),
698 sqrt(s2j_2d), sqrt(s2j_3d),
699 att[i].prop.bConstr+att[j].prop.bConstr,
701 pot1, pot2, pot3, pot);
704 /* Multiply by the number of atom pairs */
707 pot *= (double)att[i].n*(att[i].n - 1)/2;
711 pot *= (double)att[i].n*att[j].n;
713 /* We need the line density to get the energy drift of the system.
714 * The effective average r^2 is close to (rlist+sigma)^2.
716 pot *= 4*M_PI*sqr(rlist + s)/boxvol;
718 /* Add the unsigned drift to avoid cancellation of errors */
719 drift_tot += fabs(pot);
726 static real surface_frac(int cluster_size, real particle_distance, real rlist)
730 if (rlist < 0.5*particle_distance)
732 /* We have non overlapping spheres */
736 /* Half the inter-particle distance relative to rlist */
737 d = 0.5*particle_distance/rlist;
739 /* Determine the area of the surface at distance rlist to the closest
740 * particle, relative to surface of a sphere of radius rlist.
741 * The formulas below assume close to cubic cells for the pair search grid,
742 * which the pair search code tries to achieve.
743 * Note that in practice particle distances will not be delta distributed,
744 * but have some spread, often involving shorter distances,
745 * as e.g. O-H bonds in a water molecule. Thus the estimates below will
746 * usually be slightly too high and thus conservative.
748 switch (cluster_size)
751 /* One particle: trivial */
755 /* Two particles: two spheres at fractional distance 2*a */
759 /* We assume a perfect, symmetric tetrahedron geometry.
760 * The surface around a tetrahedron is too complex for a full
761 * analytical solution, so we use a Taylor expansion.
763 area_rel = (1.0 + 1/M_PI*(6*acos(1/sqrt(3))*d +
767 83.0/756.0*d*d*d*d*d*d)));
770 gmx_incons("surface_frac called with unsupported cluster_size");
774 return area_rel/cluster_size;
777 /* Returns the negative of the third derivative of a potential r^-p
778 * with a force-switch function, evaluated at the cut-off rc.
780 static real md3_force_switch(real p, real rswitch, real rc)
782 /* The switched force function is:
783 * p*r^-(p+1) + a*(r - rswitch)^2 + b*(r - rswitch)^3
786 real md3_pot, md3_sw;
788 a = -((p + 4)*rc - (p + 1)*rswitch)/(pow(rc, p+2)*pow(rc-rswitch, 2));
789 b = ((p + 3)*rc - (p + 1)*rswitch)/(pow(rc, p+2)*pow(rc-rswitch, 3));
791 md3_pot = (p + 2)*(p + 1)*p*pow(rc, p+3);
792 md3_sw = 2*a + 6*b*(rc - rswitch);
794 return md3_pot + md3_sw;
797 void calc_verlet_buffer_size(const gmx_mtop_t *mtop, real boxvol,
798 const t_inputrec *ir,
799 real reference_temperature,
800 const verletbuf_list_setup_t *list_setup,
807 real particle_distance;
808 real nb_clust_frac_pairs_not_in_list_at_cutoff;
810 verletbuf_atomtype_t *att = NULL;
813 real md1_ljd, d2_ljd, md3_ljd;
814 real md1_ljr, d2_ljr, md3_ljr;
817 real kT_fac, mass_min;
822 if (reference_temperature < 0)
824 if (EI_MD(ir->eI) && ir->etc == etcNO)
826 /* This case should be handled outside calc_verlet_buffer_size */
827 gmx_incons("calc_verlet_buffer_size called with an NVE ensemble and reference_temperature < 0");
830 /* We use the maximum temperature with multiple T-coupl groups.
831 * We could use a per particle temperature, but since particles
832 * interact, this might underestimate the buffer size.
834 reference_temperature = 0;
835 for (i = 0; i < ir->opts.ngtc; i++)
837 if (ir->opts.tau_t[i] >= 0)
839 reference_temperature = max(reference_temperature,
845 /* Resolution of the buffer size */
848 env = getenv("GMX_VERLET_BUFFER_RES");
851 sscanf(env, "%lf", &resolution);
854 /* In an atom wise pair-list there would be no pairs in the list
855 * beyond the pair-list cut-off.
856 * However, we use a pair-list of groups vs groups of atoms.
857 * For groups of 4 atoms, the parallelism of SSE instructions, only
858 * 10% of the atoms pairs are not in the list just beyond the cut-off.
859 * As this percentage increases slowly compared to the decrease of the
860 * Gaussian displacement distribution over this range, we can simply
861 * reduce the drift by this fraction.
862 * For larger groups, e.g. of 8 atoms, this fraction will be lower,
863 * so then buffer size will be on the conservative (large) side.
865 * Note that the formulas used here do not take into account
866 * cancellation of errors which could occur by missing both
867 * attractive and repulsive interactions.
869 * The only major assumption is homogeneous particle distribution.
870 * For an inhomogeneous system, such as a liquid-vapor system,
871 * the buffer will be underestimated. The actual energy drift
872 * will be higher by the factor: local/homogeneous particle density.
874 * The results of this estimate have been checked againt simulations.
875 * In most cases the real drift differs by less than a factor 2.
878 /* Worst case assumption: HCP packing of particles gives largest distance */
879 particle_distance = pow(boxvol*sqrt(2)/mtop->natoms, 1.0/3.0);
881 get_verlet_buffer_atomtypes(mtop, &att, &natt, n_nonlin_vsite);
882 assert(att != NULL && natt >= 0);
886 fprintf(debug, "particle distance assuming HCP packing: %f nm\n",
888 fprintf(debug, "energy drift atom types: %d\n", natt);
891 reppow = mtop->ffparams.reppow;
898 if (ir->vdwtype == evdwCUT)
900 real sw_range, md3_pswf;
902 switch (ir->vdw_modifier)
905 case eintmodPOTSHIFT:
906 /* -dV/dr of -r^-6 and r^-reppow */
907 md1_ljd = -6*pow(ir->rvdw, -7.0);
908 md1_ljr = reppow*pow(ir->rvdw, -(reppow+1));
909 /* The contribution of the higher derivatives is negligible */
911 case eintmodFORCESWITCH:
912 /* At the cut-off: V=V'=V''=0, so we use only V''' */
913 md3_ljd = -md3_force_switch(6.0, ir->rvdw_switch, ir->rvdw);
914 md3_ljr = md3_force_switch(reppow, ir->rvdw_switch, ir->rvdw);
916 case eintmodPOTSWITCH:
917 /* At the cut-off: V=V'=V''=0.
918 * V''' is given by the original potential times
919 * the third derivative of the switch function.
921 sw_range = ir->rvdw - ir->rvdw_switch;
922 md3_pswf = 60.0*pow(sw_range, -3.0);
924 md3_ljd = -pow(ir->rvdw, -6.0 )*md3_pswf;
925 md3_ljr = pow(ir->rvdw, -reppow)*md3_pswf;
928 gmx_incons("Unimplemented VdW modifier");
931 else if (EVDW_PME(ir->vdwtype))
933 real b, r, br, br2, br4, br6;
934 b = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
940 /* -dV/dr of g(br)*r^-6 [where g(x) = exp(-x^2)(1+x^2+x^4/2), see LJ-PME equations in manual] and r^-reppow */
941 md1_ljd = -exp(-br2)*(br6 + 3.0*br4 + 6.0*br2 + 6.0)*pow(r, -7.0);
942 md1_ljr = reppow*pow(r, -(reppow+1));
943 /* The contribution of the higher derivatives is negligible */
947 gmx_fatal(FARGS, "Energy drift calculation is only implemented for plain cut-off Lennard-Jones interactions");
950 elfac = ONE_4PI_EPS0/ir->epsilon_r;
952 /* Determine md=-dV/dr and dd=d^2V/dr^2 */
955 if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
959 if (ir->coulombtype == eelCUT)
966 eps_rf = ir->epsilon_rf/ir->epsilon_r;
969 k_rf = pow(ir->rcoulomb, -3.0)*(eps_rf - ir->epsilon_r)/(2*eps_rf + ir->epsilon_r);
973 /* epsilon_rf = infinity */
974 k_rf = 0.5*pow(ir->rcoulomb, -3.0);
980 md1_el = elfac*(pow(ir->rcoulomb, -2.0) - 2*k_rf*ir->rcoulomb);
982 d2_el = elfac*(2*pow(ir->rcoulomb, -3.0) + 2*k_rf);
984 else if (EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD)
988 b = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
991 md1_el = elfac*(b*exp(-br*br)*M_2_SQRTPI/rc + gmx_erfc(br)/(rc*rc));
992 d2_el = elfac/(rc*rc)*(2*b*(1 + br*br)*exp(-br*br)*M_2_SQRTPI + 2*gmx_erfc(br)/rc);
996 gmx_fatal(FARGS, "Energy drift calculation is only implemented for Reaction-Field and Ewald electrostatics");
999 /* Determine the variance of the atomic displacement
1000 * over nstlist-1 steps: kT_fac
1001 * For inertial dynamics (not Brownian dynamics) the mass factor
1002 * is not included in kT_fac, it is added later.
1006 /* Get the displacement distribution from the random component only.
1007 * With accurate integration the systematic (force) displacement
1008 * should be negligible (unless nstlist is extremely large, which
1009 * you wouldn't do anyhow).
1011 kT_fac = 2*BOLTZ*reference_temperature*(ir->nstlist-1)*ir->delta_t;
1012 if (ir->bd_fric > 0)
1014 /* This is directly sigma^2 of the displacement */
1015 kT_fac /= ir->bd_fric;
1017 /* Set the masses to 1 as kT_fac is the full sigma^2,
1018 * but we divide by m in ener_drift().
1020 for (i = 0; i < natt; i++)
1022 att[i].prop.mass = 1;
1029 /* Per group tau_t is not implemented yet, use the maximum */
1030 tau_t = ir->opts.tau_t[0];
1031 for (i = 1; i < ir->opts.ngtc; i++)
1033 tau_t = max(tau_t, ir->opts.tau_t[i]);
1037 /* This kT_fac needs to be divided by the mass to get sigma^2 */
1042 kT_fac = BOLTZ*reference_temperature*sqr((ir->nstlist-1)*ir->delta_t);
1045 mass_min = att[0].prop.mass;
1046 for (i = 1; i < natt; i++)
1048 mass_min = min(mass_min, att[i].prop.mass);
1053 fprintf(debug, "md1_ljd %9.2e d2_ljd %9.2e md3_ljd %9.2e\n", md1_ljd, d2_ljd, md3_ljd);
1054 fprintf(debug, "md1_ljr %9.2e d2_ljr %9.2e md3_ljr %9.2e\n", md1_ljr, d2_ljr, md3_ljr);
1055 fprintf(debug, "md1_el %9.2e d2_el %9.2e\n", md1_el, d2_el);
1056 fprintf(debug, "sqrt(kT_fac) %f\n", sqrt(kT_fac));
1057 fprintf(debug, "mass_min %f\n", mass_min);
1060 /* Search using bisection */
1062 /* The drift will be neglible at 5 times the max sigma */
1063 ib1 = (int)(5*2*sqrt(kT_fac/mass_min)/resolution) + 1;
1064 while (ib1 - ib0 > 1)
1068 rl = max(ir->rvdw, ir->rcoulomb) + rb;
1070 /* Calculate the average energy drift at the last step
1071 * of the nstlist steps at which the pair-list is used.
1073 drift = ener_drift(att, natt, &mtop->ffparams,
1075 md1_ljd, d2_ljd, md3_ljd,
1076 md1_ljr, d2_ljr, md3_ljr,
1081 /* Correct for the fact that we are using a Ni x Nj particle pair list
1082 * and not a 1 x 1 particle pair list. This reduces the drift.
1084 /* We don't have a formula for 8 (yet), use 4 which is conservative */
1085 nb_clust_frac_pairs_not_in_list_at_cutoff =
1086 surface_frac(min(list_setup->cluster_size_i, 4),
1087 particle_distance, rl)*
1088 surface_frac(min(list_setup->cluster_size_j, 4),
1089 particle_distance, rl);
1090 drift *= nb_clust_frac_pairs_not_in_list_at_cutoff;
1092 /* Convert the drift to drift per unit time per atom */
1093 drift /= ir->nstlist*ir->delta_t*mtop->natoms;
1097 fprintf(debug, "ib %3d %3d %3d rb %.3f %dx%d fac %.3f drift %.1e\n",
1099 list_setup->cluster_size_i, list_setup->cluster_size_j,
1100 nb_clust_frac_pairs_not_in_list_at_cutoff,
1104 if (fabs(drift) > ir->verletbuf_tol)
1116 *rlist = max(ir->rvdw, ir->rcoulomb) + ib1*resolution;