33a0742130c8aaf96893ad22a087d43c134bf7f9
[alexxy/gromacs.git] / src / gromacs / mdlib / ebin.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,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 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
39
40 #include <math.h>
41 #include <string.h>
42 #include "gromacs/utility/smalloc.h"
43 #include "gromacs/legacyheaders/typedefs.h"
44 #include "gromacs/utility/fatalerror.h"
45 #include "gromacs/legacyheaders/ebin.h"
46 #include "gromacs/math/utilities.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/math/units.h"
49
50 t_ebin *mk_ebin(void)
51 {
52     t_ebin *eb;
53
54     snew(eb, 1);
55
56     return eb;
57 }
58
59 int get_ebin_space(t_ebin *eb, int nener, const char *enm[], const char *unit)
60 {
61     int         index;
62     int         i, f;
63     const char *u;
64
65     index      = eb->nener;
66     eb->nener += nener;
67     srenew(eb->e, eb->nener);
68     srenew(eb->e_sim, eb->nener);
69     srenew(eb->enm, eb->nener);
70     for (i = index; (i < eb->nener); i++)
71     {
72         eb->e[i].e        = 0;
73         eb->e[i].eav      = 0;
74         eb->e[i].esum     = 0;
75         eb->e_sim[i].e    = 0;
76         eb->e_sim[i].eav  = 0;
77         eb->e_sim[i].esum = 0;
78         eb->enm[i].name   = gmx_strdup(enm[i-index]);
79         if (unit != NULL)
80         {
81             eb->enm[i].unit = gmx_strdup(unit);
82         }
83         else
84         {
85             /* Determine the unit from the longname.
86              * These units should have been defined in ifunc.c
87              * But even better would be if all interactions functions
88              * return energies and all non-interaction function
89              * entries would be removed from the ifunc array.
90              */
91             u = unit_energy;
92             for (f = 0; f < F_NRE; f++)
93             {
94                 if (strcmp(eb->enm[i].name,
95                            interaction_function[f].longname) == 0)
96                 {
97                     /* Only the terms in this list are not energies */
98                     switch (f)
99                     {
100                         case F_DISRESVIOL: u = unit_length;   break;
101                         case F_ORIRESDEV:  u = "obs";         break;
102                         case F_TEMP:       u = unit_temp_K;   break;
103                         case F_PDISPCORR:
104                         case F_PRES:       u = unit_pres_bar; break;
105                     }
106                 }
107             }
108             eb->enm[i].unit = gmx_strdup(u);
109         }
110     }
111
112     return index;
113 }
114
115 void add_ebin(t_ebin *eb, int index, int nener, real ener[], gmx_bool bSum)
116 {
117     int       i, m;
118     double    e, sum, sigma, invmm, diff;
119     t_energy *eg, *egs;
120
121     if ((index+nener > eb->nener) || (index < 0))
122     {
123         gmx_fatal(FARGS, "%s-%d: Energies out of range: index=%d nener=%d maxener=%d",
124                   __FILE__, __LINE__, index, nener, eb->nener);
125     }
126
127     eg = &(eb->e[index]);
128
129     for (i = 0; (i < nener); i++)
130     {
131         eg[i].e      = ener[i];
132     }
133
134     if (bSum)
135     {
136         egs = &(eb->e_sim[index]);
137
138         m = eb->nsum;
139
140         if (m == 0)
141         {
142             for (i = 0; (i < nener); i++)
143             {
144                 eg[i].eav    = 0;
145                 eg[i].esum   = ener[i];
146                 egs[i].esum += ener[i];
147             }
148         }
149         else
150         {
151             invmm = (1.0/(double)m)/((double)m+1.0);
152
153             for (i = 0; (i < nener); i++)
154             {
155                 /* Value for this component */
156                 e = ener[i];
157
158                 /* first update sigma, then sum */
159                 diff         = eg[i].esum - m*e;
160                 eg[i].eav   += diff*diff*invmm;
161                 eg[i].esum  += e;
162                 egs[i].esum += e;
163             }
164         }
165     }
166 }
167
168 void ebin_increase_count(t_ebin *eb, gmx_bool bSum)
169 {
170     eb->nsteps++;
171     eb->nsteps_sim++;
172
173     if (bSum)
174     {
175         eb->nsum++;
176         eb->nsum_sim++;
177     }
178 }
179
180 void reset_ebin_sums(t_ebin *eb)
181 {
182     eb->nsteps = 0;
183     eb->nsum   = 0;
184     /* The actual sums are cleared when the next frame is stored */
185 }
186
187 void pr_ebin(FILE *fp, t_ebin *eb, int index, int nener, int nperline,
188              int prmode, gmx_bool bPrHead)
189 {
190     int  i, j, i0;
191     real ee = 0;
192     int  rc;
193     char buf[30];
194
195     rc = 0;
196
197     if (index < 0)
198     {
199         gmx_fatal(FARGS, "Invalid index in pr_ebin: %d", index);
200     }
201     if (nener == -1)
202     {
203         nener = eb->nener;
204     }
205     else
206     {
207         nener = index + nener;
208     }
209     for (i = index; (i < nener) && rc >= 0; )
210     {
211         if (bPrHead)
212         {
213             i0 = i;
214             for (j = 0; (j < nperline) && (i < nener) && rc >= 0; j++, i++)
215             {
216                 if (strncmp(eb->enm[i].name, "Pres", 4) == 0)
217                 {
218                     /* Print the pressure unit to avoid confusion */
219                     sprintf(buf, "%s (%s)", eb->enm[i].name, unit_pres_bar);
220                     rc = fprintf(fp, "%15s", buf);
221                 }
222                 else
223                 {
224                     rc = fprintf(fp, "%15s", eb->enm[i].name);
225                 }
226             }
227
228             if (rc >= 0)
229             {
230                 rc = fprintf(fp, "\n");
231             }
232
233             i = i0;
234         }
235         for (j = 0; (j < nperline) && (i < nener) && rc >= 0; j++, i++)
236         {
237             switch (prmode)
238             {
239                 case eprNORMAL: ee = eb->e[i].e; break;
240                 case eprAVER:   ee = eb->e_sim[i].esum/eb->nsum_sim; break;
241                 default: gmx_fatal(FARGS, "Invalid print mode %d in pr_ebin",
242                                    prmode);
243             }
244
245             rc = fprintf(fp, "   %12.5e", ee);
246         }
247         if (rc >= 0)
248         {
249             rc = fprintf(fp, "\n");
250         }
251     }
252     if (rc < 0)
253     {
254         gmx_fatal(FARGS, "Cannot write to logfile; maybe you are out of disk space?");
255     }
256 }
257
258 #ifdef DEBUGEBIN
259 int main(int argc, char *argv[])
260 {
261 #define NE 12
262 #define NT 7
263 #define NS 5
264
265     t_ebin *eb;
266     int     i;
267     char    buf[25];
268     char   *ce[NE], *ct[NT], *cs[NS];
269     real    e[NE], t[NT], s[NS];
270     int     ie, it, is;
271
272     eb = mk_ebin();
273     for (i = 0; (i < NE); i++)
274     {
275         e[i] = i;
276         sprintf(buf, "e%d", i);
277         ce[i] = gmx_strdup(buf);
278     }
279     ie = get_ebin_space(eb, NE, ce);
280     add_ebin(eb, ie, NE, e, 0);
281     for (i = 0; (i < NS); i++)
282     {
283         s[i] = i;
284         sprintf(buf, "s%d", i);
285         cs[i] = gmx_strdup(buf);
286     }
287     is = get_ebin_space(eb, NS, cs);
288     add_ebin(eb, is, NS, s, 0);
289     for (i = 0; (i < NT); i++)
290     {
291         t[i] = i;
292         sprintf(buf, "t%d", i);
293         ct[i] = gmx_strdup(buf);
294     }
295     it = get_ebin_space(eb, NT, ct);
296     add_ebin(eb, it, NT, t, 0);
297
298     printf("Normal:\n");
299     pr_ebin(stdout, eb, 0, -1, 5, eprNORMAL, 1);
300
301     printf("Average:\n");
302     pr_ebin(stdout, eb, ie, NE, 5, eprAVER, 1);
303     pr_ebin(stdout, eb, is, NS, 3, eprAVER, 1);
304     pr_ebin(stdout, eb, it, NT, 4, eprAVER, 1);
305 }
306 #endif