2 * This file is part of the GROMACS molecular simulation package.
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.
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.
38 * Implements the HistogramSize class.
40 * \author Viveca Lindahl
41 * \author Berk Hess <hess@kth.se>
47 #include "histogramsize.h"
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"
60 #include "biasparams.h"
61 #include "pointstate.h"
66 HistogramSize::HistogramSize(const AwhBiasParams
&awhBiasParams
,
67 double histogramSizeInitial
) :
69 histogramSize_(histogramSizeInitial
),
70 inInitialStage_(awhBiasParams
.eGrowth
== eawhgrowthEXP_LINEAR
),
71 equilibrateHistogram_(awhBiasParams
.equilibrateHistogram
),
72 logScaledSampleWeight_(0),
73 maxLogScaledSampleWeight_(0),
74 havePrintedAboutCovering_(false)
78 double HistogramSize::newHistogramSizeInitialStage(const BiasParams
¶ms
,
80 bool detectedCovering
,
81 ArrayRef
<double> weightsumCovering
,
84 /* The histogram size is kept constant until the sampling region has been covered
85 and the the current sample weight is large enough and the histogram is ready. */
86 if (!detectedCovering
||
87 (logScaledSampleWeight_
< maxLogScaledSampleWeight_
) ||
88 equilibrateHistogram_
)
90 return histogramSize_
;
93 /* Reset the covering weight histogram. If we got this far we are either entering a
94 new covering stage with a new covering histogram or exiting the initial stage
96 std::fill(weightsumCovering
.begin(), weightsumCovering
.end(), 0);
98 /* The current sample weigth is now the maximum. */
99 double prevMaxLogScaledSampleWeight
= maxLogScaledSampleWeight_
;
100 maxLogScaledSampleWeight_
= logScaledSampleWeight_
;
102 /* Increase the histogram size by a constant scale factor if we can, i.e. if the sample weight
103 resulting from such a scaling is still larger than the previous maximum sample weight
104 (ensuring that the sample weights at the end of each covering stage are monotonically
105 increasing). If we cannot, exit the initial stage without changing the histogram size. */
107 /* The scale factor. The value is not very critical but should obviously be > 1 (or the exit
108 will happen very late) and probably < 5 or so (or there will be no initial stage). */
109 static const double growthFactor
= 3;
111 /* The scale factor is in most cases very close to the histogram growth factor. */
112 double scaleFactor
= growthFactor
/(1. + params
.updateWeight
*params
.localWeightScaling
/histogramSize_
);
114 bool exitInitialStage
= (logScaledSampleWeight_
- std::log(scaleFactor
) <= prevMaxLogScaledSampleWeight
);
115 double newHistogramSize
= exitInitialStage
? histogramSize_
: histogramSize_
*growthFactor
;
117 /* Update the AWH bias about the exit. */
118 inInitialStage_
= !exitInitialStage
;
120 /* Print information about coverings and if there was an exit. */
121 if (fplog
!= nullptr)
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
);
126 if (exitInitialStage
)
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
137 return newHistogramSize
;
144 * Checks if the histogram has equilibrated to the target distribution.
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.
150 * \param[in] pointStates The state of the bias points.
151 * \returns true if the histogram is equilibrated.
153 bool histogramIsEquilibrated(const std::vector
<PointState
> &pointStates
)
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
)
160 if (!pointState
.inTargetRegion())
164 totalWeight
+= pointState
.weightSumTot();
167 GMX_RELEASE_ASSERT(totalWeight
> 0, "No samples when normalizing AWH histogram.");
168 double inverseTotalWeight
= 1./totalWeight
;
170 /* Points with target weight below a certain cutoff are ignored. */
171 static const double minTargetCutoff
= 0.05;
172 double minTargetWeight
= 1./numTargetPoints
*minTargetCutoff
;
174 /* Points with error less than this tolerance pass the check.*/
175 static const double errorTolerance
= 0.2;
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
)
182 double targetWeight
= pointState
.target();
183 double sampledWeight
= pointState
.weightSumTot()*inverseTotalWeight
;
185 /* Ignore these points. */
186 if (!pointState
.inTargetRegion() || targetWeight
< minTargetWeight
)
191 if (std::abs(sampledWeight
/targetWeight
- 1) > errorTolerance
)
193 notEquilibratedWeight
+= targetWeight
;
197 equilibratedWeight
+= targetWeight
;
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;
205 return equilibratedWeight
/(equilibratedWeight
+ notEquilibratedWeight
) > minFraction
;;
210 double HistogramSize::newHistogramSize(const BiasParams
¶ms
,
213 const std::vector
<PointState
> &pointStates
,
214 ArrayRef
<double> weightsumCovering
,
217 double newHistogramSize
;
220 /* Only bother with checking equilibration if we have covered already. */
221 if (equilibrateHistogram_
&& covered
)
223 /* The histogram is equilibrated at most once. */
224 equilibrateHistogram_
= !histogramIsEquilibrated(pointStates
);
226 if (fplog
!= nullptr)
228 std::string prefix
= gmx::formatString("\nawh%d:", params
.biasIndex
+ 1);
229 if (!equilibrateHistogram_
)
231 fprintf(fplog
, "%s equilibrated histogram at t = %g ps.\n", prefix
.c_str(), t
);
233 else if (!havePrintedAboutCovering_
)
235 fprintf(fplog
, "%s covered but histogram not equilibrated at t = %g ps.\n", prefix
.c_str(), t
);
236 havePrintedAboutCovering_
= true; /* Just print once. */
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
);
246 /* If not in the initial stage, the histogram grows at a linear, possibly scaled down, rate. */
247 newHistogramSize
= histogramSize_
+ params
.updateWeight
*params
.localWeightScaling
;
250 return newHistogramSize
;
253 void HistogramSize::setHistogramSize(double histogramSize
,
254 double weightHistogramScalingFactor
)
256 GMX_ASSERT(histogramSize
> 0, "The histogram should not be empty");
257 GMX_ASSERT(weightHistogramScalingFactor
> 0, "The histogram scaling factor should be positive");
259 histogramSize_
= histogramSize
;
261 /* The weight of new samples relative to previous ones change
262 * when the histogram is rescaled. We keep the log since this number
263 * can become very large.
265 logScaledSampleWeight_
-= std::log(weightHistogramScalingFactor
);
268 void HistogramSize::restoreFromHistory(const AwhBiasStateHistory
&stateHistory
)
270 numUpdates_
= stateHistory
.numUpdates
;
271 histogramSize_
= stateHistory
.histSize
;
272 inInitialStage_
= stateHistory
.in_initial
;
273 equilibrateHistogram_
= stateHistory
.equilibrateHistogram
;
274 logScaledSampleWeight_
= stateHistory
.logScaledSampleWeight
;
275 maxLogScaledSampleWeight_
= stateHistory
.maxLogScaledSampleWeight
;
276 havePrintedAboutCovering_
= false;
279 void HistogramSize::storeState(AwhBiasStateHistory
*stateHistory
) const
281 stateHistory
->numUpdates
= numUpdates_
;
282 stateHistory
->histSize
= histogramSize_
;
283 stateHistory
->in_initial
= inInitialStage_
;
284 stateHistory
->equilibrateHistogram
= equilibrateHistogram_
;
285 stateHistory
->logScaledSampleWeight
= logScaledSampleWeight_
;
286 stateHistory
->maxLogScaledSampleWeight
= maxLogScaledSampleWeight_
;
287 /* We'll print again about covering when restoring the state */