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