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