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