Update instructions in containers.rst
[gromacs.git] / src / gromacs / mdlib / ebin.h
blobac76a2aade614fcffc7962ec951a2d61cb968ffe
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) 2013,2014,2015,2018,2019,2020, 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 /* \internal \file
38 * \brief
39 * Declares structures and interfaces to store, compute and print current and average values for thermodynamics properties.
41 * The word 'energy' is used here in wide scope and refer to any thermodynamic quantity that can benefit from
42 * averaging (e.g. temperature, pressure).
44 * \ingroup module_mdlib
46 #ifndef GMX_MDLIB_EBIN_H
47 #define GMX_MDLIB_EBIN_H
49 #include <stdio.h>
51 #include "gromacs/fileio/enxio.h"
52 #include "gromacs/trajectory/energyframe.h"
53 #include "gromacs/utility/basedefinitions.h"
55 namespace gmx
57 template<typename>
58 class ArrayRef;
61 /* \brief Running averaging structure ('energy bin') to store thermodynamic values.
63 * Collects the data on thermodynamic parameters (energy terms, temperature,
64 * pressure etc.) during the run, including their current and average values.
66 * \todo Clean this structure from unused values.
68 struct t_ebin
70 //! Number of thermodynamic terms
71 int nener;
72 //! Name and units for each term
73 gmx_enxnm_t* enm;
74 //! Number of steps used for sum (for energy history)
75 int64_t nsteps;
76 //! Number of values added to the sum so far
77 int64_t nsum;
78 //! Term values: each structure stores current, running average and sum.
79 t_energy* e;
80 //! Total number of steps saved (for energy history)
81 int64_t nsteps_sim;
82 //! Total number of values added to sum (used when printing average values at the end of the run)
83 int64_t nsum_sim;
84 //! Energy values throughout the entire simulation: structure stores current, average and sum, but only sum value is used to compute averages
85 t_energy* e_sim;
88 /* \brief Type of expected output: normal or average.
90 enum
92 eprNORMAL,
93 eprAVER,
94 eprNR
97 /*! \brief Create the structure to store thermodynamic properties*/
98 t_ebin* mk_ebin();
100 /*! \brief Destroy the \c eb structure.
102 * \param[in,out] eb Pointer to the structure to destroy.
104 void done_ebin(t_ebin* eb);
106 /*! \brief Create space for the extra thermodynamic term(s) and register its(their) name(s).
108 * The enm array must be static, because the contents are not copied, only the pointers.
110 * \param[in] eb Srtucture in which the space for the termodynamic terms shall be created..
111 * \param[in] nener Number of thermodyamic terms to allocate memory for.
112 * \param[in] enm Names of the terms.
113 * \param[in] unit Units.
115 * \returns A serial number (index) for the newly allocated terms.
117 int get_ebin_space(t_ebin* eb, int nener, const char* const enm[], const char* unit);
119 /*! \brief Add current value of the thermodynamic term(s) to the bin(s).
121 * Add nener reals (eg. energies, box-lengths, pressures) to the at position specified by \c
122 * entryIndex. If bSum is TRUE then the reals are also added to the sum and sum of squares.
124 * \param[in] eb Structure that stores the thermodynamic values.
125 * \param[in] entryIndex Internal index of the term(s) to add.
126 * \param[in] nener Number of the terms to add.
127 * \param[in] ener Value(s) of thermodynamic term(s) (nener ptc.)
128 * \param[in] bSum If the average value should be accumulated for this term(s).
130 void add_ebin(t_ebin* eb, int entryIndex, int nener, const real ener[], gmx_bool bSum);
133 /*! \brief Add values from array to the bins if the matching entry in \c shouldUse is true.
135 * Caller must ensure that \c shouldUse and \c ener to have the same
136 * size, and that \c eb has enough room for the number of true
137 * entries in \c shouldUse.
139 * \param[in] eb Structure that stores the thermodynamic values.
140 * \param[in] entryIndex Internal index of the term(s).
141 * \param[in] shouldUse Array of booleans that indicate which terms should be used.
142 * \param[in] ener Values of thermodinamic terms to add.
143 * \param[in] bSum If the average value should be accumulated for these terms.
145 void add_ebin_indexed(t_ebin* eb,
146 int entryIndex,
147 gmx::ArrayRef<bool> shouldUse,
148 gmx::ArrayRef<const real> ener,
149 gmx_bool bSum);
151 /*! \brief Increase the counters for the sums.
153 * This routine should be called after all add_ebin calls for this step.
155 * \param[in] increment How much counts should be increased
156 * \param[in] eb Structure that stores the thermodynamic values.
157 * \param[in] bSum If the sums counters should be increased as well.
159 void ebin_increase_count(int increment, t_ebin* eb, gmx_bool bSum);
162 /*! \brief Reset the average and fluctuation sums.
164 * \param[in] eb Structure that stores the thermodynamic values.
166 void reset_ebin_sums(t_ebin* eb);
169 /*! \brief Print the contents of some energy bins.
171 * We will print \c nperline entries on a text line (advisory <=
172 * 5). \c prmode may be any of the above listed enum values. \c tsteps is
173 * used only when \c eprAVER is set. If \c bPrHead than the
174 * header is printed.
176 * \c entryIndex and \c nener must be in [0,\c eb->nener), except that \c
177 * nener -1 is interpreted as \c eb->nener.
179 * \todo Callers should be refactored to pass \c eb->nener, rather than
180 * us implement and rely on this special behavior of -1.
182 * \param[in] fp I/O pointer to print to.
183 * \param[in] eb Structure that stores the thermodynamic values.
184 * \param[in] entryIndex Internal index of the term(s).
185 * \param[in] nener Number of the terms to print.
186 * \param[in] nperline Number of values per line.
187 * \param[in] prmode Print current (eprNORMAL) or average (eprAVER) values.
188 * \param[in] bPrHead If the header should be printed.
190 void pr_ebin(FILE* fp, t_ebin* eb, int entryIndex, int nener, int nperline, int prmode, gmx_bool bPrHead);
192 #endif