Update instructions in containers.rst
[gromacs.git] / src / gromacs / mdlib / ebin.cpp
blob6c830526bea1576778b38cc026603480c4ccdb44
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.
7 * Copyright (c) 2019,2020, by the GROMACS development team, led by
8 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
9 * and including many others, as listed in the AUTHORS file in the
10 * top-level source directory and at http://www.gromacs.org.
12 * GROMACS is free software; you can redistribute it and/or
13 * modify it under the terms of the GNU Lesser General Public License
14 * as published by the Free Software Foundation; either version 2.1
15 * of the License, or (at your option) any later version.
17 * GROMACS is distributed in the hope that it will be useful,
18 * but WITHOUT ANY WARRANTY; without even the implied warranty of
19 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
20 * Lesser General Public License for more details.
22 * You should have received a copy of the GNU Lesser General Public
23 * License along with GROMACS; if not, see
24 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
25 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
27 * If you want to redistribute modifications to GROMACS, please
28 * consider that scientific software is very special. Version
29 * control is crucial - bugs must be traceable. We will be happy to
30 * consider code for inclusion in the official distribution, but
31 * derived work must not be called official GROMACS. Details are found
32 * in the README & COPYING files - if they are missing, get the
33 * official version at http://www.gromacs.org.
35 * To help us fund GROMACS development, we humbly ask that you cite
36 * the research papers on the package. Check out http://www.gromacs.org.
38 /* This file is completely threadsafe - keep it that way! */
39 #include "gmxpre.h"
41 #include "ebin.h"
43 #include <cmath>
44 #include <cstring>
46 #include <algorithm>
48 #include "gromacs/math/units.h"
49 #include "gromacs/math/utilities.h"
50 #include "gromacs/math/vec.h"
51 #include "gromacs/topology/ifunc.h"
52 #include "gromacs/trajectory/energyframe.h"
53 #include "gromacs/utility/arrayref.h"
54 #include "gromacs/utility/basedefinitions.h"
55 #include "gromacs/utility/cstringutil.h"
56 #include "gromacs/utility/fatalerror.h"
57 #include "gromacs/utility/smalloc.h"
58 #include "gromacs/utility/stringutil.h"
60 t_ebin* mk_ebin()
62 t_ebin* eb;
64 snew(eb, 1);
66 return eb;
69 void done_ebin(t_ebin* eb)
71 for (int i = 0; i < eb->nener; i++)
73 sfree(eb->enm[i].name);
74 sfree(eb->enm[i].unit);
76 sfree(eb->e);
77 sfree(eb->e_sim);
78 sfree(eb->enm);
79 sfree(eb);
82 int get_ebin_space(t_ebin* eb, int nener, const char* const enm[], const char* unit)
84 int index;
85 int i, f;
86 const char* u;
88 index = eb->nener;
89 eb->nener += nener;
90 srenew(eb->e, eb->nener);
91 srenew(eb->e_sim, eb->nener);
92 srenew(eb->enm, eb->nener);
93 for (i = index; (i < eb->nener); i++)
95 eb->e[i].e = 0;
96 eb->e[i].eav = 0;
97 eb->e[i].esum = 0;
98 eb->e_sim[i].e = 0;
99 eb->e_sim[i].eav = 0;
100 eb->e_sim[i].esum = 0;
101 eb->enm[i].name = gmx_strdup(enm[i - index]);
102 if (unit != nullptr)
104 eb->enm[i].unit = gmx_strdup(unit);
106 else
108 /* Determine the unit from the longname.
109 * These units should have been defined in ifunc.c
110 * But even better would be if all interactions functions
111 * return energies and all non-interaction function
112 * entries would be removed from the ifunc array.
114 u = unit_energy;
115 for (f = 0; f < F_NRE; f++)
117 if (strcmp(eb->enm[i].name, interaction_function[f].longname) == 0)
119 /* Only the terms in this list are not energies */
120 switch (f)
122 case F_DISRESVIOL: u = unit_length; break;
123 case F_ORIRESDEV: u = "obs"; break;
124 case F_TEMP: u = unit_temp_K; break;
125 case F_PDISPCORR:
126 case F_PRES: u = unit_pres_bar; break;
130 eb->enm[i].unit = gmx_strdup(u);
134 return index;
137 // ICC 19 -O3 -msse2 generates wrong code. Lower optimization levels
138 // and other SIMD levels seem fine, however.
139 #if defined __ICC
140 # pragma intel optimization_level 2
141 #endif
142 void add_ebin(t_ebin* eb, int entryIndex, int nener, const real ener[], gmx_bool bSum)
144 int i, m;
145 double e, invmm, diff;
146 t_energy *eg, *egs;
148 if ((entryIndex + nener > eb->nener) || (entryIndex < 0))
150 gmx_fatal(FARGS, "%s-%d: Energies out of range: entryIndex=%d nener=%d maxener=%d",
151 __FILE__, __LINE__, entryIndex, nener, eb->nener);
154 eg = &(eb->e[entryIndex]);
156 for (i = 0; (i < nener); i++)
158 eg[i].e = ener[i];
161 if (bSum)
163 egs = &(eb->e_sim[entryIndex]);
165 m = eb->nsum;
167 if (m == 0)
169 for (i = 0; (i < nener); i++)
171 eg[i].eav = 0;
172 eg[i].esum = ener[i];
173 egs[i].esum += ener[i];
176 else
178 invmm = (1.0 / m) / (m + 1.0);
180 for (i = 0; (i < nener); i++)
182 /* Value for this component */
183 e = ener[i];
185 /* first update sigma, then sum */
186 diff = eg[i].esum - m * e;
187 eg[i].eav += diff * diff * invmm;
188 eg[i].esum += e;
189 egs[i].esum += e;
195 // TODO It would be faster if this function was templated on both bSum
196 // and whether eb->nsum was zero, to lift the branches out of the loop
197 // over all possible energy terms, but that is true for a lot of the
198 // ebin and mdebin functionality, so we should do it all or nothing.
199 void add_ebin_indexed(t_ebin* eb,
200 int entryIndex,
201 gmx::ArrayRef<bool> shouldUse,
202 gmx::ArrayRef<const real> ener,
203 gmx_bool bSum)
206 GMX_ASSERT(shouldUse.size() == ener.size(), "View sizes must match");
207 GMX_ASSERT(entryIndex + std::count(shouldUse.begin(), shouldUse.end(), true) <= eb->nener,
208 gmx::formatString("Energies out of range: entryIndex=%d nener=%td maxener=%d", entryIndex,
209 std::count(shouldUse.begin(), shouldUse.end(), true), eb->nener)
210 .c_str());
211 GMX_ASSERT(entryIndex >= 0, "Must have non-negative entry");
213 const int m = eb->nsum;
214 const double invmm = (m == 0) ? 0 : (1.0 / m) / (m + 1.0);
215 t_energy* energyEntry = &(eb->e[entryIndex]);
216 t_energy* simEnergyEntry = &(eb->e_sim[entryIndex]);
217 auto shouldUseIter = shouldUse.begin();
218 for (const auto& theEnergy : ener)
220 if (*shouldUseIter)
222 energyEntry->e = theEnergy;
223 if (bSum)
225 if (m == 0)
227 energyEntry->eav = 0;
228 energyEntry->esum = theEnergy;
229 simEnergyEntry->esum += theEnergy;
231 else
233 /* first update sigma, then sum */
234 double diff = energyEntry->esum - m * theEnergy;
235 energyEntry->eav += diff * diff * invmm;
236 energyEntry->esum += theEnergy;
237 simEnergyEntry->esum += theEnergy;
239 ++simEnergyEntry;
241 ++energyEntry;
243 ++shouldUseIter;
247 void ebin_increase_count(int increment, t_ebin* eb, gmx_bool bSum)
249 eb->nsteps += increment;
250 eb->nsteps_sim += increment;
252 if (bSum)
254 eb->nsum += increment;
255 eb->nsum_sim += increment;
259 void reset_ebin_sums(t_ebin* eb)
261 eb->nsteps = 0;
262 eb->nsum = 0;
263 /* The actual sums are cleared when the next frame is stored */
266 void pr_ebin(FILE* fp, t_ebin* eb, int entryIndex, int nener, int nperline, int prmode, gmx_bool bPrHead)
268 int i, j, i0;
269 int rc;
270 char buf[30];
272 rc = 0;
274 if (entryIndex < 0 || entryIndex > eb->nener)
276 gmx_fatal(FARGS, "Invalid entryIndex in pr_ebin: %d", entryIndex);
278 int start = entryIndex;
279 if (nener > eb->nener)
281 gmx_fatal(FARGS, "Invalid nener in pr_ebin: %d", nener);
283 int end = eb->nener;
284 if (nener != -1)
286 end = entryIndex + nener;
288 for (i = start; (i < end) && rc >= 0;)
290 if (bPrHead)
292 i0 = i;
293 for (j = 0; (j < nperline) && (i < end) && rc >= 0; j++, i++)
295 if (strncmp(eb->enm[i].name, "Pres", 4) == 0)
297 /* Print the pressure unit to avoid confusion */
298 sprintf(buf, "%s (%s)", eb->enm[i].name, unit_pres_bar);
299 rc = fprintf(fp, "%15s", buf);
301 else
303 rc = fprintf(fp, "%15s", eb->enm[i].name);
307 if (rc >= 0)
309 rc = fprintf(fp, "\n");
312 i = i0;
314 for (j = 0; (j < nperline) && (i < end) && rc >= 0; j++, i++)
316 switch (prmode)
318 case eprNORMAL: rc = fprintf(fp, " %12.5e", eb->e[i].e); break;
319 case eprAVER:
320 if (eb->nsum_sim > 0)
322 rc = fprintf(fp, " %12.5e", eb->e_sim[i].esum / eb->nsum_sim);
324 else
326 rc = fprintf(fp, " %-12s", "N/A");
328 break;
329 default: gmx_fatal(FARGS, "Invalid print mode %d in pr_ebin", prmode);
332 if (rc >= 0)
334 rc = fprintf(fp, "\n");
337 if (rc < 0)
339 gmx_fatal(FARGS, "Cannot write to logfile; maybe you are out of disk space?");