a7c4e38930af2502cb31158981175da4e0452444
[alexxy/gromacs.git] / src / gromacs / mdlib / nlistheuristics.c
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
5  * Copyright (c) 2001-2004, The GROMACS development team.
6  * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
7  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
8  * and including many others, as listed in the AUTHORS file in the
9  * top-level source directory and at http://www.gromacs.org.
10  *
11  * GROMACS is free software; you can redistribute it and/or
12  * modify it under the terms of the GNU Lesser General Public License
13  * as published by the Free Software Foundation; either version 2.1
14  * of the License, or (at your option) any later version.
15  *
16  * GROMACS is distributed in the hope that it will be useful,
17  * but WITHOUT ANY WARRANTY; without even the implied warranty of
18  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
19  * Lesser General Public License for more details.
20  *
21  * You should have received a copy of the GNU Lesser General Public
22  * License along with GROMACS; if not, see
23  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
24  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
25  *
26  * If you want to redistribute modifications to GROMACS, please
27  * consider that scientific software is very special. Version
28  * control is crucial - bugs must be traceable. We will be happy to
29  * consider code for inclusion in the official distribution, but
30  * derived work must not be called official GROMACS. Details are found
31  * in the README & COPYING files - if they are missing, get the
32  * official version at http://www.gromacs.org.
33  *
34  * To help us fund GROMACS development, we humbly ask that you cite
35  * the research papers on the package. Check out http://www.gromacs.org.
36  */
37 #include "gmxpre.h"
38
39 #include "config.h"
40
41 #include "gromacs/legacyheaders/typedefs.h"
42 #include "gromacs/legacyheaders/types/nlistheuristics.h"
43
44 #include "gromacs/math/vec.h"
45 #include "gromacs/utility/cstringutil.h"
46 #include "gromacs/utility/fatalerror.h"
47
48 void reset_nlistheuristics(gmx_nlheur_t *nlh, gmx_int64_t step)
49 {
50     nlh->lt_runav     = 0;
51     nlh->lt_runav2    = 0;
52     nlh->step_nscheck = step;
53 }
54
55 void init_nlistheuristics(gmx_nlheur_t *nlh,
56                           gmx_bool bGStatEveryStep, gmx_int64_t step)
57 {
58     nlh->bGStatEveryStep = bGStatEveryStep;
59     nlh->nns             = 0;
60     nlh->nabnsb          = 0;
61     nlh->s1              = 0;
62     nlh->s2              = 0;
63     nlh->ab              = 0;
64
65     reset_nlistheuristics(nlh, step);
66 }
67
68 void update_nliststatistics(gmx_nlheur_t *nlh, gmx_int64_t step)
69 {
70     gmx_int64_t     nl_lt;
71     char            sbuf[STEPSTRSIZE], sbuf2[STEPSTRSIZE];
72
73     /* Determine the neighbor list life time */
74     nl_lt = step - nlh->step_ns;
75     if (debug)
76     {
77         fprintf(debug, "%d atoms beyond ns buffer, updating neighbor list after %s steps\n", nlh->nabnsb, gmx_step_str(nl_lt, sbuf));
78     }
79     nlh->nns++;
80     nlh->s1 += nl_lt;
81     nlh->s2 += nl_lt*nl_lt;
82     nlh->ab += nlh->nabnsb;
83     if (nlh->lt_runav == 0)
84     {
85         nlh->lt_runav  = nl_lt;
86         /* Initialize the fluctuation average
87          * such that at startup we check after 0 steps.
88          */
89         nlh->lt_runav2 = sqr(nl_lt/2.0);
90     }
91     /* Running average with 0.9 gives an exp. history of 9.5 */
92     nlh->lt_runav2 = 0.9*nlh->lt_runav2 + 0.1*sqr(nlh->lt_runav - nl_lt);
93     nlh->lt_runav  = 0.9*nlh->lt_runav  + 0.1*nl_lt;
94     if (nlh->bGStatEveryStep)
95     {
96         /* Always check the nlist validity */
97         nlh->step_nscheck = step;
98     }
99     else
100     {
101         /* We check after:  <life time> - 2*sigma
102          * The factor 2 is quite conservative,
103          * but we assume that with nstlist=-1 the user
104          * prefers exact integration over performance.
105          */
106         nlh->step_nscheck = step
107             + (int)(nlh->lt_runav - 2.0*sqrt(nlh->lt_runav2)) - 1;
108     }
109     if (debug)
110     {
111         fprintf(debug, "nlist life time %s run av. %4.1f sig %3.1f check %s check with -gcom %d\n",
112                 gmx_step_str(nl_lt, sbuf), nlh->lt_runav, sqrt(nlh->lt_runav2),
113                 gmx_step_str(nlh->step_nscheck-step+1, sbuf2),
114                 (int)(nlh->lt_runav - 2.0*sqrt(nlh->lt_runav2)));
115     }
116 }
117
118 void set_nlistheuristics(gmx_nlheur_t *nlh, gmx_bool bReset, gmx_int64_t step)
119 {
120     int d;
121
122     if (bReset)
123     {
124         reset_nlistheuristics(nlh, step);
125     }
126     else
127     {
128         update_nliststatistics(nlh, step);
129     }
130
131     nlh->step_ns = step;
132     /* Initialize the cumulative coordinate scaling matrix */
133     clear_mat(nlh->scale_tot);
134     for (d = 0; d < DIM; d++)
135     {
136         nlh->scale_tot[d][d] = 1.0;
137     }
138 }