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