2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5 * Copyright (c) 2001-2004, The GROMACS development team,
6 * check out http://www.gromacs.org for more information.
7 * Copyright (c) 2012,2013, by the GROMACS development team, led by
8 * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9 * others, as listed in the AUTHORS file in the top-level source
10 * directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
44 #include <sys/types.h>
49 #include "gmx_fatal.h"
53 #include "calc_verletbuf.h"
54 #include "../mdlib/nbnxn_consts.h"
57 /* The include below sets the SIMD instruction type (precision+width)
58 * for all nbnxn SIMD search and non-bonded kernel code.
60 #ifdef GMX_NBNXN_HALF_WIDTH_SIMD
61 #define GMX_USE_HALF_WIDTH_SIMD_HERE
63 #include "gmx_simd_macros.h"
67 /* The code in this file estimates a pairlist buffer length
68 * given a target energy drift per atom per picosecond.
69 * This is done by estimating the drift given a buffer length.
70 * Ideally we would like to have a tight overestimate of the drift,
71 * but that can be difficult to achieve.
73 * Significant approximations used:
75 * Uniform particle density. UNDERESTIMATES the drift by rho_global/rho_local.
77 * Interactions don't affect particle motion. OVERESTIMATES the drift on longer
78 * time scales. This approximation probably introduces the largest errors.
80 * Only take one constraint per particle into account: OVERESTIMATES the drift.
82 * For rotating constraints assume the same functional shape for time scales
83 * where the constraints rotate significantly as the exact expression for
84 * short time scales. OVERESTIMATES the drift on long time scales.
86 * For non-linear virtual sites use the mass of the lightest constructing atom
87 * to determine the displacement. OVER/UNDERESTIMATES the drift, depending on
88 * the geometry and masses of constructing atoms.
90 * Note that the formulas for normal atoms and linear virtual sites are exact,
91 * apart from the first two approximations.
93 * Note that apart from the effect of the above approximations, the actual
94 * drift of the total energy of a system can be order of magnitude smaller
95 * due to cancellation of positive and negative drift for different pairs.
99 /* Struct for unique atom type for calculating the energy drift.
100 * The atom displacement depends on mass and constraints.
101 * The energy jump for given distance depend on LJ type and q.
105 real mass; /* mass */
106 int type; /* type (used for LJ parameters) */
108 gmx_bool bConstr; /* constrained, if TRUE, use #DOF=2 iso 3 */
109 real con_mass; /* mass of heaviest atom connected by constraints */
110 real con_len; /* constraint length to the heaviest atom */
111 } atom_nonbonded_kinetic_prop_t;
113 /* Struct for unique atom type for calculating the energy drift.
114 * The atom displacement depends on mass and constraints.
115 * The energy jump for given distance depend on LJ type and q.
119 atom_nonbonded_kinetic_prop_t prop; /* non-bonded and kinetic atom prop. */
120 int n; /* #atoms of this type in the system */
121 } verletbuf_atomtype_t;
123 void verletbuf_get_list_setup(gmx_bool bGPU,
124 verletbuf_list_setup_t *list_setup)
126 list_setup->cluster_size_i = NBNXN_CPU_CLUSTER_I_SIZE;
130 list_setup->cluster_size_j = NBNXN_GPU_CLUSTER_SIZE;
134 #ifndef GMX_NBNXN_SIMD
135 list_setup->cluster_size_j = NBNXN_CPU_CLUSTER_I_SIZE;
137 list_setup->cluster_size_j = GMX_SIMD_WIDTH_HERE;
138 #ifdef GMX_NBNXN_SIMD_2XNN
139 /* We assume the smallest cluster size to be on the safe side */
140 list_setup->cluster_size_j /= 2;
147 atom_nonbonded_kinetic_prop_equal(const atom_nonbonded_kinetic_prop_t *prop1,
148 const atom_nonbonded_kinetic_prop_t *prop2)
150 return (prop1->mass == prop2->mass &&
151 prop1->type == prop2->type &&
152 prop1->q == prop2->q &&
153 prop1->bConstr == prop2->bConstr &&
154 prop1->con_mass == prop2->con_mass &&
155 prop1->con_len == prop2->con_len);
158 static void add_at(verletbuf_atomtype_t **att_p, int *natt_p,
159 const atom_nonbonded_kinetic_prop_t *prop,
162 verletbuf_atomtype_t *att;
167 /* Ignore massless particles */
175 while (i < natt && !atom_nonbonded_kinetic_prop_equal(prop, &att[i].prop))
187 srenew(*att_p, *natt_p);
188 (*att_p)[i].prop = *prop;
189 (*att_p)[i].n = nmol;
193 static void get_vsite_masses(const gmx_moltype_t *moltype,
194 const gmx_ffparams_t *ffparams,
203 /* Check for virtual sites, determine mass from constructing atoms */
204 for (ft = 0; ft < F_NRE; ft++)
208 il = &moltype->ilist[ft];
210 for (i = 0; i < il->nr; i += 1+NRAL(ft))
213 real cam[5], inv_mass, m_aj;
214 int a1, j, aj, coeff;
216 ip = &ffparams->iparams[il->iatoms[i]];
218 a1 = il->iatoms[i+1];
222 for (j = 1; j < NRAL(ft); j++)
224 cam[j] = moltype->atoms.atom[il->iatoms[i+1+j]].m;
227 cam[j] = vsite_m[il->iatoms[i+1+j]];
231 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.",
233 interaction_function[ft].longname,
234 il->iatoms[i+1+j]+1);
243 vsite_m[a1] = (cam[1]*cam[2])/(cam[2]*sqr(1-ip->vsite.a) + cam[1]*sqr(ip->vsite.a));
247 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));
252 for (j = 0; j < 3*ip->vsiten.n; j += 3)
254 aj = il->iatoms[i+j+2];
255 coeff = ip[il->iatoms[i+j]].vsiten.a;
256 if (moltype->atoms.atom[aj].ptype == eptVSite)
262 m_aj = moltype->atoms.atom[aj].m;
266 gmx_incons("The mass of a vsiten constructing atom is <= 0");
268 inv_mass += coeff*coeff/m_aj;
270 vsite_m[a1] = 1/inv_mass;
273 /* Use the mass of the lightest constructing atom.
274 * This is an approximation.
275 * If the distance of the virtual site to the
276 * constructing atom is less than all distances
277 * between constructing atoms, this is a safe
278 * over-estimate of the displacement of the vsite.
279 * This condition holds for all H mass replacement
280 * vsite constructions, except for SP2/3 groups.
281 * In SP3 groups one H will have a F_VSITE3
282 * construction, so even there the total drift
283 * estimate shouldn't be far off.
286 vsite_m[a1] = cam[1];
287 for (j = 2; j < NRAL(ft); j++)
289 vsite_m[a1] = min(vsite_m[a1], cam[j]);
296 fprintf(debug, "atom %4d %-20s mass %6.3f\n",
297 a1, interaction_function[ft].longname, vsite_m[a1]);
304 static void get_verlet_buffer_atomtypes(const gmx_mtop_t *mtop,
305 verletbuf_atomtype_t **att_p,
309 verletbuf_atomtype_t *att;
311 int mb, nmol, ft, i, a1, a2, a3, a;
312 const t_atoms *atoms;
315 atom_nonbonded_kinetic_prop_t *prop;
317 int n_nonlin_vsite_mol;
322 if (n_nonlin_vsite != NULL)
327 for (mb = 0; mb < mtop->nmolblock; mb++)
329 nmol = mtop->molblock[mb].nmol;
331 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
333 /* Check for constraints, as they affect the kinetic energy.
334 * For virtual sites we need the masses and geometry of
335 * the constructing atoms to determine their velocity distribution.
337 snew(prop, atoms->nr);
338 snew(vsite_m, atoms->nr);
340 for (ft = F_CONSTR; ft <= F_CONSTRNC; ft++)
342 il = &mtop->moltype[mtop->molblock[mb].type].ilist[ft];
344 for (i = 0; i < il->nr; i += 1+NRAL(ft))
346 ip = &mtop->ffparams.iparams[il->iatoms[i]];
347 a1 = il->iatoms[i+1];
348 a2 = il->iatoms[i+2];
349 if (atoms->atom[a2].m > prop[a1].con_mass)
351 prop[a1].con_mass = atoms->atom[a2].m;
352 prop[a1].con_len = ip->constr.dA;
354 if (atoms->atom[a1].m > prop[a2].con_mass)
356 prop[a2].con_mass = atoms->atom[a1].m;
357 prop[a2].con_len = ip->constr.dA;
362 il = &mtop->moltype[mtop->molblock[mb].type].ilist[F_SETTLE];
364 for (i = 0; i < il->nr; i += 1+NRAL(F_SETTLE))
366 ip = &mtop->ffparams.iparams[il->iatoms[i]];
367 a1 = il->iatoms[i+1];
368 a2 = il->iatoms[i+2];
369 a3 = il->iatoms[i+3];
370 /* Usually the mass of a1 (usually oxygen) is larger than a2/a3.
371 * If this is not the case, we overestimate the displacement,
372 * which leads to a larger buffer (ok since this is an exotic case).
374 prop[a1].con_mass = atoms->atom[a2].m;
375 prop[a1].con_len = ip->settle.doh;
377 prop[a2].con_mass = atoms->atom[a1].m;
378 prop[a2].con_len = ip->settle.doh;
380 prop[a3].con_mass = atoms->atom[a1].m;
381 prop[a3].con_len = ip->settle.doh;
384 get_vsite_masses(&mtop->moltype[mtop->molblock[mb].type],
387 &n_nonlin_vsite_mol);
388 if (n_nonlin_vsite != NULL)
390 *n_nonlin_vsite += nmol*n_nonlin_vsite_mol;
393 for (a = 0; a < atoms->nr; a++)
395 if (atoms->atom[a].ptype == eptVSite)
397 prop[a].mass = vsite_m[a];
401 prop[a].mass = atoms->atom[a].m;
403 prop[a].type = atoms->atom[a].type;
404 prop[a].q = atoms->atom[a].q;
405 /* We consider an atom constrained, #DOF=2, when it is
406 * connected with constraints to (at least one) atom with
407 * a mass of more than 0.4x its own mass. This is not a critical
408 * parameter, since with roughly equal masses the unconstrained
409 * and constrained displacement will not differ much (and both
410 * overestimate the displacement).
412 prop[a].bConstr = (prop[a].con_mass > 0.4*prop[a].mass);
414 add_at(&att, &natt, &prop[a], nmol);
423 for (a = 0; a < natt; a++)
425 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",
426 a, att[a].prop.mass, att[a].prop.type, att[a].prop.q,
427 att[a].prop.bConstr, att[a].prop.con_mass, att[a].prop.con_len,
436 /* This function computes two components of the estimate of the variance
437 * in the displacement of one atom in a system of two constrained atoms.
438 * Returns in sigma2_2d the variance due to rotation of the constrained
439 * atom around the atom to which it constrained.
440 * Returns in sigma2_3d the variance due to displacement of the COM
441 * of the whole system of the two constrained atoms.
443 * Note that we only take a single constraint (the one to the heaviest atom)
444 * into account. If an atom has multiple constraints, this will result in
445 * an overestimate of the displacement, which gives a larger drift and buffer.
447 static void constrained_atom_sigma2(real kT_fac,
448 const atom_nonbonded_kinetic_prop_t *prop,
457 /* Here we decompose the motion of a constrained atom into two
458 * components: rotation around the COM and translation of the COM.
461 /* Determine the variance for the displacement of the rotational mode */
462 sigma2_rot = kT_fac/(prop->mass*(prop->mass + prop->con_mass)/prop->con_mass);
464 /* The distance from the atom to the COM, i.e. the rotational arm */
465 com_dist = prop->con_len*prop->con_mass/(prop->mass + prop->con_mass);
467 /* The variance relative to the arm */
468 sigma2_rel = sigma2_rot/(com_dist*com_dist);
469 /* At 6 the scaling formula has slope 0,
470 * so we keep sigma2_2d constant after that.
474 /* A constrained atom rotates around the atom it is constrained to.
475 * This results in a smaller linear displacement than for a free atom.
476 * For a perfectly circular displacement, this lowers the displacement
477 * by: 1/arcsin(arc_length)
478 * and arcsin(x) = 1 + x^2/6 + ...
479 * For sigma2_rel<<1 the displacement distribution is erfc
480 * (exact formula is provided below). For larger sigma, it is clear
481 * that the displacement can't be larger than 2*com_dist.
482 * It turns out that the distribution becomes nearly uniform.
483 * For intermediate sigma2_rel, scaling down sigma with the third
484 * order expansion of arcsin with argument sigma_rel turns out
485 * to give a very good approximation of the distribution and variance.
486 * Even for larger values, the variance is only slightly overestimated.
487 * Note that the most relevant displacements are in the long tail.
488 * This rotation approximation always overestimates the tail (which
489 * runs to infinity, whereas it should be <= 2*com_dist).
490 * Thus we always overestimate the drift and the buffer size.
492 scale = 1/(1 + sigma2_rel/6);
493 *sigma2_2d = sigma2_rot*scale*scale;
497 /* sigma_2d is set to the maximum given by the scaling above.
498 * For large sigma2 the real displacement distribution is close
499 * to uniform over -2*con_len to 2*com_dist.
500 * Our erfc with sigma_2d=sqrt(1.5)*com_dist (which means the sigma
501 * of the erfc output distribution is con_dist) overestimates
502 * the variance and additionally has a long tail. This means
503 * we have a (safe) overestimation of the drift.
505 *sigma2_2d = 1.5*com_dist*com_dist;
508 /* The constrained atom also moves (in 3D) with the COM of both atoms */
509 *sigma2_3d = kT_fac/(prop->mass + prop->con_mass);
512 static void get_atom_sigma2(real kT_fac,
513 const atom_nonbonded_kinetic_prop_t *prop,
519 /* Complicated constraint calculation in a separate function */
520 constrained_atom_sigma2(kT_fac, prop, sigma2_2d, sigma2_3d);
524 /* Unconstrained atom: trivial */
526 *sigma2_3d = kT_fac/prop->mass;
530 static void approx_2dof(real s2, real x, real *shift, real *scale)
532 /* A particle with 1 DOF constrained has 2 DOFs instead of 3.
533 * This code is also used for particles with multiple constraints,
534 * in which case we overestimate the displacement.
535 * The 2DOF distribution is sqrt(pi/2)*erfc(r/(sqrt(2)*s))/(2*s).
536 * We approximate this with scale*Gaussian(s,r+shift),
537 * by matching the distribution value and derivative at x.
538 * This is a tight overestimate for all r>=0 at any s and x.
542 ex = exp(-x*x/(2*s2));
543 er = gmx_erfc(x/sqrt(2*s2));
545 *shift = -x + sqrt(2*s2/M_PI)*ex/er;
546 *scale = 0.5*M_PI*exp(ex*ex/(M_PI*er*er))*er;
549 static real ener_drift(const verletbuf_atomtype_t *att, int natt,
550 const gmx_ffparams_t *ffp,
552 real md_ljd, real md_ljr, real md_el, real dd_el,
554 real rlist, real boxvol)
556 double drift_tot, pot1, pot2, pot;
558 real s2i_2d, s2i_3d, s2j_2d, s2j_3d, s2, s;
562 double c_exp, c_erfc;
566 /* Loop over the different atom type pairs */
567 for (i = 0; i < natt; i++)
569 get_atom_sigma2(kT_fac, &att[i].prop, &s2i_2d, &s2i_3d);
570 ti = att[i].prop.type;
572 for (j = i; j < natt; j++)
574 get_atom_sigma2(kT_fac, &att[j].prop, &s2j_2d, &s2j_3d);
575 tj = att[j].prop.type;
577 /* Add up the up to four independent variances */
578 s2 = s2i_2d + s2i_3d + s2j_2d + s2j_3d;
580 /* Note that attractive and repulsive potentials for individual
581 * pairs will partially cancel.
583 /* -dV/dr at the cut-off for LJ + Coulomb */
585 md_ljd*ffp->iparams[ti*ffp->atnr+tj].lj.c6 +
586 md_ljr*ffp->iparams[ti*ffp->atnr+tj].lj.c12 +
587 md_el*att[i].prop.q*att[j].prop.q;
589 /* d2V/dr2 at the cut-off for Coulomb, we neglect LJ */
590 dd = dd_el*att[i].prop.q*att[j].prop.q;
594 /* For constraints: adapt r and scaling for the Gaussian */
595 if (att[i].prop.bConstr)
599 approx_2dof(s2i_2d, r_buffer*s2i_2d/s2, &sh, &sc);
603 if (att[j].prop.bConstr)
607 approx_2dof(s2j_2d, r_buffer*s2j_2d/s2, &sh, &sc);
612 /* Exact contribution of an atom pair with Gaussian displacement
613 * with sigma s to the energy drift for a potential with
614 * derivative -md and second derivative dd at the cut-off.
615 * The only catch is that for potentials that change sign
616 * near the cut-off there could be an unlucky compensation
617 * of positive and negative energy drift.
618 * Such potentials are extremely rare though.
620 * Note that pot has unit energy*length, as the linear
621 * atom density still needs to be put in.
623 c_exp = exp(-rsh*rsh/(2*s2))/sqrt(2*M_PI);
624 c_erfc = 0.5*gmx_erfc(rsh/(sqrt(2*s2)));
628 md/2*((rsh*rsh + s2)*c_erfc - rsh*s*c_exp);
630 dd/6*(s*(rsh*rsh + 2*s2)*c_exp - rsh*(rsh*rsh + 3*s2)*c_erfc);
635 fprintf(debug, "n %d %d d s %.3f %.3f %.3f %.3f con %d md %8.1e dd %8.1e pot1 %8.1e pot2 %8.1e pot %8.1e\n",
637 sqrt(s2i_2d), sqrt(s2i_3d),
638 sqrt(s2j_2d), sqrt(s2j_3d),
639 att[i].prop.bConstr+att[j].prop.bConstr,
640 md, dd, pot1, pot2, pot);
643 /* Multiply by the number of atom pairs */
646 pot *= (double)att[i].n*(att[i].n - 1)/2;
650 pot *= (double)att[i].n*att[j].n;
652 /* We need the line density to get the energy drift of the system.
653 * The effective average r^2 is close to (rlist+sigma)^2.
655 pot *= 4*M_PI*sqr(rlist + s)/boxvol;
657 /* Add the unsigned drift to avoid cancellation of errors */
658 drift_tot += fabs(pot);
665 static real surface_frac(int cluster_size, real particle_distance, real rlist)
669 if (rlist < 0.5*particle_distance)
671 /* We have non overlapping spheres */
675 /* Half the inter-particle distance relative to rlist */
676 d = 0.5*particle_distance/rlist;
678 /* Determine the area of the surface at distance rlist to the closest
679 * particle, relative to surface of a sphere of radius rlist.
680 * The formulas below assume close to cubic cells for the pair search grid,
681 * which the pair search code tries to achieve.
682 * Note that in practice particle distances will not be delta distributed,
683 * but have some spread, often involving shorter distances,
684 * as e.g. O-H bonds in a water molecule. Thus the estimates below will
685 * usually be slightly too high and thus conservative.
687 switch (cluster_size)
690 /* One particle: trivial */
694 /* Two particles: two spheres at fractional distance 2*a */
698 /* We assume a perfect, symmetric tetrahedron geometry.
699 * The surface around a tetrahedron is too complex for a full
700 * analytical solution, so we use a Taylor expansion.
702 area_rel = (1.0 + 1/M_PI*(6*acos(1/sqrt(3))*d +
706 83.0/756.0*d*d*d*d*d*d)));
709 gmx_incons("surface_frac called with unsupported cluster_size");
713 return area_rel/cluster_size;
716 void calc_verlet_buffer_size(const gmx_mtop_t *mtop, real boxvol,
717 const t_inputrec *ir, real drift_target,
718 const verletbuf_list_setup_t *list_setup,
725 real particle_distance;
726 real nb_clust_frac_pairs_not_in_list_at_cutoff;
728 verletbuf_atomtype_t *att = NULL;
731 real md_ljd, md_ljr, md_el, dd_el;
733 real kT_fac, mass_min;
738 /* Resolution of the buffer size */
741 env = getenv("GMX_VERLET_BUFFER_RES");
744 sscanf(env, "%lf", &resolution);
747 /* In an atom wise pair-list there would be no pairs in the list
748 * beyond the pair-list cut-off.
749 * However, we use a pair-list of groups vs groups of atoms.
750 * For groups of 4 atoms, the parallelism of SSE instructions, only
751 * 10% of the atoms pairs are not in the list just beyond the cut-off.
752 * As this percentage increases slowly compared to the decrease of the
753 * Gaussian displacement distribution over this range, we can simply
754 * reduce the drift by this fraction.
755 * For larger groups, e.g. of 8 atoms, this fraction will be lower,
756 * so then buffer size will be on the conservative (large) side.
758 * Note that the formulas used here do not take into account
759 * cancellation of errors which could occur by missing both
760 * attractive and repulsive interactions.
762 * The only major assumption is homogeneous particle distribution.
763 * For an inhomogeneous system, such as a liquid-vapor system,
764 * the buffer will be underestimated. The actual energy drift
765 * will be higher by the factor: local/homogeneous particle density.
767 * The results of this estimate have been checked againt simulations.
768 * In most cases the real drift differs by less than a factor 2.
771 /* Worst case assumption: HCP packing of particles gives largest distance */
772 particle_distance = pow(boxvol*sqrt(2)/mtop->natoms, 1.0/3.0);
774 get_verlet_buffer_atomtypes(mtop, &att, &natt, n_nonlin_vsite);
775 assert(att != NULL && natt >= 0);
779 fprintf(debug, "particle distance assuming HCP packing: %f nm\n",
781 fprintf(debug, "energy drift atom types: %d\n", natt);
784 reppow = mtop->ffparams.reppow;
787 if (ir->vdwtype == evdwCUT)
789 /* -dV/dr of -r^-6 and r^-repporw */
790 md_ljd = -6*pow(ir->rvdw, -7.0);
791 md_ljr = reppow*pow(ir->rvdw, -(reppow+1));
792 /* The contribution of the second derivative is negligible */
796 gmx_fatal(FARGS, "Energy drift calculation is only implemented for plain cut-off Lennard-Jones interactions");
799 elfac = ONE_4PI_EPS0/ir->epsilon_r;
801 /* Determine md=-dV/dr and dd=d^2V/dr^2 */
804 if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
808 if (ir->coulombtype == eelCUT)
815 eps_rf = ir->epsilon_rf/ir->epsilon_r;
818 k_rf = pow(ir->rcoulomb, -3.0)*(eps_rf - ir->epsilon_r)/(2*eps_rf + ir->epsilon_r);
822 /* epsilon_rf = infinity */
823 k_rf = 0.5*pow(ir->rcoulomb, -3.0);
829 md_el = elfac*(pow(ir->rcoulomb, -2.0) - 2*k_rf*ir->rcoulomb);
831 dd_el = elfac*(2*pow(ir->rcoulomb, -3.0) + 2*k_rf);
833 else if (EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD)
837 b = calc_ewaldcoeff(ir->rcoulomb, ir->ewald_rtol);
840 md_el = elfac*(b*exp(-br*br)*M_2_SQRTPI/rc + gmx_erfc(br)/(rc*rc));
841 dd_el = elfac/(rc*rc)*(2*b*(1 + br*br)*exp(-br*br)*M_2_SQRTPI + 2*gmx_erfc(br)/rc);
845 gmx_fatal(FARGS, "Energy drift calculation is only implemented for Reaction-Field and Ewald electrostatics");
848 /* Determine the variance of the atomic displacement
849 * over nstlist-1 steps: kT_fac
850 * For inertial dynamics (not Brownian dynamics) the mass factor
851 * is not included in kT_fac, it is added later.
855 /* Get the displacement distribution from the random component only.
856 * With accurate integration the systematic (force) displacement
857 * should be negligible (unless nstlist is extremely large, which
858 * you wouldn't do anyhow).
860 kT_fac = 2*BOLTZ*ir->opts.ref_t[0]*(ir->nstlist-1)*ir->delta_t;
863 /* This is directly sigma^2 of the displacement */
864 kT_fac /= ir->bd_fric;
866 /* Set the masses to 1 as kT_fac is the full sigma^2,
867 * but we divide by m in ener_drift().
869 for (i = 0; i < natt; i++)
871 att[i].prop.mass = 1;
878 /* Per group tau_t is not implemented yet, use the maximum */
879 tau_t = ir->opts.tau_t[0];
880 for (i = 1; i < ir->opts.ngtc; i++)
882 tau_t = max(tau_t, ir->opts.tau_t[i]);
886 /* This kT_fac needs to be divided by the mass to get sigma^2 */
891 kT_fac = BOLTZ*ir->opts.ref_t[0]*sqr((ir->nstlist-1)*ir->delta_t);
894 mass_min = att[0].prop.mass;
895 for (i = 1; i < natt; i++)
897 mass_min = min(mass_min, att[i].prop.mass);
902 fprintf(debug, "md_ljd %e md_ljr %e\n", md_ljd, md_ljr);
903 fprintf(debug, "md_el %e dd_el %e\n", md_el, dd_el);
904 fprintf(debug, "sqrt(kT_fac) %f\n", sqrt(kT_fac));
905 fprintf(debug, "mass_min %f\n", mass_min);
908 /* Search using bisection */
910 /* The drift will be neglible at 5 times the max sigma */
911 ib1 = (int)(5*2*sqrt(kT_fac/mass_min)/resolution) + 1;
912 while (ib1 - ib0 > 1)
916 rl = max(ir->rvdw, ir->rcoulomb) + rb;
918 /* Calculate the average energy drift at the last step
919 * of the nstlist steps at which the pair-list is used.
921 drift = ener_drift(att, natt, &mtop->ffparams,
923 md_ljd, md_ljr, md_el, dd_el, rb,
926 /* Correct for the fact that we are using a Ni x Nj particle pair list
927 * and not a 1 x 1 particle pair list. This reduces the drift.
929 /* We don't have a formula for 8 (yet), use 4 which is conservative */
930 nb_clust_frac_pairs_not_in_list_at_cutoff =
931 surface_frac(min(list_setup->cluster_size_i, 4),
932 particle_distance, rl)*
933 surface_frac(min(list_setup->cluster_size_j, 4),
934 particle_distance, rl);
935 drift *= nb_clust_frac_pairs_not_in_list_at_cutoff;
937 /* Convert the drift to drift per unit time per atom */
938 drift /= ir->nstlist*ir->delta_t*mtop->natoms;
942 fprintf(debug, "ib %3d %3d %3d rb %.3f %dx%d fac %.3f drift %f\n",
944 list_setup->cluster_size_i, list_setup->cluster_size_j,
945 nb_clust_frac_pairs_not_in_list_at_cutoff,
949 if (fabs(drift) > drift_target)
961 *rlist = max(ir->rvdw, ir->rcoulomb) + ib1*resolution;