Fixes issue with vsiten and verlet buffers.
[alexxy/gromacs.git] / src / kernel / calc_verletbuf.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team,
6  * check out http://www.gromacs.org for more information.
7  * Copyright (c) 2012,2013, by the GROMACS development team, led by
8  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
9  * others, as listed in the AUTHORS file in the top-level source
10  * directory and at http://www.gromacs.org.
11  *
12  * GROMACS is free software; you can redistribute it and/or
13  * modify it under the terms of the GNU Lesser General Public License
14  * as published by the Free Software Foundation; either version 2.1
15  * of the License, or (at your option) any later version.
16  *
17  * GROMACS is distributed in the hope that it will be useful,
18  * but WITHOUT ANY WARRANTY; without even the implied warranty of
19  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
20  * Lesser General Public License for more details.
21  *
22  * You should have received a copy of the GNU Lesser General Public
23  * License along with GROMACS; if not, see
24  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
26  *
27  * If you want to redistribute modifications to GROMACS, please
28  * consider that scientific software is very special. Version
29  * control is crucial - bugs must be traceable. We will be happy to
30  * consider code for inclusion in the official distribution, but
31  * derived work must not be called official GROMACS. Details are found
32  * in the README & COPYING files - if they are missing, get the
33  * official version at http://www.gromacs.org.
34  *
35  * To help us fund GROMACS development, we humbly ask that you cite
36  * the research papers on the package. Check out http://www.gromacs.org.
37  */
38 #ifdef HAVE_CONFIG_H
39 #include <config.h>
40 #endif
41
42 #include "assert.h"
43
44 #include <sys/types.h>
45 #include <math.h>
46 #include "typedefs.h"
47 #include "physics.h"
48 #include "smalloc.h"
49 #include "gmx_fatal.h"
50 #include "macros.h"
51 #include "vec.h"
52 #include "coulomb.h"
53 #include "calc_verletbuf.h"
54 #include "../mdlib/nbnxn_consts.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 "gmx_simd_macros.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_WIDTH_HERE;
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], 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         sfree(vsite_m);
420         sfree(prop);
421     }
422
423     if (gmx_debug_at)
424     {
425         for (a = 0; a < natt; a++)
426         {
427             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",
428                     a, att[a].prop.mass, att[a].prop.type, att[a].prop.q,
429                     att[a].prop.bConstr, att[a].prop.con_mass, att[a].prop.con_len, 
430                     att[a].n);
431         }
432     }
433
434     *att_p  = att;
435     *natt_p = natt;
436 }
437
438 /* This function computes two components of the estimate of the variance
439  * in the displacement of one atom in a system of two constrained atoms.
440  * Returns in sigma2_2d the variance due to rotation of the constrained
441  * atom around the atom to which it constrained.
442  * Returns in sigma2_3d the variance due to displacement of the COM
443  * of the whole system of the two constrained atoms.
444  *
445  * Note that we only take a single constraint (the one to the heaviest atom)
446  * into account. If an atom has multiple constraints, this will result in
447  * an overestimate of the displacement, which gives a larger drift and buffer.
448  */
449 static void constrained_atom_sigma2(real kT_fac,
450                                     const atom_nonbonded_kinetic_prop_t *prop,
451                                     real *sigma2_2d,
452                                     real *sigma2_3d)
453 {
454     real sigma2_rot;
455     real com_dist;
456     real sigma2_rel;
457     real scale;
458
459     /* Here we decompose the motion of a constrained atom into two
460      * components: rotation around the COM and translation of the COM.
461      */
462
463     /* Determine the variance for the displacement of the rotational mode */
464     sigma2_rot = kT_fac/(prop->mass*(prop->mass + prop->con_mass)/prop->con_mass);
465
466     /* The distance from the atom to the COM, i.e. the rotational arm */
467     com_dist = prop->con_len*prop->con_mass/(prop->mass + prop->con_mass);
468
469     /* The variance relative to the arm */
470     sigma2_rel = sigma2_rot/(com_dist*com_dist);
471     /* At 6 the scaling formula has slope 0,
472      * so we keep sigma2_2d constant after that.
473      */
474     if (sigma2_rel < 6)
475     {
476         /* A constrained atom rotates around the atom it is constrained to.
477          * This results in a smaller linear displacement than for a free atom.
478          * For a perfectly circular displacement, this lowers the displacement
479          * by: 1/arcsin(arc_length)
480          * and arcsin(x) = 1 + x^2/6 + ...
481          * For sigma2_rel<<1 the displacement distribution is erfc
482          * (exact formula is provided below). For larger sigma, it is clear
483          * that the displacement can't be larger than 2*com_dist.
484          * It turns out that the distribution becomes nearly uniform.
485          * For intermediate sigma2_rel, scaling down sigma with the third
486          * order expansion of arcsin with argument sigma_rel turns out
487          * to give a very good approximation of the distribution and variance.
488          * Even for larger values, the variance is only slightly overestimated.
489          * Note that the most relevant displacements are in the long tail.
490          * This rotation approximation always overestimates the tail (which
491          * runs to infinity, whereas it should be <= 2*com_dist).
492          * Thus we always overestimate the drift and the buffer size.
493          */
494         scale      = 1/(1 + sigma2_rel/6);
495         *sigma2_2d = sigma2_rot*scale*scale;
496     }
497     else
498     {
499         /* sigma_2d is set to the maximum given by the scaling above.
500          * For large sigma2 the real displacement distribution is close
501          * to uniform over -2*con_len to 2*com_dist.
502          * Our erfc with sigma_2d=sqrt(1.5)*com_dist (which means the sigma
503          * of the erfc output distribution is con_dist) overestimates
504          * the variance and additionally has a long tail. This means
505          * we have a (safe) overestimation of the drift.
506          */
507         *sigma2_2d = 1.5*com_dist*com_dist;
508     }
509
510     /* The constrained atom also moves (in 3D) with the COM of both atoms */
511     *sigma2_3d = kT_fac/(prop->mass + prop->con_mass);
512 }
513
514 static void get_atom_sigma2(real kT_fac,
515                             const atom_nonbonded_kinetic_prop_t *prop,
516                             real *sigma2_2d,
517                             real *sigma2_3d)
518 {
519     if (prop->bConstr)
520     {
521         /* Complicated constraint calculation in a separate function */
522         constrained_atom_sigma2(kT_fac, prop, sigma2_2d, sigma2_3d);
523     }
524     else
525     {
526         /* Unconstrained atom: trivial */
527         *sigma2_2d = 0;
528         *sigma2_3d = kT_fac/prop->mass;
529     }
530 }
531
532 static void approx_2dof(real s2, real x, real *shift, real *scale)
533 {
534     /* A particle with 1 DOF constrained has 2 DOFs instead of 3.
535      * This code is also used for particles with multiple constraints,
536      * in which case we overestimate the displacement.
537      * The 2DOF distribution is sqrt(pi/2)*erfc(r/(sqrt(2)*s))/(2*s).
538      * We approximate this with scale*Gaussian(s,r+shift),
539      * by matching the distribution value and derivative at x.
540      * This is a tight overestimate for all r>=0 at any s and x.
541      */
542     real ex, er;
543
544     ex = exp(-x*x/(2*s2));
545     er = gmx_erfc(x/sqrt(2*s2));
546
547     *shift = -x + sqrt(2*s2/M_PI)*ex/er;
548     *scale = 0.5*M_PI*exp(ex*ex/(M_PI*er*er))*er;
549 }
550
551 static real ener_drift(const verletbuf_atomtype_t *att, int natt,
552                        const gmx_ffparams_t *ffp,
553                        real kT_fac,
554                        real md_ljd, real md_ljr, real md_el, real dd_el,
555                        real r_buffer,
556                        real rlist, real boxvol)
557 {
558     double drift_tot, pot1, pot2, pot;
559     int    i, j;
560     real   s2i_2d, s2i_3d, s2j_2d, s2j_3d, s2, s;
561     int    ti, tj;
562     real   md, dd;
563     real   sc_fac, rsh;
564     double c_exp, c_erfc;
565
566     drift_tot = 0;
567
568     /* Loop over the different atom type pairs */
569     for (i = 0; i < natt; i++)
570     {
571         get_atom_sigma2(kT_fac, &att[i].prop, &s2i_2d, &s2i_3d);
572         ti = att[i].prop.type;
573
574         for (j = i; j < natt; j++)
575         {
576             get_atom_sigma2(kT_fac, &att[j].prop, &s2j_2d, &s2j_3d);
577             tj = att[j].prop.type;
578
579             /* Add up the up to four independent variances */
580             s2 = s2i_2d + s2i_3d + s2j_2d + s2j_3d; 
581
582             /* Note that attractive and repulsive potentials for individual
583              * pairs will partially cancel.
584              */
585             /* -dV/dr at the cut-off for LJ + Coulomb */
586             md =
587                 md_ljd*ffp->iparams[ti*ffp->atnr+tj].lj.c6 +
588                 md_ljr*ffp->iparams[ti*ffp->atnr+tj].lj.c12 +
589                 md_el*att[i].prop.q*att[j].prop.q;
590
591             /* d2V/dr2 at the cut-off for Coulomb, we neglect LJ */
592             dd = dd_el*att[i].prop.q*att[j].prop.q;
593
594             rsh    = r_buffer;
595             sc_fac = 1.0;
596             /* For constraints: adapt r and scaling for the Gaussian */
597             if (att[i].prop.bConstr)
598             {
599                 real sh, sc;
600
601                 approx_2dof(s2i_2d, r_buffer*s2i_2d/s2, &sh, &sc);
602                 rsh    += sh;
603                 sc_fac *= sc;
604             }
605             if (att[j].prop.bConstr)
606             {
607                 real sh, sc;
608
609                 approx_2dof(s2j_2d, r_buffer*s2j_2d/s2, &sh, &sc);
610                 rsh    += sh;
611                 sc_fac *= sc;
612             }
613
614             /* Exact contribution of an atom pair with Gaussian displacement
615              * with sigma s to the energy drift for a potential with
616              * derivative -md and second derivative dd at the cut-off.
617              * The only catch is that for potentials that change sign
618              * near the cut-off there could be an unlucky compensation
619              * of positive and negative energy drift.
620              * Such potentials are extremely rare though.
621              *
622              * Note that pot has unit energy*length, as the linear
623              * atom density still needs to be put in.
624              */
625             c_exp  = exp(-rsh*rsh/(2*s2))/sqrt(2*M_PI);
626             c_erfc = 0.5*gmx_erfc(rsh/(sqrt(2*s2)));
627             s      = sqrt(s2);
628
629             pot1 = sc_fac*
630                 md/2*((rsh*rsh + s2)*c_erfc - rsh*s*c_exp);
631             pot2 = sc_fac*
632                 dd/6*(s*(rsh*rsh + 2*s2)*c_exp - rsh*(rsh*rsh + 3*s2)*c_erfc);
633             pot = pot1 + pot2;
634
635             if (gmx_debug_at)
636             {
637                 fprintf(debug, "n %d %d d s %.3f %.3f %.3f %.3f con %d md %8.1e dd %8.1e pot1 %8.1e pot2 %8.1e pot %8.1e\n",
638                         att[i].n, att[j].n,
639                         sqrt(s2i_2d), sqrt(s2i_3d),
640                         sqrt(s2j_2d), sqrt(s2j_3d),
641                         att[i].prop.bConstr+att[j].prop.bConstr,
642                         md, dd, pot1, pot2, pot);
643             }
644
645             /* Multiply by the number of atom pairs */
646             if (j == i)
647             {
648                 pot *= (double)att[i].n*(att[i].n - 1)/2;
649             }
650             else
651             {
652                 pot *= (double)att[i].n*att[j].n;
653             }
654             /* We need the line density to get the energy drift of the system.
655              * The effective average r^2 is close to (rlist+sigma)^2.
656              */
657             pot *= 4*M_PI*sqr(rlist + s)/boxvol;
658
659             /* Add the unsigned drift to avoid cancellation of errors */
660             drift_tot += fabs(pot);
661         }
662     }
663
664     return drift_tot;
665 }
666
667 static real surface_frac(int cluster_size, real particle_distance, real rlist)
668 {
669     real d, area_rel;
670
671     if (rlist < 0.5*particle_distance)
672     {
673         /* We have non overlapping spheres */
674         return 1.0;
675     }
676
677     /* Half the inter-particle distance relative to rlist */
678     d = 0.5*particle_distance/rlist;
679
680     /* Determine the area of the surface at distance rlist to the closest
681      * particle, relative to surface of a sphere of radius rlist.
682      * The formulas below assume close to cubic cells for the pair search grid,
683      * which the pair search code tries to achieve.
684      * Note that in practice particle distances will not be delta distributed,
685      * but have some spread, often involving shorter distances,
686      * as e.g. O-H bonds in a water molecule. Thus the estimates below will
687      * usually be slightly too high and thus conservative.
688      */
689     switch (cluster_size)
690     {
691         case 1:
692             /* One particle: trivial */
693             area_rel = 1.0;
694             break;
695         case 2:
696             /* Two particles: two spheres at fractional distance 2*a */
697             area_rel = 1.0 + d;
698             break;
699         case 4:
700             /* We assume a perfect, symmetric tetrahedron geometry.
701              * The surface around a tetrahedron is too complex for a full
702              * analytical solution, so we use a Taylor expansion.
703              */
704             area_rel = (1.0 + 1/M_PI*(6*acos(1/sqrt(3))*d +
705                                       sqrt(3)*d*d*(1.0 +
706                                                    5.0/18.0*d*d +
707                                                    7.0/45.0*d*d*d*d +
708                                                    83.0/756.0*d*d*d*d*d*d)));
709             break;
710         default:
711             gmx_incons("surface_frac called with unsupported cluster_size");
712             area_rel = 1.0;
713     }
714
715     return area_rel/cluster_size;
716 }
717
718 void calc_verlet_buffer_size(const gmx_mtop_t *mtop, real boxvol,
719                              const t_inputrec *ir, real drift_target,
720                              const verletbuf_list_setup_t *list_setup,
721                              int *n_nonlin_vsite,
722                              real *rlist)
723 {
724     double                resolution;
725     char                 *env;
726
727     real                  particle_distance;
728     real                  nb_clust_frac_pairs_not_in_list_at_cutoff;
729
730     verletbuf_atomtype_t *att  = NULL;
731     int                   natt = -1, i;
732     double                reppow;
733     real                  md_ljd, md_ljr, md_el, dd_el;
734     real                  elfac;
735     real                  kT_fac, mass_min;
736     int                   ib0, ib1, ib;
737     real                  rb, rl;
738     real                  drift;
739
740     /* Resolution of the buffer size */
741     resolution = 0.001;
742
743     env = getenv("GMX_VERLET_BUFFER_RES");
744     if (env != NULL)
745     {
746         sscanf(env, "%lf", &resolution);
747     }
748
749     /* In an atom wise pair-list there would be no pairs in the list
750      * beyond the pair-list cut-off.
751      * However, we use a pair-list of groups vs groups of atoms.
752      * For groups of 4 atoms, the parallelism of SSE instructions, only
753      * 10% of the atoms pairs are not in the list just beyond the cut-off.
754      * As this percentage increases slowly compared to the decrease of the
755      * Gaussian displacement distribution over this range, we can simply
756      * reduce the drift by this fraction.
757      * For larger groups, e.g. of 8 atoms, this fraction will be lower,
758      * so then buffer size will be on the conservative (large) side.
759      *
760      * Note that the formulas used here do not take into account
761      * cancellation of errors which could occur by missing both
762      * attractive and repulsive interactions.
763      *
764      * The only major assumption is homogeneous particle distribution.
765      * For an inhomogeneous system, such as a liquid-vapor system,
766      * the buffer will be underestimated. The actual energy drift
767      * will be higher by the factor: local/homogeneous particle density.
768      *
769      * The results of this estimate have been checked againt simulations.
770      * In most cases the real drift differs by less than a factor 2.
771      */
772
773     /* Worst case assumption: HCP packing of particles gives largest distance */
774     particle_distance = pow(boxvol*sqrt(2)/mtop->natoms, 1.0/3.0);
775
776     get_verlet_buffer_atomtypes(mtop, &att, &natt, n_nonlin_vsite);
777     assert(att != NULL && natt >= 0);
778
779     if (debug)
780     {
781         fprintf(debug, "particle distance assuming HCP packing: %f nm\n",
782                 particle_distance);
783         fprintf(debug, "energy drift atom types: %d\n", natt);
784     }
785
786     reppow = mtop->ffparams.reppow;
787     md_ljd = 0;
788     md_ljr = 0;
789     if (ir->vdwtype == evdwCUT)
790     {
791         /* -dV/dr of -r^-6 and r^-repporw */
792         md_ljd = -6*pow(ir->rvdw, -7.0);
793         md_ljr = reppow*pow(ir->rvdw, -(reppow+1));
794         /* The contribution of the second derivative is negligible */
795     }
796     else
797     {
798         gmx_fatal(FARGS, "Energy drift calculation is only implemented for plain cut-off Lennard-Jones interactions");
799     }
800
801     elfac = ONE_4PI_EPS0/ir->epsilon_r;
802
803     /* Determine md=-dV/dr and dd=d^2V/dr^2 */
804     md_el = 0;
805     dd_el = 0;
806     if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
807     {
808         real eps_rf, k_rf;
809
810         if (ir->coulombtype == eelCUT)
811         {
812             eps_rf = 1;
813             k_rf   = 0;
814         }
815         else
816         {
817             eps_rf = ir->epsilon_rf/ir->epsilon_r;
818             if (eps_rf != 0)
819             {
820                 k_rf = pow(ir->rcoulomb, -3.0)*(eps_rf - ir->epsilon_r)/(2*eps_rf + ir->epsilon_r);
821             }
822             else
823             {
824                 /* epsilon_rf = infinity */
825                 k_rf = 0.5*pow(ir->rcoulomb, -3.0);
826             }
827         }
828
829         if (eps_rf > 0)
830         {
831             md_el = elfac*(pow(ir->rcoulomb, -2.0) - 2*k_rf*ir->rcoulomb);
832         }
833         dd_el = elfac*(2*pow(ir->rcoulomb, -3.0) + 2*k_rf);
834     }
835     else if (EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD)
836     {
837         real b, rc, br;
838
839         b     = calc_ewaldcoeff(ir->rcoulomb, ir->ewald_rtol);
840         rc    = ir->rcoulomb;
841         br    = b*rc;
842         md_el = elfac*(b*exp(-br*br)*M_2_SQRTPI/rc + gmx_erfc(br)/(rc*rc));
843         dd_el = elfac/(rc*rc)*(2*b*(1 + br*br)*exp(-br*br)*M_2_SQRTPI + 2*gmx_erfc(br)/rc);
844     }
845     else
846     {
847         gmx_fatal(FARGS, "Energy drift calculation is only implemented for Reaction-Field and Ewald electrostatics");
848     }
849
850     /* Determine the variance of the atomic displacement
851      * over nstlist-1 steps: kT_fac
852      * For inertial dynamics (not Brownian dynamics) the mass factor
853      * is not included in kT_fac, it is added later.
854      */
855     if (ir->eI == eiBD)
856     {
857         /* Get the displacement distribution from the random component only.
858          * With accurate integration the systematic (force) displacement
859          * should be negligible (unless nstlist is extremely large, which
860          * you wouldn't do anyhow).
861          */
862         kT_fac = 2*BOLTZ*ir->opts.ref_t[0]*(ir->nstlist-1)*ir->delta_t;
863         if (ir->bd_fric > 0)
864         {
865             /* This is directly sigma^2 of the displacement */
866             kT_fac /= ir->bd_fric;
867
868             /* Set the masses to 1 as kT_fac is the full sigma^2,
869              * but we divide by m in ener_drift().
870              */
871             for (i = 0; i < natt; i++)
872             {
873                 att[i].prop.mass = 1;
874             }
875         }
876         else
877         {
878             real tau_t;
879
880             /* Per group tau_t is not implemented yet, use the maximum */
881             tau_t = ir->opts.tau_t[0];
882             for (i = 1; i < ir->opts.ngtc; i++)
883             {
884                 tau_t = max(tau_t, ir->opts.tau_t[i]);
885             }
886
887             kT_fac *= tau_t;
888             /* This kT_fac needs to be divided by the mass to get sigma^2 */
889         }
890     }
891     else
892     {
893         kT_fac = BOLTZ*ir->opts.ref_t[0]*sqr((ir->nstlist-1)*ir->delta_t);
894     }
895
896     mass_min = att[0].prop.mass;
897     for (i = 1; i < natt; i++)
898     {
899         mass_min = min(mass_min, att[i].prop.mass);
900     }
901
902     if (debug)
903     {
904         fprintf(debug, "md_ljd %e md_ljr %e\n", md_ljd, md_ljr);
905         fprintf(debug, "md_el %e dd_el %e\n", md_el, dd_el);
906         fprintf(debug, "sqrt(kT_fac) %f\n", sqrt(kT_fac));
907         fprintf(debug, "mass_min %f\n", mass_min);
908     }
909
910     /* Search using bisection */
911     ib0 = -1;
912     /* The drift will be neglible at 5 times the max sigma */
913     ib1 = (int)(5*2*sqrt(kT_fac/mass_min)/resolution) + 1;
914     while (ib1 - ib0 > 1)
915     {
916         ib = (ib0 + ib1)/2;
917         rb = ib*resolution;
918         rl = max(ir->rvdw, ir->rcoulomb) + rb;
919
920         /* Calculate the average energy drift at the last step
921          * of the nstlist steps at which the pair-list is used.
922          */
923         drift = ener_drift(att, natt, &mtop->ffparams,
924                            kT_fac,
925                            md_ljd, md_ljr, md_el, dd_el, rb,
926                            rl, boxvol);
927
928         /* Correct for the fact that we are using a Ni x Nj particle pair list
929          * and not a 1 x 1 particle pair list. This reduces the drift.
930          */
931         /* We don't have a formula for 8 (yet), use 4 which is conservative */
932         nb_clust_frac_pairs_not_in_list_at_cutoff =
933             surface_frac(min(list_setup->cluster_size_i, 4),
934                          particle_distance, rl)*
935             surface_frac(min(list_setup->cluster_size_j, 4),
936                          particle_distance, rl);
937         drift *= nb_clust_frac_pairs_not_in_list_at_cutoff;
938
939         /* Convert the drift to drift per unit time per atom */
940         drift /= ir->nstlist*ir->delta_t*mtop->natoms;
941
942         if (debug)
943         {
944             fprintf(debug, "ib %3d %3d %3d rb %.3f %dx%d fac %.3f drift %f\n",
945                     ib0, ib, ib1, rb,
946                     list_setup->cluster_size_i, list_setup->cluster_size_j,
947                     nb_clust_frac_pairs_not_in_list_at_cutoff,
948                     drift);
949         }
950
951         if (fabs(drift) > drift_target)
952         {
953             ib0 = ib;
954         }
955         else
956         {
957             ib1 = ib;
958         }
959     }
960
961     sfree(att);
962
963     *rlist = max(ir->rvdw, ir->rcoulomb) + ib1*resolution;
964 }