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