OCL: Make variables const
[gromacs.git] / src / gromacs / mdlib / ebin.cpp
blob45633bbf7eda11dfcc499318c1737bfbc323cd1e
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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,2017,2018, 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.
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.
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.
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.
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.
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.
37 /* This file is completely threadsafe - keep it that way! */
38 #include "gmxpre.h"
40 #include "ebin.h"
42 #include <cmath>
43 #include <cstring>
45 #include "gromacs/math/units.h"
46 #include "gromacs/math/utilities.h"
47 #include "gromacs/math/vec.h"
48 #include "gromacs/topology/ifunc.h"
49 #include "gromacs/trajectory/energyframe.h"
50 #include "gromacs/utility/cstringutil.h"
51 #include "gromacs/utility/fatalerror.h"
52 #include "gromacs/utility/smalloc.h"
54 t_ebin *mk_ebin()
56 t_ebin *eb;
58 snew(eb, 1);
60 return eb;
63 void done_ebin(t_ebin *eb)
65 for (int i = 0; i < eb->nener; i++)
67 sfree(eb->enm[i].name);
68 sfree(eb->enm[i].unit);
70 sfree(eb->e);
71 sfree(eb->e_sim);
72 sfree(eb->enm);
75 int get_ebin_space(t_ebin *eb, int nener, const char *const enm[], const char *unit)
77 int index;
78 int i, f;
79 const char *u;
81 index = eb->nener;
82 eb->nener += nener;
83 srenew(eb->e, eb->nener);
84 srenew(eb->e_sim, eb->nener);
85 srenew(eb->enm, eb->nener);
86 for (i = index; (i < eb->nener); i++)
88 eb->e[i].e = 0;
89 eb->e[i].eav = 0;
90 eb->e[i].esum = 0;
91 eb->e_sim[i].e = 0;
92 eb->e_sim[i].eav = 0;
93 eb->e_sim[i].esum = 0;
94 eb->enm[i].name = gmx_strdup(enm[i-index]);
95 if (unit != nullptr)
97 eb->enm[i].unit = gmx_strdup(unit);
99 else
101 /* Determine the unit from the longname.
102 * These units should have been defined in ifunc.c
103 * But even better would be if all interactions functions
104 * return energies and all non-interaction function
105 * entries would be removed from the ifunc array.
107 u = unit_energy;
108 for (f = 0; f < F_NRE; f++)
110 if (strcmp(eb->enm[i].name,
111 interaction_function[f].longname) == 0)
113 /* Only the terms in this list are not energies */
114 switch (f)
116 case F_DISRESVIOL: u = unit_length; break;
117 case F_ORIRESDEV: u = "obs"; break;
118 case F_TEMP: u = unit_temp_K; break;
119 case F_PDISPCORR:
120 case F_PRES: u = unit_pres_bar; break;
124 eb->enm[i].unit = gmx_strdup(u);
128 return index;
131 void add_ebin(t_ebin *eb, int index, int nener, const real ener[], gmx_bool bSum)
133 int i, m;
134 double e, invmm, diff;
135 t_energy *eg, *egs;
137 if ((index+nener > eb->nener) || (index < 0))
139 gmx_fatal(FARGS, "%s-%d: Energies out of range: index=%d nener=%d maxener=%d",
140 __FILE__, __LINE__, index, nener, eb->nener);
143 eg = &(eb->e[index]);
145 for (i = 0; (i < nener); i++)
147 eg[i].e = ener[i];
150 if (bSum)
152 egs = &(eb->e_sim[index]);
154 m = eb->nsum;
156 if (m == 0)
158 for (i = 0; (i < nener); i++)
160 eg[i].eav = 0;
161 eg[i].esum = ener[i];
162 egs[i].esum += ener[i];
165 else
167 invmm = (1.0/m)/(m+1.0);
169 for (i = 0; (i < nener); i++)
171 /* Value for this component */
172 e = ener[i];
174 /* first update sigma, then sum */
175 diff = eg[i].esum - m*e;
176 eg[i].eav += diff*diff*invmm;
177 eg[i].esum += e;
178 egs[i].esum += e;
184 void ebin_increase_count(t_ebin *eb, gmx_bool bSum)
186 eb->nsteps++;
187 eb->nsteps_sim++;
189 if (bSum)
191 eb->nsum++;
192 eb->nsum_sim++;
196 void reset_ebin_sums(t_ebin *eb)
198 eb->nsteps = 0;
199 eb->nsum = 0;
200 /* The actual sums are cleared when the next frame is stored */
203 void pr_ebin(FILE *fp, t_ebin *eb, int index, int nener, int nperline,
204 int prmode, gmx_bool bPrHead)
206 int i, j, i0;
207 int rc;
208 char buf[30];
210 rc = 0;
212 if (index < 0 || index > eb->nener)
214 gmx_fatal(FARGS, "Invalid index in pr_ebin: %d", index);
216 int start = index;
217 if (nener > eb->nener)
219 gmx_fatal(FARGS, "Invalid nener in pr_ebin: %d", nener);
221 int end = eb->nener;
222 if (nener != -1)
224 end = index + nener;
226 for (i = start; (i < end) && rc >= 0; )
228 if (bPrHead)
230 i0 = i;
231 for (j = 0; (j < nperline) && (i < end) && rc >= 0; j++, i++)
233 if (strncmp(eb->enm[i].name, "Pres", 4) == 0)
235 /* Print the pressure unit to avoid confusion */
236 sprintf(buf, "%s (%s)", eb->enm[i].name, unit_pres_bar);
237 rc = fprintf(fp, "%15s", buf);
239 else
241 rc = fprintf(fp, "%15s", eb->enm[i].name);
245 if (rc >= 0)
247 rc = fprintf(fp, "\n");
250 i = i0;
252 for (j = 0; (j < nperline) && (i < end) && rc >= 0; j++, i++)
254 switch (prmode)
256 case eprNORMAL:
257 rc = fprintf(fp, " %12.5e", eb->e[i].e);
258 break;
259 case eprAVER:
260 if (eb->nsum_sim > 0)
262 rc = fprintf(fp, " %12.5e", eb->e_sim[i].esum/eb->nsum_sim);
264 else
266 rc = fprintf(fp, " %-12s", "N/A");
268 break;
269 default: gmx_fatal(FARGS, "Invalid print mode %d in pr_ebin",
270 prmode);
274 if (rc >= 0)
276 rc = fprintf(fp, "\n");
279 if (rc < 0)
281 gmx_fatal(FARGS, "Cannot write to logfile; maybe you are out of disk space?");