Add C++ version of t_ilist
[gromacs.git] / src / gromacs / mdlib / calc_verletbuf.cpp
blob80a7f19f7b2f1ab69da819694b0fb28a1d08a042
1 /*
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.
35 #include "gmxpre.h"
37 #include "calc_verletbuf.h"
39 #include <cassert>
40 #include <cmath>
41 #include <cstdlib>
43 #include <algorithm>
45 #include "gromacs/ewald/ewald-utils.h"
46 #include "gromacs/math/functions.h"
47 #include "gromacs/math/units.h"
48 #include "gromacs/math/vec.h"
49 #include "gromacs/mdlib/nb_verlet.h"
50 #include "gromacs/mdlib/nbnxn_simd.h"
51 #include "gromacs/mdlib/nbnxn_util.h"
52 #include "gromacs/mdtypes/inputrec.h"
53 #include "gromacs/mdtypes/md_enums.h"
54 #include "gromacs/topology/ifunc.h"
55 #include "gromacs/topology/topology.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/smalloc.h"
58 #include "gromacs/utility/strconvert.h"
60 /* The code in this file estimates a pairlist buffer length
61 * given a target energy drift per atom per picosecond.
62 * This is done by estimating the drift given a buffer length.
63 * Ideally we would like to have a tight overestimate of the drift,
64 * but that can be difficult to achieve.
66 * Significant approximations used:
68 * Uniform particle density. UNDERESTIMATES the drift by rho_global/rho_local.
70 * Interactions don't affect particle motion. OVERESTIMATES the drift on longer
71 * time scales. This approximation probably introduces the largest errors.
73 * Only take one constraint per particle into account: OVERESTIMATES the drift.
75 * For rotating constraints assume the same functional shape for time scales
76 * where the constraints rotate significantly as the exact expression for
77 * short time scales. OVERESTIMATES the drift on long time scales.
79 * For non-linear virtual sites use the mass of the lightest constructing atom
80 * to determine the displacement. OVER/UNDERESTIMATES the drift, depending on
81 * the geometry and masses of constructing atoms.
83 * Note that the formulas for normal atoms and linear virtual sites are exact,
84 * apart from the first two approximations.
86 * Note that apart from the effect of the above approximations, the actual
87 * drift of the total energy of a system can be orders of magnitude smaller
88 * due to cancellation of positive and negative drift for different pairs.
92 /* Struct for unique atom type for calculating the energy drift.
93 * The atom displacement depends on mass and constraints.
94 * The energy jump for given distance depend on LJ type and q.
96 typedef struct
98 atom_nonbonded_kinetic_prop_t prop; /* non-bonded and kinetic atom prop. */
99 int n; /* #atoms of this type in the system */
100 } verletbuf_atomtype_t;
102 // Struct for derivatives of a non-bonded interaction potential
103 typedef struct
105 real md1; // -V' at the cutoff
106 real d2; // V'' at the cutoff
107 real md3; // -V''' at the cutoff
108 } pot_derivatives_t;
110 VerletbufListSetup verletbufGetListSetup(int nbnxnKernelType)
112 /* Note that the current buffer estimation code only handles clusters
113 * of size 1, 2 or 4, so for 4x8 or 8x8 we use the estimate for 4x4.
115 VerletbufListSetup listSetup;
117 listSetup.cluster_size_i = nbnxn_kernel_to_cluster_i_size(nbnxnKernelType);
118 listSetup.cluster_size_j = nbnxn_kernel_to_cluster_j_size(nbnxnKernelType);
120 if (nbnxnKernelType == nbnxnk8x8x8_GPU ||
121 nbnxnKernelType == nbnxnk8x8x8_PlainC)
123 /* The GPU kernels (except for OpenCL) split the j-clusters in two halves */
124 listSetup.cluster_size_j /= 2;
127 return listSetup;
130 VerletbufListSetup verletbufGetSafeListSetup(ListSetupType listType)
132 /* When calling this function we often don't know which kernel type we
133 * are going to use. We choose the kernel type with the smallest possible
134 * i- and j-cluster sizes, so we potentially overestimate, but never
135 * underestimate, the buffer drift.
137 int nbnxnKernelType;
139 if (listType == ListSetupType::Gpu)
141 nbnxnKernelType = nbnxnk8x8x8_GPU;
143 else if (GMX_SIMD && listType == ListSetupType::CpuSimdWhenSupported)
145 #ifdef GMX_NBNXN_SIMD_2XNN
146 /* We use the smallest cluster size to be on the safe side */
147 nbnxnKernelType = nbnxnk4xN_SIMD_2xNN;
148 #else
149 nbnxnKernelType = nbnxnk4xN_SIMD_4xN;
150 #endif
152 else
154 nbnxnKernelType = nbnxnk4x4_PlainC;
157 return verletbufGetListSetup(nbnxnKernelType);
160 static gmx_bool
161 atom_nonbonded_kinetic_prop_equal(const atom_nonbonded_kinetic_prop_t *prop1,
162 const atom_nonbonded_kinetic_prop_t *prop2)
164 return (prop1->mass == prop2->mass &&
165 prop1->type == prop2->type &&
166 prop1->q == prop2->q &&
167 prop1->bConstr == prop2->bConstr &&
168 prop1->con_mass == prop2->con_mass &&
169 prop1->con_len == prop2->con_len);
172 static void add_at(verletbuf_atomtype_t **att_p, int *natt_p,
173 const atom_nonbonded_kinetic_prop_t *prop,
174 int nmol)
176 verletbuf_atomtype_t *att;
177 int natt, i;
179 if (prop->mass == 0)
181 /* Ignore massless particles */
182 return;
185 att = *att_p;
186 natt = *natt_p;
188 i = 0;
189 while (i < natt && !atom_nonbonded_kinetic_prop_equal(prop, &att[i].prop))
191 i++;
194 if (i < natt)
196 att[i].n += nmol;
198 else
200 (*natt_p)++;
201 srenew(*att_p, *natt_p);
202 (*att_p)[i].prop = *prop;
203 (*att_p)[i].n = nmol;
207 /* Returns the mass of atom atomIndex or 1 when setMassesToOne=true */
208 static real getMass(const t_atoms &atoms,
209 int atomIndex,
210 bool setMassesToOne)
212 if (!setMassesToOne)
214 return atoms.atom[atomIndex].m;
216 else
218 return 1;
222 static void get_vsite_masses(const gmx_moltype_t *moltype,
223 const gmx_ffparams_t *ffparams,
224 bool setMassesToOne,
225 real *vsite_m,
226 int *n_nonlin_vsite)
228 int ft, i;
230 *n_nonlin_vsite = 0;
232 /* Check for virtual sites, determine mass from constructing atoms */
233 for (ft = 0; ft < F_NRE; ft++)
235 if (IS_VSITE(ft))
237 const InteractionList *il = &moltype->ilist[ft];
239 for (i = 0; i < il->size(); i += 1+NRAL(ft))
241 const t_iparams *ip;
242 real inv_mass, coeff, m_aj;
243 int a1;
245 ip = &ffparams->iparams[il->iatoms[i]];
247 a1 = il->iatoms[i+1];
249 if (ft != F_VSITEN)
251 /* Only vsiten can have more than four
252 constructing atoms, so NRAL(ft) <= 5 */
253 int j;
254 real *cam;
255 const int maxj = NRAL(ft);
257 snew(cam, maxj);
258 assert(maxj <= 5);
259 for (j = 1; j < maxj; j++)
261 int aj = il->iatoms[i + 1 + j];
262 cam[j] = getMass(moltype->atoms, aj, setMassesToOne);
263 if (cam[j] == 0)
265 cam[j] = vsite_m[aj];
267 if (cam[j] == 0)
269 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.",
270 *moltype->name,
271 interaction_function[ft].longname,
272 aj + 1);
276 switch (ft)
278 case F_VSITE2:
279 /* Exact */
280 vsite_m[a1] = (cam[1]*cam[2])/(cam[2]*gmx::square(1-ip->vsite.a) + cam[1]*gmx::square(ip->vsite.a));
281 break;
282 case F_VSITE3:
283 /* Exact */
284 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));
285 break;
286 case F_VSITEN:
287 gmx_incons("Invalid vsite type");
288 default:
289 /* Use the mass of the lightest constructing atom.
290 * This is an approximation.
291 * If the distance of the virtual site to the
292 * constructing atom is less than all distances
293 * between constructing atoms, this is a safe
294 * over-estimate of the displacement of the vsite.
295 * This condition holds for all H mass replacement
296 * vsite constructions, except for SP2/3 groups.
297 * In SP3 groups one H will have a F_VSITE3
298 * construction, so even there the total drift
299 * estimate shouldn't be far off.
301 vsite_m[a1] = cam[1];
302 for (j = 2; j < maxj; j++)
304 vsite_m[a1] = std::min(vsite_m[a1], cam[j]);
306 (*n_nonlin_vsite)++;
307 break;
309 sfree(cam);
311 else
313 int j;
315 /* Exact */
316 inv_mass = 0;
317 for (j = 0; j < 3*ffparams->iparams[il->iatoms[i]].vsiten.n; j += 3)
319 int aj = il->iatoms[i + j + 2];
320 coeff = ffparams->iparams[il->iatoms[i+j]].vsiten.a;
321 if (moltype->atoms.atom[aj].ptype == eptVSite)
323 m_aj = vsite_m[aj];
325 else
327 m_aj = moltype->atoms.atom[aj].m;
329 if (m_aj <= 0)
331 gmx_incons("The mass of a vsiten constructing atom is <= 0");
333 inv_mass += coeff*coeff/m_aj;
335 vsite_m[a1] = 1/inv_mass;
336 /* Correct for loop increment of i */
337 i += j - 1 - NRAL(ft);
339 if (gmx_debug_at)
341 fprintf(debug, "atom %4d %-20s mass %6.3f\n",
342 a1, interaction_function[ft].longname, vsite_m[a1]);
349 static void get_verlet_buffer_atomtypes(const gmx_mtop_t *mtop,
350 bool setMassesToOne,
351 verletbuf_atomtype_t **att_p,
352 int *natt_p,
353 int *n_nonlin_vsite)
355 verletbuf_atomtype_t *att;
356 int natt;
357 int ft, i, a1, a2, a3, a;
358 const t_iparams *ip;
359 atom_nonbonded_kinetic_prop_t *prop;
360 real *vsite_m;
361 int n_nonlin_vsite_mol;
363 att = nullptr;
364 natt = 0;
366 if (n_nonlin_vsite != nullptr)
368 *n_nonlin_vsite = 0;
371 for (const gmx_molblock_t &molblock : mtop->molblock)
373 int nmol = molblock.nmol;
374 const gmx_moltype_t &moltype = mtop->moltype[molblock.type];
375 const t_atoms *atoms = &moltype.atoms;
377 /* Check for constraints, as they affect the kinetic energy.
378 * For virtual sites we need the masses and geometry of
379 * the constructing atoms to determine their velocity distribution.
381 snew(prop, atoms->nr);
382 snew(vsite_m, atoms->nr);
384 for (ft = F_CONSTR; ft <= F_CONSTRNC; ft++)
386 const InteractionList *il = &moltype.ilist[ft];
388 for (i = 0; i < il->size(); i += 1+NRAL(ft))
390 ip = &mtop->ffparams.iparams[il->iatoms[i]];
391 a1 = il->iatoms[i+1];
392 a2 = il->iatoms[i+2];
393 real mass1 = getMass(*atoms, a1, setMassesToOne);
394 real mass2 = getMass(*atoms, a2, setMassesToOne);
395 if (mass2 > prop[a1].con_mass)
397 prop[a1].con_mass = mass2;
398 prop[a1].con_len = ip->constr.dA;
400 if (mass1 > prop[a2].con_mass)
402 prop[a2].con_mass = mass1;
403 prop[a2].con_len = ip->constr.dA;
408 const InteractionList *il = &moltype.ilist[F_SETTLE];
410 for (i = 0; i < il->size(); i += 1+NRAL(F_SETTLE))
412 ip = &mtop->ffparams.iparams[il->iatoms[i]];
413 a1 = il->iatoms[i+1];
414 a2 = il->iatoms[i+2];
415 a3 = il->iatoms[i+3];
416 /* Usually the mass of a1 (usually oxygen) is larger than a2/a3.
417 * If this is not the case, we overestimate the displacement,
418 * which leads to a larger buffer (ok since this is an exotic case).
420 prop[a1].con_mass = getMass(*atoms, a2, setMassesToOne);
421 prop[a1].con_len = ip->settle.doh;
423 prop[a2].con_mass = getMass(*atoms, a1, setMassesToOne);
424 prop[a2].con_len = ip->settle.doh;
426 prop[a3].con_mass = getMass(*atoms, a1, setMassesToOne);
427 prop[a3].con_len = ip->settle.doh;
430 get_vsite_masses(&moltype,
431 &mtop->ffparams,
432 setMassesToOne,
433 vsite_m,
434 &n_nonlin_vsite_mol);
435 if (n_nonlin_vsite != nullptr)
437 *n_nonlin_vsite += nmol*n_nonlin_vsite_mol;
440 for (a = 0; a < atoms->nr; a++)
442 if (atoms->atom[a].ptype == eptVSite)
444 prop[a].mass = vsite_m[a];
446 else
448 prop[a].mass = getMass(*atoms, a, setMassesToOne);
450 prop[a].type = atoms->atom[a].type;
451 prop[a].q = atoms->atom[a].q;
452 /* We consider an atom constrained, #DOF=2, when it is
453 * connected with constraints to (at least one) atom with
454 * a mass of more than 0.4x its own mass. This is not a critical
455 * parameter, since with roughly equal masses the unconstrained
456 * and constrained displacement will not differ much (and both
457 * overestimate the displacement).
459 prop[a].bConstr = (prop[a].con_mass > 0.4*prop[a].mass);
461 add_at(&att, &natt, &prop[a], nmol);
464 sfree(vsite_m);
465 sfree(prop);
468 if (gmx_debug_at)
470 for (a = 0; a < natt; a++)
472 fprintf(debug, "type %d: m %5.2f t %d q %6.3f con %s con_m %5.3f con_l %5.3f n %d\n",
473 a, att[a].prop.mass, att[a].prop.type, att[a].prop.q,
474 gmx::boolToString(att[a].prop.bConstr), att[a].prop.con_mass, att[a].prop.con_len,
475 att[a].n);
479 *att_p = att;
480 *natt_p = natt;
483 /* This function computes two components of the estimate of the variance
484 * in the displacement of one atom in a system of two constrained atoms.
485 * Returns in sigma2_2d the variance due to rotation of the constrained
486 * atom around the atom to which it constrained.
487 * Returns in sigma2_3d the variance due to displacement of the COM
488 * of the whole system of the two constrained atoms.
490 * Note that we only take a single constraint (the one to the heaviest atom)
491 * into account. If an atom has multiple constraints, this will result in
492 * an overestimate of the displacement, which gives a larger drift and buffer.
494 void constrained_atom_sigma2(real kT_fac,
495 const atom_nonbonded_kinetic_prop_t *prop,
496 real *sigma2_2d,
497 real *sigma2_3d)
499 /* Here we decompose the motion of a constrained atom into two
500 * components: rotation around the COM and translation of the COM.
503 /* Determine the variance of the arc length for the two rotational DOFs */
504 real massFraction = prop->con_mass/(prop->mass + prop->con_mass);
505 real sigma2_rot = kT_fac*massFraction/prop->mass;
507 /* The distance from the atom to the COM, i.e. the rotational arm */
508 real comDistance = prop->con_len*massFraction;
510 /* The variance relative to the arm */
511 real sigma2_rel = sigma2_rot/gmx::square(comDistance);
513 /* For sigma2_rel << 1 we don't notice the rotational effect and
514 * we have a normal, Gaussian displacement distribution.
515 * For larger sigma2_rel the displacement is much less, in fact it can
516 * not exceed 2*comDistance. We can calculate MSD/arm^2 as:
517 * integral_x=0-inf distance2(x) x/sigma2_rel exp(-x^2/(2 sigma2_rel)) dx
518 * where x is angular displacement and distance2(x) is the distance^2
519 * between points at angle 0 and x:
520 * distance2(x) = (sin(x) - sin(0))^2 + (cos(x) - cos(0))^2
521 * The limiting value of this MSD is 2, which is also the value for
522 * a uniform rotation distribution that would be reached at long time.
523 * The maximum is 2.5695 at sigma2_rel = 4.5119.
524 * We approximate this integral with a rational polynomial with
525 * coefficients from a Taylor expansion. This approximation is an
526 * overestimate for all values of sigma2_rel. Its maximum value
527 * of 2.6491 is reached at sigma2_rel = sqrt(45/2) = 4.7434.
528 * We keep the approximation constant after that.
529 * We use this approximate MSD as the variance for a Gaussian distribution.
531 * NOTE: For any sensible buffer tolerance this will result in a (large)
532 * overestimate of the buffer size, since the Gaussian has a long tail,
533 * whereas the actual distribution can not reach values larger than 2.
535 /* Coeffients obtained from a Taylor expansion */
536 const real a = 1.0/3.0;
537 const real b = 2.0/45.0;
539 /* Our approximation is constant after sigma2_rel = 1/sqrt(b) */
540 sigma2_rel = std::min(sigma2_rel, 1/std::sqrt(b));
542 /* Compute the approximate sigma^2 for 2D motion due to the rotation */
543 *sigma2_2d = gmx::square(comDistance)*
544 sigma2_rel/(1 + a*sigma2_rel + b*gmx::square(sigma2_rel));
546 /* The constrained atom also moves (in 3D) with the COM of both atoms */
547 *sigma2_3d = kT_fac/(prop->mass + prop->con_mass);
550 static void get_atom_sigma2(real kT_fac,
551 const atom_nonbonded_kinetic_prop_t *prop,
552 real *sigma2_2d,
553 real *sigma2_3d)
555 if (prop->bConstr)
557 /* Complicated constraint calculation in a separate function */
558 constrained_atom_sigma2(kT_fac, prop, sigma2_2d, sigma2_3d);
560 else
562 /* Unconstrained atom: trivial */
563 *sigma2_2d = 0;
564 *sigma2_3d = kT_fac/prop->mass;
568 static void approx_2dof(real s2, real x, real *shift, real *scale)
570 /* A particle with 1 DOF constrained has 2 DOFs instead of 3.
571 * This code is also used for particles with multiple constraints,
572 * in which case we overestimate the displacement.
573 * The 2DOF distribution is sqrt(pi/2)*erfc(r/(sqrt(2)*s))/(2*s).
574 * We approximate this with scale*Gaussian(s,r+shift),
575 * by matching the distribution value and derivative at x.
576 * This is a tight overestimate for all r>=0 at any s and x.
578 real ex, er;
580 ex = std::exp(-x*x/(2*s2));
581 er = std::erfc(x/std::sqrt(2*s2));
583 *shift = -x + std::sqrt(2*s2/M_PI)*ex/er;
584 *scale = 0.5*M_PI*std::exp(ex*ex/(M_PI*er*er))*er;
587 // Returns an (over)estimate of the energy drift for a single atom pair,
588 // given the kinetic properties, displacement variances and list buffer.
589 static real energyDriftAtomPair(bool isConstrained_i,
590 bool isConstrained_j,
591 real s2, real s2i_2d, real s2j_2d,
592 real r_buffer,
593 const pot_derivatives_t *der)
595 // For relatively small arguments erfc() is so small that if will be 0.0
596 // when stored in a float. We set an argument limit of 8 (Erfc(8)=1e-29),
597 // such that we can divide by erfc and have some space left for arithmetic.
598 const real erfc_arg_max = 8.0;
600 real rsh = r_buffer;
601 real sc_fac = 1.0;
603 real c_exp, c_erfc;
605 if (rsh*rsh > 2*s2*erfc_arg_max*erfc_arg_max)
607 // Below we calculate c_erfc = 0.5*erfc(rsh/sqrt(2*s2))
608 // When rsh/sqrt(2*s2) increases, this erfc will be the first
609 // result that underflows and becomes 0.0. To avoid this,
610 // we set c_exp=0 and c_erfc=0 for large arguments.
611 // This also avoids NaN in approx_2dof().
612 // In any relevant case this has no effect on the results,
613 // since c_exp < 6e-29, so the displacement is completely
614 // negligible for such atom pairs (and an overestimate).
615 // In nearly all use cases, there will be other atom pairs
616 // that contribute much more to the total, so zeroing
617 // this particular contribution has no effect at all.
618 c_exp = 0;
619 c_erfc = 0;
621 else
623 /* For constraints: adapt r and scaling for the Gaussian */
624 if (isConstrained_i)
626 real sh, sc;
628 approx_2dof(s2i_2d, r_buffer*s2i_2d/s2, &sh, &sc);
629 rsh += sh;
630 sc_fac *= sc;
632 if (isConstrained_j)
634 real sh, sc;
636 approx_2dof(s2j_2d, r_buffer*s2j_2d/s2, &sh, &sc);
637 rsh += sh;
638 sc_fac *= sc;
641 /* Exact contribution of an atom pair with Gaussian displacement
642 * with sigma s to the energy drift for a potential with
643 * derivative -md and second derivative dd at the cut-off.
644 * The only catch is that for potentials that change sign
645 * near the cut-off there could be an unlucky compensation
646 * of positive and negative energy drift.
647 * Such potentials are extremely rare though.
649 * Note that pot has unit energy*length, as the linear
650 * atom density still needs to be put in.
652 c_exp = std::exp(-rsh*rsh/(2*s2))/std::sqrt(2*M_PI);
653 c_erfc = 0.5*std::erfc(rsh/(std::sqrt(2*s2)));
655 real s = std::sqrt(s2);
656 real rsh2 = rsh*rsh;
658 real pot1 = sc_fac*
659 der->md1/2*((rsh2 + s2)*c_erfc - rsh*s*c_exp);
660 real pot2 = sc_fac*
661 der->d2/6*(s*(rsh2 + 2*s2)*c_exp - rsh*(rsh2 + 3*s2)*c_erfc);
662 real pot3 = sc_fac*
663 der->md3/24*((rsh2*rsh2 + 6*rsh2*s2 + 3*s2*s2)*c_erfc - rsh*s*(rsh2 + 5*s2)*c_exp);
665 return pot1 + pot2 + pot3;
668 static real energyDrift(const verletbuf_atomtype_t *att, int natt,
669 const gmx_ffparams_t *ffp,
670 real kT_fac,
671 const pot_derivatives_t *ljDisp,
672 const pot_derivatives_t *ljRep,
673 const pot_derivatives_t *elec,
674 real rlj, real rcoulomb,
675 real rlist, real boxvol)
677 double drift_tot = 0;
679 if (kT_fac == 0)
681 /* No atom displacements: no drift, avoid division by 0 */
682 return drift_tot;
685 // Here add up the contribution of all atom pairs in the system to
686 // (estimated) energy drift by looping over all atom type pairs.
687 for (int i = 0; i < natt; i++)
689 // Get the thermal displacement variance for the i-atom type
690 const atom_nonbonded_kinetic_prop_t *prop_i = &att[i].prop;
691 real s2i_2d, s2i_3d;
692 get_atom_sigma2(kT_fac, prop_i, &s2i_2d, &s2i_3d);
694 for (int j = i; j < natt; j++)
696 // Get the thermal displacement variance for the j-atom type
697 const atom_nonbonded_kinetic_prop_t *prop_j = &att[j].prop;
698 real s2j_2d, s2j_3d;
699 get_atom_sigma2(kT_fac, prop_j, &s2j_2d, &s2j_3d);
701 /* Add up the up to four independent variances */
702 real s2 = s2i_2d + s2i_3d + s2j_2d + s2j_3d;
704 // Set -V', V'' and -V''' at the cut-off for LJ */
705 real c6 = ffp->iparams[prop_i->type*ffp->atnr + prop_j->type].lj.c6;
706 real c12 = ffp->iparams[prop_i->type*ffp->atnr + prop_j->type].lj.c12;
707 pot_derivatives_t lj;
708 lj.md1 = c6*ljDisp->md1 + c12*ljRep->md1;
709 lj.d2 = c6*ljDisp->d2 + c12*ljRep->d2;
710 lj.md3 = c6*ljDisp->md3 + c12*ljRep->md3;
712 real pot_lj = energyDriftAtomPair(prop_i->bConstr, prop_j->bConstr,
713 s2, s2i_2d, s2j_2d,
714 rlist - rlj,
715 &lj);
717 // Set -V' and V'' at the cut-off for Coulomb
718 pot_derivatives_t elec_qq;
719 elec_qq.md1 = elec->md1*prop_i->q*prop_j->q;
720 elec_qq.d2 = elec->d2 *prop_i->q*prop_j->q;
721 elec_qq.md3 = 0;
723 real pot_q = energyDriftAtomPair(prop_i->bConstr, prop_j->bConstr,
724 s2, s2i_2d, s2j_2d,
725 rlist - rcoulomb,
726 &elec_qq);
728 // Note that attractive and repulsive potentials for individual
729 // pairs can partially cancel.
730 real pot = pot_lj + pot_q;
732 /* Multiply by the number of atom pairs */
733 if (j == i)
735 pot *= static_cast<double>(att[i].n)*(att[i].n - 1)/2;
737 else
739 pot *= static_cast<double>(att[i].n)*att[j].n;
741 /* We need the line density to get the energy drift of the system.
742 * The effective average r^2 is close to (rlist+sigma)^2.
744 pot *= 4*M_PI*gmx::square(rlist + std::sqrt(s2))/boxvol;
746 /* Add the unsigned drift to avoid cancellation of errors */
747 drift_tot += std::abs(pot);
751 return drift_tot;
754 static real surface_frac(int cluster_size, real particle_distance, real rlist)
756 real d, area_rel;
758 if (rlist < 0.5*particle_distance)
760 /* We have non overlapping spheres */
761 return 1.0;
764 /* Half the inter-particle distance relative to rlist */
765 d = 0.5*particle_distance/rlist;
767 /* Determine the area of the surface at distance rlist to the closest
768 * particle, relative to surface of a sphere of radius rlist.
769 * The formulas below assume close to cubic cells for the pair search grid,
770 * which the pair search code tries to achieve.
771 * Note that in practice particle distances will not be delta distributed,
772 * but have some spread, often involving shorter distances,
773 * as e.g. O-H bonds in a water molecule. Thus the estimates below will
774 * usually be slightly too high and thus conservative.
776 switch (cluster_size)
778 case 1:
779 /* One particle: trivial */
780 area_rel = 1.0;
781 break;
782 case 2:
783 /* Two particles: two spheres at fractional distance 2*a */
784 area_rel = 1.0 + d;
785 break;
786 case 4:
787 /* We assume a perfect, symmetric tetrahedron geometry.
788 * The surface around a tetrahedron is too complex for a full
789 * analytical solution, so we use a Taylor expansion.
791 area_rel = (1.0 + 1/M_PI*(6*std::acos(1/std::sqrt(3))*d +
792 std::sqrt(3)*d*d*(1.0 +
793 5.0/18.0*d*d +
794 7.0/45.0*d*d*d*d +
795 83.0/756.0*d*d*d*d*d*d)));
796 break;
797 default:
798 gmx_incons("surface_frac called with unsupported cluster_size");
801 return area_rel/cluster_size;
804 /* Returns the negative of the third derivative of a potential r^-p
805 * with a force-switch function, evaluated at the cut-off rc.
807 static real md3_force_switch(real p, real rswitch, real rc)
809 /* The switched force function is:
810 * p*r^-(p+1) + a*(r - rswitch)^2 + b*(r - rswitch)^3
812 real a, b;
813 real md3_pot, md3_sw;
815 a = -((p + 4)*rc - (p + 1)*rswitch)/(pow(rc, p+2)*gmx::square(rc-rswitch));
816 b = ((p + 3)*rc - (p + 1)*rswitch)/(pow(rc, p+2)*gmx::power3(rc-rswitch));
818 md3_pot = (p + 2)*(p + 1)*p*pow(rc, p+3);
819 md3_sw = 2*a + 6*b*(rc - rswitch);
821 return md3_pot + md3_sw;
824 /* Returns the maximum reference temperature over all coupled groups */
825 static real maxReferenceTemperature(const t_inputrec &ir)
827 if (EI_MD(ir.eI) && ir.etc == etcNO)
829 /* This case should be handled outside calc_verlet_buffer_size */
830 gmx_incons("calc_verlet_buffer_size called with an NVE ensemble and reference_temperature < 0");
833 real maxTemperature = 0;
834 for (int i = 0; i < ir.opts.ngtc; i++)
836 if (ir.opts.tau_t[i] >= 0)
838 maxTemperature = std::max(maxTemperature, ir.opts.ref_t[i]);
842 return maxTemperature;
845 /* Returns the variance of the atomic displacement over timePeriod.
847 * Note: When not using BD with a non-mass dependendent friction coefficient,
848 * the return value still needs to be divided by the particle mass.
850 static real displacementVariance(const t_inputrec &ir,
851 real temperature,
852 real timePeriod)
854 real kT_fac;
856 if (ir.eI == eiBD)
858 /* Get the displacement distribution from the random component only.
859 * With accurate integration the systematic (force) displacement
860 * should be negligible (unless nstlist is extremely large, which
861 * you wouldn't do anyhow).
863 kT_fac = 2*BOLTZ*temperature*timePeriod;
864 if (ir.bd_fric > 0)
866 /* This is directly sigma^2 of the displacement */
867 kT_fac /= ir.bd_fric;
869 else
871 /* Per group tau_t is not implemented yet, use the maximum */
872 real tau_t = ir.opts.tau_t[0];
873 for (int i = 1; i < ir.opts.ngtc; i++)
875 tau_t = std::max(tau_t, ir.opts.tau_t[i]);
878 kT_fac *= tau_t;
879 /* This kT_fac needs to be divided by the mass to get sigma^2 */
882 else
884 kT_fac = BOLTZ*temperature*gmx::square(timePeriod);
887 return kT_fac;
890 /* Returns the largest sigma of the Gaussian displacement over all particle
891 * types. This ignores constraints, so is an overestimate.
893 static real maxSigma(real kT_fac,
894 int natt,
895 const verletbuf_atomtype_t *att)
897 assert(att);
898 real smallestMass = att[0].prop.mass;
899 for (int i = 1; i < natt; i++)
901 smallestMass = std::min(smallestMass, att[i].prop.mass);
904 return 2*std::sqrt(kT_fac/smallestMass);
907 void calc_verlet_buffer_size(const gmx_mtop_t *mtop, real boxvol,
908 const t_inputrec *ir,
909 int nstlist,
910 int list_lifetime,
911 real reference_temperature,
912 const VerletbufListSetup *list_setup,
913 int *n_nonlin_vsite,
914 real *rlist)
916 double resolution;
917 char *env;
919 real particle_distance;
920 real nb_clust_frac_pairs_not_in_list_at_cutoff;
922 verletbuf_atomtype_t *att = nullptr;
923 int natt = -1;
924 real elfac;
925 int ib0, ib1, ib;
926 real rb, rl;
927 real drift;
929 if (!EI_DYNAMICS(ir->eI))
931 gmx_incons("Can only determine the Verlet buffer size for integrators that perform dynamics");
933 if (ir->verletbuf_tol <= 0)
935 gmx_incons("The Verlet buffer tolerance needs to be larger than zero");
938 if (reference_temperature < 0)
940 /* We use the maximum temperature with multiple T-coupl groups.
941 * We could use a per particle temperature, but since particles
942 * interact, this might underestimate the buffer size.
944 reference_temperature = maxReferenceTemperature(*ir);
947 /* Resolution of the buffer size */
948 resolution = 0.001;
950 env = getenv("GMX_VERLET_BUFFER_RES");
951 if (env != nullptr)
953 sscanf(env, "%lf", &resolution);
956 /* In an atom wise pair-list there would be no pairs in the list
957 * beyond the pair-list cut-off.
958 * However, we use a pair-list of groups vs groups of atoms.
959 * For groups of 4 atoms, the parallelism of SSE instructions, only
960 * 10% of the atoms pairs are not in the list just beyond the cut-off.
961 * As this percentage increases slowly compared to the decrease of the
962 * Gaussian displacement distribution over this range, we can simply
963 * reduce the drift by this fraction.
964 * For larger groups, e.g. of 8 atoms, this fraction will be lower,
965 * so then buffer size will be on the conservative (large) side.
967 * Note that the formulas used here do not take into account
968 * cancellation of errors which could occur by missing both
969 * attractive and repulsive interactions.
971 * The only major assumption is homogeneous particle distribution.
972 * For an inhomogeneous system, such as a liquid-vapor system,
973 * the buffer will be underestimated. The actual energy drift
974 * will be higher by the factor: local/homogeneous particle density.
976 * The results of this estimate have been checked againt simulations.
977 * In most cases the real drift differs by less than a factor 2.
980 /* Worst case assumption: HCP packing of particles gives largest distance */
981 particle_distance = std::cbrt(boxvol*std::sqrt(2)/mtop->natoms);
983 /* TODO: Obtain masses through (future) integrator functionality
984 * to avoid scattering the code with (or forgetting) checks.
986 const bool setMassesToOne = (ir->eI == eiBD && ir->bd_fric > 0);
987 get_verlet_buffer_atomtypes(mtop, setMassesToOne, &att, &natt, n_nonlin_vsite);
988 assert(att != nullptr && natt >= 0);
990 if (debug)
992 fprintf(debug, "particle distance assuming HCP packing: %f nm\n",
993 particle_distance);
994 fprintf(debug, "energy drift atom types: %d\n", natt);
997 pot_derivatives_t ljDisp = { 0, 0, 0 };
998 pot_derivatives_t ljRep = { 0, 0, 0 };
999 real repPow = mtop->ffparams.reppow;
1001 if (ir->vdwtype == evdwCUT)
1003 real sw_range, md3_pswf;
1005 switch (ir->vdw_modifier)
1007 case eintmodNONE:
1008 case eintmodPOTSHIFT:
1009 /* -dV/dr of -r^-6 and r^-reppow */
1010 ljDisp.md1 = -6*std::pow(ir->rvdw, -7.0);
1011 ljRep.md1 = repPow*std::pow(ir->rvdw, -(repPow + 1));
1012 /* The contribution of the higher derivatives is negligible */
1013 break;
1014 case eintmodFORCESWITCH:
1015 /* At the cut-off: V=V'=V''=0, so we use only V''' */
1016 ljDisp.md3 = -md3_force_switch(6.0, ir->rvdw_switch, ir->rvdw);
1017 ljRep.md3 = md3_force_switch(repPow, ir->rvdw_switch, ir->rvdw);
1018 break;
1019 case eintmodPOTSWITCH:
1020 /* At the cut-off: V=V'=V''=0.
1021 * V''' is given by the original potential times
1022 * the third derivative of the switch function.
1024 sw_range = ir->rvdw - ir->rvdw_switch;
1025 md3_pswf = 60.0/gmx::power3(sw_range);
1027 ljDisp.md3 = -std::pow(ir->rvdw, -6.0 )*md3_pswf;
1028 ljRep.md3 = std::pow(ir->rvdw, -repPow)*md3_pswf;
1029 break;
1030 default:
1031 gmx_incons("Unimplemented VdW modifier");
1034 else if (EVDW_PME(ir->vdwtype))
1036 real b = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
1037 real r = ir->rvdw;
1038 real br = b*r;
1039 real br2 = br*br;
1040 real br4 = br2*br2;
1041 real br6 = br4*br2;
1042 // -dV/dr of g(br)*r^-6 [where g(x) = exp(-x^2)(1+x^2+x^4/2),
1043 // see LJ-PME equations in manual] and r^-reppow
1044 ljDisp.md1 = -std::exp(-br2)*(br6 + 3.0*br4 + 6.0*br2 + 6.0)*std::pow(r, -7.0);
1045 ljRep.md1 = repPow*pow(r, -(repPow + 1));
1046 // The contribution of the higher derivatives is negligible
1048 else
1050 gmx_fatal(FARGS, "Energy drift calculation is only implemented for plain cut-off Lennard-Jones interactions");
1053 elfac = ONE_4PI_EPS0/ir->epsilon_r;
1055 // Determine the 1st and 2nd derivative for the electostatics
1056 pot_derivatives_t elec = { 0, 0, 0 };
1058 if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
1060 real eps_rf, k_rf;
1062 if (ir->coulombtype == eelCUT)
1064 eps_rf = 1;
1065 k_rf = 0;
1067 else
1069 eps_rf = ir->epsilon_rf/ir->epsilon_r;
1070 if (eps_rf != 0)
1072 k_rf = (eps_rf - ir->epsilon_r)/( gmx::power3(ir->rcoulomb) * (2*eps_rf + ir->epsilon_r) );
1074 else
1076 /* epsilon_rf = infinity */
1077 k_rf = 0.5/gmx::power3(ir->rcoulomb);
1081 if (eps_rf > 0)
1083 elec.md1 = elfac*(1.0/gmx::square(ir->rcoulomb) - 2*k_rf*ir->rcoulomb);
1085 elec.d2 = elfac*(2.0/gmx::power3(ir->rcoulomb) + 2*k_rf);
1087 else if (EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD)
1089 real b, rc, br;
1091 b = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
1092 rc = ir->rcoulomb;
1093 br = b*rc;
1094 elec.md1 = elfac*(b*std::exp(-br*br)*M_2_SQRTPI/rc + std::erfc(br)/(rc*rc));
1095 elec.d2 = elfac/(rc*rc)*(2*b*(1 + br*br)*std::exp(-br*br)*M_2_SQRTPI + 2*std::erfc(br)/rc);
1097 else
1099 gmx_fatal(FARGS, "Energy drift calculation is only implemented for Reaction-Field and Ewald electrostatics");
1102 /* Determine the variance of the atomic displacement
1103 * over list_lifetime steps: kT_fac
1104 * For inertial dynamics (not Brownian dynamics) the mass factor
1105 * is not included in kT_fac, it is added later.
1107 const real kT_fac = displacementVariance(*ir, reference_temperature,
1108 list_lifetime*ir->delta_t);
1110 if (debug)
1112 fprintf(debug, "Derivatives of non-bonded potentials at the cut-off:\n");
1113 fprintf(debug, "LJ disp. -V' %9.2e V'' %9.2e -V''' %9.2e\n", ljDisp.md1, ljDisp.d2, ljDisp.md3);
1114 fprintf(debug, "LJ rep. -V' %9.2e V'' %9.2e -V''' %9.2e\n", ljRep.md1, ljRep.d2, ljRep.md3);
1115 fprintf(debug, "Electro. -V' %9.2e V'' %9.2e\n", elec.md1, elec.d2);
1116 fprintf(debug, "sqrt(kT_fac) %f\n", std::sqrt(kT_fac));
1119 /* Search using bisection */
1120 ib0 = -1;
1121 /* The drift will be neglible at 5 times the max sigma */
1122 ib1 = static_cast<int>(5*maxSigma(kT_fac, natt, att)/resolution) + 1;
1123 while (ib1 - ib0 > 1)
1125 ib = (ib0 + ib1)/2;
1126 rb = ib*resolution;
1127 rl = std::max(ir->rvdw, ir->rcoulomb) + rb;
1129 /* Calculate the average energy drift at the last step
1130 * of the nstlist steps at which the pair-list is used.
1132 drift = energyDrift(att, natt, &mtop->ffparams,
1133 kT_fac,
1134 &ljDisp, &ljRep, &elec,
1135 ir->rvdw, ir->rcoulomb,
1136 rl, boxvol);
1138 /* Correct for the fact that we are using a Ni x Nj particle pair list
1139 * and not a 1 x 1 particle pair list. This reduces the drift.
1141 /* We don't have a formula for 8 (yet), use 4 which is conservative */
1142 nb_clust_frac_pairs_not_in_list_at_cutoff =
1143 surface_frac(std::min(list_setup->cluster_size_i, 4),
1144 particle_distance, rl)*
1145 surface_frac(std::min(list_setup->cluster_size_j, 4),
1146 particle_distance, rl);
1147 drift *= nb_clust_frac_pairs_not_in_list_at_cutoff;
1149 /* Convert the drift to drift per unit time per atom */
1150 drift /= nstlist*ir->delta_t*mtop->natoms;
1152 if (debug)
1154 fprintf(debug, "ib %3d %3d %3d rb %.3f %dx%d fac %.3f drift %.1e\n",
1155 ib0, ib, ib1, rb,
1156 list_setup->cluster_size_i, list_setup->cluster_size_j,
1157 nb_clust_frac_pairs_not_in_list_at_cutoff,
1158 drift);
1161 if (std::abs(drift) > ir->verletbuf_tol)
1163 ib0 = ib;
1165 else
1167 ib1 = ib;
1171 sfree(att);
1173 *rlist = std::max(ir->rvdw, ir->rcoulomb) + ib1*resolution;
1176 /* Returns the pairlist buffer size for use as a minimum buffer size
1178 * Note that this is a rather crude estimate. It is ok for a buffer
1179 * set for good energy conservation or RF electrostatics. But it is
1180 * too small with PME and the buffer set with the default tolerance.
1182 static real minCellSizeFromPairlistBuffer(const t_inputrec &ir)
1184 return ir.rlist - std::max(ir.rvdw, ir.rcoulomb);
1187 real minCellSizeForAtomDisplacement(const gmx_mtop_t &mtop,
1188 const t_inputrec &ir,
1189 real chanceRequested)
1191 if (!EI_DYNAMICS(ir.eI) || (EI_MD(ir.eI) && ir.etc == etcNO))
1193 return minCellSizeFromPairlistBuffer(ir);
1196 /* We use the maximum temperature with multiple T-coupl groups.
1197 * We could use a per particle temperature, but since particles
1198 * interact, this might underestimate the displacements.
1200 const real temperature = maxReferenceTemperature(ir);
1202 const bool setMassesToOne = (ir.eI == eiBD && ir.bd_fric > 0);
1204 verletbuf_atomtype_t *att = nullptr;
1205 int natt = -1;
1206 get_verlet_buffer_atomtypes(&mtop, setMassesToOne, &att, &natt, nullptr);
1208 const real kT_fac = displacementVariance(ir, temperature,
1209 ir.nstlist*ir.delta_t);
1211 /* Resolution of the cell size */
1212 real resolution = 0.001;
1214 /* Search using bisection, avoid 0 and start at 1 */
1215 int ib0 = 0;
1216 /* The chance will be neglible at 10 times the max sigma */
1217 int ib1 = int(10*maxSigma(kT_fac, natt, att)/resolution) + 1;
1218 real cellSize = 0;
1219 while (ib1 - ib0 > 1)
1221 int ib = (ib0 + ib1)/2;
1222 cellSize = ib*resolution;
1224 /* We assumes atom are distributed uniformly over the cell width.
1225 * Once an atom has moved by more than the cellSize (as passed
1226 * as the buffer argument to energyDriftAtomPair() below),
1227 * the chance of crossing the boundary of the neighbor cell
1228 * thus increases as 1/cellSize with the additional displacement
1229 * on to of cellSize. We thus create a linear interaction with
1230 * derivative = -1/cellSize. Using this in the energyDriftAtomPair
1231 * function will return the chance of crossing the next boundary.
1233 const pot_derivatives_t boundaryInteraction = { 1/cellSize, 0, 0 };
1235 real chance = 0;
1236 for (int i = 0; i < natt; i++)
1238 const atom_nonbonded_kinetic_prop_t &propAtom = att[i].prop;
1239 real s2_2d;
1240 real s2_3d;
1241 get_atom_sigma2(kT_fac, &propAtom, &s2_2d, &s2_3d);
1243 real chancePerAtom = energyDriftAtomPair(propAtom.bConstr, false,
1244 s2_2d + s2_3d, s2_2d, 0,
1245 cellSize,
1246 &boundaryInteraction);
1248 if (propAtom.bConstr)
1250 /* energyDriftAtomPair() uses an unlimited Gaussian displacement
1251 * distribution for constrained atoms, whereas they can
1252 * actually not move more than the COM of the two constrained
1253 * atoms plus twice the distance from the COM.
1254 * Use this maximum, limited displacement when this results in
1255 * a smaller chance (note that this is still an overestimate).
1257 real massFraction = propAtom.con_mass/(propAtom.mass + propAtom.con_mass);
1258 real comDistance = propAtom.con_len*massFraction;
1260 real chanceWithMaxDistance =
1261 energyDriftAtomPair(false, false,
1262 s2_3d, 0, 0,
1263 cellSize - 2*comDistance,
1264 &boundaryInteraction);
1265 chancePerAtom = std::min(chancePerAtom, chanceWithMaxDistance);
1268 /* Take into account the line density of the boundary */
1269 chancePerAtom /= cellSize;
1271 chance += att[i].n*chancePerAtom;
1274 /* Note: chance is for every nstlist steps */
1275 if (chance > chanceRequested*ir.nstlist)
1277 ib0 = ib;
1279 else
1281 ib1 = ib;
1285 sfree(att);
1287 return cellSize;