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