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