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