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