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