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