f0a48401fc798d746bae2b7a19b8916f42663746
[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 <cassert>
40 #include <cmath>
41 #include <cstdlib>
42
43 #include <algorithm>
44
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"
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                         default:
290                             /* Use the mass of the lightest constructing atom.
291                              * This is an approximation.
292                              * If the distance of the virtual site to the
293                              * constructing atom is less than all distances
294                              * between constructing atoms, this is a safe
295                              * over-estimate of the displacement of the vsite.
296                              * This condition holds for all H mass replacement
297                              * vsite constructions, except for SP2/3 groups.
298                              * In SP3 groups one H will have a F_VSITE3
299                              * construction, so even there the total drift
300                              * estimate shouldn't be far off.
301                              */
302                             vsite_m[a1] = cam[1];
303                             for (j = 2; j < maxj; j++)
304                             {
305                                 vsite_m[a1] = std::min(vsite_m[a1], cam[j]);
306                             }
307                             (*n_nonlin_vsite)++;
308                             break;
309                     }
310                     sfree(cam);
311                 }
312                 else
313                 {
314                     int j;
315
316                     /* Exact */
317                     inv_mass = 0;
318                     for (j = 0; j < 3*ffparams->iparams[il->iatoms[i]].vsiten.n; j += 3)
319                     {
320                         int aj = il->iatoms[i + j + 2];
321                         coeff  = ffparams->iparams[il->iatoms[i+j]].vsiten.a;
322                         if (moltype->atoms.atom[aj].ptype == eptVSite)
323                         {
324                             m_aj = vsite_m[aj];
325                         }
326                         else
327                         {
328                             m_aj = moltype->atoms.atom[aj].m;
329                         }
330                         if (m_aj <= 0)
331                         {
332                             gmx_incons("The mass of a vsiten constructing atom is <= 0");
333                         }
334                         inv_mass += coeff*coeff/m_aj;
335                     }
336                     vsite_m[a1] = 1/inv_mass;
337                     /* Correct for loop increment of i */
338                     i += j - 1 - NRAL(ft);
339                 }
340                 if (gmx_debug_at)
341                 {
342                     fprintf(debug, "atom %4d %-20s mass %6.3f\n",
343                             a1, interaction_function[ft].longname, vsite_m[a1]);
344                 }
345             }
346         }
347     }
348 }
349
350 static void get_verlet_buffer_atomtypes(const gmx_mtop_t      *mtop,
351                                         bool                   setMassesToOne,
352                                         verletbuf_atomtype_t **att_p,
353                                         int                   *natt_p,
354                                         int                   *n_nonlin_vsite)
355 {
356     verletbuf_atomtype_t          *att;
357     int                            natt;
358     int                            ft, i, a1, a2, a3, a;
359     const t_ilist                 *il;
360     const t_iparams               *ip;
361     atom_nonbonded_kinetic_prop_t *prop;
362     real                          *vsite_m;
363     int                            n_nonlin_vsite_mol;
364
365     att  = nullptr;
366     natt = 0;
367
368     if (n_nonlin_vsite != nullptr)
369     {
370         *n_nonlin_vsite = 0;
371     }
372
373     for (const gmx_molblock_t &molblock : mtop->molblock)
374     {
375         int                  nmol    = molblock.nmol;
376         const gmx_moltype_t &moltype = mtop->moltype[molblock.type];
377         const t_atoms       *atoms   = &moltype.atoms;
378
379         /* Check for constraints, as they affect the kinetic energy.
380          * For virtual sites we need the masses and geometry of
381          * the constructing atoms to determine their velocity distribution.
382          */
383         snew(prop, atoms->nr);
384         snew(vsite_m, atoms->nr);
385
386         for (ft = F_CONSTR; ft <= F_CONSTRNC; ft++)
387         {
388             il = &moltype.ilist[ft];
389
390             for (i = 0; i < il->nr; i += 1+NRAL(ft))
391             {
392                 ip         = &mtop->ffparams.iparams[il->iatoms[i]];
393                 a1         = il->iatoms[i+1];
394                 a2         = il->iatoms[i+2];
395                 real mass1 = getMass(*atoms, a1, setMassesToOne);
396                 real mass2 = getMass(*atoms, a2, setMassesToOne);
397                 if (mass2 > prop[a1].con_mass)
398                 {
399                     prop[a1].con_mass = mass2;
400                     prop[a1].con_len  = ip->constr.dA;
401                 }
402                 if (mass1 > prop[a2].con_mass)
403                 {
404                     prop[a2].con_mass = mass1;
405                     prop[a2].con_len  = ip->constr.dA;
406                 }
407             }
408         }
409
410         il = &moltype.ilist[F_SETTLE];
411
412         for (i = 0; i < il->nr; i += 1+NRAL(F_SETTLE))
413         {
414             ip         = &mtop->ffparams.iparams[il->iatoms[i]];
415             a1         = il->iatoms[i+1];
416             a2         = il->iatoms[i+2];
417             a3         = il->iatoms[i+3];
418             /* Usually the mass of a1 (usually oxygen) is larger than a2/a3.
419              * If this is not the case, we overestimate the displacement,
420              * which leads to a larger buffer (ok since this is an exotic case).
421              */
422             prop[a1].con_mass = getMass(*atoms, a2, setMassesToOne);
423             prop[a1].con_len  = ip->settle.doh;
424
425             prop[a2].con_mass = getMass(*atoms, a1, setMassesToOne);
426             prop[a2].con_len  = ip->settle.doh;
427
428             prop[a3].con_mass = getMass(*atoms, a1, setMassesToOne);
429             prop[a3].con_len  = ip->settle.doh;
430         }
431
432         get_vsite_masses(&moltype,
433                          &mtop->ffparams,
434                          setMassesToOne,
435                          vsite_m,
436                          &n_nonlin_vsite_mol);
437         if (n_nonlin_vsite != nullptr)
438         {
439             *n_nonlin_vsite += nmol*n_nonlin_vsite_mol;
440         }
441
442         for (a = 0; a < atoms->nr; a++)
443         {
444             if (atoms->atom[a].ptype == eptVSite)
445             {
446                 prop[a].mass = vsite_m[a];
447             }
448             else
449             {
450                 prop[a].mass = getMass(*atoms, a, setMassesToOne);
451             }
452             prop[a].type     = atoms->atom[a].type;
453             prop[a].q        = atoms->atom[a].q;
454             /* We consider an atom constrained, #DOF=2, when it is
455              * connected with constraints to (at least one) atom with
456              * a mass of more than 0.4x its own mass. This is not a critical
457              * parameter, since with roughly equal masses the unconstrained
458              * and constrained displacement will not differ much (and both
459              * overestimate the displacement).
460              */
461             prop[a].bConstr = (prop[a].con_mass > 0.4*prop[a].mass);
462
463             add_at(&att, &natt, &prop[a], nmol);
464         }
465
466         sfree(vsite_m);
467         sfree(prop);
468     }
469
470     if (gmx_debug_at)
471     {
472         for (a = 0; a < natt; a++)
473         {
474             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",
475                     a, att[a].prop.mass, att[a].prop.type, att[a].prop.q,
476                     gmx::boolToString(att[a].prop.bConstr), att[a].prop.con_mass, att[a].prop.con_len,
477                     att[a].n);
478         }
479     }
480
481     *att_p  = att;
482     *natt_p = natt;
483 }
484
485 /* This function computes two components of the estimate of the variance
486  * in the displacement of one atom in a system of two constrained atoms.
487  * Returns in sigma2_2d the variance due to rotation of the constrained
488  * atom around the atom to which it constrained.
489  * Returns in sigma2_3d the variance due to displacement of the COM
490  * of the whole system of the two constrained atoms.
491  *
492  * Note that we only take a single constraint (the one to the heaviest atom)
493  * into account. If an atom has multiple constraints, this will result in
494  * an overestimate of the displacement, which gives a larger drift and buffer.
495  */
496 void constrained_atom_sigma2(real                                 kT_fac,
497                              const atom_nonbonded_kinetic_prop_t *prop,
498                              real                                *sigma2_2d,
499                              real                                *sigma2_3d)
500 {
501     /* Here we decompose the motion of a constrained atom into two
502      * components: rotation around the COM and translation of the COM.
503      */
504
505     /* Determine the variance of the arc length for the two rotational DOFs */
506     real massFraction = prop->con_mass/(prop->mass + prop->con_mass);
507     real sigma2_rot   = kT_fac*massFraction/prop->mass;
508
509     /* The distance from the atom to the COM, i.e. the rotational arm */
510     real comDistance  = prop->con_len*massFraction;
511
512     /* The variance relative to the arm */
513     real sigma2_rel   = sigma2_rot/gmx::square(comDistance);
514
515     /* For sigma2_rel << 1 we don't notice the rotational effect and
516      * we have a normal, Gaussian displacement distribution.
517      * For larger sigma2_rel the displacement is much less, in fact it can
518      * not exceed 2*comDistance. We can calculate MSD/arm^2 as:
519      *   integral_x=0-inf distance2(x) x/sigma2_rel exp(-x^2/(2 sigma2_rel)) dx
520      * where x is angular displacement and distance2(x) is the distance^2
521      * between points at angle 0 and x:
522      *   distance2(x) = (sin(x) - sin(0))^2 + (cos(x) - cos(0))^2
523      * The limiting value of this MSD is 2, which is also the value for
524      * a uniform rotation distribution that would be reached at long time.
525      * The maximum is 2.5695 at sigma2_rel = 4.5119.
526      * We approximate this integral with a rational polynomial with
527      * coefficients from a Taylor expansion. This approximation is an
528      * overestimate for all values of sigma2_rel. Its maximum value
529      * of 2.6491 is reached at sigma2_rel = sqrt(45/2) = 4.7434.
530      * We keep the approximation constant after that.
531      * We use this approximate MSD as the variance for a Gaussian distribution.
532      *
533      * NOTE: For any sensible buffer tolerance this will result in a (large)
534      * overestimate of the buffer size, since the Gaussian has a long tail,
535      * whereas the actual distribution can not reach values larger than 2.
536      */
537     /* Coeffients obtained from a Taylor expansion */
538     const real a = 1.0/3.0;
539     const real b = 2.0/45.0;
540
541     /* Our approximation is constant after sigma2_rel = 1/sqrt(b) */
542     sigma2_rel   = std::min(sigma2_rel, 1/std::sqrt(b));
543
544     /* Compute the approximate sigma^2 for 2D motion due to the rotation */
545     *sigma2_2d   = gmx::square(comDistance)*
546         sigma2_rel/(1 + a*sigma2_rel + b*gmx::square(sigma2_rel));
547
548     /* The constrained atom also moves (in 3D) with the COM of both atoms */
549     *sigma2_3d   = kT_fac/(prop->mass + prop->con_mass);
550 }
551
552 static void get_atom_sigma2(real                                 kT_fac,
553                             const atom_nonbonded_kinetic_prop_t *prop,
554                             real                                *sigma2_2d,
555                             real                                *sigma2_3d)
556 {
557     if (prop->bConstr)
558     {
559         /* Complicated constraint calculation in a separate function */
560         constrained_atom_sigma2(kT_fac, prop, sigma2_2d, sigma2_3d);
561     }
562     else
563     {
564         /* Unconstrained atom: trivial */
565         *sigma2_2d = 0;
566         *sigma2_3d = kT_fac/prop->mass;
567     }
568 }
569
570 static void approx_2dof(real s2, real x, real *shift, real *scale)
571 {
572     /* A particle with 1 DOF constrained has 2 DOFs instead of 3.
573      * This code is also used for particles with multiple constraints,
574      * in which case we overestimate the displacement.
575      * The 2DOF distribution is sqrt(pi/2)*erfc(r/(sqrt(2)*s))/(2*s).
576      * We approximate this with scale*Gaussian(s,r+shift),
577      * by matching the distribution value and derivative at x.
578      * This is a tight overestimate for all r>=0 at any s and x.
579      */
580     real ex, er;
581
582     ex = std::exp(-x*x/(2*s2));
583     er = std::erfc(x/std::sqrt(2*s2));
584
585     *shift = -x + std::sqrt(2*s2/M_PI)*ex/er;
586     *scale = 0.5*M_PI*std::exp(ex*ex/(M_PI*er*er))*er;
587 }
588
589 // Returns an (over)estimate of the energy drift for a single atom pair,
590 // given the kinetic properties, displacement variances and list buffer.
591 static real energyDriftAtomPair(bool isConstrained_i,
592                                 bool isConstrained_j,
593                                 real s2, real s2i_2d, real s2j_2d,
594                                 real r_buffer,
595                                 const pot_derivatives_t *der)
596 {
597     // For relatively small arguments erfc() is so small that if will be 0.0
598     // when stored in a float. We set an argument limit of 8 (Erfc(8)=1e-29),
599     // such that we can divide by erfc and have some space left for arithmetic.
600     const real erfc_arg_max = 8.0;
601
602     real       rsh    = r_buffer;
603     real       sc_fac = 1.0;
604
605     real       c_exp, c_erfc;
606
607     if (rsh*rsh > 2*s2*erfc_arg_max*erfc_arg_max)
608     {
609         // Below we calculate c_erfc = 0.5*erfc(rsh/sqrt(2*s2))
610         // When rsh/sqrt(2*s2) increases, this erfc will be the first
611         // result that underflows and becomes 0.0. To avoid this,
612         // we set c_exp=0 and c_erfc=0 for large arguments.
613         // This also avoids NaN in approx_2dof().
614         // In any relevant case this has no effect on the results,
615         // since c_exp < 6e-29, so the displacement is completely
616         // negligible for such atom pairs (and an overestimate).
617         // In nearly all use cases, there will be other atom pairs
618         // that contribute much more to the total, so zeroing
619         // this particular contribution has no effect at all.
620         c_exp  = 0;
621         c_erfc = 0;
622     }
623     else
624     {
625         /* For constraints: adapt r and scaling for the Gaussian */
626         if (isConstrained_i)
627         {
628             real sh, sc;
629
630             approx_2dof(s2i_2d, r_buffer*s2i_2d/s2, &sh, &sc);
631             rsh    += sh;
632             sc_fac *= sc;
633         }
634         if (isConstrained_j)
635         {
636             real sh, sc;
637
638             approx_2dof(s2j_2d, r_buffer*s2j_2d/s2, &sh, &sc);
639             rsh    += sh;
640             sc_fac *= sc;
641         }
642
643         /* Exact contribution of an atom pair with Gaussian displacement
644          * with sigma s to the energy drift for a potential with
645          * derivative -md and second derivative dd at the cut-off.
646          * The only catch is that for potentials that change sign
647          * near the cut-off there could be an unlucky compensation
648          * of positive and negative energy drift.
649          * Such potentials are extremely rare though.
650          *
651          * Note that pot has unit energy*length, as the linear
652          * atom density still needs to be put in.
653          */
654         c_exp  = std::exp(-rsh*rsh/(2*s2))/std::sqrt(2*M_PI);
655         c_erfc = 0.5*std::erfc(rsh/(std::sqrt(2*s2)));
656     }
657     real s    = std::sqrt(s2);
658     real rsh2 = rsh*rsh;
659
660     real pot1 = sc_fac*
661         der->md1/2*((rsh2 + s2)*c_erfc - rsh*s*c_exp);
662     real pot2 = sc_fac*
663         der->d2/6*(s*(rsh2 + 2*s2)*c_exp - rsh*(rsh2 + 3*s2)*c_erfc);
664     real pot3 = sc_fac*
665         der->md3/24*((rsh2*rsh2 + 6*rsh2*s2 + 3*s2*s2)*c_erfc - rsh*s*(rsh2 + 5*s2)*c_exp);
666
667     return pot1 + pot2 + pot3;
668 }
669
670 static real energyDrift(const verletbuf_atomtype_t *att, int natt,
671                         const gmx_ffparams_t *ffp,
672                         real kT_fac,
673                         const pot_derivatives_t *ljDisp,
674                         const pot_derivatives_t *ljRep,
675                         const pot_derivatives_t *elec,
676                         real rlj, real rcoulomb,
677                         real rlist, real boxvol)
678 {
679     double drift_tot = 0;
680
681     if (kT_fac == 0)
682     {
683         /* No atom displacements: no drift, avoid division by 0 */
684         return drift_tot;
685     }
686
687     // Here add up the contribution of all atom pairs in the system to
688     // (estimated) energy drift by looping over all atom type pairs.
689     for (int i = 0; i < natt; i++)
690     {
691         // Get the thermal displacement variance for the i-atom type
692         const atom_nonbonded_kinetic_prop_t *prop_i = &att[i].prop;
693         real                                 s2i_2d, s2i_3d;
694         get_atom_sigma2(kT_fac, prop_i, &s2i_2d, &s2i_3d);
695
696         for (int j = i; j < natt; j++)
697         {
698             // Get the thermal displacement variance for the j-atom type
699             const atom_nonbonded_kinetic_prop_t *prop_j = &att[j].prop;
700             real                                 s2j_2d, s2j_3d;
701             get_atom_sigma2(kT_fac, prop_j, &s2j_2d, &s2j_3d);
702
703             /* Add up the up to four independent variances */
704             real s2 = s2i_2d + s2i_3d + s2j_2d + s2j_3d;
705
706             // Set -V', V'' and -V''' at the cut-off for LJ */
707             real              c6  = ffp->iparams[prop_i->type*ffp->atnr + prop_j->type].lj.c6;
708             real              c12 = ffp->iparams[prop_i->type*ffp->atnr + prop_j->type].lj.c12;
709             pot_derivatives_t lj;
710             lj.md1 = c6*ljDisp->md1 + c12*ljRep->md1;
711             lj.d2  = c6*ljDisp->d2  + c12*ljRep->d2;
712             lj.md3 = c6*ljDisp->md3 + c12*ljRep->md3;
713
714             real pot_lj = energyDriftAtomPair(prop_i->bConstr, prop_j->bConstr,
715                                               s2, s2i_2d, s2j_2d,
716                                               rlist - rlj,
717                                               &lj);
718
719             // Set -V' and V'' at the cut-off for Coulomb
720             pot_derivatives_t elec_qq;
721             elec_qq.md1 = elec->md1*prop_i->q*prop_j->q;
722             elec_qq.d2  = elec->d2 *prop_i->q*prop_j->q;
723             elec_qq.md3 = 0;
724
725             real pot_q  = energyDriftAtomPair(prop_i->bConstr, prop_j->bConstr,
726                                               s2, s2i_2d, s2j_2d,
727                                               rlist - rcoulomb,
728                                               &elec_qq);
729
730             // Note that attractive and repulsive potentials for individual
731             // pairs can partially cancel.
732             real pot = pot_lj + pot_q;
733
734             /* Multiply by the number of atom pairs */
735             if (j == i)
736             {
737                 pot *= static_cast<double>(att[i].n)*(att[i].n - 1)/2;
738             }
739             else
740             {
741                 pot *= static_cast<double>(att[i].n)*att[j].n;
742             }
743             /* We need the line density to get the energy drift of the system.
744              * The effective average r^2 is close to (rlist+sigma)^2.
745              */
746             pot *= 4*M_PI*gmx::square(rlist + std::sqrt(s2))/boxvol;
747
748             /* Add the unsigned drift to avoid cancellation of errors */
749             drift_tot += std::abs(pot);
750         }
751     }
752
753     return drift_tot;
754 }
755
756 static real surface_frac(int cluster_size, real particle_distance, real rlist)
757 {
758     real d, area_rel;
759
760     if (rlist < 0.5*particle_distance)
761     {
762         /* We have non overlapping spheres */
763         return 1.0;
764     }
765
766     /* Half the inter-particle distance relative to rlist */
767     d = 0.5*particle_distance/rlist;
768
769     /* Determine the area of the surface at distance rlist to the closest
770      * particle, relative to surface of a sphere of radius rlist.
771      * The formulas below assume close to cubic cells for the pair search grid,
772      * which the pair search code tries to achieve.
773      * Note that in practice particle distances will not be delta distributed,
774      * but have some spread, often involving shorter distances,
775      * as e.g. O-H bonds in a water molecule. Thus the estimates below will
776      * usually be slightly too high and thus conservative.
777      */
778     switch (cluster_size)
779     {
780         case 1:
781             /* One particle: trivial */
782             area_rel = 1.0;
783             break;
784         case 2:
785             /* Two particles: two spheres at fractional distance 2*a */
786             area_rel = 1.0 + d;
787             break;
788         case 4:
789             /* We assume a perfect, symmetric tetrahedron geometry.
790              * The surface around a tetrahedron is too complex for a full
791              * analytical solution, so we use a Taylor expansion.
792              */
793             area_rel = (1.0 + 1/M_PI*(6*std::acos(1/std::sqrt(3))*d +
794                                       std::sqrt(3)*d*d*(1.0 +
795                                                         5.0/18.0*d*d +
796                                                         7.0/45.0*d*d*d*d +
797                                                         83.0/756.0*d*d*d*d*d*d)));
798             break;
799         default:
800             gmx_incons("surface_frac called with unsupported cluster_size");
801     }
802
803     return area_rel/cluster_size;
804 }
805
806 /* Returns the negative of the third derivative of a potential r^-p
807  * with a force-switch function, evaluated at the cut-off rc.
808  */
809 static real md3_force_switch(real p, real rswitch, real rc)
810 {
811     /* The switched force function is:
812      * p*r^-(p+1) + a*(r - rswitch)^2 + b*(r - rswitch)^3
813      */
814     real a, b;
815     real md3_pot, md3_sw;
816
817     a = -((p + 4)*rc - (p + 1)*rswitch)/(pow(rc, p+2)*gmx::square(rc-rswitch));
818     b =  ((p + 3)*rc - (p + 1)*rswitch)/(pow(rc, p+2)*gmx::power3(rc-rswitch));
819
820     md3_pot = (p + 2)*(p + 1)*p*pow(rc, p+3);
821     md3_sw  = 2*a + 6*b*(rc - rswitch);
822
823     return md3_pot + md3_sw;
824 }
825
826 /* Returns the maximum reference temperature over all coupled groups */
827 static real maxReferenceTemperature(const t_inputrec &ir)
828 {
829     if (EI_MD(ir.eI) && ir.etc == etcNO)
830     {
831         /* This case should be handled outside calc_verlet_buffer_size */
832         gmx_incons("calc_verlet_buffer_size called with an NVE ensemble and reference_temperature < 0");
833     }
834
835     real maxTemperature = 0;
836     for (int i = 0; i < ir.opts.ngtc; i++)
837     {
838         if (ir.opts.tau_t[i] >= 0)
839         {
840             maxTemperature = std::max(maxTemperature, ir.opts.ref_t[i]);
841         }
842     }
843
844     return maxTemperature;
845 }
846
847 /* Returns the variance of the atomic displacement over timePeriod.
848  *
849  * Note: When not using BD with a non-mass dependendent friction coefficient,
850  *       the return value still needs to be divided by the particle mass.
851  */
852 static real displacementVariance(const t_inputrec &ir,
853                                  real              temperature,
854                                  real              timePeriod)
855 {
856     real kT_fac;
857
858     if (ir.eI == eiBD)
859     {
860         /* Get the displacement distribution from the random component only.
861          * With accurate integration the systematic (force) displacement
862          * should be negligible (unless nstlist is extremely large, which
863          * you wouldn't do anyhow).
864          */
865         kT_fac = 2*BOLTZ*temperature*timePeriod;
866         if (ir.bd_fric > 0)
867         {
868             /* This is directly sigma^2 of the displacement */
869             kT_fac /= ir.bd_fric;
870         }
871         else
872         {
873             /* Per group tau_t is not implemented yet, use the maximum */
874             real tau_t = ir.opts.tau_t[0];
875             for (int i = 1; i < ir.opts.ngtc; i++)
876             {
877                 tau_t = std::max(tau_t, ir.opts.tau_t[i]);
878             }
879
880             kT_fac *= tau_t;
881             /* This kT_fac needs to be divided by the mass to get sigma^2 */
882         }
883     }
884     else
885     {
886         kT_fac = BOLTZ*temperature*gmx::square(timePeriod);
887     }
888
889     return kT_fac;
890 }
891
892 /* Returns the largest sigma of the Gaussian displacement over all particle
893  * types. This ignores constraints, so is an overestimate.
894  */
895 static real maxSigma(real                        kT_fac,
896                      int                         natt,
897                      const verletbuf_atomtype_t *att)
898 {
899     assert(att);
900     real smallestMass = att[0].prop.mass;
901     for (int i = 1; i < natt; i++)
902     {
903         smallestMass = std::min(smallestMass, att[i].prop.mass);
904     }
905
906     return 2*std::sqrt(kT_fac/smallestMass);
907 }
908
909 void calc_verlet_buffer_size(const gmx_mtop_t *mtop, real boxvol,
910                              const t_inputrec *ir,
911                              int               nstlist,
912                              int               list_lifetime,
913                              real reference_temperature,
914                              const VerletbufListSetup *list_setup,
915                              int *n_nonlin_vsite,
916                              real *rlist)
917 {
918     double                resolution;
919     char                 *env;
920
921     real                  particle_distance;
922     real                  nb_clust_frac_pairs_not_in_list_at_cutoff;
923
924     verletbuf_atomtype_t *att  = nullptr;
925     int                   natt = -1;
926     real                  elfac;
927     int                   ib0, ib1, ib;
928     real                  rb, rl;
929     real                  drift;
930
931     if (!EI_DYNAMICS(ir->eI))
932     {
933         gmx_incons("Can only determine the Verlet buffer size for integrators that perform dynamics");
934     }
935     if (ir->verletbuf_tol <= 0)
936     {
937         gmx_incons("The Verlet buffer tolerance needs to be larger than zero");
938     }
939
940     if (reference_temperature < 0)
941     {
942         /* We use the maximum temperature with multiple T-coupl groups.
943          * We could use a per particle temperature, but since particles
944          * interact, this might underestimate the buffer size.
945          */
946         reference_temperature = maxReferenceTemperature(*ir);
947     }
948
949     /* Resolution of the buffer size */
950     resolution = 0.001;
951
952     env = getenv("GMX_VERLET_BUFFER_RES");
953     if (env != nullptr)
954     {
955         sscanf(env, "%lf", &resolution);
956     }
957
958     /* In an atom wise pair-list there would be no pairs in the list
959      * beyond the pair-list cut-off.
960      * However, we use a pair-list of groups vs groups of atoms.
961      * For groups of 4 atoms, the parallelism of SSE instructions, only
962      * 10% of the atoms pairs are not in the list just beyond the cut-off.
963      * As this percentage increases slowly compared to the decrease of the
964      * Gaussian displacement distribution over this range, we can simply
965      * reduce the drift by this fraction.
966      * For larger groups, e.g. of 8 atoms, this fraction will be lower,
967      * so then buffer size will be on the conservative (large) side.
968      *
969      * Note that the formulas used here do not take into account
970      * cancellation of errors which could occur by missing both
971      * attractive and repulsive interactions.
972      *
973      * The only major assumption is homogeneous particle distribution.
974      * For an inhomogeneous system, such as a liquid-vapor system,
975      * the buffer will be underestimated. The actual energy drift
976      * will be higher by the factor: local/homogeneous particle density.
977      *
978      * The results of this estimate have been checked againt simulations.
979      * In most cases the real drift differs by less than a factor 2.
980      */
981
982     /* Worst case assumption: HCP packing of particles gives largest distance */
983     particle_distance = std::cbrt(boxvol*std::sqrt(2)/mtop->natoms);
984
985     /* TODO: Obtain masses through (future) integrator functionality
986      *       to avoid scattering the code with (or forgetting) checks.
987      */
988     const bool setMassesToOne = (ir->eI == eiBD && ir->bd_fric > 0);
989     get_verlet_buffer_atomtypes(mtop, setMassesToOne, &att, &natt, n_nonlin_vsite);
990     assert(att != nullptr && natt >= 0);
991
992     if (debug)
993     {
994         fprintf(debug, "particle distance assuming HCP packing: %f nm\n",
995                 particle_distance);
996         fprintf(debug, "energy drift atom types: %d\n", natt);
997     }
998
999     pot_derivatives_t ljDisp = { 0, 0, 0 };
1000     pot_derivatives_t ljRep  = { 0, 0, 0 };
1001     real              repPow = mtop->ffparams.reppow;
1002
1003     if (ir->vdwtype == evdwCUT)
1004     {
1005         real sw_range, md3_pswf;
1006
1007         switch (ir->vdw_modifier)
1008         {
1009             case eintmodNONE:
1010             case eintmodPOTSHIFT:
1011                 /* -dV/dr of -r^-6 and r^-reppow */
1012                 ljDisp.md1 =     -6*std::pow(ir->rvdw, -7.0);
1013                 ljRep.md1  = repPow*std::pow(ir->rvdw, -(repPow + 1));
1014                 /* The contribution of the higher derivatives is negligible */
1015                 break;
1016             case eintmodFORCESWITCH:
1017                 /* At the cut-off: V=V'=V''=0, so we use only V''' */
1018                 ljDisp.md3 = -md3_force_switch(6.0,    ir->rvdw_switch, ir->rvdw);
1019                 ljRep.md3  =  md3_force_switch(repPow, ir->rvdw_switch, ir->rvdw);
1020                 break;
1021             case eintmodPOTSWITCH:
1022                 /* At the cut-off: V=V'=V''=0.
1023                  * V''' is given by the original potential times
1024                  * the third derivative of the switch function.
1025                  */
1026                 sw_range   = ir->rvdw - ir->rvdw_switch;
1027                 md3_pswf   = 60.0/gmx::power3(sw_range);
1028
1029                 ljDisp.md3 = -std::pow(ir->rvdw, -6.0   )*md3_pswf;
1030                 ljRep.md3  =  std::pow(ir->rvdw, -repPow)*md3_pswf;
1031                 break;
1032             default:
1033                 gmx_incons("Unimplemented VdW modifier");
1034         }
1035     }
1036     else if (EVDW_PME(ir->vdwtype))
1037     {
1038         real b     = calc_ewaldcoeff_lj(ir->rvdw, ir->ewald_rtol_lj);
1039         real r     = ir->rvdw;
1040         real br    = b*r;
1041         real br2   = br*br;
1042         real br4   = br2*br2;
1043         real br6   = br4*br2;
1044         // -dV/dr of g(br)*r^-6 [where g(x) = exp(-x^2)(1+x^2+x^4/2),
1045         // see LJ-PME equations in manual] and r^-reppow
1046         ljDisp.md1 = -std::exp(-br2)*(br6 + 3.0*br4 + 6.0*br2 + 6.0)*std::pow(r, -7.0);
1047         ljRep.md1  = repPow*pow(r, -(repPow + 1));
1048         // The contribution of the higher derivatives is negligible
1049     }
1050     else
1051     {
1052         gmx_fatal(FARGS, "Energy drift calculation is only implemented for plain cut-off Lennard-Jones interactions");
1053     }
1054
1055     elfac = ONE_4PI_EPS0/ir->epsilon_r;
1056
1057     // Determine the 1st and 2nd derivative for the electostatics
1058     pot_derivatives_t elec = { 0, 0, 0 };
1059
1060     if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
1061     {
1062         real eps_rf, k_rf;
1063
1064         if (ir->coulombtype == eelCUT)
1065         {
1066             eps_rf = 1;
1067             k_rf   = 0;
1068         }
1069         else
1070         {
1071             eps_rf = ir->epsilon_rf/ir->epsilon_r;
1072             if (eps_rf != 0)
1073             {
1074                 k_rf = (eps_rf - ir->epsilon_r)/( gmx::power3(ir->rcoulomb) * (2*eps_rf + ir->epsilon_r) );
1075             }
1076             else
1077             {
1078                 /* epsilon_rf = infinity */
1079                 k_rf = 0.5/gmx::power3(ir->rcoulomb);
1080             }
1081         }
1082
1083         if (eps_rf > 0)
1084         {
1085             elec.md1 = elfac*(1.0/gmx::square(ir->rcoulomb) - 2*k_rf*ir->rcoulomb);
1086         }
1087         elec.d2      = elfac*(2.0/gmx::power3(ir->rcoulomb) + 2*k_rf);
1088     }
1089     else if (EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD)
1090     {
1091         real b, rc, br;
1092
1093         b        = calc_ewaldcoeff_q(ir->rcoulomb, ir->ewald_rtol);
1094         rc       = ir->rcoulomb;
1095         br       = b*rc;
1096         elec.md1 = elfac*(b*std::exp(-br*br)*M_2_SQRTPI/rc + std::erfc(br)/(rc*rc));
1097         elec.d2  = elfac/(rc*rc)*(2*b*(1 + br*br)*std::exp(-br*br)*M_2_SQRTPI + 2*std::erfc(br)/rc);
1098     }
1099     else
1100     {
1101         gmx_fatal(FARGS, "Energy drift calculation is only implemented for Reaction-Field and Ewald electrostatics");
1102     }
1103
1104     /* Determine the variance of the atomic displacement
1105      * over list_lifetime steps: kT_fac
1106      * For inertial dynamics (not Brownian dynamics) the mass factor
1107      * is not included in kT_fac, it is added later.
1108      */
1109     const real kT_fac = displacementVariance(*ir, reference_temperature,
1110                                              list_lifetime*ir->delta_t);
1111
1112     if (debug)
1113     {
1114         fprintf(debug, "Derivatives of non-bonded potentials at the cut-off:\n");
1115         fprintf(debug, "LJ disp. -V' %9.2e V'' %9.2e -V''' %9.2e\n", ljDisp.md1, ljDisp.d2, ljDisp.md3);
1116         fprintf(debug, "LJ rep.  -V' %9.2e V'' %9.2e -V''' %9.2e\n", ljRep.md1, ljRep.d2, ljRep.md3);
1117         fprintf(debug, "Electro. -V' %9.2e V'' %9.2e\n", elec.md1, elec.d2);
1118         fprintf(debug, "sqrt(kT_fac) %f\n", std::sqrt(kT_fac));
1119     }
1120
1121     /* Search using bisection */
1122     ib0 = -1;
1123     /* The drift will be neglible at 5 times the max sigma */
1124     ib1 = static_cast<int>(5*maxSigma(kT_fac, natt, att)/resolution) + 1;
1125     while (ib1 - ib0 > 1)
1126     {
1127         ib = (ib0 + ib1)/2;
1128         rb = ib*resolution;
1129         rl = std::max(ir->rvdw, ir->rcoulomb) + rb;
1130
1131         /* Calculate the average energy drift at the last step
1132          * of the nstlist steps at which the pair-list is used.
1133          */
1134         drift = energyDrift(att, natt, &mtop->ffparams,
1135                             kT_fac,
1136                             &ljDisp, &ljRep, &elec,
1137                             ir->rvdw, ir->rcoulomb,
1138                             rl, boxvol);
1139
1140         /* Correct for the fact that we are using a Ni x Nj particle pair list
1141          * and not a 1 x 1 particle pair list. This reduces the drift.
1142          */
1143         /* We don't have a formula for 8 (yet), use 4 which is conservative */
1144         nb_clust_frac_pairs_not_in_list_at_cutoff =
1145             surface_frac(std::min(list_setup->cluster_size_i, 4),
1146                          particle_distance, rl)*
1147             surface_frac(std::min(list_setup->cluster_size_j, 4),
1148                          particle_distance, rl);
1149         drift *= nb_clust_frac_pairs_not_in_list_at_cutoff;
1150
1151         /* Convert the drift to drift per unit time per atom */
1152         drift /= nstlist*ir->delta_t*mtop->natoms;
1153
1154         if (debug)
1155         {
1156             fprintf(debug, "ib %3d %3d %3d rb %.3f %dx%d fac %.3f drift %.1e\n",
1157                     ib0, ib, ib1, rb,
1158                     list_setup->cluster_size_i, list_setup->cluster_size_j,
1159                     nb_clust_frac_pairs_not_in_list_at_cutoff,
1160                     drift);
1161         }
1162
1163         if (std::abs(drift) > ir->verletbuf_tol)
1164         {
1165             ib0 = ib;
1166         }
1167         else
1168         {
1169             ib1 = ib;
1170         }
1171     }
1172
1173     sfree(att);
1174
1175     *rlist = std::max(ir->rvdw, ir->rcoulomb) + ib1*resolution;
1176 }
1177
1178 /* Returns the pairlist buffer size for use as a minimum buffer size
1179  *
1180  * Note that this is a rather crude estimate. It is ok for a buffer
1181  * set for good energy conservation or RF electrostatics. But it is
1182  * too small with PME and the buffer set with the default tolerance.
1183  */
1184 static real minCellSizeFromPairlistBuffer(const t_inputrec &ir)
1185 {
1186     return ir.rlist - std::max(ir.rvdw, ir.rcoulomb);
1187 }
1188
1189 real minCellSizeForAtomDisplacement(const gmx_mtop_t &mtop,
1190                                     const t_inputrec &ir,
1191                                     real              chanceRequested)
1192 {
1193     if (!EI_DYNAMICS(ir.eI) || (EI_MD(ir.eI) && ir.etc == etcNO))
1194     {
1195         return minCellSizeFromPairlistBuffer(ir);
1196     }
1197
1198     /* We use the maximum temperature with multiple T-coupl groups.
1199      * We could use a per particle temperature, but since particles
1200      * interact, this might underestimate the displacements.
1201      */
1202     const real            temperature = maxReferenceTemperature(ir);
1203
1204     const bool            setMassesToOne = (ir.eI == eiBD && ir.bd_fric > 0);
1205
1206     verletbuf_atomtype_t *att  = nullptr;
1207     int                   natt = -1;
1208     get_verlet_buffer_atomtypes(&mtop, setMassesToOne, &att, &natt, nullptr);
1209
1210     const real kT_fac = displacementVariance(ir, temperature,
1211                                              ir.nstlist*ir.delta_t);
1212
1213     /* Resolution of the cell size */
1214     real resolution = 0.001;
1215
1216     /* Search using bisection, avoid 0 and start at 1 */
1217     int  ib0      = 0;
1218     /* The chance will be neglible at 10 times the max sigma */
1219     int  ib1      = int(10*maxSigma(kT_fac, natt, att)/resolution) + 1;
1220     real cellSize = 0;
1221     while (ib1 - ib0 > 1)
1222     {
1223         int ib = (ib0 + ib1)/2;
1224         cellSize = ib*resolution;
1225
1226         /* We assumes atom are distributed uniformly over the cell width.
1227          * Once an atom has moved by more than the cellSize (as passed
1228          * as the buffer argument to energyDriftAtomPair() below),
1229          * the chance of crossing the boundary of the neighbor cell
1230          * thus increases as 1/cellSize with the additional displacement
1231          * on to of cellSize. We thus create a linear interaction with
1232          * derivative = -1/cellSize. Using this in the energyDriftAtomPair
1233          * function will return the chance of crossing the next boundary.
1234          */
1235         const pot_derivatives_t boundaryInteraction = { 1/cellSize, 0, 0 };
1236
1237         real                    chance = 0;
1238         for (int i = 0; i < natt; i++)
1239         {
1240             const atom_nonbonded_kinetic_prop_t &propAtom = att[i].prop;
1241             real s2_2d;
1242             real s2_3d;
1243             get_atom_sigma2(kT_fac, &propAtom, &s2_2d, &s2_3d);
1244
1245             real chancePerAtom = energyDriftAtomPair(propAtom.bConstr, false,
1246                                                      s2_2d + s2_3d, s2_2d, 0,
1247                                                      cellSize,
1248                                                      &boundaryInteraction);
1249
1250             if (propAtom.bConstr)
1251             {
1252                 /* energyDriftAtomPair() uses an unlimited Gaussian displacement
1253                  * distribution for constrained atoms, whereas they can
1254                  * actually not move more than the COM of the two constrained
1255                  * atoms plus twice the distance from the COM.
1256                  * Use this maximum, limited displacement when this results in
1257                  * a smaller chance (note that this is still an overestimate).
1258                  */
1259                 real massFraction = propAtom.con_mass/(propAtom.mass + propAtom.con_mass);
1260                 real comDistance  = propAtom.con_len*massFraction;
1261
1262                 real chanceWithMaxDistance =
1263                     energyDriftAtomPair(false, false,
1264                                         s2_3d, 0, 0,
1265                                         cellSize - 2*comDistance,
1266                                         &boundaryInteraction);
1267                 chancePerAtom = std::min(chancePerAtom, chanceWithMaxDistance);
1268             }
1269
1270             /* Take into account the line density of the boundary */
1271             chancePerAtom /= cellSize;
1272
1273             chance        += att[i].n*chancePerAtom;
1274         }
1275
1276         /* Note: chance is for every nstlist steps */
1277         if (chance > chanceRequested*ir.nstlist)
1278         {
1279             ib0 = ib;
1280         }
1281         else
1282         {
1283             ib1 = ib;
1284         }
1285     }
1286
1287     sfree(att);
1288
1289     return cellSize;
1290 }