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