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