Merge "Merge branch release-4-6 into master"
[alexxy/gromacs.git] / src / gromacs / gmxpreprocess / calc_verletbuf.c
1 /*  -*- mode: c; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4; c-file-style: "stroustrup"; -*-
2  *
3  *
4  *                This source code is part of
5  *
6  *                 G   R   O   M   A   C   S
7  *
8  *          GROningen MAchine for Chemical Simulations
9  *
10  *                        VERSION 3.2.03
11  * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
12  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
13  * Copyright (c) 2001-2004, The GROMACS development team,
14  * check out http://www.gromacs.org for more information.
15
16  * This program is free software; you can redistribute it and/or
17  * modify it under the terms of the GNU General Public License
18  * as published by the Free Software Foundation; either version 2
19  * of the License, or (at your option) any later version.
20  *
21  * If you want to redistribute modifications, please consider that
22  * scientific software is very special. Version control is crucial -
23  * bugs must be traceable. We will be happy to consider code for
24  * inclusion in the official distribution, but derived work must not
25  * be called official GROMACS. Details are found in the README & COPYING
26  * files - if they are missing, get the official version at www.gromacs.org.
27  *
28  * To help us fund GROMACS development, we humbly ask that you cite
29  * the papers on the package - you can find them in the top README file.
30  *
31  * For more info, check our website at http://www.gromacs.org
32  *
33  * And Hey:
34  * Gallium Rubidium Oxygen Manganese Argon Carbon Silicon
35  */
36 #ifdef HAVE_CONFIG_H
37 #include <config.h>
38 #endif
39
40 #include <assert.h>
41
42 #include <sys/types.h>
43 #include <math.h>
44 #include "typedefs.h"
45 #include "physics.h"
46 #include "smalloc.h"
47 #include "gmx_fatal.h"
48 #include "macros.h"
49 #include "vec.h"
50 #include "coulomb.h"
51 #include "calc_verletbuf.h"
52 #include "../mdlib/nbnxn_consts.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 "gmx_simd_macros.h"
62 #endif
63
64 /* Struct for unique atom type for calculating the energy drift.
65  * The atom displacement depends on mass and constraints.
66  * The energy jump for given distance depend on LJ type and q.
67  */
68 typedef struct
69 {
70     real     mass; /* mass */
71     int      type; /* type (used for LJ parameters) */
72     real     q;    /* charge */
73     int      con;  /* constrained: 0, else 1, if 1, use #DOF=2 iso 3 */
74     int      n;    /* total #atoms of this type in the system */
75 } verletbuf_atomtype_t;
76
77
78 void verletbuf_get_list_setup(gmx_bool                bGPU,
79                               verletbuf_list_setup_t *list_setup)
80 {
81     list_setup->cluster_size_i     = NBNXN_CPU_CLUSTER_I_SIZE;
82
83     if (bGPU)
84     {
85         list_setup->cluster_size_j = NBNXN_GPU_CLUSTER_SIZE;
86     }
87     else
88     {
89 #ifndef GMX_NBNXN_SIMD
90         list_setup->cluster_size_j = NBNXN_CPU_CLUSTER_I_SIZE;
91 #else
92         list_setup->cluster_size_j = GMX_SIMD_WIDTH_HERE;
93 #ifdef GMX_NBNXN_SIMD_2XNN
94         /* We assume the smallest cluster size to be on the safe side */
95         list_setup->cluster_size_j /= 2;
96 #endif
97 #endif
98     }
99 }
100
101 static void add_at(verletbuf_atomtype_t **att_p, int *natt_p,
102                    real mass, int type, real q, int con, int nmol)
103 {
104     verletbuf_atomtype_t *att;
105     int                   natt, i;
106
107     if (mass == 0)
108     {
109         /* Ignore massless particles */
110         return;
111     }
112
113     att  = *att_p;
114     natt = *natt_p;
115
116     i = 0;
117     while (i < natt &&
118            !(mass == att[i].mass &&
119              type == att[i].type &&
120              q    == att[i].q &&
121              con  == att[i].con))
122     {
123         i++;
124     }
125
126     if (i < natt)
127     {
128         att[i].n += nmol;
129     }
130     else
131     {
132         (*natt_p)++;
133         srenew(*att_p, *natt_p);
134         (*att_p)[i].mass = mass;
135         (*att_p)[i].type = type;
136         (*att_p)[i].q    = q;
137         (*att_p)[i].con  = con;
138         (*att_p)[i].n    = nmol;
139     }
140 }
141
142 static void get_verlet_buffer_atomtypes(const gmx_mtop_t      *mtop,
143                                         verletbuf_atomtype_t **att_p,
144                                         int                   *natt_p,
145                                         int                   *n_nonlin_vsite)
146 {
147     verletbuf_atomtype_t *att;
148     int                   natt;
149     int                   mb, nmol, ft, i, j, a1, a2, a3, a;
150     const t_atoms        *atoms;
151     const t_ilist        *il;
152     const t_atom         *at;
153     const t_iparams      *ip;
154     real                 *con_m, *vsite_m, cam[5];
155
156     att  = NULL;
157     natt = 0;
158
159     if (n_nonlin_vsite != NULL)
160     {
161         *n_nonlin_vsite = 0;
162     }
163
164     for (mb = 0; mb < mtop->nmolblock; mb++)
165     {
166         nmol = mtop->molblock[mb].nmol;
167
168         atoms = &mtop->moltype[mtop->molblock[mb].type].atoms;
169
170         /* Check for constraints, as they affect the kinetic energy */
171         snew(con_m, atoms->nr);
172         snew(vsite_m, atoms->nr);
173
174         for (ft = F_CONSTR; ft <= F_CONSTRNC; ft++)
175         {
176             il = &mtop->moltype[mtop->molblock[mb].type].ilist[ft];
177
178             for (i = 0; i < il->nr; i += 1+NRAL(ft))
179             {
180                 a1         = il->iatoms[i+1];
181                 a2         = il->iatoms[i+2];
182                 con_m[a1] += atoms->atom[a2].m;
183                 con_m[a2] += atoms->atom[a1].m;
184             }
185         }
186
187         il = &mtop->moltype[mtop->molblock[mb].type].ilist[F_SETTLE];
188
189         for (i = 0; i < il->nr; i += 1+NRAL(F_SETTLE))
190         {
191             a1         = il->iatoms[i+1];
192             a2         = il->iatoms[i+2];
193             a3         = il->iatoms[i+3];
194             con_m[a1] += atoms->atom[a2].m + atoms->atom[a3].m;
195             con_m[a2] += atoms->atom[a1].m + atoms->atom[a3].m;
196             con_m[a3] += atoms->atom[a1].m + atoms->atom[a2].m;
197         }
198
199         /* Check for virtual sites, determine mass from constructing atoms */
200         for (ft = 0; ft < F_NRE; ft++)
201         {
202             if (IS_VSITE(ft))
203             {
204                 il = &mtop->moltype[mtop->molblock[mb].type].ilist[ft];
205
206                 for (i = 0; i < il->nr; i += 1+NRAL(ft))
207                 {
208                     ip = &mtop->ffparams.iparams[il->iatoms[i]];
209
210                     a1 = il->iatoms[i+1];
211
212                     for (j = 1; j < NRAL(ft); j++)
213                     {
214                         cam[j] = atoms->atom[il->iatoms[i+1+j]].m;
215                         if (cam[j] == 0)
216                         {
217                             cam[j] = vsite_m[il->iatoms[i+1+j]];
218                         }
219                         if (cam[j] == 0)
220                         {
221                             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.",
222                                       *mtop->moltype[mtop->molblock[mb].type].name,
223                                       interaction_function[ft].longname,
224                                       il->iatoms[i+1+j]+1);
225                         }
226                     }
227
228                     switch (ft)
229                     {
230                         case F_VSITE2:
231                             /* Exact except for ignoring constraints */
232                             vsite_m[a1] = (cam[2]*sqr(1-ip->vsite.a) + cam[1]*sqr(ip->vsite.a))/(cam[1]*cam[2]);
233                             break;
234                         case F_VSITE3:
235                             /* Exact except for ignoring constraints */
236                             vsite_m[a1] = (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))/(cam[1]*cam[2]*cam[3]);
237                             break;
238                         default:
239                             /* Use the mass of the lightest constructing atom.
240                              * This is an approximation.
241                              * If the distance of the virtual site to the
242                              * constructing atom is less than all distances
243                              * between constructing atoms, this is a safe
244                              * over-estimate of the displacement of the vsite.
245                              * This condition holds for all H mass replacement
246                              * replacement vsite constructions, except for SP2/3
247                              * groups. In SP3 groups one H will have a F_VSITE3
248                              * construction, so even there the total drift
249                              * estimation shouldn't be far off.
250                              */
251                             assert(j >= 1);
252                             vsite_m[a1] = cam[1];
253                             for (j = 2; j < NRAL(ft); j++)
254                             {
255                                 vsite_m[a1] = min(vsite_m[a1], cam[j]);
256                             }
257                             if (n_nonlin_vsite != NULL)
258                             {
259                                 *n_nonlin_vsite += nmol;
260                             }
261                             break;
262                     }
263                 }
264             }
265         }
266
267         for (a = 0; a < atoms->nr; a++)
268         {
269             at = &atoms->atom[a];
270             /* We consider an atom constrained, #DOF=2, when it is
271              * connected with constraints to one or more atoms with
272              * total mass larger than 1.5 that of the atom itself.
273              */
274             add_at(&att, &natt,
275                    at->m, at->type, at->q, con_m[a] > 1.5*at->m, nmol);
276         }
277
278         sfree(vsite_m);
279         sfree(con_m);
280     }
281
282     if (gmx_debug_at)
283     {
284         for (a = 0; a < natt; a++)
285         {
286             fprintf(debug, "type %d: m %5.2f t %d q %6.3f con %d n %d\n",
287                     a, att[a].mass, att[a].type, att[a].q, att[a].con, att[a].n);
288         }
289     }
290
291     *att_p  = att;
292     *natt_p = natt;
293 }
294
295 static void approx_2dof(real s2, real x,
296                         real *shift, real *scale)
297 {
298     /* A particle with 1 DOF constrained has 2 DOFs instead of 3.
299      * This code is also used for particles with multiple constraints,
300      * in which case we overestimate the displacement.
301      * The 2DOF distribution is sqrt(pi/2)*erfc(r/(sqrt(2)*s))/(2*s).
302      * We approximate this with scale*Gaussian(s,r+shift),
303      * by matching the distribution value and derivative at x.
304      * This is a tight overestimate for all r>=0 at any s and x.
305      */
306     real ex, er;
307
308     ex = exp(-x*x/(2*s2));
309     er = gmx_erfc(x/sqrt(2*s2));
310
311     *shift = -x + sqrt(2*s2/M_PI)*ex/er;
312     *scale = 0.5*M_PI*exp(ex*ex/(M_PI*er*er))*er;
313 }
314
315 static real ener_drift(const verletbuf_atomtype_t *att, int natt,
316                        const gmx_ffparams_t *ffp,
317                        real kT_fac,
318                        real md_ljd, real md_ljr, real md_el, real dd_el,
319                        real r_buffer,
320                        real rlist, real boxvol)
321 {
322     double drift_tot, pot1, pot2, pot;
323     int    i, j;
324     real   s2i, s2j, s2, s;
325     int    ti, tj;
326     real   md, dd;
327     real   sc_fac, rsh;
328     double c_exp, c_erfc;
329
330     drift_tot = 0;
331
332     /* Loop over the different atom type pairs */
333     for (i = 0; i < natt; i++)
334     {
335         s2i = kT_fac/att[i].mass;
336         ti  = att[i].type;
337
338         for (j = i; j < natt; j++)
339         {
340             s2j = kT_fac/att[j].mass;
341             tj  = att[j].type;
342
343             /* Note that attractive and repulsive potentials for individual
344              * pairs will partially cancel.
345              */
346             /* -dV/dr at the cut-off for LJ + Coulomb */
347             md =
348                 md_ljd*ffp->iparams[ti*ffp->atnr+tj].lj.c6 +
349                 md_ljr*ffp->iparams[ti*ffp->atnr+tj].lj.c12 +
350                 md_el*att[i].q*att[j].q;
351
352             /* d2V/dr2 at the cut-off for Coulomb, we neglect LJ */
353             dd = dd_el*att[i].q*att[j].q;
354
355             s2  = s2i + s2j;
356
357             rsh    = r_buffer;
358             sc_fac = 1.0;
359             /* For constraints: adapt r and scaling for the Gaussian */
360             if (att[i].con)
361             {
362                 real sh, sc;
363                 approx_2dof(s2i, r_buffer*s2i/s2, &sh, &sc);
364                 rsh    += sh;
365                 sc_fac *= sc;
366             }
367             if (att[j].con)
368             {
369                 real sh, sc;
370                 approx_2dof(s2j, r_buffer*s2j/s2, &sh, &sc);
371                 rsh    += sh;
372                 sc_fac *= sc;
373             }
374
375             /* Exact contribution of an atom pair with Gaussian displacement
376              * with sigma s to the energy drift for a potential with
377              * derivative -md and second derivative dd at the cut-off.
378              * The only catch is that for potentials that change sign
379              * near the cut-off there could be an unlucky compensation
380              * of positive and negative energy drift.
381              * Such potentials are extremely rare though.
382              *
383              * Note that pot has unit energy*length, as the linear
384              * atom density still needs to be put in.
385              */
386             c_exp  = exp(-rsh*rsh/(2*s2))/sqrt(2*M_PI);
387             c_erfc = 0.5*gmx_erfc(rsh/(sqrt(2*s2)));
388             s      = sqrt(s2);
389
390             pot1 = sc_fac*
391                 md/2*((rsh*rsh + s2)*c_erfc - rsh*s*c_exp);
392             pot2 = sc_fac*
393                 dd/6*(s*(rsh*rsh + 2*s2)*c_exp - rsh*(rsh*rsh + 3*s2)*c_erfc);
394             pot = pot1 + pot2;
395
396             if (gmx_debug_at)
397             {
398                 fprintf(debug, "n %d %d d s %.3f %.3f con %d md %8.1e dd %8.1e pot1 %8.1e pot2 %8.1e pot %8.1e\n",
399                         att[i].n, att[j].n, sqrt(s2i), sqrt(s2j),
400                         att[i].con+att[j].con,
401                         md, dd, pot1, pot2, pot);
402             }
403
404             /* Multiply by the number of atom pairs */
405             if (j == i)
406             {
407                 pot *= (double)att[i].n*(att[i].n - 1)/2;
408             }
409             else
410             {
411                 pot *= (double)att[i].n*att[j].n;
412             }
413             /* We need the line density to get the energy drift of the system.
414              * The effective average r^2 is close to (rlist+sigma)^2.
415              */
416             pot *= 4*M_PI*sqr(rlist + s)/boxvol;
417
418             /* Add the unsigned drift to avoid cancellation of errors */
419             drift_tot += fabs(pot);
420         }
421     }
422
423     return drift_tot;
424 }
425
426 static real surface_frac(int cluster_size, real particle_distance, real rlist)
427 {
428     real d, area_rel;
429
430     if (rlist < 0.5*particle_distance)
431     {
432         /* We have non overlapping spheres */
433         return 1.0;
434     }
435
436     /* Half the inter-particle distance relative to rlist */
437     d = 0.5*particle_distance/rlist;
438
439     /* Determine the area of the surface at distance rlist to the closest
440      * particle, relative to surface of a sphere of radius rlist.
441      * The formulas below assume close to cubic cells for the pair search grid,
442      * which the pair search code tries to achieve.
443      * Note that in practice particle distances will not be delta distributed,
444      * but have some spread, often involving shorter distances,
445      * as e.g. O-H bonds in a water molecule. Thus the estimates below will
446      * usually be slightly too high and thus conservative.
447      */
448     switch (cluster_size)
449     {
450         case 1:
451             /* One particle: trivial */
452             area_rel = 1.0;
453             break;
454         case 2:
455             /* Two particles: two spheres at fractional distance 2*a */
456             area_rel = 1.0 + d;
457             break;
458         case 4:
459             /* We assume a perfect, symmetric tetrahedron geometry.
460              * The surface around a tetrahedron is too complex for a full
461              * analytical solution, so we use a Taylor expansion.
462              */
463             area_rel = (1.0 + 1/M_PI*(6*acos(1/sqrt(3))*d +
464                                       sqrt(3)*d*d*(1.0 +
465                                                    5.0/18.0*d*d +
466                                                    7.0/45.0*d*d*d*d +
467                                                    83.0/756.0*d*d*d*d*d*d)));
468             break;
469         default:
470             gmx_incons("surface_frac called with unsupported cluster_size");
471             area_rel = 1.0;
472     }
473
474     return area_rel/cluster_size;
475 }
476
477 void calc_verlet_buffer_size(const gmx_mtop_t *mtop, real boxvol,
478                              const t_inputrec *ir, real drift_target,
479                              const verletbuf_list_setup_t *list_setup,
480                              int *n_nonlin_vsite,
481                              real *rlist)
482 {
483     double                resolution;
484     char                 *env;
485
486     real                  particle_distance;
487     real                  nb_clust_frac_pairs_not_in_list_at_cutoff;
488
489     verletbuf_atomtype_t *att  = NULL;
490     int                   natt = -1, i;
491     double                reppow;
492     real                  md_ljd, md_ljr, md_el, dd_el;
493     real                  elfac;
494     real                  kT_fac, mass_min;
495     int                   ib0, ib1, ib;
496     real                  rb, rl;
497     real                  drift;
498
499     /* Resolution of the buffer size */
500     resolution = 0.001;
501
502     env = getenv("GMX_VERLET_BUFFER_RES");
503     if (env != NULL)
504     {
505         sscanf(env, "%lf", &resolution);
506     }
507
508     /* In an atom wise pair-list there would be no pairs in the list
509      * beyond the pair-list cut-off.
510      * However, we use a pair-list of groups vs groups of atoms.
511      * For groups of 4 atoms, the parallelism of SSE instructions, only
512      * 10% of the atoms pairs are not in the list just beyond the cut-off.
513      * As this percentage increases slowly compared to the decrease of the
514      * Gaussian displacement distribution over this range, we can simply
515      * reduce the drift by this fraction.
516      * For larger groups, e.g. of 8 atoms, this fraction will be lower,
517      * so then buffer size will be on the conservative (large) side.
518      *
519      * Note that the formulas used here do not take into account
520      * cancellation of errors which could occur by missing both
521      * attractive and repulsive interactions.
522      *
523      * The only major assumption is homogeneous particle distribution.
524      * For an inhomogeneous system, such as a liquid-vapor system,
525      * the buffer will be underestimated. The actual energy drift
526      * will be higher by the factor: local/homogeneous particle density.
527      *
528      * The results of this estimate have been checked againt simulations.
529      * In most cases the real drift differs by less than a factor 2.
530      */
531
532     /* Worst case assumption: HCP packing of particles gives largest distance */
533     particle_distance = pow(boxvol*sqrt(2)/mtop->natoms, 1.0/3.0);
534
535     get_verlet_buffer_atomtypes(mtop, &att, &natt, n_nonlin_vsite);
536     assert(att != NULL && natt >= 0);
537
538     if (debug)
539     {
540         fprintf(debug, "particle distance assuming HCP packing: %f nm\n",
541                 particle_distance);
542         fprintf(debug, "energy drift atom types: %d\n", natt);
543     }
544
545     reppow = mtop->ffparams.reppow;
546     md_ljd = 0;
547     md_ljr = 0;
548     if (ir->vdwtype == evdwCUT)
549     {
550         /* -dV/dr of -r^-6 and r^-repporw */
551         md_ljd = -6*pow(ir->rvdw, -7.0);
552         md_ljr = reppow*pow(ir->rvdw, -(reppow+1));
553         /* The contribution of the second derivative is negligible */
554     }
555     else
556     {
557         gmx_fatal(FARGS, "Energy drift calculation is only implemented for plain cut-off Lennard-Jones interactions");
558     }
559
560     elfac = ONE_4PI_EPS0/ir->epsilon_r;
561
562     /* Determine md=-dV/dr and dd=d^2V/dr^2 */
563     md_el = 0;
564     dd_el = 0;
565     if (ir->coulombtype == eelCUT || EEL_RF(ir->coulombtype))
566     {
567         real eps_rf, k_rf;
568
569         if (ir->coulombtype == eelCUT)
570         {
571             eps_rf = 1;
572             k_rf   = 0;
573         }
574         else
575         {
576             eps_rf = ir->epsilon_rf/ir->epsilon_r;
577             if (eps_rf != 0)
578             {
579                 k_rf = pow(ir->rcoulomb, -3.0)*(eps_rf - ir->epsilon_r)/(2*eps_rf + ir->epsilon_r);
580             }
581             else
582             {
583                 /* epsilon_rf = infinity */
584                 k_rf = 0.5*pow(ir->rcoulomb, -3.0);
585             }
586         }
587
588         if (eps_rf > 0)
589         {
590             md_el = elfac*(pow(ir->rcoulomb, -2.0) - 2*k_rf*ir->rcoulomb);
591         }
592         dd_el = elfac*(2*pow(ir->rcoulomb, -3.0) + 2*k_rf);
593     }
594     else if (EEL_PME(ir->coulombtype) || ir->coulombtype == eelEWALD)
595     {
596         real b, rc, br;
597
598         b     = calc_ewaldcoeff(ir->rcoulomb, ir->ewald_rtol);
599         rc    = ir->rcoulomb;
600         br    = b*rc;
601         md_el = elfac*(b*exp(-br*br)*M_2_SQRTPI/rc + gmx_erfc(br)/(rc*rc));
602         dd_el = elfac/(rc*rc)*(2*b*(1 + br*br)*exp(-br*br)*M_2_SQRTPI + 2*gmx_erfc(br)/rc);
603     }
604     else
605     {
606         gmx_fatal(FARGS, "Energy drift calculation is only implemented for Reaction-Field and Ewald electrostatics");
607     }
608
609     /* Determine the variance of the atomic displacement
610      * over nstlist-1 steps: kT_fac
611      * For inertial dynamics (not Brownian dynamics) the mass factor
612      * is not included in kT_fac, it is added later.
613      */
614     if (ir->eI == eiBD)
615     {
616         /* Get the displacement distribution from the random component only.
617          * With accurate integration the systematic (force) displacement
618          * should be negligible (unless nstlist is extremely large, which
619          * you wouldn't do anyhow).
620          */
621         kT_fac = 2*BOLTZ*ir->opts.ref_t[0]*(ir->nstlist-1)*ir->delta_t;
622         if (ir->bd_fric > 0)
623         {
624             /* This is directly sigma^2 of the displacement */
625             kT_fac /= ir->bd_fric;
626
627             /* Set the masses to 1 as kT_fac is the full sigma^2,
628              * but we divide by m in ener_drift().
629              */
630             for (i = 0; i < natt; i++)
631             {
632                 att[i].mass = 1;
633             }
634         }
635         else
636         {
637             real tau_t;
638
639             /* Per group tau_t is not implemented yet, use the maximum */
640             tau_t = ir->opts.tau_t[0];
641             for (i = 1; i < ir->opts.ngtc; i++)
642             {
643                 tau_t = max(tau_t, ir->opts.tau_t[i]);
644             }
645
646             kT_fac *= tau_t;
647             /* This kT_fac needs to be divided by the mass to get sigma^2 */
648         }
649     }
650     else
651     {
652         kT_fac = BOLTZ*ir->opts.ref_t[0]*sqr((ir->nstlist-1)*ir->delta_t);
653     }
654
655     mass_min = att[0].mass;
656     for (i = 1; i < natt; i++)
657     {
658         mass_min = min(mass_min, att[i].mass);
659     }
660
661     if (debug)
662     {
663         fprintf(debug, "md_ljd %e md_ljr %e\n", md_ljd, md_ljr);
664         fprintf(debug, "md_el %e dd_el %e\n", md_el, dd_el);
665         fprintf(debug, "sqrt(kT_fac) %f\n", sqrt(kT_fac));
666         fprintf(debug, "mass_min %f\n", mass_min);
667     }
668
669     /* Search using bisection */
670     ib0 = -1;
671     /* The drift will be neglible at 5 times the max sigma */
672     ib1 = (int)(5*2*sqrt(kT_fac/mass_min)/resolution) + 1;
673     while (ib1 - ib0 > 1)
674     {
675         ib = (ib0 + ib1)/2;
676         rb = ib*resolution;
677         rl = max(ir->rvdw, ir->rcoulomb) + rb;
678
679         /* Calculate the average energy drift at the last step
680          * of the nstlist steps at which the pair-list is used.
681          */
682         drift = ener_drift(att, natt, &mtop->ffparams,
683                            kT_fac,
684                            md_ljd, md_ljr, md_el, dd_el, rb,
685                            rl, boxvol);
686
687         /* Correct for the fact that we are using a Ni x Nj particle pair list
688          * and not a 1 x 1 particle pair list. This reduces the drift.
689          */
690         /* We don't have a formula for 8 (yet), use 4 which is conservative */
691         nb_clust_frac_pairs_not_in_list_at_cutoff =
692             surface_frac(min(list_setup->cluster_size_i, 4),
693                          particle_distance, rl)*
694             surface_frac(min(list_setup->cluster_size_j, 4),
695                          particle_distance, rl);
696         drift *= nb_clust_frac_pairs_not_in_list_at_cutoff;
697
698         /* Convert the drift to drift per unit time per atom */
699         drift /= ir->nstlist*ir->delta_t*mtop->natoms;
700
701         if (debug)
702         {
703             fprintf(debug, "ib %3d %3d %3d rb %.3f %dx%d fac %.3f drift %f\n",
704                     ib0, ib, ib1, rb,
705                     list_setup->cluster_size_i, list_setup->cluster_size_j,
706                     nb_clust_frac_pairs_not_in_list_at_cutoff,
707                     drift);
708         }
709
710         if (fabs(drift) > drift_target)
711         {
712             ib0 = ib;
713         }
714         else
715         {
716             ib1 = ib;
717         }
718     }
719
720     sfree(att);
721
722     *rlist = max(ir->rvdw, ir->rcoulomb) + ib1*resolution;
723 }