Merge remote-tracking branch 'gerrit/release-4-6'
[alexxy/gromacs.git] / src / gromacs / mdlib / nlistheuristics.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.0
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 "typedefs.h"
41 #include "gmx_fatal.h"
42 #include "vec.h"
43
44 void reset_nlistheuristics(gmx_nlheur_t *nlh,gmx_large_int_t step)
45 {
46     nlh->lt_runav  = 0;
47     nlh->lt_runav2 = 0;
48     nlh->step_nscheck = step;
49 }
50
51 void init_nlistheuristics(gmx_nlheur_t *nlh,
52                           gmx_bool bGStatEveryStep,gmx_large_int_t step)
53 {
54     nlh->bGStatEveryStep = bGStatEveryStep;
55     nlh->nns       = 0;
56     nlh->nabnsb    = 0;
57     nlh->s1        = 0;
58     nlh->s2        = 0;
59     nlh->ab        = 0;
60
61     reset_nlistheuristics(nlh,step);
62 }
63
64 void update_nliststatistics(gmx_nlheur_t *nlh,gmx_large_int_t step)
65 {
66     gmx_large_int_t nl_lt;
67     char sbuf[STEPSTRSIZE],sbuf2[STEPSTRSIZE];
68
69     /* Determine the neighbor list life time */
70     nl_lt = step - nlh->step_ns;
71     if (debug)
72     {
73         fprintf(debug,"%d atoms beyond ns buffer, updating neighbor list after %s steps\n",nlh->nabnsb,gmx_step_str(nl_lt,sbuf));
74     }
75     nlh->nns++;
76     nlh->s1 += nl_lt;
77     nlh->s2 += nl_lt*nl_lt;
78     nlh->ab += nlh->nabnsb;
79     if (nlh->lt_runav == 0)
80     {
81         nlh->lt_runav  = nl_lt;
82         /* Initialize the fluctuation average
83          * such that at startup we check after 0 steps.
84          */
85         nlh->lt_runav2 = sqr(nl_lt/2.0);
86     }
87     /* Running average with 0.9 gives an exp. history of 9.5 */
88     nlh->lt_runav2 = 0.9*nlh->lt_runav2 + 0.1*sqr(nlh->lt_runav - nl_lt);
89     nlh->lt_runav  = 0.9*nlh->lt_runav  + 0.1*nl_lt;
90     if (nlh->bGStatEveryStep)
91     {
92         /* Always check the nlist validity */
93         nlh->step_nscheck = step;
94     }
95     else
96     {
97         /* We check after:  <life time> - 2*sigma
98          * The factor 2 is quite conservative,
99          * but we assume that with nstlist=-1 the user
100          * prefers exact integration over performance.
101          */
102         nlh->step_nscheck = step
103                   + (int)(nlh->lt_runav - 2.0*sqrt(nlh->lt_runav2)) - 1;
104     }
105     if (debug)
106     {
107         fprintf(debug,"nlist life time %s run av. %4.1f sig %3.1f check %s check with -gcom %d\n",
108                 gmx_step_str(nl_lt,sbuf),nlh->lt_runav,sqrt(nlh->lt_runav2),
109                 gmx_step_str(nlh->step_nscheck-step+1,sbuf2),
110                 (int)(nlh->lt_runav - 2.0*sqrt(nlh->lt_runav2)));
111     }
112 }
113
114 void set_nlistheuristics(gmx_nlheur_t *nlh,gmx_bool bReset,gmx_large_int_t step)
115 {
116     int d;
117
118     if (bReset)
119     {
120         reset_nlistheuristics(nlh,step);
121     }
122     else
123     {
124         update_nliststatistics(nlh,step);
125     }
126
127     nlh->step_ns = step;
128     /* Initialize the cumulative coordinate scaling matrix */
129     clear_mat(nlh->scale_tot);
130     for(d=0; d<DIM; d++)
131     {
132         nlh->scale_tot[d][d] = 1.0;
133     }
134 }