Clang-11 fixes for applied forces
[alexxy/gromacs.git] / src / gromacs / applied_forces / awh / biaswriter.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2015,2016,2017,2018,2019,2020,2021, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  *
38  * \brief
39  * This file contains the BiasWriter class that prepares and writes data of a Bias to an energy file.
40  *
41  * \author Viveca Lindahl
42  * \author Berk Hess <hess@kth.se>
43  * \ingroup module_awh
44  */
45
46 #ifndef GMX_AWH_BIASWRITER_H
47 #define GMX_AWH_BIASWRITER_H
48
49 #include <map>
50 #include <vector>
51
52 #include "gromacs/fileio/enxio.h"
53 #include "gromacs/utility/arrayref.h"
54 #include "gromacs/utility/basedefinitions.h"
55
56 struct gmx_multisim_t;
57 struct t_enxsubblock;
58
59 namespace gmx
60 {
61 class Bias;
62
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.
67  */
68
69 //! Enum with the AWH variables to write.
70 enum class AwhOutputEntryType
71 {
72     MetaData,               //!< Meta data.
73     CoordValue,             //!< Coordinate value.
74     Pmf,                    //!< The pmf.
75     Bias,                   //!< The bias.
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.
81 };
82
83 //! Enum with the types of metadata to write.
84 enum class AwhOutputMetaData
85 {
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.
90 };
91
92 //! Enum with different ways of normalizing the output.
93 enum class Normalization
94 {
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.
99 };
100
101 /*! \internal \brief AWH output data block that can be written to an energy file block.
102  */
103 class AwhEnergyBlock
104 {
105 public:
106     /*! \brief Constructor
107      *
108      * \param[in] numPoints           Number of points in block.
109      * \param[in] normalizationType   Type of normalization.
110      * \param[in] normalizationValue  Value to normalization with.
111      */
112     AwhEnergyBlock(int numPoints, Normalization normalizationType, float normalizationValue);
113
114     /*! \brief Returns an ArrarRef to the data in the block.
115      */
116     gmx::ArrayRef<float> data() { return data_; }
117
118     const Normalization normalizationType;  /**< How to normalize the output data */
119     const float         normalizationValue; /**< The normalization value */
120 private:
121     std::vector<float> data_; /**< The data, always float which is enough since this is statistical data */
122 };
123
124 /*! \internal \brief Class organizing the output data storing and writing of an AWH bias.
125  */
126 class BiasWriter
127 {
128 public:
129     /*! \brief Constructor.
130      *
131      * \param[in] bias  The AWH bias.
132      */
133     BiasWriter(const Bias& bias);
134
135     /*! \brief Returns the number of data blocks.
136      *
137      * \returns the number of data blocks.
138      */
139     int numBlocks() const { return block_.size(); }
140
141     /*! \brief Collect AWH bias data in blocks and write to energy subblocks.
142      *
143      * \param[in]     bias      The AWH Bias.
144      * \param[in,out] subblock  Energy subblocks to write to.
145      * \returns the number of blocks written.
146      */
147     int writeToEnergySubblocks(const Bias& bias, t_enxsubblock* subblock);
148
149 private:
150     /*! \brief Query if the writer has a block for the given variable.
151      *
152      * \param[in] outputType  Output entry type.
153      */
154     bool hasVarBlock(AwhOutputEntryType outputType) const
155     {
156         return outputTypeToBlock_.find(outputType)->second >= 0;
157     }
158
159     /*! \brief* Find the first block containing the given variable.
160      *
161      * \param[in] outputType  Output entry type.
162      * \returns the first block index for the variable, or -1 there is no block.
163      */
164     int getVarStartBlock(AwhOutputEntryType outputType) const
165     {
166         return outputTypeToBlock_.find(outputType)->second;
167     }
168
169     /*! \brief Transfer AWH point data to writer data blocks.
170      *
171      * \param[in] metaDataIndex  Meta data block index.
172      * \param[in] metaDataType   The type of meta data to write.
173      * \param[in] bias           The AWH Bias.
174      */
175     void transferMetaDataToWriter(gmx::index metaDataIndex, AwhOutputMetaData metaDataType, const Bias& bias);
176
177     /*! \brief Transfer AWH point data to writer data blocks.
178      *
179      * \param[in] outputType  Output entry type.
180      * \param[in] pointIndex  The point index.
181      * \param[in] bias        The AWH Bias.
182      * \param[in] pmf         PMF values.
183      */
184     void transferPointDataToWriter(AwhOutputEntryType         outputType,
185                                    int                        pointIndex,
186                                    const Bias&                bias,
187                                    gmx::ArrayRef<const float> pmf);
188
189     /*! \brief
190      * Prepare the bias output data.
191      *
192      * \param[in] bias  The AWH Bias.
193      */
194     void prepareBiasOutput(const Bias& bias);
195
196     std::vector<AwhEnergyBlock>       block_; /**< The data blocks */
197     std::map<AwhOutputEntryType, int> outputTypeToBlock_; /**< Start block index for each variable, -1 when variable should not be written */
198 };
199
200 } // namespace gmx
201
202 #endif /* GMX_AWH_BIASWRITER_H */