2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017,2018, by the GROMACS development team, led by
5 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
6 * and including many others, as listed in the AUTHORS file in the
7 * top-level source directory and at http://www.gromacs.org.
9 * GROMACS is free software; you can redistribute it and/or
10 * modify it under the terms of the GNU Lesser General Public License
11 * as published by the Free Software Foundation; either version 2.1
12 * of the License, or (at your option) any later version.
14 * GROMACS is distributed in the hope that it will be useful,
15 * but WITHOUT ANY WARRANTY; without even the implied warranty of
16 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
17 * Lesser General Public License for more details.
19 * You should have received a copy of the GNU Lesser General Public
20 * License along with GROMACS; if not, see
21 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
22 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
24 * If you want to redistribute modifications to GROMACS, please
25 * consider that scientific software is very special. Version
26 * control is crucial - bugs must be traceable. We will be happy to
27 * consider code for inclusion in the official distribution, but
28 * derived work must not be called official GROMACS. Details are found
29 * in the README & COPYING files - if they are missing, get the
30 * official version at http://www.gromacs.org.
32 * To help us fund GROMACS development, we humbly ask that you cite
33 * the research papers on the package. Check out http://www.gromacs.org.
39 * This file contains the BiasWriter class that prepares and writes data of a Bias to an energy file.
41 * \author Viveca Lindahl
42 * \author Berk Hess <hess@kth.se>
46 #ifndef GMX_AWH_BIASWRITER_H
47 #define GMX_AWH_BIASWRITER_H
52 #include "gromacs/fileio/enxio.h"
53 #include "gromacs/utility/arrayref.h"
54 #include "gromacs/utility/basedefinitions.h"
56 struct gmx_multisim_t
;
63 /* TODO: the post-simulations AWH reader and this AWH writer are totally
64 * disconnected although they read/write the same data. I'm not sure how
65 * to handle that or if it should be left as it is until the writing is done
66 * in a differen format (i.e. TNG) than the current energy file.
69 //! Enum with the AWH variables to write.
70 enum class AwhOutputEntryType
72 MetaData
, //!< Meta data.
73 CoordValue
, //!< Coordinate value.
76 Visits
, //!< The number of visits.
77 Weights
, //!< The weights.
78 Target
, //!< The target distribition.
79 ForceCorrelationVolume
, //!< The volume of the force correlation tensor.
80 FrictionTensor
//!< The full friction tensor.
83 //! Enum with the types of metadata to write.
84 enum class AwhOutputMetaData
86 NumBlock
, //!< The number of blocks.
87 TargetError
, //!< The target error.
88 ScaledSampleWeight
, //!< The logarithm of the sample weight relative to a sample weight of 1 at the initial time.
89 Count
//!< The number of enum values, not including Count.
92 //! Enum with different ways of normalizing the output.
93 enum class Normalization
95 None
, //!< No normalization.
96 Coordinate
, //!< Scale using the internal/user input coordinate scaling factor.
97 FreeEnergy
, //!< Normalize free energy values by subtracting the minimum value.
98 Distribution
//!< Normalize the distribution to 1.
101 /*! \internal \brief AWH output data block that can be written to an energy file block.
106 /*! \brief Constructor
108 * \param[in] numPoints Number of points in block.
109 * \param[in] normalizationType Type of normalization.
110 * \param[in] normalizationValue Value to normalization with.
112 AwhEnergyBlock(int numPoints
,
113 Normalization normalizationType
,
114 float normalizationValue
);
116 /*! \brief Returns an ArrarRef to the data in the block.
118 gmx::ArrayRef
<float> data()
123 const Normalization normalizationType
; /**< How to normalize the output data */
124 const float normalizationValue
; /**< The normalization value */
126 std::vector
<float> data_
; /**< The data, always float which is enough since this is statistical data */
129 /*! \internal \brief Class organizing the output data storing and writing of an AWH bias.
134 /*! \brief Constructor.
136 * \param[in] bias The AWH bias.
138 BiasWriter(const Bias
&bias
);
140 /*! \brief Returns the number of data blocks.
142 * \returns the number of data blocks.
144 int numBlocks() const
146 return block_
.size();
149 /*! \brief Collect AWH bias data in blocks and write to energy subblocks.
151 * \param[in] bias The AWH Bias.
152 * \param[in,out] subblock Energy subblocks to write to.
153 * \returns the number of blocks written.
155 int writeToEnergySubblocks(const Bias
&bias
, t_enxsubblock
*subblock
);
158 /*! \brief Query if the writer has a block for the given variable.
160 * \param[in] outputType Output entry type.
162 bool hasVarBlock(AwhOutputEntryType outputType
) const
164 return outputTypeToBlock_
.find(outputType
)->second
>= 0;
167 /*! \brief* Find the first block containing the given variable.
169 * \param[in] outputType Output entry type.
170 * \returns the first block index for the variable, or -1 there is no block.
172 int getVarStartBlock(AwhOutputEntryType outputType
) const
174 return outputTypeToBlock_
.find(outputType
)->second
;
177 /*! \brief Transfer AWH point data to writer data blocks.
179 * \param[in] metaDataIndex Meta data block index.
180 * \param[in] metaDataType The type of meta data to write.
181 * \param[in] bias The AWH Bias.
183 void transferMetaDataToWriter(gmx::index metaDataIndex
,
184 AwhOutputMetaData metaDataType
,
187 /*! \brief Transfer AWH point data to writer data blocks.
189 * \param[in] outputType Output entry type.
190 * \param[in] pointIndex The point index.
191 * \param[in] bias The AWH Bias.
192 * \param[in] pmf PMF values.
194 void transferPointDataToWriter(AwhOutputEntryType outputType
,
197 gmx::ArrayRef
<const float> pmf
);
200 * Prepare the bias output data.
202 * \param[in] bias The AWH Bias.
204 void prepareBiasOutput(const Bias
&bias
);
207 std::vector
<AwhEnergyBlock
> block_
; /**< The data blocks */
208 std::map
<AwhOutputEntryType
, int> outputTypeToBlock_
; /**< Start block index for each variable, -1 when variable should not be written */
213 #endif /* GMX_AWH_BIASWRITER_H */