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