2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2012,2013,2014,2015,2016,2017,2018, 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"
44 #include "gromacs/ewald/ewald-utils.h"
45 #include "gromacs/math/functions.h"
46 #include "gromacs/math/units.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/mdlib/nb_verlet.h"
49 #include "gromacs/mdlib/nbnxn_simd.h"
50 #include "gromacs/mdlib/nbnxn_util.h"
51 #include "gromacs/mdtypes/inputrec.h"
52 #include "gromacs/mdtypes/md_enums.h"
53 #include "gromacs/topology/ifunc.h"
54 #include "gromacs/topology/topology.h"
55 #include "gromacs/utility/fatalerror.h"
56 #include "gromacs/utility/strconvert.h"
58 /* The code in this file estimates a pairlist buffer length
59 * given a target energy drift per atom per picosecond.
60 * This is done by estimating the drift given a buffer length.
61 * Ideally we would like to have a tight overestimate of the drift,
62 * but that can be difficult to achieve.
64 * Significant approximations used:
66 * Uniform particle density. UNDERESTIMATES the drift by rho_global/rho_local.
68 * Interactions don't affect particle motion. OVERESTIMATES the drift on longer
69 * time scales. This approximation probably introduces the largest errors.
71 * Only take one constraint per particle into account: OVERESTIMATES the drift.
73 * For rotating constraints assume the same functional shape for time scales
74 * where the constraints rotate significantly as the exact expression for
75 * short time scales. OVERESTIMATES the drift on long time scales.
77 * For non-linear virtual sites use the mass of the lightest constructing atom
78 * to determine the displacement. OVER/UNDERESTIMATES the drift, depending on
79 * the geometry and masses of constructing atoms.
81 * Note that the formulas for normal atoms and linear virtual sites are exact,
82 * apart from the first two approximations.
84 * Note that apart from the effect of the above approximations, the actual
85 * drift of the total energy of a system can be orders of magnitude smaller
86 * due to cancellation of positive and negative drift for different pairs.
90 /* Struct for unique atom type for calculating the energy drift.
91 * The atom displacement depends on mass and constraints.
92 * The energy jump for given distance depend on LJ type and q.
94 struct VerletbufAtomtype
96 atom_nonbonded_kinetic_prop_t prop; /* non-bonded and kinetic atom prop. */
97 int n; /* #atoms of this type in the system */
100 // Struct for derivatives of a non-bonded interaction potential
101 struct pot_derivatives_t
103 real md1; // -V' at the cutoff
104 real d2; // V'' at the cutoff
105 real md3; // -V''' at the cutoff
108 VerletbufListSetup verletbufGetListSetup(int nbnxnKernelType)
110 /* Note that the current buffer estimation code only handles clusters
111 * of size 1, 2 or 4, so for 4x8 or 8x8 we use the estimate for 4x4.
113 VerletbufListSetup listSetup;
115 listSetup.cluster_size_i = nbnxn_kernel_to_cluster_i_size(nbnxnKernelType);
116 listSetup.cluster_size_j = nbnxn_kernel_to_cluster_j_size(nbnxnKernelType);
118 if (nbnxnKernelType == nbnxnk8x8x8_GPU ||
119 nbnxnKernelType == nbnxnk8x8x8_PlainC)
121 /* The GPU kernels (except for OpenCL) split the j-clusters in two halves */
122 listSetup.cluster_size_j /= 2;
128 VerletbufListSetup verletbufGetSafeListSetup(ListSetupType listType)
130 /* When calling this function we often don't know which kernel type we
131 * are going to use. We choose the kernel type with the smallest possible
132 * i- and j-cluster sizes, so we potentially overestimate, but never
133 * underestimate, the buffer drift.
137 if (listType == ListSetupType::Gpu)
139 nbnxnKernelType = nbnxnk8x8x8_GPU;
141 else if (GMX_SIMD && listType == ListSetupType::CpuSimdWhenSupported)
143 #ifdef GMX_NBNXN_SIMD_2XNN
144 /* We use the smallest cluster size to be on the safe side */
145 nbnxnKernelType = nbnxnk4xN_SIMD_2xNN;
147 nbnxnKernelType = nbnxnk4xN_SIMD_4xN;
152 nbnxnKernelType = nbnxnk4x4_PlainC;
155 return verletbufGetListSetup(nbnxnKernelType);
158 // Returns whether prop1 and prop2 are identical
160 atom_nonbonded_kinetic_prop_equal(const atom_nonbonded_kinetic_prop_t &prop1,
161 const atom_nonbonded_kinetic_prop_t &prop2)
163 return (prop1.mass == prop2.mass &&
164 prop1.type == prop2.type &&
165 prop1.q == prop2.q &&
166 prop1.bConstr == prop2.bConstr &&
167 prop1.con_mass == prop2.con_mass &&
168 prop1.con_len == prop2.con_len);
171 static void addAtomtype(std::vector<VerletbufAtomtype> *att,
172 const atom_nonbonded_kinetic_prop_t &prop,
177 /* Ignore massless particles */
182 while (i < att->size() &&
183 !atom_nonbonded_kinetic_prop_equal(prop, (*att)[i].prop))
194 att->push_back({ prop, nmol });
198 /* Returns the mass of atom atomIndex or 1 when setMassesToOne=true */
199 static real getMass(const t_atoms &atoms,
205 return atoms.atom[atomIndex].m;
213 // Set the masses of a vsites in vsite_m and the non-linear vsite count in n_nonlin_vsite
214 static void get_vsite_masses(const gmx_moltype_t &moltype,
215 const gmx_ffparams_t &ffparams,
217 gmx::ArrayRef<real> vsite_m,
220 GMX_RELEASE_ASSERT(n_nonlin_vsite, "Expect a valid pointer");
224 /* Check for virtual sites, determine mass from constructing atoms */
225 for (const auto &ilist : extractILists(moltype.ilist, IF_VSITE))
227 for (size_t i = 0; i < ilist.iatoms.size(); i += ilistStride(ilist))
229 const t_iparams &ip = ffparams.iparams[ilist.iatoms[i]];
230 const int a1 = ilist.iatoms[i + 1];
232 if (ilist.functionType != F_VSITEN)
234 /* Only vsiten can have more than four
235 constructing atoms, so NRAL(ft) <= 5 */
236 const int maxj = NRAL(ilist.functionType);
237 std::vector<real> cam(maxj, 0);
238 GMX_ASSERT(maxj <= 5, "This code expect at most 5 atoms in a vsite");
239 for (int j = 1; j < maxj; j++)
241 const int aj = ilist.iatoms[i + 1 + j];
242 cam[j] = getMass(moltype.atoms, aj, setMassesToOne);
245 cam[j] = vsite_m[aj];
247 /* A vsite should be constructed from normal atoms or
248 * vsites of lower complexity, which we have processed
249 * in a previous iteration.
251 GMX_ASSERT(cam[j] != 0, "We should have a non-zero mass");
254 switch (ilist.functionType)
258 vsite_m[a1] = (cam[1]*cam[2])/(cam[2]*gmx::square(1 - ip.vsite.a) + cam[1]*gmx::square(ip.vsite.a));
262 vsite_m[a1] = (cam[1]*cam[2]*cam[3])/(cam[2]*cam[3]*gmx::square(1 - ip.vsite.a - ip.vsite.b) + cam[1]*cam[3]*gmx::square(ip.vsite.a) + cam[1]*cam[2]*gmx::square(ip.vsite.b));
265 GMX_RELEASE_ASSERT(false, "VsiteN should not end up in this code path");
268 /* Use the mass of the lightest constructing atom.
269 * This is an approximation.
270 * If the distance of the virtual site to the
271 * constructing atom is less than all distances
272 * between constructing atoms, this is a safe
273 * over-estimate of the displacement of the vsite.
274 * This condition holds for all H mass replacement
275 * vsite constructions, except for SP2/3 groups.
276 * In SP3 groups one H will have a F_VSITE3
277 * construction, so even there the total drift
278 * estimate shouldn't be far off.
280 vsite_m[a1] = cam[1];
281 for (int j = 2; j < maxj; j++)
283 vsite_m[a1] = std::min(vsite_m[a1], cam[j]);
293 int numConstructingAtoms = ffparams.iparams[ilist.iatoms[i]].vsiten.n;
294 for (int j = 0; j < 3*numConstructingAtoms; j += 3)
296 int aj = ilist.iatoms[i + j + 2];
297 real coeff = ffparams.iparams[ilist.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 the loop increment of i for processes more than 1 entry */
315 i += (numConstructingAtoms - 1)*ilistStride(ilist);
319 fprintf(debug, "atom %4d %-20s mass %6.3f\n",
320 a1, interaction_function[ilist.functionType].longname, vsite_m[a1]);
326 static std::vector<VerletbufAtomtype>
327 get_verlet_buffer_atomtypes(const gmx_mtop_t *mtop,
331 std::vector<VerletbufAtomtype> att;
332 int ft, i, a1, a2, a3, a;
334 int n_nonlin_vsite_mol;
336 if (n_nonlin_vsite != nullptr)
341 for (const gmx_molblock_t &molblock : mtop->molblock)
343 int nmol = molblock.nmol;
344 const gmx_moltype_t &moltype = mtop->moltype[molblock.type];
345 const t_atoms *atoms = &moltype.atoms;
347 /* Check for constraints, as they affect the kinetic energy.
348 * For virtual sites we need the masses and geometry of
349 * the constructing atoms to determine their velocity distribution.
350 * Thus we need a list of properties for all atoms which
351 * we partially fill when looping over constraints.
353 std::vector<atom_nonbonded_kinetic_prop_t> prop(atoms->nr);
355 for (ft = F_CONSTR; ft <= F_CONSTRNC; ft++)
357 const InteractionList &il = moltype.ilist[ft];
359 for (i = 0; i < il.size(); i += 1+NRAL(ft))
361 ip = &mtop->ffparams.iparams[il.iatoms[i]];
364 real mass1 = getMass(*atoms, a1, setMassesToOne);
365 real mass2 = getMass(*atoms, a2, setMassesToOne);
366 if (mass2 > prop[a1].con_mass)
368 prop[a1].con_mass = mass2;
369 prop[a1].con_len = ip->constr.dA;
371 if (mass1 > prop[a2].con_mass)
373 prop[a2].con_mass = mass1;
374 prop[a2].con_len = ip->constr.dA;
379 const InteractionList &il = moltype.ilist[F_SETTLE];
381 for (i = 0; i < il.size(); i += 1+NRAL(F_SETTLE))
383 ip = &mtop->ffparams.iparams[il.iatoms[i]];
387 /* Usually the mass of a1 (usually oxygen) is larger than a2/a3.
388 * If this is not the case, we overestimate the displacement,
389 * which leads to a larger buffer (ok since this is an exotic case).
391 prop[a1].con_mass = getMass(*atoms, a2, setMassesToOne);
392 prop[a1].con_len = ip->settle.doh;
394 prop[a2].con_mass = getMass(*atoms, a1, setMassesToOne);
395 prop[a2].con_len = ip->settle.doh;
397 prop[a3].con_mass = getMass(*atoms, a1, setMassesToOne);
398 prop[a3].con_len = ip->settle.doh;
401 std::vector<real> vsite_m(atoms->nr);
402 get_vsite_masses(moltype,
406 &n_nonlin_vsite_mol);
407 if (n_nonlin_vsite != nullptr)
409 *n_nonlin_vsite += nmol*n_nonlin_vsite_mol;
412 for (a = 0; a < atoms->nr; a++)
414 if (atoms->atom[a].ptype == eptVSite)
416 prop[a].mass = vsite_m[a];
420 prop[a].mass = getMass(*atoms, a, setMassesToOne);
422 prop[a].type = atoms->atom[a].type;
423 prop[a].q = atoms->atom[a].q;
424 /* We consider an atom constrained, #DOF=2, when it is
425 * connected with constraints to (at least one) atom with
426 * a mass of more than 0.4x its own mass. This is not a critical
427 * parameter, since with roughly equal masses the unconstrained
428 * and constrained displacement will not differ much (and both
429 * overestimate the displacement).
431 prop[a].bConstr = (prop[a].con_mass > 0.4*prop[a].mass);
433 addAtomtype(&att, prop[a], nmol);
439 for (size_t a = 0; a < att.size(); a++)
441 fprintf(debug, "type %zu: m %5.2f t %d q %6.3f con %s con_m %5.3f con_l %5.3f n %d\n",
442 a, att[a].prop.mass, att[a].prop.type, att[a].prop.q,
443 gmx::boolToString(att[a].prop.bConstr), att[a].prop.con_mass, att[a].prop.con_len,
451 /* This function computes two components of the estimate of the variance
452 * in the displacement of one atom in a system of two constrained atoms.
453 * Returns in sigma2_2d the variance due to rotation of the constrained
454 * atom around the atom to which it constrained.
455 * Returns in sigma2_3d the variance due to displacement of the COM
456 * of the whole system of the two constrained atoms.
458 * Note that we only take a single constraint (the one to the heaviest atom)
459 * into account. If an atom has multiple constraints, this will result in
460 * an overestimate of the displacement, which gives a larger drift and buffer.
462 void constrained_atom_sigma2(real kT_fac,
463 const atom_nonbonded_kinetic_prop_t *prop,
467 /* Here we decompose the motion of a constrained atom into two
468 * components: rotation around the COM and translation of the COM.
471 /* Determine the variance of the arc length for the two rotational DOFs */
472 real massFraction = prop->con_mass/(prop->mass + prop->con_mass);
473 real sigma2_rot = kT_fac*massFraction/prop->mass;
475 /* The distance from the atom to the COM, i.e. the rotational arm */
476 real comDistance = prop->con_len*massFraction;
478 /* The variance relative to the arm */
479 real sigma2_rel = sigma2_rot/gmx::square(comDistance);
481 /* For sigma2_rel << 1 we don't notice the rotational effect and
482 * we have a normal, Gaussian displacement distribution.
483 * For larger sigma2_rel the displacement is much less, in fact it can
484 * not exceed 2*comDistance. We can calculate MSD/arm^2 as:
485 * integral_x=0-inf distance2(x) x/sigma2_rel exp(-x^2/(2 sigma2_rel)) dx
486 * where x is angular displacement and distance2(x) is the distance^2
487 * between points at angle 0 and x:
488 * distance2(x) = (sin(x) - sin(0))^2 + (cos(x) - cos(0))^2
489 * The limiting value of this MSD is 2, which is also the value for
490 * a uniform rotation distribution that would be reached at long time.
491 * The maximum is 2.5695 at sigma2_rel = 4.5119.
492 * We approximate this integral with a rational polynomial with
493 * coefficients from a Taylor expansion. This approximation is an
494 * overestimate for all values of sigma2_rel. Its maximum value
495 * of 2.6491 is reached at sigma2_rel = sqrt(45/2) = 4.7434.
496 * We keep the approximation constant after that.
497 * We use this approximate MSD as the variance for a Gaussian distribution.
499 * NOTE: For any sensible buffer tolerance this will result in a (large)
500 * overestimate of the buffer size, since the Gaussian has a long tail,
501 * whereas the actual distribution can not reach values larger than 2.
503 /* Coeffients obtained from a Taylor expansion */
504 const real a = 1.0/3.0;
505 const real b = 2.0/45.0;
507 /* Our approximation is constant after sigma2_rel = 1/sqrt(b) */
508 sigma2_rel = std::min(sigma2_rel, 1/std::sqrt(b));
510 /* Compute the approximate sigma^2 for 2D motion due to the rotation */
511 *sigma2_2d = gmx::square(comDistance)*
512 sigma2_rel/(1 + a*sigma2_rel + b*gmx::square(sigma2_rel));
514 /* The constrained atom also moves (in 3D) with the COM of both atoms */
515 *sigma2_3d = kT_fac/(prop->mass + prop->con_mass);
518 static void get_atom_sigma2(real kT_fac,
519 const atom_nonbonded_kinetic_prop_t *prop,
525 /* Complicated constraint calculation in a separate function */
526 constrained_atom_sigma2(kT_fac, prop, sigma2_2d, sigma2_3d);
530 /* Unconstrained atom: trivial */
532 *sigma2_3d = kT_fac/prop->mass;
536 static void approx_2dof(real s2, real x, real *shift, real *scale)
538 /* A particle with 1 DOF constrained has 2 DOFs instead of 3.
539 * This code is also used for particles with multiple constraints,
540 * in which case we overestimate the displacement.
541 * The 2DOF distribution is sqrt(pi/2)*erfc(r/(sqrt(2)*s))/(2*s).
542 * We approximate this with scale*Gaussian(s,r+shift),
543 * by matching the distribution value and derivative at x.
544 * This is a tight overestimate for all r>=0 at any s and x.
548 ex = std::exp(-x*x/(2*s2));
549 er = std::erfc(x/std::sqrt(2*s2));
551 *shift = -x + std::sqrt(2*s2/M_PI)*ex/er;
552 *scale = 0.5*M_PI*std::exp(ex*ex/(M_PI*er*er))*er;
555 // Returns an (over)estimate of the energy drift for a single atom pair,
556 // given the kinetic properties, displacement variances and list buffer.
557 static real energyDriftAtomPair(bool isConstrained_i,
558 bool isConstrained_j,
559 real s2, real s2i_2d, real s2j_2d,
561 const pot_derivatives_t *der)
563 // For relatively small arguments erfc() is so small that if will be 0.0
564 // when stored in a float. We set an argument limit of 8 (Erfc(8)=1e-29),
565 // such that we can divide by erfc and have some space left for arithmetic.
566 const real erfc_arg_max = 8.0;
573 if (rsh*rsh > 2*s2*erfc_arg_max*erfc_arg_max)
575 // Below we calculate c_erfc = 0.5*erfc(rsh/sqrt(2*s2))
576 // When rsh/sqrt(2*s2) increases, this erfc will be the first
577 // result that underflows and becomes 0.0. To avoid this,
578 // we set c_exp=0 and c_erfc=0 for large arguments.
579 // This also avoids NaN in approx_2dof().
580 // In any relevant case this has no effect on the results,
581 // since c_exp < 6e-29, so the displacement is completely
582 // negligible for such atom pairs (and an overestimate).
583 // In nearly all use cases, there will be other atom pairs
584 // that contribute much more to the total, so zeroing
585 // this particular contribution has no effect at all.
591 /* For constraints: adapt r and scaling for the Gaussian */
596 approx_2dof(s2i_2d, r_buffer*s2i_2d/s2, &sh, &sc);
604 approx_2dof(s2j_2d, r_buffer*s2j_2d/s2, &sh, &sc);
609 /* Exact contribution of an atom pair with Gaussian displacement
610 * with sigma s to the energy drift for a potential with
611 * derivative -md and second derivative dd at the cut-off.
612 * The only catch is that for potentials that change sign
613 * near the cut-off there could be an unlucky compensation
614 * of positive and negative energy drift.
615 * Such potentials are extremely rare though.
617 * Note that pot has unit energy*length, as the linear
618 * atom density still needs to be put in.
620 c_exp = std::exp(-rsh*rsh/(2*s2))/std::sqrt(2*M_PI);
621 c_erfc = 0.5*std::erfc(rsh/(std::sqrt(2*s2)));
623 real s = std::sqrt(s2);
627 der->md1/2*((rsh2 + s2)*c_erfc - rsh*s*c_exp);
629 der->d2/6*(s*(rsh2 + 2*s2)*c_exp - rsh*(rsh2 + 3*s2)*c_erfc);
631 der->md3/24*((rsh2*rsh2 + 6*rsh2*s2 + 3*s2*s2)*c_erfc - rsh*s*(rsh2 + 5*s2)*c_exp);
633 return pot1 + pot2 + pot3;
636 // Computes and returns an estimate of the energy drift for the whole system
637 static real energyDrift(gmx::ArrayRef<const VerletbufAtomtype> att,
638 const gmx_ffparams_t *ffp,
640 const pot_derivatives_t *ljDisp,
641 const pot_derivatives_t *ljRep,
642 const pot_derivatives_t *elec,
643 real rlj, real rcoulomb,
644 real rlist, real boxvol)
646 double drift_tot = 0;
650 /* No atom displacements: no drift, avoid division by 0 */
654 // Here add up the contribution of all atom pairs in the system to
655 // (estimated) energy drift by looping over all atom type pairs.
656 for (int i = 0; i < att.size(); i++)
658 // Get the thermal displacement variance for the i-atom type
659 const atom_nonbonded_kinetic_prop_t *prop_i = &att[i].prop;
661 get_atom_sigma2(kT_fac, prop_i, &s2i_2d, &s2i_3d);
663 for (int j = i; j < att.size(); j++)
665 // Get the thermal displacement variance for the j-atom type
666 const atom_nonbonded_kinetic_prop_t *prop_j = &att[j].prop;
668 get_atom_sigma2(kT_fac, prop_j, &s2j_2d, &s2j_3d);
670 /* Add up the up to four independent variances */
671 real s2 = s2i_2d + s2i_3d + s2j_2d + s2j_3d;
673 // Set -V', V'' and -V''' at the cut-off for LJ */
674 real c6 = ffp->iparams[prop_i->type*ffp->atnr + prop_j->type].lj.c6;
675 real c12 = ffp->iparams[prop_i->type*ffp->atnr + prop_j->type].lj.c12;
676 pot_derivatives_t lj;
677 lj.md1 = c6*ljDisp->md1 + c12*ljRep->md1;
678 lj.d2 = c6*ljDisp->d2 + c12*ljRep->d2;
679 lj.md3 = c6*ljDisp->md3 + c12*ljRep->md3;
681 real pot_lj = energyDriftAtomPair(prop_i->bConstr, prop_j->bConstr,
686 // Set -V' and V'' at the cut-off for Coulomb
687 pot_derivatives_t elec_qq;
688 elec_qq.md1 = elec->md1*prop_i->q*prop_j->q;
689 elec_qq.d2 = elec->d2 *prop_i->q*prop_j->q;
692 real pot_q = energyDriftAtomPair(prop_i->bConstr, prop_j->bConstr,
697 // Note that attractive and repulsive potentials for individual
698 // pairs can partially cancel.
699 real pot = pot_lj + pot_q;
701 /* Multiply by the number of atom pairs */
704 pot *= static_cast<double>(att[i].n)*(att[i].n - 1)/2;
708 pot *= static_cast<double>(att[i].n)*att[j].n;
710 /* We need the line density to get the energy drift of the system.
711 * The effective average r^2 is close to (rlist+sigma)^2.
713 pot *= 4*M_PI*gmx::square(rlist + std::sqrt(s2))/boxvol;
715 /* Add the unsigned drift to avoid cancellation of errors */
716 drift_tot += std::abs(pot);
723 // Returns the chance that a particle in a cluster is at distance rlist
724 // when the cluster is at distance rlist
725 static real surface_frac(int cluster_size, real particle_distance, real rlist)
729 if (rlist < 0.5*particle_distance)
731 /* We have non overlapping spheres */
735 /* Half the inter-particle distance relative to rlist */
736 d = 0.5*particle_distance/rlist;
738 /* Determine the area of the surface at distance rlist to the closest
739 * particle, relative to surface of a sphere of radius rlist.
740 * The formulas below assume close to cubic cells for the pair search grid,
741 * which the pair search code tries to achieve.
742 * Note that in practice particle distances will not be delta distributed,
743 * but have some spread, often involving shorter distances,
744 * as e.g. O-H bonds in a water molecule. Thus the estimates below will
745 * usually be slightly too high and thus conservative.
747 switch (cluster_size)
750 /* One particle: trivial */
754 /* Two particles: two spheres at fractional distance 2*a */
758 /* We assume a perfect, symmetric tetrahedron geometry.
759 * The surface around a tetrahedron is too complex for a full
760 * analytical solution, so we use a Taylor expansion.
762 area_rel = (1.0 + 1/M_PI*(6*std::acos(1/std::sqrt(3))*d +
763 std::sqrt(3)*d*d*(1.0 +
766 83.0/756.0*d*d*d*d*d*d)));
769 gmx_incons("surface_frac called with unsupported cluster_size");
772 return area_rel/cluster_size;
775 /* Returns the negative of the third derivative of a potential r^-p
776 * with a force-switch function, evaluated at the cut-off rc.
778 static real md3_force_switch(real p, real rswitch, real rc)
780 /* The switched force function is:
781 * p*r^-(p+1) + a*(r - rswitch)^2 + b*(r - rswitch)^3
784 real md3_pot, md3_sw;
786 a = -((p + 4)*rc - (p + 1)*rswitch)/(pow(rc, p+2)*gmx::square(rc-rswitch));
787 b = ((p + 3)*rc - (p + 1)*rswitch)/(pow(rc, p+2)*gmx::power3(rc-rswitch));
789 md3_pot = (p + 2)*(p + 1)*p*pow(rc, p+3);
790 md3_sw = 2*a + 6*b*(rc - rswitch);
792 return md3_pot + md3_sw;
795 /* Returns the variance of the atomic displacement over timePeriod.
797 * Note: When not using BD with a non-mass dependendent friction coefficient,
798 * the return value still needs to be divided by the particle mass.
800 static real displacementVariance(const t_inputrec &ir,
808 /* Get the displacement distribution from the random component only.
809 * With accurate integration the systematic (force) displacement
810 * should be negligible (unless nstlist is extremely large, which
811 * you wouldn't do anyhow).
813 kT_fac = 2*BOLTZ*temperature*timePeriod;
816 /* This is directly sigma^2 of the displacement */
817 kT_fac /= ir.bd_fric;
821 /* Per group tau_t is not implemented yet, use the maximum */
822 real tau_t = ir.opts.tau_t[0];
823 for (int i = 1; i < ir.opts.ngtc; i++)
825 tau_t = std::max(tau_t, ir.opts.tau_t[i]);
829 /* This kT_fac needs to be divided by the mass to get sigma^2 */
834 kT_fac = BOLTZ*temperature*gmx::square(timePeriod);
840 /* Returns the largest sigma of the Gaussian displacement over all particle
841 * types. This ignores constraints, so is an overestimate.
843 static real maxSigma(real kT_fac,
844 gmx::ArrayRef<const VerletbufAtomtype> att)
846 GMX_ASSERT(!att.empty(), "We should have at least one type");
847 real smallestMass = att[0].prop.mass;
848 for (int i = 1; i < att.size(); i++)
850 smallestMass = std::min(smallestMass, att[i].prop.mass);
853 return 2*std::sqrt(kT_fac/smallestMass);
856 void calc_verlet_buffer_size(const gmx_mtop_t *mtop, real boxvol,
857 const t_inputrec *ir,
860 real reference_temperature,
861 const VerletbufListSetup *list_setup,
868 real particle_distance;
869 real nb_clust_frac_pairs_not_in_list_at_cutoff;
876 if (!EI_DYNAMICS(ir->eI))
878 gmx_incons("Can only determine the Verlet buffer size for integrators that perform dynamics");
880 if (ir->verletbuf_tol <= 0)
882 gmx_incons("The Verlet buffer tolerance needs to be larger than zero");
885 if (reference_temperature < 0)
887 /* We use the maximum temperature with multiple T-coupl groups.
888 * We could use a per particle temperature, but since particles
889 * interact, this might underestimate the buffer size.
891 reference_temperature = maxReferenceTemperature(*ir);
893 GMX_RELEASE_ASSERT(reference_temperature >= 0, "Without T-coupling we should not end up here");
896 /* Resolution of the buffer size */
899 env = getenv("GMX_VERLET_BUFFER_RES");
902 sscanf(env, "%lf", &resolution);
905 /* In an atom wise pair-list there would be no pairs in the list
906 * beyond the pair-list cut-off.
907 * However, we use a pair-list of groups vs groups of atoms.
908 * For groups of 4 atoms, the parallelism of SSE instructions, only
909 * 10% of the atoms pairs are not in the list just beyond the cut-off.
910 * As this percentage increases slowly compared to the decrease of the
911 * Gaussian displacement distribution over this range, we can simply
912 * reduce the drift by this fraction.
913 * For larger groups, e.g. of 8 atoms, this fraction will be lower,
914 * so then buffer size will be on the conservative (large) side.
916 * Note that the formulas used here do not take into account
917 * cancellation of errors which could occur by missing both
918 * attractive and repulsive interactions.
920 * The only major assumption is homogeneous particle distribution.
921 * For an inhomogeneous system, such as a liquid-vapor system,
922 * the buffer will be underestimated. The actual energy drift
923 * will be higher by the factor: local/homogeneous particle density.
925 * The results of this estimate have been checked againt simulations.
926 * In most cases the real drift differs by less than a factor 2.
929 /* Worst case assumption: HCP packing of particles gives largest distance */
930 particle_distance = std::cbrt(boxvol*std::sqrt(2)/mtop->natoms);
932 /* TODO: Obtain masses through (future) integrator functionality
933 * to avoid scattering the code with (or forgetting) checks.
935 const bool setMassesToOne = (ir->eI == eiBD && ir->bd_fric > 0);
937 get_verlet_buffer_atomtypes(mtop, setMassesToOne, n_nonlin_vsite);
938 GMX_ASSERT(!att.empty(), "We expect at least one type");
942 fprintf(debug, "particle distance assuming HCP packing: %f nm\n",
944 fprintf(debug, "energy drift atom types: %zu\n", att.size());
947 pot_derivatives_t ljDisp = { 0, 0, 0 };
948 pot_derivatives_t ljRep = { 0, 0, 0 };
949 real repPow = mtop->ffparams.reppow;
951 if (ir->vdwtype == evdwCUT)
953 real sw_range, md3_pswf;
955 switch (ir->vdw_modifier)
958 case eintmodPOTSHIFT:
959 /* -dV/dr of -r^-6 and r^-reppow */
960 ljDisp.md1 = -6*std::pow(ir->rvdw, -7.0);
961 ljRep.md1 = repPow*std::pow(ir->rvdw, -(repPow + 1));
962 /* The contribution of the higher derivatives is negligible */
964 case eintmodFORCESWITCH:
965 /* At the cut-off: V=V'=V''=0, so we use only V''' */
966 ljDisp.md3 = -md3_force_switch(6.0, ir->rvdw_switch, ir->rvdw);
967 ljRep.md3 = md3_force_switch(repPow, ir->rvdw_switch, ir->rvdw);
969 case eintmodPOTSWITCH:
970 /* At the cut-off: V=V'=V''=0.
971 * V''' is given by the original potential times
972 * the third derivative of the switch function.
974 sw_range = ir->rvdw - ir->rvdw_switch;
975 md3_pswf = 60.0/gmx::power3(sw_range);
977 ljDisp.md3 = -std::pow(ir->rvdw, -6.0 )*md3_pswf;
978 ljRep.md3 = std::pow(ir->rvdw, -repPow)*md3_pswf;
981 gmx_incons("Unimplemented VdW modifier");
984 else if (EVDW_PME(ir->vdwtype))
986 real b = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
992 // -dV/dr of g(br)*r^-6 [where g(x) = exp(-x^2)(1+x^2+x^4/2),
993 // see LJ-PME equations in manual] and r^-reppow
994 ljDisp.md1 = -std::exp(-br2)*(br6 + 3.0*br4 + 6.0*br2 + 6.0)*std::pow(r, -7.0);
995 ljRep.md1 = repPow*pow(r, -(repPow + 1));
996 // The contribution of the higher derivatives is negligible
1000 gmx_fatal(FARGS, "Energy drift calculation is only implemented for plain cut-off Lennard-Jones interactions");
1003 elfac = ONE_4PI_EPS0/ir->epsilon_r;
1005 // Determine the 1st and 2nd derivative for the electostatics
1006 pot_derivatives_t elec = { 0, 0, 0 };
1008 if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
1012 if (ir->coulombtype == eelCUT)
1019 eps_rf = ir->epsilon_rf/ir->epsilon_r;
1022 k_rf = (eps_rf - ir->epsilon_r)/( gmx::power3(ir->rcoulomb) * (2*eps_rf + ir->epsilon_r) );
1026 /* epsilon_rf = infinity */
1027 k_rf = 0.5/gmx::power3(ir->rcoulomb);
1033 elec.md1 = elfac*(1.0/gmx::square(ir->rcoulomb) - 2*k_rf*ir->rcoulomb);
1035 elec.d2 = elfac*(2.0/gmx::power3(ir->rcoulomb) + 2*k_rf);
1037 else if (EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD)
1041 b = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
1044 elec.md1 = elfac*(b*std::exp(-br*br)*M_2_SQRTPI/rc + std::erfc(br)/(rc*rc));
1045 elec.d2 = elfac/(rc*rc)*(2*b*(1 + br*br)*std::exp(-br*br)*M_2_SQRTPI + 2*std::erfc(br)/rc);
1049 gmx_fatal(FARGS, "Energy drift calculation is only implemented for Reaction-Field and Ewald electrostatics");
1052 /* Determine the variance of the atomic displacement
1053 * over list_lifetime steps: kT_fac
1054 * For inertial dynamics (not Brownian dynamics) the mass factor
1055 * is not included in kT_fac, it is added later.
1057 const real kT_fac = displacementVariance(*ir, reference_temperature,
1058 list_lifetime*ir->delta_t);
1062 fprintf(debug, "Derivatives of non-bonded potentials at the cut-off:\n");
1063 fprintf(debug, "LJ disp. -V' %9.2e V'' %9.2e -V''' %9.2e\n", ljDisp.md1, ljDisp.d2, ljDisp.md3);
1064 fprintf(debug, "LJ rep. -V' %9.2e V'' %9.2e -V''' %9.2e\n", ljRep.md1, ljRep.d2, ljRep.md3);
1065 fprintf(debug, "Electro. -V' %9.2e V'' %9.2e\n", elec.md1, elec.d2);
1066 fprintf(debug, "sqrt(kT_fac) %f\n", std::sqrt(kT_fac));
1069 /* Search using bisection */
1071 /* The drift will be neglible at 5 times the max sigma */
1072 ib1 = static_cast<int>(5*maxSigma(kT_fac, att)/resolution) + 1;
1073 while (ib1 - ib0 > 1)
1077 rl = std::max(ir->rvdw, ir->rcoulomb) + rb;
1079 /* Calculate the average energy drift at the last step
1080 * of the nstlist steps at which the pair-list is used.
1082 drift = energyDrift(att, &mtop->ffparams,
1084 &ljDisp, &ljRep, &elec,
1085 ir->rvdw, ir->rcoulomb,
1088 /* Correct for the fact that we are using a Ni x Nj particle pair list
1089 * and not a 1 x 1 particle pair list. This reduces the drift.
1091 /* We don't have a formula for 8 (yet), use 4 which is conservative */
1092 nb_clust_frac_pairs_not_in_list_at_cutoff =
1093 surface_frac(std::min(list_setup->cluster_size_i, 4),
1094 particle_distance, rl)*
1095 surface_frac(std::min(list_setup->cluster_size_j, 4),
1096 particle_distance, rl);
1097 drift *= nb_clust_frac_pairs_not_in_list_at_cutoff;
1099 /* Convert the drift to drift per unit time per atom */
1100 drift /= nstlist*ir->delta_t*mtop->natoms;
1104 fprintf(debug, "ib %3d %3d %3d rb %.3f %dx%d fac %.3f drift %.1e\n",
1106 list_setup->cluster_size_i, list_setup->cluster_size_j,
1107 nb_clust_frac_pairs_not_in_list_at_cutoff,
1111 if (std::abs(drift) > ir->verletbuf_tol)
1121 *rlist = std::max(ir->rvdw, ir->rcoulomb) + ib1*resolution;
1124 /* Returns the pairlist buffer size for use as a minimum buffer size
1126 * Note that this is a rather crude estimate. It is ok for a buffer
1127 * set for good energy conservation or RF electrostatics. But it is
1128 * too small with PME and the buffer set with the default tolerance.
1130 static real minCellSizeFromPairlistBuffer(const t_inputrec &ir)
1132 return ir.rlist - std::max(ir.rvdw, ir.rcoulomb);
1135 real minCellSizeForAtomDisplacement(const gmx_mtop_t &mtop,
1136 const t_inputrec &ir,
1137 real chanceRequested)
1139 if (!EI_DYNAMICS(ir.eI) || (EI_MD(ir.eI) && ir.etc == etcNO))
1141 return minCellSizeFromPairlistBuffer(ir);
1144 /* We use the maximum temperature with multiple T-coupl groups.
1145 * We could use a per particle temperature, but since particles
1146 * interact, this might underestimate the displacements.
1148 const real temperature = maxReferenceTemperature(ir);
1150 const bool setMassesToOne = (ir.eI == eiBD && ir.bd_fric > 0);
1152 const auto atomtypes = get_verlet_buffer_atomtypes(&mtop, setMassesToOne, nullptr);
1154 const real kT_fac = displacementVariance(ir, temperature,
1155 ir.nstlist*ir.delta_t);
1157 /* Resolution of the cell size */
1158 real resolution = 0.001;
1160 /* Search using bisection, avoid 0 and start at 1 */
1162 /* The chance will be neglible at 10 times the max sigma */
1163 int ib1 = int(10*maxSigma(kT_fac, atomtypes)/resolution) + 1;
1165 while (ib1 - ib0 > 1)
1167 int ib = (ib0 + ib1)/2;
1168 cellSize = ib*resolution;
1170 /* We assumes atom are distributed uniformly over the cell width.
1171 * Once an atom has moved by more than the cellSize (as passed
1172 * as the buffer argument to energyDriftAtomPair() below),
1173 * the chance of crossing the boundary of the neighbor cell
1174 * thus increases as 1/cellSize with the additional displacement
1175 * on to of cellSize. We thus create a linear interaction with
1176 * derivative = -1/cellSize. Using this in the energyDriftAtomPair
1177 * function will return the chance of crossing the next boundary.
1179 const pot_derivatives_t boundaryInteraction = { 1/cellSize, 0, 0 };
1182 for (const VerletbufAtomtype &att : atomtypes)
1184 const atom_nonbonded_kinetic_prop_t &propAtom = att.prop;
1187 get_atom_sigma2(kT_fac, &propAtom, &s2_2d, &s2_3d);
1189 real chancePerAtom = energyDriftAtomPair(propAtom.bConstr, false,
1190 s2_2d + s2_3d, s2_2d, 0,
1192 &boundaryInteraction);
1194 if (propAtom.bConstr)
1196 /* energyDriftAtomPair() uses an unlimited Gaussian displacement
1197 * distribution for constrained atoms, whereas they can
1198 * actually not move more than the COM of the two constrained
1199 * atoms plus twice the distance from the COM.
1200 * Use this maximum, limited displacement when this results in
1201 * a smaller chance (note that this is still an overestimate).
1203 real massFraction = propAtom.con_mass/(propAtom.mass + propAtom.con_mass);
1204 real comDistance = propAtom.con_len*massFraction;
1206 real chanceWithMaxDistance =
1207 energyDriftAtomPair(false, false,
1209 cellSize - 2*comDistance,
1210 &boundaryInteraction);
1211 chancePerAtom = std::min(chancePerAtom, chanceWithMaxDistance);
1214 /* Take into account the line density of the boundary */
1215 chancePerAtom /= cellSize;
1217 chance += att.n*chancePerAtom;
1220 /* Note: chance is for every nstlist steps */
1221 if (chance > chanceRequested*ir.nstlist)