1 /* -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
4 * This source code is part of
8 * GROningen MAchine for Chemical Simulations
11 * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12 * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13 * Copyright (c) 2001-2004, The GROMACS development team,
14 * check out http://www.gromacs.org for more information.
16 * This program is free software; you can redistribute it and/or
17 * modify it under the terms of the GNU General Public License
18 * as published by the Free Software Foundation; either version 2
19 * of the License, or (at your option) any later version.
21 * If you want to redistribute modifications, please consider that
22 * scientific software is very special. Version control is crucial -
23 * bugs must be traceable. We will be happy to consider code for
24 * inclusion in the official distribution, but derived work must not
25 * be called official GROMACS. Details are found in the README & COPYING
26 * files - if they are missing, get the official version at www.gromacs.org.
28 * To help us fund GROMACS development, we humbly ask that you cite
29 * the papers on the package - you can find them in the top README file.
31 * For more info, check our website at http://www.gromacs.org
34 * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
42 #include <sys/types.h>
47 #include "gmx_fatal.h"
51 #include "calc_verletbuf.h"
52 #include "../mdlib/nbnxn_consts.h"
55 /* The include below sets the SIMD instruction type (precision+width)
56 * for all nbnxn SIMD search and non-bonded kernel code.
58 #ifdef GMX_NBNXN_HALF_WIDTH_SIMD
59 #define GMX_USE_HALF_WIDTH_SIMD_HERE
61 #include "gmx_simd_macros.h"
64 /* Struct for unique atom type for calculating the energy drift.
65 * The atom displacement depends on mass and constraints.
66 * The energy jump for given distance depend on LJ type and q.
71 int type; /* type (used for LJ parameters) */
73 int con; /* constrained: 0, else 1, if 1, use #DOF=2 iso 3 */
74 int n; /* total #atoms of this type in the system */
75 } verletbuf_atomtype_t;
78 void verletbuf_get_list_setup(gmx_bool bGPU,
79 verletbuf_list_setup_t *list_setup)
81 list_setup->cluster_size_i = NBNXN_CPU_CLUSTER_I_SIZE;
85 list_setup->cluster_size_j = NBNXN_GPU_CLUSTER_SIZE;
89 #ifndef GMX_NBNXN_SIMD
90 list_setup->cluster_size_j = NBNXN_CPU_CLUSTER_I_SIZE;
92 list_setup->cluster_size_j = GMX_SIMD_WIDTH_HERE;
93 #ifdef GMX_NBNXN_SIMD_2XNN
94 /* We assume the smallest cluster size to be on the safe side */
95 list_setup->cluster_size_j /= 2;
101 static void add_at(verletbuf_atomtype_t **att_p, int *natt_p,
102 real mass, int type, real q, int con, int nmol)
104 verletbuf_atomtype_t *att;
109 /* Ignore massless particles */
118 !(mass == att[i].mass &&
119 type == att[i].type &&
133 srenew(*att_p, *natt_p);
134 (*att_p)[i].mass = mass;
135 (*att_p)[i].type = type;
137 (*att_p)[i].con = con;
138 (*att_p)[i].n = nmol;
142 static void get_verlet_buffer_atomtypes(const gmx_mtop_t *mtop,
143 verletbuf_atomtype_t **att_p,
147 verletbuf_atomtype_t *att;
149 int mb, nmol, ft, i, j, a1, a2, a3, a;
150 const t_atoms *atoms;
154 real *con_m, *vsite_m, cam[5];
159 if (n_nonlin_vsite != NULL)
164 for (mb = 0; mb < mtop->nmolblock; mb++)
166 nmol = mtop->molblock[mb].nmol;
168 atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
170 /* Check for constraints, as they affect the kinetic energy */
171 snew(con_m, atoms->nr);
172 snew(vsite_m, atoms->nr);
174 for (ft = F_CONSTR; ft <= F_CONSTRNC; ft++)
176 il = &mtop->moltype[mtop->molblock[mb].type].ilist[ft];
178 for (i = 0; i < il->nr; i += 1+NRAL(ft))
180 a1 = il->iatoms[i+1];
181 a2 = il->iatoms[i+2];
182 con_m[a1] += atoms->atom[a2].m;
183 con_m[a2] += atoms->atom[a1].m;
187 il = &mtop->moltype[mtop->molblock[mb].type].ilist[F_SETTLE];
189 for (i = 0; i < il->nr; i += 1+NRAL(F_SETTLE))
191 a1 = il->iatoms[i+1];
192 a2 = il->iatoms[i+2];
193 a3 = il->iatoms[i+3];
194 con_m[a1] += atoms->atom[a2].m + atoms->atom[a3].m;
195 con_m[a2] += atoms->atom[a1].m + atoms->atom[a3].m;
196 con_m[a3] += atoms->atom[a1].m + atoms->atom[a2].m;
199 /* Check for virtual sites, determine mass from constructing atoms */
200 for (ft = 0; ft < F_NRE; ft++)
204 il = &mtop->moltype[mtop->molblock[mb].type].ilist[ft];
206 for (i = 0; i < il->nr; i += 1+NRAL(ft))
208 ip = &mtop->ffparams.iparams[il->iatoms[i]];
210 a1 = il->iatoms[i+1];
212 for (j = 1; j < NRAL(ft); j++)
214 cam[j] = atoms->atom[il->iatoms[i+1+j]].m;
217 cam[j] = vsite_m[il->iatoms[i+1+j]];
221 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.",
222 *mtop->moltype[mtop->molblock[mb].type].name,
223 interaction_function[ft].longname,
224 il->iatoms[i+1+j]+1);
231 /* Exact except for ignoring constraints */
232 vsite_m[a1] = (cam[2]*sqr(1-ip->vsite.a) + cam[1]*sqr(ip->vsite.a))/(cam[1]*cam[2]);
235 /* Exact except for ignoring constraints */
236 vsite_m[a1] = (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))/(cam[1]*cam[2]*cam[3]);
239 /* Use the mass of the lightest constructing atom.
240 * This is an approximation.
241 * If the distance of the virtual site to the
242 * constructing atom is less than all distances
243 * between constructing atoms, this is a safe
244 * over-estimate of the displacement of the vsite.
245 * This condition holds for all H mass replacement
246 * replacement vsite constructions, except for SP2/3
247 * groups. In SP3 groups one H will have a F_VSITE3
248 * construction, so even there the total drift
249 * estimation shouldn't be far off.
252 vsite_m[a1] = cam[1];
253 for (j = 2; j < NRAL(ft); j++)
255 vsite_m[a1] = min(vsite_m[a1], cam[j]);
257 if (n_nonlin_vsite != NULL)
259 *n_nonlin_vsite += nmol;
267 for (a = 0; a < atoms->nr; a++)
269 at = &atoms->atom[a];
270 /* We consider an atom constrained, #DOF=2, when it is
271 * connected with constraints to one or more atoms with
272 * total mass larger than 1.5 that of the atom itself.
275 at->m, at->type, at->q, con_m[a] > 1.5*at->m, nmol);
284 for (a = 0; a < natt; a++)
286 fprintf(debug, "type %d: m %5.2f t %d q %6.3f con %d n %d\n",
287 a, att[a].mass, att[a].type, att[a].q, att[a].con, att[a].n);
295 static void approx_2dof(real s2, real x,
296 real *shift, real *scale)
298 /* A particle with 1 DOF constrained has 2 DOFs instead of 3.
299 * This code is also used for particles with multiple constraints,
300 * in which case we overestimate the displacement.
301 * The 2DOF distribution is sqrt(pi/2)*erfc(r/(sqrt(2)*s))/(2*s).
302 * We approximate this with scale*Gaussian(s,r+shift),
303 * by matching the distribution value and derivative at x.
304 * This is a tight overestimate for all r>=0 at any s and x.
308 ex = exp(-x*x/(2*s2));
309 er = gmx_erfc(x/sqrt(2*s2));
311 *shift = -x + sqrt(2*s2/M_PI)*ex/er;
312 *scale = 0.5*M_PI*exp(ex*ex/(M_PI*er*er))*er;
315 static real ener_drift(const verletbuf_atomtype_t *att, int natt,
316 const gmx_ffparams_t *ffp,
318 real md_ljd, real md_ljr, real md_el, real dd_el,
320 real rlist, real boxvol)
322 double drift_tot, pot1, pot2, pot;
324 real s2i, s2j, s2, s;
328 double c_exp, c_erfc;
332 /* Loop over the different atom type pairs */
333 for (i = 0; i < natt; i++)
335 s2i = kT_fac/att[i].mass;
338 for (j = i; j < natt; j++)
340 s2j = kT_fac/att[j].mass;
343 /* Note that attractive and repulsive potentials for individual
344 * pairs will partially cancel.
346 /* -dV/dr at the cut-off for LJ + Coulomb */
348 md_ljd*ffp->iparams[ti*ffp->atnr+tj].lj.c6 +
349 md_ljr*ffp->iparams[ti*ffp->atnr+tj].lj.c12 +
350 md_el*att[i].q*att[j].q;
352 /* d2V/dr2 at the cut-off for Coulomb, we neglect LJ */
353 dd = dd_el*att[i].q*att[j].q;
359 /* For constraints: adapt r and scaling for the Gaussian */
363 approx_2dof(s2i, r_buffer*s2i/s2, &sh, &sc);
370 approx_2dof(s2j, r_buffer*s2j/s2, &sh, &sc);
375 /* Exact contribution of an atom pair with Gaussian displacement
376 * with sigma s to the energy drift for a potential with
377 * derivative -md and second derivative dd at the cut-off.
378 * The only catch is that for potentials that change sign
379 * near the cut-off there could be an unlucky compensation
380 * of positive and negative energy drift.
381 * Such potentials are extremely rare though.
383 * Note that pot has unit energy*length, as the linear
384 * atom density still needs to be put in.
386 c_exp = exp(-rsh*rsh/(2*s2))/sqrt(2*M_PI);
387 c_erfc = 0.5*gmx_erfc(rsh/(sqrt(2*s2)));
391 md/2*((rsh*rsh + s2)*c_erfc - rsh*s*c_exp);
393 dd/6*(s*(rsh*rsh + 2*s2)*c_exp - rsh*(rsh*rsh + 3*s2)*c_erfc);
398 fprintf(debug, "n %d %d d s %.3f %.3f con %d md %8.1e dd %8.1e pot1 %8.1e pot2 %8.1e pot %8.1e\n",
399 att[i].n, att[j].n, sqrt(s2i), sqrt(s2j),
400 att[i].con+att[j].con,
401 md, dd, pot1, pot2, pot);
404 /* Multiply by the number of atom pairs */
407 pot *= (double)att[i].n*(att[i].n - 1)/2;
411 pot *= (double)att[i].n*att[j].n;
413 /* We need the line density to get the energy drift of the system.
414 * The effective average r^2 is close to (rlist+sigma)^2.
416 pot *= 4*M_PI*sqr(rlist + s)/boxvol;
418 /* Add the unsigned drift to avoid cancellation of errors */
419 drift_tot += fabs(pot);
426 static real surface_frac(int cluster_size, real particle_distance, real rlist)
430 if (rlist < 0.5*particle_distance)
432 /* We have non overlapping spheres */
436 /* Half the inter-particle distance relative to rlist */
437 d = 0.5*particle_distance/rlist;
439 /* Determine the area of the surface at distance rlist to the closest
440 * particle, relative to surface of a sphere of radius rlist.
441 * The formulas below assume close to cubic cells for the pair search grid,
442 * which the pair search code tries to achieve.
443 * Note that in practice particle distances will not be delta distributed,
444 * but have some spread, often involving shorter distances,
445 * as e.g. O-H bonds in a water molecule. Thus the estimates below will
446 * usually be slightly too high and thus conservative.
448 switch (cluster_size)
451 /* One particle: trivial */
455 /* Two particles: two spheres at fractional distance 2*a */
459 /* We assume a perfect, symmetric tetrahedron geometry.
460 * The surface around a tetrahedron is too complex for a full
461 * analytical solution, so we use a Taylor expansion.
463 area_rel = (1.0 + 1/M_PI*(6*acos(1/sqrt(3))*d +
467 83.0/756.0*d*d*d*d*d*d)));
470 gmx_incons("surface_frac called with unsupported cluster_size");
474 return area_rel/cluster_size;
477 void calc_verlet_buffer_size(const gmx_mtop_t *mtop, real boxvol,
478 const t_inputrec *ir, real drift_target,
479 const verletbuf_list_setup_t *list_setup,
486 real particle_distance;
487 real nb_clust_frac_pairs_not_in_list_at_cutoff;
489 verletbuf_atomtype_t *att = NULL;
492 real md_ljd, md_ljr, md_el, dd_el;
494 real kT_fac, mass_min;
499 /* Resolution of the buffer size */
502 env = getenv("GMX_VERLET_BUFFER_RES");
505 sscanf(env, "%lf", &resolution);
508 /* In an atom wise pair-list there would be no pairs in the list
509 * beyond the pair-list cut-off.
510 * However, we use a pair-list of groups vs groups of atoms.
511 * For groups of 4 atoms, the parallelism of SSE instructions, only
512 * 10% of the atoms pairs are not in the list just beyond the cut-off.
513 * As this percentage increases slowly compared to the decrease of the
514 * Gaussian displacement distribution over this range, we can simply
515 * reduce the drift by this fraction.
516 * For larger groups, e.g. of 8 atoms, this fraction will be lower,
517 * so then buffer size will be on the conservative (large) side.
519 * Note that the formulas used here do not take into account
520 * cancellation of errors which could occur by missing both
521 * attractive and repulsive interactions.
523 * The only major assumption is homogeneous particle distribution.
524 * For an inhomogeneous system, such as a liquid-vapor system,
525 * the buffer will be underestimated. The actual energy drift
526 * will be higher by the factor: local/homogeneous particle density.
528 * The results of this estimate have been checked againt simulations.
529 * In most cases the real drift differs by less than a factor 2.
532 /* Worst case assumption: HCP packing of particles gives largest distance */
533 particle_distance = pow(boxvol*sqrt(2)/mtop->natoms, 1.0/3.0);
535 get_verlet_buffer_atomtypes(mtop, &att, &natt, n_nonlin_vsite);
536 assert(att != NULL && natt >= 0);
540 fprintf(debug, "particle distance assuming HCP packing: %f nm\n",
542 fprintf(debug, "energy drift atom types: %d\n", natt);
545 reppow = mtop->ffparams.reppow;
548 if (ir->vdwtype == evdwCUT)
550 /* -dV/dr of -r^-6 and r^-repporw */
551 md_ljd = -6*pow(ir->rvdw, -7.0);
552 md_ljr = reppow*pow(ir->rvdw, -(reppow+1));
553 /* The contribution of the second derivative is negligible */
557 gmx_fatal(FARGS, "Energy drift calculation is only implemented for plain cut-off Lennard-Jones interactions");
560 elfac = ONE_4PI_EPS0/ir->epsilon_r;
562 /* Determine md=-dV/dr and dd=d^2V/dr^2 */
565 if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
569 if (ir->coulombtype == eelCUT)
576 eps_rf = ir->epsilon_rf/ir->epsilon_r;
579 k_rf = pow(ir->rcoulomb, -3.0)*(eps_rf - ir->epsilon_r)/(2*eps_rf + ir->epsilon_r);
583 /* epsilon_rf = infinity */
584 k_rf = 0.5*pow(ir->rcoulomb, -3.0);
590 md_el = elfac*(pow(ir->rcoulomb, -2.0) - 2*k_rf*ir->rcoulomb);
592 dd_el = elfac*(2*pow(ir->rcoulomb, -3.0) + 2*k_rf);
594 else if (EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD)
598 b = calc_ewaldcoeff(ir->rcoulomb, ir->ewald_rtol);
601 md_el = elfac*(b*exp(-br*br)*M_2_SQRTPI/rc + gmx_erfc(br)/(rc*rc));
602 dd_el = elfac/(rc*rc)*(2*b*(1 + br*br)*exp(-br*br)*M_2_SQRTPI + 2*gmx_erfc(br)/rc);
606 gmx_fatal(FARGS, "Energy drift calculation is only implemented for Reaction-Field and Ewald electrostatics");
609 /* Determine the variance of the atomic displacement
610 * over nstlist-1 steps: kT_fac
611 * For inertial dynamics (not Brownian dynamics) the mass factor
612 * is not included in kT_fac, it is added later.
616 /* Get the displacement distribution from the random component only.
617 * With accurate integration the systematic (force) displacement
618 * should be negligible (unless nstlist is extremely large, which
619 * you wouldn't do anyhow).
621 kT_fac = 2*BOLTZ*ir->opts.ref_t[0]*(ir->nstlist-1)*ir->delta_t;
624 /* This is directly sigma^2 of the displacement */
625 kT_fac /= ir->bd_fric;
627 /* Set the masses to 1 as kT_fac is the full sigma^2,
628 * but we divide by m in ener_drift().
630 for (i = 0; i < natt; i++)
639 /* Per group tau_t is not implemented yet, use the maximum */
640 tau_t = ir->opts.tau_t[0];
641 for (i = 1; i < ir->opts.ngtc; i++)
643 tau_t = max(tau_t, ir->opts.tau_t[i]);
647 /* This kT_fac needs to be divided by the mass to get sigma^2 */
652 kT_fac = BOLTZ*ir->opts.ref_t[0]*sqr((ir->nstlist-1)*ir->delta_t);
655 mass_min = att[0].mass;
656 for (i = 1; i < natt; i++)
658 mass_min = min(mass_min, att[i].mass);
663 fprintf(debug, "md_ljd %e md_ljr %e\n", md_ljd, md_ljr);
664 fprintf(debug, "md_el %e dd_el %e\n", md_el, dd_el);
665 fprintf(debug, "sqrt(kT_fac) %f\n", sqrt(kT_fac));
666 fprintf(debug, "mass_min %f\n", mass_min);
669 /* Search using bisection */
671 /* The drift will be neglible at 5 times the max sigma */
672 ib1 = (int)(5*2*sqrt(kT_fac/mass_min)/resolution) + 1;
673 while (ib1 - ib0 > 1)
677 rl = max(ir->rvdw, ir->rcoulomb) + rb;
679 /* Calculate the average energy drift at the last step
680 * of the nstlist steps at which the pair-list is used.
682 drift = ener_drift(att, natt, &mtop->ffparams,
684 md_ljd, md_ljr, md_el, dd_el, rb,
687 /* Correct for the fact that we are using a Ni x Nj particle pair list
688 * and not a 1 x 1 particle pair list. This reduces the drift.
690 /* We don't have a formula for 8 (yet), use 4 which is conservative */
691 nb_clust_frac_pairs_not_in_list_at_cutoff =
692 surface_frac(min(list_setup->cluster_size_i, 4),
693 particle_distance, rl)*
694 surface_frac(min(list_setup->cluster_size_j, 4),
695 particle_distance, rl);
696 drift *= nb_clust_frac_pairs_not_in_list_at_cutoff;
698 /* Convert the drift to drift per unit time per atom */
699 drift /= ir->nstlist*ir->delta_t*mtop->natoms;
703 fprintf(debug, "ib %3d %3d %3d rb %.3f %dx%d fac %.3f drift %f\n",
705 list_setup->cluster_size_i, list_setup->cluster_size_j,
706 nb_clust_frac_pairs_not_in_list_at_cutoff,
710 if (fabs(drift) > drift_target)
722 *rlist = max(ir->rvdw, ir->rcoulomb) + ib1*resolution;