Remove gmx custom fixed int (e.g. gmx_int64_t) types
[alexxy/gromacs.git] / src / gromacs / awh / histogramsize.h
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
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  * Declares the HistogramSize class.
40  *
41  * The data members of this class keep track of global size and update related
42  * properties of the bias histogram and the evolution of the histogram size.
43  * Initially histogramSize_ (and thus the convergence rate) is controlled
44  * heuristically to get good initial estimates,  i.e. increase the robustness
45  * of the method.
46  *
47  * \author Viveca Lindahl
48  * \author Berk Hess <hess@kth.se>
49  * \ingroup module_awh
50  */
51
52 #ifndef GMX_AWH_HISTOGRAMSIZE_H
53 #define GMX_AWH_HISTOGRAMSIZE_H
54
55 #include <cstdio>
56
57 #include <vector>
58
59 #include "gromacs/math/vectypes.h"
60 #include "gromacs/utility/arrayref.h"
61
62 namespace gmx
63 {
64
65 struct AwhBiasStateHistory;
66 struct AwhBiasParams;
67 class BiasParams;
68 class PointState;
69
70 /*! \internal
71  * \brief Tracks global size related properties of the bias histogram.
72  *
73  * Tracks the number of updates and the histogram size.
74  * Also keep track of the stage (initial/final of the AWH method
75  * and printing warnings about covering.
76  *
77  * \note Histogram sizes are floating-point values, since the histogram uses weighted
78  *        entries and we can assign a floating-point scaling factor when changing it.
79  */
80 class HistogramSize
81 {
82     public:
83         /*! \brief Constructor.
84          *
85          * \param[in] awhBiasParams         The Bias parameters from inputrec.
86          * \param[in] histogramSizeInitial  The initial histogram size.
87          */
88         HistogramSize(const AwhBiasParams &awhBiasParams,
89                       double               histogramSizeInitial);
90
91     private:
92         /*! \brief
93          * Returns the new size of the reference weight histogram in the initial stage.
94          *
95          * This function also takes care resetting the histogram used for covering checks
96          * and for exiting the initial stage.
97          *
98          * \param[in]     params             The bias parameters.
99          * \param[in]     t                  Time.
100          * \param[in]     detectedCovering   True if we detected that the sampling interval has been sufficiently covered.
101          * \param[in,out] weightsumCovering  The weight sum for checking covering.
102          * \param[in,out] fplog              Log file.
103          * \returns the new histogram size.
104          */
105         double newHistogramSizeInitialStage(const BiasParams &params,
106                                             double            t,
107                                             bool              detectedCovering,
108                                             ArrayRef<double>  weightsumCovering,
109                                             FILE             *fplog);
110
111     public:
112         /*! \brief
113          * Return the new reference weight histogram size for the current update.
114          *
115          * This function also takes care of checking for covering in the initial stage.
116          *
117          * \param[in]     params             The bias parameters.
118          * \param[in]     t                  Time.
119          * \param[in]     covered            True if the sampling interval has been covered enough.
120          * \param[in]     pointStates        The state of the grid points.
121          * \param[in,out] weightsumCovering  The weight sum for checking covering.
122          * \param[in,out] fplog              Log file.
123          * \returns the new histogram size.
124          */
125         double newHistogramSize(const BiasParams              &params,
126                                 double                         t,
127                                 bool                           covered,
128                                 const std::vector<PointState> &pointStates,
129                                 ArrayRef<double>               weightsumCovering,
130                                 FILE                          *fplog);
131
132         /*! \brief Restores the histogram size from history.
133          *
134          * \param[in] stateHistory  The AWH bias state history.
135          */
136         void restoreFromHistory(const AwhBiasStateHistory &stateHistory);
137
138         /*! \brief Store the histogram size state in a history struct.
139          *
140          * \param[in,out] stateHistory  The AWH bias state history.
141          */
142         void storeState(AwhBiasStateHistory *stateHistory) const;
143
144         /*! \brief Returns the number of updates since the start of the simulation.
145          */
146         int numUpdates() const
147         {
148             return numUpdates_;
149         }
150
151         /*! \brief Increments the number of updates by 1.
152          */
153         void incrementNumUpdates()
154         {
155             numUpdates_ += 1;
156         }
157
158         /*! \brief Returns the histogram size.
159          */
160         double histogramSize() const
161         {
162             return histogramSize_;
163         }
164
165         /*! \brief Sets the histogram size.
166          *
167          * \param[in] histogramSize                 The new histogram size.
168          * \param[in] weightHistogramScalingFactor  The factor to scale the weight by.
169          */
170         void setHistogramSize(double histogramSize,
171                               double weightHistogramScalingFactor);
172
173         /*! \brief Returns true if we are in the initial stage of the AWH method.
174          */
175         bool inInitialStage() const
176         {
177             return inInitialStage_;
178         }
179
180         /*! \brief Returns The log of the current sample weight, scaled because of the histogram rescaling.
181          */
182         double logScaledSampleWeight() const
183         {
184             return logScaledSampleWeight_;
185         }
186
187     private:
188         int64_t numUpdates_; /**< The number of updates performed since the start of the simulation. */
189
190         /* The histogram size sets the update size and so controls the convergence rate of the free energy and bias. */
191         double      histogramSize_; /**< Size of reference weight histogram. */
192
193         /* Values that control the evolution of the histogram size. */
194         bool      inInitialStage_;           /**< True if in the intial stage. */
195         bool      equilibrateHistogram_;     /**< True if samples are kept from accumulating until the sampled distribution is close enough to the target. */
196         double    logScaledSampleWeight_;    /**< The log of the current sample weight, scaled because of the histogram rescaling. */
197         double    maxLogScaledSampleWeight_; /**< Maximum sample weight obtained for previous (smaller) histogram sizes. */
198
199         /* Bool to avoid printing multiple, not so useful, messages to log */
200         bool      havePrintedAboutCovering_; /**< True if we have printed about covering to the log while equilibrateHistogram==true */
201 };
202
203 }      // namespace gmx
204
205 #endif /* GMX_AWH_HISTOGRAMSIZE_H */