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