cf886d8079ae4611b324f9c482a2b4e6bca8a99f
[alexxy/gromacs.git] / src / gromacs / mdlib / calc_verletbuf.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35 #include "gmxpre.h"
36
37 #include "calc_verletbuf.h"
38
39 #include <assert.h>
40 #include <stdlib.h>
41
42 #include <cmath>
43
44 #include <algorithm>
45
46 #include "gromacs/ewald/ewald-utils.h"
47 #include "gromacs/math/functions.h"
48 #include "gromacs/math/units.h"
49 #include "gromacs/math/vec.h"
50 #include "gromacs/mdlib/nb_verlet.h"
51 #include "gromacs/mdlib/nbnxn_simd.h"
52 #include "gromacs/mdlib/nbnxn_util.h"
53 #include "gromacs/mdtypes/inputrec.h"
54 #include "gromacs/mdtypes/md_enums.h"
55 #include "gromacs/topology/ifunc.h"
56 #include "gromacs/topology/topology.h"
57 #include "gromacs/utility/fatalerror.h"
58 #include "gromacs/utility/smalloc.h"
59
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.
65  *
66  * Significant approximations used:
67  *
68  * Uniform particle density. UNDERESTIMATES the drift by rho_global/rho_local.
69  *
70  * Interactions don't affect particle motion. OVERESTIMATES the drift on longer
71  * time scales. This approximation probably introduces the largest errors.
72  *
73  * Only take one constraint per particle into account: OVERESTIMATES the drift.
74  *
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.
78  *
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.
82  *
83  * Note that the formulas for normal atoms and linear virtual sites are exact,
84  * apart from the first two approximations.
85  *
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.
89  */
90
91
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.
95  */
96 typedef struct
97 {
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;
101
102 // Struct for derivatives of a non-bonded interaction potential
103 typedef struct
104 {
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;
109
110 VerletbufListSetup verletbufGetListSetup(int nbnxnKernelType)
111 {
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.
114      */
115     VerletbufListSetup listSetup;
116
117     listSetup.cluster_size_i = nbnxn_kernel_to_cluster_i_size(nbnxnKernelType);
118     listSetup.cluster_size_j = nbnxn_kernel_to_cluster_j_size(nbnxnKernelType);
119
120     if (nbnxnKernelType == nbnxnk8x8x8_GPU ||
121         nbnxnKernelType == nbnxnk8x8x8_PlainC)
122     {
123         /* The GPU kernels (except for OpenCL) split the j-clusters in two halves */
124         listSetup.cluster_size_j /= 2;
125     }
126
127     return listSetup;
128 }
129
130 VerletbufListSetup verletbufGetSafeListSetup(ListSetupType listType)
131 {
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.
136      */
137     int nbnxnKernelType;
138
139     if (listType == ListSetupType::Gpu)
140     {
141         nbnxnKernelType = nbnxnk8x8x8_GPU;
142     }
143     else if (GMX_SIMD && listType == ListSetupType::CpuSimdWhenSupported)
144     {
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
151     }
152     else
153     {
154         nbnxnKernelType = nbnxnk4x4_PlainC;
155     }
156
157     return verletbufGetListSetup(nbnxnKernelType);
158 }
159
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)
163 {
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);
170 }
171
172 static void add_at(verletbuf_atomtype_t **att_p, int *natt_p,
173                    const atom_nonbonded_kinetic_prop_t *prop,
174                    int nmol)
175 {
176     verletbuf_atomtype_t   *att;
177     int                     natt, i;
178
179     if (prop->mass == 0)
180     {
181         /* Ignore massless particles */
182         return;
183     }
184
185     att  = *att_p;
186     natt = *natt_p;
187
188     i = 0;
189     while (i < natt && !atom_nonbonded_kinetic_prop_equal(prop, &att[i].prop))
190     {
191         i++;
192     }
193
194     if (i < natt)
195     {
196         att[i].n += nmol;
197     }
198     else
199     {
200         (*natt_p)++;
201         srenew(*att_p, *natt_p);
202         (*att_p)[i].prop = *prop;
203         (*att_p)[i].n    = nmol;
204     }
205 }
206
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)
211 {
212     if (!setMassesToOne)
213     {
214         return atoms.atom[atomIndex].m;
215     }
216     else
217     {
218         return 1;
219     }
220 }
221
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)
227 {
228     int            ft, i;
229     const t_ilist *il;
230
231     *n_nonlin_vsite = 0;
232
233     /* Check for virtual sites, determine mass from constructing atoms */
234     for (ft = 0; ft < F_NRE; ft++)
235     {
236         if (IS_VSITE(ft))
237         {
238             il = &moltype->ilist[ft];
239
240             for (i = 0; i < il->nr; i += 1+NRAL(ft))
241             {
242                 const t_iparams *ip;
243                 real             inv_mass, coeff, m_aj;
244                 int              a1;
245
246                 ip = &ffparams->iparams[il->iatoms[i]];
247
248                 a1 = il->iatoms[i+1];
249
250                 if (ft != F_VSITEN)
251                 {
252                     /* Only vsiten can have more than four
253                        constructing atoms, so NRAL(ft) <= 5 */
254                     int        j;
255                     real      *cam;
256                     const int  maxj = NRAL(ft);
257
258                     snew(cam, maxj);
259                     assert(maxj <= 5);
260                     for (j = 1; j < maxj; j++)
261                     {
262                         int aj = il->iatoms[i + 1 + j];
263                         cam[j] = getMass(moltype->atoms, aj, setMassesToOne);
264                         if (cam[j] == 0)
265                         {
266                             cam[j] = vsite_m[aj];
267                         }
268                         if (cam[j] == 0)
269                         {
270                             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.",
271                                       *moltype->name,
272                                       interaction_function[ft].longname,
273                                       aj + 1);
274                         }
275                     }
276
277                     switch (ft)
278                     {
279                         case F_VSITE2:
280                             /* Exact */
281                             vsite_m[a1] = (cam[1]*cam[2])/(cam[2]*gmx::square(1-ip->vsite.a) + cam[1]*gmx::square(ip->vsite.a));
282                             break;
283                         case F_VSITE3:
284                             /* Exact */
285                             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));
286                             break;
287                         case F_VSITEN:
288                             gmx_incons("Invalid vsite type");
289                             break;
290                         default:
291                             /* Use the mass of the lightest constructing atom.
292                              * This is an approximation.
293                              * If the distance of the virtual site to the
294                              * constructing atom is less than all distances
295                              * between constructing atoms, this is a safe
296                              * over-estimate of the displacement of the vsite.
297                              * This condition holds for all H mass replacement
298                              * vsite constructions, except for SP2/3 groups.
299                              * In SP3 groups one H will have a F_VSITE3
300                              * construction, so even there the total drift
301                              * estimate shouldn't be far off.
302                              */
303                             vsite_m[a1] = cam[1];
304                             for (j = 2; j < maxj; j++)
305                             {
306                                 vsite_m[a1] = std::min(vsite_m[a1], cam[j]);
307                             }
308                             (*n_nonlin_vsite)++;
309                             break;
310                     }
311                     sfree(cam);
312                 }
313                 else
314                 {
315                     int j;
316
317                     /* Exact */
318                     inv_mass = 0;
319                     for (j = 0; j < 3*ffparams->iparams[il->iatoms[i]].vsiten.n; j += 3)
320                     {
321                         int aj = il->iatoms[i + j + 2];
322                         coeff  = ffparams->iparams[il->iatoms[i+j]].vsiten.a;
323                         if (moltype->atoms.atom[aj].ptype == eptVSite)
324                         {
325                             m_aj = vsite_m[aj];
326                         }
327                         else
328                         {
329                             m_aj = moltype->atoms.atom[aj].m;
330                         }
331                         if (m_aj <= 0)
332                         {
333                             gmx_incons("The mass of a vsiten constructing atom is <= 0");
334                         }
335                         inv_mass += coeff*coeff/m_aj;
336                     }
337                     vsite_m[a1] = 1/inv_mass;
338                     /* Correct for loop increment of i */
339                     i += j - 1 - NRAL(ft);
340                 }
341                 if (gmx_debug_at)
342                 {
343                     fprintf(debug, "atom %4d %-20s mass %6.3f\n",
344                             a1, interaction_function[ft].longname, vsite_m[a1]);
345                 }
346             }
347         }
348     }
349 }
350
351 static void get_verlet_buffer_atomtypes(const gmx_mtop_t      *mtop,
352                                         bool                   setMassesToOne,
353                                         verletbuf_atomtype_t **att_p,
354                                         int                   *natt_p,
355                                         int                   *n_nonlin_vsite)
356 {
357     verletbuf_atomtype_t          *att;
358     int                            natt;
359     int                            mb, nmol, ft, i, a1, a2, a3, a;
360     const t_atoms                 *atoms;
361     const t_ilist                 *il;
362     const t_iparams               *ip;
363     atom_nonbonded_kinetic_prop_t *prop;
364     real                          *vsite_m;
365     int                            n_nonlin_vsite_mol;
366
367     att  = nullptr;
368     natt = 0;
369
370     if (n_nonlin_vsite != nullptr)
371     {
372         *n_nonlin_vsite = 0;
373     }
374
375     for (mb = 0; mb < mtop->nmolblock; mb++)
376     {
377         nmol = mtop->molblock[mb].nmol;
378
379         atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
380
381         /* Check for constraints, as they affect the kinetic energy.
382          * For virtual sites we need the masses and geometry of
383          * the constructing atoms to determine their velocity distribution.
384          */
385         snew(prop, atoms->nr);
386         snew(vsite_m, atoms->nr);
387
388         for (ft = F_CONSTR; ft <= F_CONSTRNC; ft++)
389         {
390             il = &mtop->moltype[mtop->molblock[mb].type].ilist[ft];
391
392             for (i = 0; i < il->nr; i += 1+NRAL(ft))
393             {
394                 ip         = &mtop->ffparams.iparams[il->iatoms[i]];
395                 a1         = il->iatoms[i+1];
396                 a2         = il->iatoms[i+2];
397                 real mass1 = getMass(*atoms, a1, setMassesToOne);
398                 real mass2 = getMass(*atoms, a2, setMassesToOne);
399                 if (mass2 > prop[a1].con_mass)
400                 {
401                     prop[a1].con_mass = mass2;
402                     prop[a1].con_len  = ip->constr.dA;
403                 }
404                 if (mass1 > prop[a2].con_mass)
405                 {
406                     prop[a2].con_mass = mass1;
407                     prop[a2].con_len  = ip->constr.dA;
408                 }
409             }
410         }
411
412         il = &mtop->moltype[mtop->molblock[mb].type].ilist[F_SETTLE];
413
414         for (i = 0; i < il->nr; i += 1+NRAL(F_SETTLE))
415         {
416             ip         = &mtop->ffparams.iparams[il->iatoms[i]];
417             a1         = il->iatoms[i+1];
418             a2         = il->iatoms[i+2];
419             a3         = il->iatoms[i+3];
420             /* Usually the mass of a1 (usually oxygen) is larger than a2/a3.
421              * If this is not the case, we overestimate the displacement,
422              * which leads to a larger buffer (ok since this is an exotic case).
423              */
424             prop[a1].con_mass = getMass(*atoms, a2, setMassesToOne);
425             prop[a1].con_len  = ip->settle.doh;
426
427             prop[a2].con_mass = getMass(*atoms, a1, setMassesToOne);
428             prop[a2].con_len  = ip->settle.doh;
429
430             prop[a3].con_mass = getMass(*atoms, a1, setMassesToOne);
431             prop[a3].con_len  = ip->settle.doh;
432         }
433
434         get_vsite_masses(&mtop->moltype[mtop->molblock[mb].type],
435                          &mtop->ffparams,
436                          setMassesToOne,
437                          vsite_m,
438                          &n_nonlin_vsite_mol);
439         if (n_nonlin_vsite != nullptr)
440         {
441             *n_nonlin_vsite += nmol*n_nonlin_vsite_mol;
442         }
443
444         for (a = 0; a < atoms->nr; a++)
445         {
446             if (atoms->atom[a].ptype == eptVSite)
447             {
448                 prop[a].mass = vsite_m[a];
449             }
450             else
451             {
452                 prop[a].mass = getMass(*atoms, a, setMassesToOne);
453             }
454             prop[a].type     = atoms->atom[a].type;
455             prop[a].q        = atoms->atom[a].q;
456             /* We consider an atom constrained, #DOF=2, when it is
457              * connected with constraints to (at least one) atom with
458              * a mass of more than 0.4x its own mass. This is not a critical
459              * parameter, since with roughly equal masses the unconstrained
460              * and constrained displacement will not differ much (and both
461              * overestimate the displacement).
462              */
463             prop[a].bConstr = (prop[a].con_mass > 0.4*prop[a].mass);
464
465             add_at(&att, &natt, &prop[a], nmol);
466         }
467
468         sfree(vsite_m);
469         sfree(prop);
470     }
471
472     if (gmx_debug_at)
473     {
474         for (a = 0; a < natt; a++)
475         {
476             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",
477                     a, att[a].prop.mass, att[a].prop.type, att[a].prop.q,
478                     att[a].prop.bConstr, att[a].prop.con_mass, att[a].prop.con_len,
479                     att[a].n);
480         }
481     }
482
483     *att_p  = att;
484     *natt_p = natt;
485 }
486
487 /* This function computes two components of the estimate of the variance
488  * in the displacement of one atom in a system of two constrained atoms.
489  * Returns in sigma2_2d the variance due to rotation of the constrained
490  * atom around the atom to which it constrained.
491  * Returns in sigma2_3d the variance due to displacement of the COM
492  * of the whole system of the two constrained atoms.
493  *
494  * Note that we only take a single constraint (the one to the heaviest atom)
495  * into account. If an atom has multiple constraints, this will result in
496  * an overestimate of the displacement, which gives a larger drift and buffer.
497  */
498 void constrained_atom_sigma2(real                                 kT_fac,
499                              const atom_nonbonded_kinetic_prop_t *prop,
500                              real                                *sigma2_2d,
501                              real                                *sigma2_3d)
502 {
503     /* Here we decompose the motion of a constrained atom into two
504      * components: rotation around the COM and translation of the COM.
505      */
506
507     /* Determine the variance of the arc length for the two rotational DOFs */
508     real massFraction = prop->con_mass/(prop->mass + prop->con_mass);
509     real sigma2_rot   = kT_fac*massFraction/prop->mass;
510
511     /* The distance from the atom to the COM, i.e. the rotational arm */
512     real comDistance  = prop->con_len*massFraction;
513
514     /* The variance relative to the arm */
515     real sigma2_rel   = sigma2_rot/gmx::square(comDistance);
516
517     /* For sigma2_rel << 1 we don't notice the rotational effect and
518      * we have a normal, Gaussian displacement distribution.
519      * For larger sigma2_rel the displacement is much less, in fact it can
520      * not exceed 2*comDistance. We can calculate MSD/arm^2 as:
521      *   integral_x=0-inf distance2(x) x/sigma2_rel exp(-x^2/(2 sigma2_rel)) dx
522      * where x is angular displacement and distance2(x) is the distance^2
523      * between points at angle 0 and x:
524      *   distance2(x) = (sin(x) - sin(0))^2 + (cos(x) - cos(0))^2
525      * The limiting value of this MSD is 2, which is also the value for
526      * a uniform rotation distribution that would be reached at long time.
527      * The maximum is 2.5695 at sigma2_rel = 4.5119.
528      * We approximate this integral with a rational polynomial with
529      * coefficients from a Taylor expansion. This approximation is an
530      * overestimate for all values of sigma2_rel. Its maximum value
531      * of 2.6491 is reached at sigma2_rel = sqrt(45/2) = 4.7434.
532      * We keep the approximation constant after that.
533      * We use this approximate MSD as the variance for a Gaussian distribution.
534      *
535      * NOTE: For any sensible buffer tolerance this will result in a (large)
536      * overestimate of the buffer size, since the Gaussian has a long tail,
537      * whereas the actual distribution can not reach values larger than 2.
538      */
539     /* Coeffients obtained from a Taylor expansion */
540     const real a = 1.0/3.0;
541     const real b = 2.0/45.0;
542
543     /* Our approximation is constant after sigma2_rel = 1/sqrt(b) */
544     sigma2_rel   = std::min(sigma2_rel, 1/std::sqrt(b));
545
546     /* Compute the approximate sigma^2 for 2D motion due to the rotation */
547     *sigma2_2d   = gmx::square(comDistance)*
548         sigma2_rel/(1 + a*sigma2_rel + b*gmx::square(sigma2_rel));
549
550     /* The constrained atom also moves (in 3D) with the COM of both atoms */
551     *sigma2_3d   = kT_fac/(prop->mass + prop->con_mass);
552 }
553
554 static void get_atom_sigma2(real                                 kT_fac,
555                             const atom_nonbonded_kinetic_prop_t *prop,
556                             real                                *sigma2_2d,
557                             real                                *sigma2_3d)
558 {
559     if (prop->bConstr)
560     {
561         /* Complicated constraint calculation in a separate function */
562         constrained_atom_sigma2(kT_fac, prop, sigma2_2d, sigma2_3d);
563     }
564     else
565     {
566         /* Unconstrained atom: trivial */
567         *sigma2_2d = 0;
568         *sigma2_3d = kT_fac/prop->mass;
569     }
570 }
571
572 static void approx_2dof(real s2, real x, real *shift, real *scale)
573 {
574     /* A particle with 1 DOF constrained has 2 DOFs instead of 3.
575      * This code is also used for particles with multiple constraints,
576      * in which case we overestimate the displacement.
577      * The 2DOF distribution is sqrt(pi/2)*erfc(r/(sqrt(2)*s))/(2*s).
578      * We approximate this with scale*Gaussian(s,r+shift),
579      * by matching the distribution value and derivative at x.
580      * This is a tight overestimate for all r>=0 at any s and x.
581      */
582     real ex, er;
583
584     ex = std::exp(-x*x/(2*s2));
585     er = std::erfc(x/std::sqrt(2*s2));
586
587     *shift = -x + std::sqrt(2*s2/M_PI)*ex/er;
588     *scale = 0.5*M_PI*std::exp(ex*ex/(M_PI*er*er))*er;
589 }
590
591 // Returns an (over)estimate of the energy drift for a single atom pair,
592 // given the kinetic properties, displacement variances and list buffer.
593 static real energyDriftAtomPair(bool isConstrained_i,
594                                 bool isConstrained_j,
595                                 real s2, real s2i_2d, real s2j_2d,
596                                 real r_buffer,
597                                 const pot_derivatives_t *der)
598 {
599     // For relatively small arguments erfc() is so small that if will be 0.0
600     // when stored in a float. We set an argument limit of 8 (Erfc(8)=1e-29),
601     // such that we can divide by erfc and have some space left for arithmetic.
602     const real erfc_arg_max = 8.0;
603
604     real       rsh    = r_buffer;
605     real       sc_fac = 1.0;
606
607     real       c_exp, c_erfc;
608
609     if (rsh*rsh > 2*s2*erfc_arg_max*erfc_arg_max)
610     {
611         // Below we calculate c_erfc = 0.5*erfc(rsh/sqrt(2*s2))
612         // When rsh/sqrt(2*s2) increases, this erfc will be the first
613         // result that underflows and becomes 0.0. To avoid this,
614         // we set c_exp=0 and c_erfc=0 for large arguments.
615         // This also avoids NaN in approx_2dof().
616         // In any relevant case this has no effect on the results,
617         // since c_exp < 6e-29, so the displacement is completely
618         // negligible for such atom pairs (and an overestimate).
619         // In nearly all use cases, there will be other atom pairs
620         // that contribute much more to the total, so zeroing
621         // this particular contribution has no effect at all.
622         c_exp  = 0;
623         c_erfc = 0;
624     }
625     else
626     {
627         /* For constraints: adapt r and scaling for the Gaussian */
628         if (isConstrained_i)
629         {
630             real sh, sc;
631
632             approx_2dof(s2i_2d, r_buffer*s2i_2d/s2, &sh, &sc);
633             rsh    += sh;
634             sc_fac *= sc;
635         }
636         if (isConstrained_j)
637         {
638             real sh, sc;
639
640             approx_2dof(s2j_2d, r_buffer*s2j_2d/s2, &sh, &sc);
641             rsh    += sh;
642             sc_fac *= sc;
643         }
644
645         /* Exact contribution of an atom pair with Gaussian displacement
646          * with sigma s to the energy drift for a potential with
647          * derivative -md and second derivative dd at the cut-off.
648          * The only catch is that for potentials that change sign
649          * near the cut-off there could be an unlucky compensation
650          * of positive and negative energy drift.
651          * Such potentials are extremely rare though.
652          *
653          * Note that pot has unit energy*length, as the linear
654          * atom density still needs to be put in.
655          */
656         c_exp  = std::exp(-rsh*rsh/(2*s2))/std::sqrt(2*M_PI);
657         c_erfc = 0.5*std::erfc(rsh/(std::sqrt(2*s2)));
658     }
659     real s    = std::sqrt(s2);
660     real rsh2 = rsh*rsh;
661
662     real pot1 = sc_fac*
663         der->md1/2*((rsh2 + s2)*c_erfc - rsh*s*c_exp);
664     real pot2 = sc_fac*
665         der->d2/6*(s*(rsh2 + 2*s2)*c_exp - rsh*(rsh2 + 3*s2)*c_erfc);
666     real pot3 = sc_fac*
667         der->md3/24*((rsh2*rsh2 + 6*rsh2*s2 + 3*s2*s2)*c_erfc - rsh*s*(rsh2 + 5*s2)*c_exp);
668
669     return pot1 + pot2 + pot3;
670 }
671
672 static real energyDrift(const verletbuf_atomtype_t *att, int natt,
673                         const gmx_ffparams_t *ffp,
674                         real kT_fac,
675                         const pot_derivatives_t *ljDisp,
676                         const pot_derivatives_t *ljRep,
677                         const pot_derivatives_t *elec,
678                         real rlj, real rcoulomb,
679                         real rlist, real boxvol)
680 {
681     double drift_tot = 0;
682
683     if (kT_fac == 0)
684     {
685         /* No atom displacements: no drift, avoid division by 0 */
686         return drift_tot;
687     }
688
689     // Here add up the contribution of all atom pairs in the system to
690     // (estimated) energy drift by looping over all atom type pairs.
691     for (int i = 0; i < natt; i++)
692     {
693         // Get the thermal displacement variance for the i-atom type
694         const atom_nonbonded_kinetic_prop_t *prop_i = &att[i].prop;
695         real                                 s2i_2d, s2i_3d;
696         get_atom_sigma2(kT_fac, prop_i, &s2i_2d, &s2i_3d);
697
698         for (int j = i; j < natt; j++)
699         {
700             // Get the thermal displacement variance for the j-atom type
701             const atom_nonbonded_kinetic_prop_t *prop_j = &att[j].prop;
702             real                                 s2j_2d, s2j_3d;
703             get_atom_sigma2(kT_fac, prop_j, &s2j_2d, &s2j_3d);
704
705             /* Add up the up to four independent variances */
706             real s2 = s2i_2d + s2i_3d + s2j_2d + s2j_3d;
707
708             // Set -V', V'' and -V''' at the cut-off for LJ */
709             real              c6  = ffp->iparams[prop_i->type*ffp->atnr + prop_j->type].lj.c6;
710             real              c12 = ffp->iparams[prop_i->type*ffp->atnr + prop_j->type].lj.c12;
711             pot_derivatives_t lj;
712             lj.md1 = c6*ljDisp->md1 + c12*ljRep->md1;
713             lj.d2  = c6*ljDisp->d2  + c12*ljRep->d2;
714             lj.md3 = c6*ljDisp->md3 + c12*ljRep->md3;
715
716             real pot_lj = energyDriftAtomPair(prop_i->bConstr, prop_j->bConstr,
717                                               s2, s2i_2d, s2j_2d,
718                                               rlist - rlj,
719                                               &lj);
720
721             // Set -V' and V'' at the cut-off for Coulomb
722             pot_derivatives_t elec_qq;
723             elec_qq.md1 = elec->md1*prop_i->q*prop_j->q;
724             elec_qq.d2  = elec->d2 *prop_i->q*prop_j->q;
725             elec_qq.md3 = 0;
726
727             real pot_q  = energyDriftAtomPair(prop_i->bConstr, prop_j->bConstr,
728                                               s2, s2i_2d, s2j_2d,
729                                               rlist - rcoulomb,
730                                               &elec_qq);
731
732             // Note that attractive and repulsive potentials for individual
733             // pairs can partially cancel.
734             real pot = pot_lj + pot_q;
735
736             /* Multiply by the number of atom pairs */
737             if (j == i)
738             {
739                 pot *= (double)att[i].n*(att[i].n - 1)/2;
740             }
741             else
742             {
743                 pot *= (double)att[i].n*att[j].n;
744             }
745             /* We need the line density to get the energy drift of the system.
746              * The effective average r^2 is close to (rlist+sigma)^2.
747              */
748             pot *= 4*M_PI*gmx::square(rlist + std::sqrt(s2))/boxvol;
749
750             /* Add the unsigned drift to avoid cancellation of errors */
751             drift_tot += std::abs(pot);
752         }
753     }
754
755     return drift_tot;
756 }
757
758 static real surface_frac(int cluster_size, real particle_distance, real rlist)
759 {
760     real d, area_rel;
761
762     if (rlist < 0.5*particle_distance)
763     {
764         /* We have non overlapping spheres */
765         return 1.0;
766     }
767
768     /* Half the inter-particle distance relative to rlist */
769     d = 0.5*particle_distance/rlist;
770
771     /* Determine the area of the surface at distance rlist to the closest
772      * particle, relative to surface of a sphere of radius rlist.
773      * The formulas below assume close to cubic cells for the pair search grid,
774      * which the pair search code tries to achieve.
775      * Note that in practice particle distances will not be delta distributed,
776      * but have some spread, often involving shorter distances,
777      * as e.g. O-H bonds in a water molecule. Thus the estimates below will
778      * usually be slightly too high and thus conservative.
779      */
780     switch (cluster_size)
781     {
782         case 1:
783             /* One particle: trivial */
784             area_rel = 1.0;
785             break;
786         case 2:
787             /* Two particles: two spheres at fractional distance 2*a */
788             area_rel = 1.0 + d;
789             break;
790         case 4:
791             /* We assume a perfect, symmetric tetrahedron geometry.
792              * The surface around a tetrahedron is too complex for a full
793              * analytical solution, so we use a Taylor expansion.
794              */
795             area_rel = (1.0 + 1/M_PI*(6*std::acos(1/std::sqrt(3))*d +
796                                       std::sqrt(3)*d*d*(1.0 +
797                                                         5.0/18.0*d*d +
798                                                         7.0/45.0*d*d*d*d +
799                                                         83.0/756.0*d*d*d*d*d*d)));
800             break;
801         default:
802             gmx_incons("surface_frac called with unsupported cluster_size");
803             area_rel = 1.0;
804     }
805
806     return area_rel/cluster_size;
807 }
808
809 /* Returns the negative of the third derivative of a potential r^-p
810  * with a force-switch function, evaluated at the cut-off rc.
811  */
812 static real md3_force_switch(real p, real rswitch, real rc)
813 {
814     /* The switched force function is:
815      * p*r^-(p+1) + a*(r - rswitch)^2 + b*(r - rswitch)^3
816      */
817     real a, b;
818     real md3_pot, md3_sw;
819
820     a = -((p + 4)*rc - (p + 1)*rswitch)/(pow(rc, p+2)*gmx::square(rc-rswitch));
821     b =  ((p + 3)*rc - (p + 1)*rswitch)/(pow(rc, p+2)*gmx::power3(rc-rswitch));
822
823     md3_pot = (p + 2)*(p + 1)*p*pow(rc, p+3);
824     md3_sw  = 2*a + 6*b*(rc - rswitch);
825
826     return md3_pot + md3_sw;
827 }
828
829 /* Returns the maximum reference temperature over all coupled groups */
830 static real maxReferenceTemperature(const t_inputrec &ir)
831 {
832     if (EI_MD(ir.eI) && ir.etc == etcNO)
833     {
834         /* This case should be handled outside calc_verlet_buffer_size */
835         gmx_incons("calc_verlet_buffer_size called with an NVE ensemble and reference_temperature < 0");
836     }
837
838     real maxTemperature = 0;
839     for (int i = 0; i < ir.opts.ngtc; i++)
840     {
841         if (ir.opts.tau_t[i] >= 0)
842         {
843             maxTemperature = std::max(maxTemperature, ir.opts.ref_t[i]);
844         }
845     }
846
847     return maxTemperature;
848 }
849
850 /* Returns the variance of the atomic displacement over timePeriod.
851  *
852  * Note: When not using BD with a non-mass dependendent friction coefficient,
853  *       the return value still needs to be divided by the particle mass.
854  */
855 static real displacementVariance(const t_inputrec &ir,
856                                  real              temperature,
857                                  real              timePeriod)
858 {
859     real kT_fac;
860
861     if (ir.eI == eiBD)
862     {
863         /* Get the displacement distribution from the random component only.
864          * With accurate integration the systematic (force) displacement
865          * should be negligible (unless nstlist is extremely large, which
866          * you wouldn't do anyhow).
867          */
868         kT_fac = 2*BOLTZ*temperature*timePeriod;
869         if (ir.bd_fric > 0)
870         {
871             /* This is directly sigma^2 of the displacement */
872             kT_fac /= ir.bd_fric;
873         }
874         else
875         {
876             /* Per group tau_t is not implemented yet, use the maximum */
877             real tau_t = ir.opts.tau_t[0];
878             for (int i = 1; i < ir.opts.ngtc; i++)
879             {
880                 tau_t = std::max(tau_t, ir.opts.tau_t[i]);
881             }
882
883             kT_fac *= tau_t;
884             /* This kT_fac needs to be divided by the mass to get sigma^2 */
885         }
886     }
887     else
888     {
889         kT_fac = BOLTZ*temperature*gmx::square(timePeriod);
890     }
891
892     return kT_fac;
893 }
894
895 /* Returns the largest sigma of the Gaussian displacement over all particle
896  * types. This ignores constraints, so is an overestimate.
897  */
898 static real maxSigma(real                        kT_fac,
899                      int                         natt,
900                      const verletbuf_atomtype_t *att)
901 {
902     assert(att);
903     real smallestMass = att[0].prop.mass;
904     for (int i = 1; i < natt; i++)
905     {
906         smallestMass = std::min(smallestMass, att[i].prop.mass);
907     }
908
909     return 2*std::sqrt(kT_fac/smallestMass);
910 }
911
912 void calc_verlet_buffer_size(const gmx_mtop_t *mtop, real boxvol,
913                              const t_inputrec *ir,
914                              int               nstlist,
915                              int               list_lifetime,
916                              real reference_temperature,
917                              const VerletbufListSetup *list_setup,
918                              int *n_nonlin_vsite,
919                              real *rlist)
920 {
921     double                resolution;
922     char                 *env;
923
924     real                  particle_distance;
925     real                  nb_clust_frac_pairs_not_in_list_at_cutoff;
926
927     verletbuf_atomtype_t *att  = nullptr;
928     int                   natt = -1;
929     real                  elfac;
930     int                   ib0, ib1, ib;
931     real                  rb, rl;
932     real                  drift;
933
934     if (!EI_DYNAMICS(ir->eI))
935     {
936         gmx_incons("Can only determine the Verlet buffer size for integrators that perform dynamics");
937     }
938     if (ir->verletbuf_tol <= 0)
939     {
940         gmx_incons("The Verlet buffer tolerance needs to be larger than zero");
941     }
942
943     if (reference_temperature < 0)
944     {
945         /* We use the maximum temperature with multiple T-coupl groups.
946          * We could use a per particle temperature, but since particles
947          * interact, this might underestimate the buffer size.
948          */
949         reference_temperature = maxReferenceTemperature(*ir);
950     }
951
952     /* Resolution of the buffer size */
953     resolution = 0.001;
954
955     env = getenv("GMX_VERLET_BUFFER_RES");
956     if (env != nullptr)
957     {
958         sscanf(env, "%lf", &resolution);
959     }
960
961     /* In an atom wise pair-list there would be no pairs in the list
962      * beyond the pair-list cut-off.
963      * However, we use a pair-list of groups vs groups of atoms.
964      * For groups of 4 atoms, the parallelism of SSE instructions, only
965      * 10% of the atoms pairs are not in the list just beyond the cut-off.
966      * As this percentage increases slowly compared to the decrease of the
967      * Gaussian displacement distribution over this range, we can simply
968      * reduce the drift by this fraction.
969      * For larger groups, e.g. of 8 atoms, this fraction will be lower,
970      * so then buffer size will be on the conservative (large) side.
971      *
972      * Note that the formulas used here do not take into account
973      * cancellation of errors which could occur by missing both
974      * attractive and repulsive interactions.
975      *
976      * The only major assumption is homogeneous particle distribution.
977      * For an inhomogeneous system, such as a liquid-vapor system,
978      * the buffer will be underestimated. The actual energy drift
979      * will be higher by the factor: local/homogeneous particle density.
980      *
981      * The results of this estimate have been checked againt simulations.
982      * In most cases the real drift differs by less than a factor 2.
983      */
984
985     /* Worst case assumption: HCP packing of particles gives largest distance */
986     particle_distance = std::cbrt(boxvol*std::sqrt(2)/mtop->natoms);
987
988     /* TODO: Obtain masses through (future) integrator functionality
989      *       to avoid scattering the code with (or forgetting) checks.
990      */
991     const bool setMassesToOne = (ir->eI == eiBD && ir->bd_fric > 0);
992     get_verlet_buffer_atomtypes(mtop, setMassesToOne, &att, &natt, n_nonlin_vsite);
993     assert(att != NULL && natt >= 0);
994
995     if (debug)
996     {
997         fprintf(debug, "particle distance assuming HCP packing: %f nm\n",
998                 particle_distance);
999         fprintf(debug, "energy drift atom types: %d\n", natt);
1000     }
1001
1002     pot_derivatives_t ljDisp = { 0, 0, 0 };
1003     pot_derivatives_t ljRep  = { 0, 0, 0 };
1004     real              repPow = mtop->ffparams.reppow;
1005
1006     if (ir->vdwtype == evdwCUT)
1007     {
1008         real sw_range, md3_pswf;
1009
1010         switch (ir->vdw_modifier)
1011         {
1012             case eintmodNONE:
1013             case eintmodPOTSHIFT:
1014                 /* -dV/dr of -r^-6 and r^-reppow */
1015                 ljDisp.md1 =     -6*std::pow(ir->rvdw, -7.0);
1016                 ljRep.md1  = repPow*std::pow(ir->rvdw, -(repPow + 1));
1017                 /* The contribution of the higher derivatives is negligible */
1018                 break;
1019             case eintmodFORCESWITCH:
1020                 /* At the cut-off: V=V'=V''=0, so we use only V''' */
1021                 ljDisp.md3 = -md3_force_switch(6.0,    ir->rvdw_switch, ir->rvdw);
1022                 ljRep.md3  =  md3_force_switch(repPow, ir->rvdw_switch, ir->rvdw);
1023                 break;
1024             case eintmodPOTSWITCH:
1025                 /* At the cut-off: V=V'=V''=0.
1026                  * V''' is given by the original potential times
1027                  * the third derivative of the switch function.
1028                  */
1029                 sw_range   = ir->rvdw - ir->rvdw_switch;
1030                 md3_pswf   = 60.0/gmx::power3(sw_range);
1031
1032                 ljDisp.md3 = -std::pow(ir->rvdw, -6.0   )*md3_pswf;
1033                 ljRep.md3  =  std::pow(ir->rvdw, -repPow)*md3_pswf;
1034                 break;
1035             default:
1036                 gmx_incons("Unimplemented VdW modifier");
1037         }
1038     }
1039     else if (EVDW_PME(ir->vdwtype))
1040     {
1041         real b     = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
1042         real r     = ir->rvdw;
1043         real br    = b*r;
1044         real br2   = br*br;
1045         real br4   = br2*br2;
1046         real br6   = br4*br2;
1047         // -dV/dr of g(br)*r^-6 [where g(x) = exp(-x^2)(1+x^2+x^4/2),
1048         // see LJ-PME equations in manual] and r^-reppow
1049         ljDisp.md1 = -std::exp(-br2)*(br6 + 3.0*br4 + 6.0*br2 + 6.0)*std::pow(r, -7.0);
1050         ljRep.md1  = repPow*pow(r, -(repPow + 1));
1051         // The contribution of the higher derivatives is negligible
1052     }
1053     else
1054     {
1055         gmx_fatal(FARGS, "Energy drift calculation is only implemented for plain cut-off Lennard-Jones interactions");
1056     }
1057
1058     elfac = ONE_4PI_EPS0/ir->epsilon_r;
1059
1060     // Determine the 1st and 2nd derivative for the electostatics
1061     pot_derivatives_t elec = { 0, 0, 0 };
1062
1063     if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
1064     {
1065         real eps_rf, k_rf;
1066
1067         if (ir->coulombtype == eelCUT)
1068         {
1069             eps_rf = 1;
1070             k_rf   = 0;
1071         }
1072         else
1073         {
1074             eps_rf = ir->epsilon_rf/ir->epsilon_r;
1075             if (eps_rf != 0)
1076             {
1077                 k_rf = (eps_rf - ir->epsilon_r)/( gmx::power3(ir->rcoulomb) * (2*eps_rf + ir->epsilon_r) );
1078             }
1079             else
1080             {
1081                 /* epsilon_rf = infinity */
1082                 k_rf = 0.5/gmx::power3(ir->rcoulomb);
1083             }
1084         }
1085
1086         if (eps_rf > 0)
1087         {
1088             elec.md1 = elfac*(1.0/gmx::square(ir->rcoulomb) - 2*k_rf*ir->rcoulomb);
1089         }
1090         elec.d2      = elfac*(2.0/gmx::power3(ir->rcoulomb) + 2*k_rf);
1091     }
1092     else if (EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD)
1093     {
1094         real b, rc, br;
1095
1096         b        = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
1097         rc       = ir->rcoulomb;
1098         br       = b*rc;
1099         elec.md1 = elfac*(b*std::exp(-br*br)*M_2_SQRTPI/rc + std::erfc(br)/(rc*rc));
1100         elec.d2  = elfac/(rc*rc)*(2*b*(1 + br*br)*std::exp(-br*br)*M_2_SQRTPI + 2*std::erfc(br)/rc);
1101     }
1102     else
1103     {
1104         gmx_fatal(FARGS, "Energy drift calculation is only implemented for Reaction-Field and Ewald electrostatics");
1105     }
1106
1107     /* Determine the variance of the atomic displacement
1108      * over list_lifetime steps: kT_fac
1109      * For inertial dynamics (not Brownian dynamics) the mass factor
1110      * is not included in kT_fac, it is added later.
1111      */
1112     const real kT_fac = displacementVariance(*ir, reference_temperature,
1113                                              list_lifetime*ir->delta_t);
1114
1115     if (debug)
1116     {
1117         fprintf(debug, "Derivatives of non-bonded potentials at the cut-off:\n");
1118         fprintf(debug, "LJ disp. -V' %9.2e V'' %9.2e -V''' %9.2e\n", ljDisp.md1, ljDisp.d2, ljDisp.md3);
1119         fprintf(debug, "LJ rep.  -V' %9.2e V'' %9.2e -V''' %9.2e\n", ljRep.md1, ljRep.d2, ljRep.md3);
1120         fprintf(debug, "Electro. -V' %9.2e V'' %9.2e\n", elec.md1, elec.d2);
1121         fprintf(debug, "sqrt(kT_fac) %f\n", std::sqrt(kT_fac));
1122     }
1123
1124     /* Search using bisection */
1125     ib0 = -1;
1126     /* The drift will be neglible at 5 times the max sigma */
1127     ib1 = (int)(5*maxSigma(kT_fac, natt, att)/resolution) + 1;
1128     while (ib1 - ib0 > 1)
1129     {
1130         ib = (ib0 + ib1)/2;
1131         rb = ib*resolution;
1132         rl = std::max(ir->rvdw, ir->rcoulomb) + rb;
1133
1134         /* Calculate the average energy drift at the last step
1135          * of the nstlist steps at which the pair-list is used.
1136          */
1137         drift = energyDrift(att, natt, &mtop->ffparams,
1138                             kT_fac,
1139                             &ljDisp, &ljRep, &elec,
1140                             ir->rvdw, ir->rcoulomb,
1141                             rl, boxvol);
1142
1143         /* Correct for the fact that we are using a Ni x Nj particle pair list
1144          * and not a 1 x 1 particle pair list. This reduces the drift.
1145          */
1146         /* We don't have a formula for 8 (yet), use 4 which is conservative */
1147         nb_clust_frac_pairs_not_in_list_at_cutoff =
1148             surface_frac(std::min(list_setup->cluster_size_i, 4),
1149                          particle_distance, rl)*
1150             surface_frac(std::min(list_setup->cluster_size_j, 4),
1151                          particle_distance, rl);
1152         drift *= nb_clust_frac_pairs_not_in_list_at_cutoff;
1153
1154         /* Convert the drift to drift per unit time per atom */
1155         drift /= nstlist*ir->delta_t*mtop->natoms;
1156
1157         if (debug)
1158         {
1159             fprintf(debug, "ib %3d %3d %3d rb %.3f %dx%d fac %.3f drift %.1e\n",
1160                     ib0, ib, ib1, rb,
1161                     list_setup->cluster_size_i, list_setup->cluster_size_j,
1162                     nb_clust_frac_pairs_not_in_list_at_cutoff,
1163                     drift);
1164         }
1165
1166         if (std::abs(drift) > ir->verletbuf_tol)
1167         {
1168             ib0 = ib;
1169         }
1170         else
1171         {
1172             ib1 = ib;
1173         }
1174     }
1175
1176     sfree(att);
1177
1178     *rlist = std::max(ir->rvdw, ir->rcoulomb) + ib1*resolution;
1179 }