Modernize AWH enumerations
[alexxy/gromacs.git] / src / gromacs / applied_forces / awh / histogramsize.cpp
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  * \brief
38  * Implements the HistogramSize class.
39  *
40  * \author Viveca Lindahl
41  * \author Berk Hess <hess@kth.se>
42  * \ingroup module_awh
43  */
44
45 #include "gmxpre.h"
46
47 #include "histogramsize.h"
48
49 #include <cstdio>
50 #include <cstdlib>
51 #include <cstring>
52
53 #include <algorithm>
54
55 #include "gromacs/mdtypes/awh_history.h"
56 #include "gromacs/mdtypes/awh_params.h"
57 #include "gromacs/utility/arrayref.h"
58 #include "gromacs/utility/gmxassert.h"
59 #include "gromacs/utility/stringutil.h"
60
61 #include "biasparams.h"
62 #include "pointstate.h"
63
64 namespace gmx
65 {
66
67 HistogramSize::HistogramSize(const AwhBiasParams& awhBiasParams, double histogramSizeInitial) :
68     numUpdates_(0),
69     histogramSize_(histogramSizeInitial),
70     inInitialStage_(awhBiasParams.eGrowth == AwhHistogramGrowthType::ExponentialLinear),
71     equilibrateHistogram_(awhBiasParams.equilibrateHistogram),
72     logScaledSampleWeight_(0),
73     maxLogScaledSampleWeight_(0),
74     havePrintedAboutCovering_(false)
75 {
76 }
77
78 double HistogramSize::newHistogramSizeInitialStage(const BiasParams& params,
79                                                    double            t,
80                                                    bool              detectedCovering,
81                                                    ArrayRef<double>  weightsumCovering,
82                                                    FILE*             fplog)
83 {
84     /* The histogram size is kept constant until the sampling region has been covered
85        and the current sample weight is large enough and the histogram is ready. */
86     if (!detectedCovering || (logScaledSampleWeight_ < maxLogScaledSampleWeight_) || equilibrateHistogram_)
87     {
88         return histogramSize_;
89     }
90
91     /* Reset the covering weight histogram. If we got this far we are either entering a
92        new covering stage with a new covering histogram or exiting the initial stage
93        altogether. */
94     std::fill(weightsumCovering.begin(), weightsumCovering.end(), 0);
95
96     /*  The current sample weigth is now the maximum. */
97     double prevMaxLogScaledSampleWeight = maxLogScaledSampleWeight_;
98     maxLogScaledSampleWeight_           = logScaledSampleWeight_;
99
100     /* Increase the histogram size by a constant scale factor if we can, i.e. if the sample weight
101        resulting from such a scaling is still larger than the previous maximum sample weight
102        (ensuring that the sample weights at the end of each covering stage are monotonically
103        increasing). If we cannot, exit the initial stage without changing the histogram size. */
104
105     /* The scale factor. The value is not very critical but should obviously be > 1 (or the exit
106        will happen very late) and probably < 5 or so (or there will be no initial stage). */
107     static const double growthFactor = 3;
108
109     /* The scale factor is in most cases very close to the histogram growth factor. */
110     double scaleFactor =
111             growthFactor / (1. + params.updateWeight * params.localWeightScaling / histogramSize_);
112
113     bool exitInitialStage =
114             (logScaledSampleWeight_ - std::log(scaleFactor) <= prevMaxLogScaledSampleWeight);
115     double newHistogramSize = exitInitialStage ? histogramSize_ : histogramSize_ * growthFactor;
116
117     /* Update the AWH bias about the exit. */
118     inInitialStage_ = !exitInitialStage;
119
120     /* Print information about coverings and if there was an exit. */
121     if (fplog != nullptr)
122     {
123         std::string prefix = gmx::formatString("\nawh%d:", params.biasIndex + 1);
124         fprintf(fplog, "%s covering at t = %g ps. Decreased the update size.\n", prefix.c_str(), t);
125
126         if (exitInitialStage)
127         {
128             fprintf(fplog, "%s out of the initial stage at t = %g.\n", prefix.c_str(), t);
129             /* It would be nice to have a way of estimating a minimum time until exit but it
130                is difficult because the exit time is determined by how long it takes to cover
131                relative to the time it takes to "regaining" enough sample weight. The latter
132                is easy to calculate, but how the former depends on the histogram size
133                is not known. */
134         }
135         fflush(fplog);
136     }
137     return newHistogramSize;
138 }
139
140 namespace
141 {
142
143 /*! \brief
144  * Checks if the histogram has equilibrated to the target distribution.
145  *
146  * The histogram is considered equilibrated if, for a minimum fraction of
147  * the target region, the relative error of the sampled weight relative
148  * to the target is less than a tolerance value.
149  *
150  * \param[in] pointStates  The state of the bias points.
151  * \returns true if the histogram is equilibrated.
152  */
153 bool histogramIsEquilibrated(const std::vector<PointState>& pointStates)
154 {
155     /* Get the total weight of the total weight histogram; needed for normalization. */
156     double totalWeight     = 0;
157     int    numTargetPoints = 0;
158     for (auto& pointState : pointStates)
159     {
160         if (!pointState.inTargetRegion())
161         {
162             continue;
163         }
164         totalWeight += pointState.weightSumTot();
165         numTargetPoints++;
166     }
167     GMX_RELEASE_ASSERT(totalWeight > 0, "No samples when normalizing AWH histogram.");
168     double inverseTotalWeight = 1. / totalWeight;
169
170     /* Points with target weight below a certain cutoff are ignored. */
171     static const double minTargetCutoff = 0.05;
172     double              minTargetWeight = 1. / numTargetPoints * minTargetCutoff;
173
174     /* Points with error less than this tolerance pass the check.*/
175     static const double errorTolerance = 0.2;
176
177     /* Sum up weight of points that do or don't pass the check. */
178     double equilibratedWeight    = 0;
179     double notEquilibratedWeight = 0;
180     for (auto& pointState : pointStates)
181     {
182         double targetWeight  = pointState.target();
183         double sampledWeight = pointState.weightSumTot() * inverseTotalWeight;
184
185         /* Ignore these points. */
186         if (!pointState.inTargetRegion() || targetWeight < minTargetWeight)
187         {
188             continue;
189         }
190
191         if (std::abs(sampledWeight / targetWeight - 1) > errorTolerance)
192         {
193             notEquilibratedWeight += targetWeight;
194         }
195         else
196         {
197             equilibratedWeight += targetWeight;
198         }
199     }
200
201     /* It is enough if sampling in at least a fraction of the target region follows the target
202        distribution. Boundaries will in general fail and this should be ignored (to some extent). */
203     static const double minFraction = 0.8;
204
205     return equilibratedWeight / (equilibratedWeight + notEquilibratedWeight) > minFraction;
206 }
207
208 } // namespace
209
210 double HistogramSize::newHistogramSize(const BiasParams&              params,
211                                        double                         t,
212                                        bool                           covered,
213                                        const std::vector<PointState>& pointStates,
214                                        ArrayRef<double>               weightsumCovering,
215                                        FILE*                          fplog)
216 {
217     double newHistogramSize;
218     if (inInitialStage_)
219     {
220         /* Only bother with checking equilibration if we have covered already. */
221         if (equilibrateHistogram_ && covered)
222         {
223             /* The histogram is equilibrated at most once. */
224             equilibrateHistogram_ = !histogramIsEquilibrated(pointStates);
225
226             if (fplog != nullptr)
227             {
228                 std::string prefix = gmx::formatString("\nawh%d:", params.biasIndex + 1);
229                 if (!equilibrateHistogram_)
230                 {
231                     fprintf(fplog, "%s equilibrated histogram at t = %g ps.\n", prefix.c_str(), t);
232                 }
233                 else if (!havePrintedAboutCovering_)
234                 {
235                     fprintf(fplog,
236                             "%s covered but histogram not equilibrated at t = %g ps.\n",
237                             prefix.c_str(),
238                             t);
239                     havePrintedAboutCovering_ = true; /* Just print once. */
240                 }
241             }
242         }
243
244         /* In the initial stage, the histogram grows dynamically as a function of the number of coverings. */
245         newHistogramSize = newHistogramSizeInitialStage(params, t, covered, weightsumCovering, fplog);
246     }
247     else
248     {
249         /* If not in the initial stage, the histogram grows at a linear, possibly scaled down, rate. */
250         newHistogramSize = histogramSize_ + params.updateWeight * params.localWeightScaling;
251     }
252
253     return newHistogramSize;
254 }
255
256 void HistogramSize::setHistogramSize(double histogramSize, double weightHistogramScalingFactor)
257 {
258     GMX_ASSERT(histogramSize > 0, "The histogram should not be empty");
259     GMX_ASSERT(weightHistogramScalingFactor > 0, "The histogram scaling factor should be positive");
260
261     histogramSize_ = histogramSize;
262
263     /* The weight of new samples relative to previous ones change
264      * when the histogram is rescaled. We keep the log since this number
265      * can become very large.
266      */
267     logScaledSampleWeight_ -= std::log(weightHistogramScalingFactor);
268 };
269
270 void HistogramSize::restoreFromHistory(const AwhBiasStateHistory& stateHistory)
271 {
272     numUpdates_               = stateHistory.numUpdates;
273     histogramSize_            = stateHistory.histSize;
274     inInitialStage_           = stateHistory.in_initial;
275     equilibrateHistogram_     = stateHistory.equilibrateHistogram;
276     logScaledSampleWeight_    = stateHistory.logScaledSampleWeight;
277     maxLogScaledSampleWeight_ = stateHistory.maxLogScaledSampleWeight;
278     havePrintedAboutCovering_ = false;
279 }
280
281 void HistogramSize::storeState(AwhBiasStateHistory* stateHistory) const
282 {
283     stateHistory->numUpdates               = numUpdates_;
284     stateHistory->histSize                 = histogramSize_;
285     stateHistory->in_initial               = inInitialStage_;
286     stateHistory->equilibrateHistogram     = equilibrateHistogram_;
287     stateHistory->logScaledSampleWeight    = logScaledSampleWeight_;
288     stateHistory->maxLogScaledSampleWeight = maxLogScaledSampleWeight_;
289     /* We'll print again about covering when restoring the state */
290 }
291
292 } // namespace gmx