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 Bias class.
40 * \author Viveca Lindahl
41 * \author Berk Hess <hess@kth.se>
59 #include "gromacs/fileio/gmxfio.h"
60 #include "gromacs/gmxlib/network.h"
61 #include "gromacs/math/functions.h"
62 #include "gromacs/math/utilities.h"
63 #include "gromacs/mdtypes/awh_history.h"
64 #include "gromacs/mdtypes/awh_params.h"
65 #include "gromacs/mdtypes/commrec.h"
66 #include "gromacs/utility/exceptions.h"
67 #include "gromacs/utility/gmxassert.h"
68 #include "gromacs/utility/stringutil.h"
70 #include "correlationgrid.h"
71 #include "correlationhistory.h"
72 #include "pointstate.h"
77 void Bias::warnForHistogramAnomalies(double t, int64_t step, FILE* fplog)
79 const int maxNumWarningsInCheck = 1; /* The maximum number of warnings to print per check */
80 const int maxNumWarningsInRun = 10; /* The maximum number of warnings to print in a run */
82 if (fplog == nullptr || numWarningsIssued_ >= maxNumWarningsInRun || state_.inInitialStage()
83 || !params_.isCheckHistogramForAnomaliesStep(step))
89 state_.warnForHistogramAnomalies(grid_, biasIndex(), t, fplog, maxNumWarningsInCheck);
91 if (numWarningsIssued_ >= maxNumWarningsInRun)
93 fprintf(fplog, "\nawh%d: suppressing future AWH warnings.\n", biasIndex() + 1);
97 void Bias::doSkippedUpdatesForAllPoints()
99 if (params_.skipUpdates())
101 state_.doSkippedUpdatesForAllPoints(params_);
105 gmx::ArrayRef<const double> Bias::calcForceAndUpdateBias(const awh_dvec coordValue,
106 double* awhPotential,
107 double* potentialJump,
108 const t_commrec* commRecord,
109 const gmx_multisim_t* ms,
117 GMX_THROW(InvalidInputError(
118 "The step number is negative which is not supported by the AWH code."));
121 state_.setCoordValue(grid_, coordValue);
123 std::vector<double, AlignedAllocator<double>>& probWeightNeighbor = alignedTempWorkSpace_;
125 /* If the convolved force is needed or this is a sampling step,
126 * the bias in the current neighborhood needs to be up-to-date
127 * and the probablity weights need to be calculated.
129 const bool sampleCoord = params_.isSampleCoordStep(step);
130 const bool moveUmbrella = (sampleCoord || step == 0);
131 double convolvedBias = 0;
132 if (params_.convolveForce || moveUmbrella || sampleCoord)
134 if (params_.skipUpdates())
136 state_.doSkippedUpdatesInNeighborhood(params_, grid_);
140 state_.updateProbabilityWeightsAndConvolvedBias(dimParams_, grid_, &probWeightNeighbor);
144 updateForceCorrelationGrid(probWeightNeighbor, t);
146 state_.sampleCoordAndPmf(grid_, probWeightNeighbor, convolvedBias);
150 const CoordState& coordState = state_.coordState();
152 /* Set the bias force and get the potential contribution from this bias.
153 * The potential jump occurs at different times depending on how
154 * the force is applied (and how the potential is normalized).
155 * For the convolved force it happens when the bias is updated,
156 * for the umbrella when the umbrella is moved.
160 if (params_.convolveForce)
162 state_.calcConvolvedForce(dimParams_, grid_, probWeightNeighbor, tempForce_, biasForce_);
164 potential = -convolvedBias * params_.invBeta;
169 GMX_RELEASE_ASSERT(state_.points()[coordState.umbrellaGridpoint()].inTargetRegion(),
170 "AWH bias grid point for the umbrella reference value is outside of the "
172 potential = state_.calcUmbrellaForceAndPotential(
173 dimParams_, grid_, coordState.umbrellaGridpoint(), biasForce_);
175 /* Moving the umbrella results in a force correction and
176 * a new potential. The umbrella center is sampled as often as
177 * the coordinate so we know the probability weights needed
178 * for moving the umbrella are up-to-date.
182 double newPotential = state_.moveUmbrella(dimParams_, grid_, probWeightNeighbor,
183 biasForce_, step, seed, params_.biasIndex);
184 *potentialJump = newPotential - potential;
188 /* Update the free energy estimates and bias and other history dependent method parameters */
189 if (params_.isUpdateFreeEnergyStep(step))
191 state_.updateFreeEnergyAndAddSamplesToHistogram(dimParams_, grid_, params_, commRecord, ms,
192 t, step, fplog, &updateList_);
194 if (params_.convolveForce)
196 /* The update results in a potential jump, so we need the new convolved potential. */
197 double newPotential = -calcConvolvedBias(coordState.coordValue()) * params_.invBeta;
198 *potentialJump = newPotential - potential;
202 /* Return the potential. */
203 *awhPotential = potential;
205 /* Check the sampled histograms and potentially warn user if something is suspicious */
206 warnForHistogramAnomalies(t, step, fplog);
212 * Count the total number of samples / sample weight over all grid points.
214 * \param[in] pointState The state of the points in a bias.
215 * \returns the total sample count.
217 static int64_t countSamples(const std::vector<PointState>& pointState)
219 double numSamples = 0;
220 for (const PointState& point : pointState)
222 numSamples += point.weightSumTot();
225 return gmx::roundToInt64(numSamples);
229 * Check if the state (loaded from checkpoint) and the run are consistent.
231 * When the state and the run setup are inconsistent, an exception is thrown.
233 * \param[in] params The parameters of the bias.
234 * \param[in] state The state of the bias.
236 static void ensureStateAndRunConsistency(const BiasParams& params, const BiasState& state)
238 int64_t numSamples = countSamples(state.points());
239 int64_t numUpdatesFromSamples =
240 numSamples / (params.numSamplesUpdateFreeEnergy_ * params.numSharedUpdate);
241 int64_t numUpdatesExpected = state.histogramSize().numUpdates();
242 if (numUpdatesFromSamples != numUpdatesExpected)
244 std::string mesg = gmx::formatString(
245 "The number of AWH updates in the checkpoint file (%" PRId64
246 ") does not match the total number of AWH samples divided by the number of samples "
247 "per update for %d sharing AWH bias(es) (%" PRId64 "/%d=%" PRId64 ")",
248 numUpdatesExpected, params.numSharedUpdate, numSamples,
249 params.numSamplesUpdateFreeEnergy_ * params.numSharedUpdate, numUpdatesFromSamples);
250 mesg += " Maybe you changed AWH parameters.";
251 /* Unfortunately we currently do not store the number of simulations
252 * sharing the bias or the state to checkpoint. But we can hint at
253 * a mismatch in the number of sharing simulations.
255 if (numUpdatesFromSamples % state.histogramSize().numUpdates() == 0)
257 mesg += gmx::formatString(
258 " Or the run you continued from used %" PRId64
259 " sharing simulations, whereas you now specified %d sharing simulations.",
260 numUpdatesFromSamples / state.histogramSize().numUpdates(), params.numSharedUpdate);
262 GMX_THROW(InvalidInputError(mesg));
266 void Bias::restoreStateFromHistory(const AwhBiasHistory* biasHistory, const t_commrec* cr)
268 GMX_RELEASE_ASSERT(thisRankDoesIO_ == MASTER(cr),
269 "The master rank should do I/O, the other ranks should not");
273 GMX_RELEASE_ASSERT(biasHistory != nullptr,
274 "On the master rank we need a valid history object to restore from");
275 state_.restoreFromHistory(*biasHistory, grid_);
277 /* Ensure that the state is consistent with our current run setup,
278 * since the user can have changed AWH parameters or the number
279 * of simulations sharing the bias.
281 ensureStateAndRunConsistency(params_, state_);
283 if (forceCorrelationGrid_ != nullptr)
285 forceCorrelationGrid_->restoreStateFromHistory(biasHistory->forceCorrelationGrid);
291 state_.broadcast(cr);
295 void Bias::initHistoryFromState(AwhBiasHistory* biasHistory) const
297 GMX_RELEASE_ASSERT(biasHistory != nullptr, "Need a valid biasHistory");
299 state_.initHistoryFromState(biasHistory);
301 if (forceCorrelationGrid_ != nullptr)
303 biasHistory->forceCorrelationGrid = initCorrelationGridHistoryFromState(forceCorrelationGrid());
307 void Bias::updateHistory(AwhBiasHistory* biasHistory) const
309 GMX_RELEASE_ASSERT(biasHistory != nullptr, "Need a valid biasHistory");
311 state_.updateHistory(biasHistory, grid_);
313 if (forceCorrelationGrid_ != nullptr)
315 updateCorrelationGridHistory(&biasHistory->forceCorrelationGrid, forceCorrelationGrid());
319 Bias::Bias(int biasIndexInCollection,
320 const AwhParams& awhParams,
321 const AwhBiasParams& awhBiasParams,
322 const std::vector<DimParams>& dimParamsInit,
325 int numSharingSimulations,
326 const std::string& biasInitFilename,
327 ThisRankWillDoIO thisRankWillDoIO,
328 BiasParams::DisableUpdateSkips disableUpdateSkips) :
329 dimParams_(dimParamsInit),
330 grid_(dimParamsInit, awhBiasParams.dimParams),
337 numSharingSimulations,
339 biasIndexInCollection),
340 state_(awhBiasParams, params_.initialHistogramSize, dimParams_, grid_),
341 thisRankDoesIO_(thisRankWillDoIO == ThisRankWillDoIO::Yes),
345 numWarningsIssued_(0)
347 /* For a global update updateList covers all points, so reserve that */
348 updateList_.reserve(grid_.numPoints());
350 state_.initGridPointState(awhBiasParams, dimParams_, grid_, params_, biasInitFilename, awhParams.numBias);
354 /* Set up the force correlation object. */
356 /* We let the correlation init function set its parameters
357 * to something useful for now.
359 double blockLength = 0;
360 /* Construct the force correlation object. */
361 forceCorrelationGrid_ = std::make_unique<CorrelationGrid>(
362 state_.points().size(), ndim(), blockLength,
363 CorrelationGrid::BlockLengthMeasure::Time, awhParams.nstSampleCoord * mdTimeStep);
365 writer_ = std::make_unique<BiasWriter>(*this);
369 void Bias::printInitializationToLog(FILE* fplog) const
371 if (fplog != nullptr && forceCorrelationGrid_ != nullptr)
373 std::string prefix = gmx::formatString("\nawh%d:", params_.biasIndex + 1);
376 "%s initial force correlation block length = %g %s"
377 "%s force correlation number of blocks = %d",
378 prefix.c_str(), forceCorrelationGrid().getBlockLength(),
379 forceCorrelationGrid().blockLengthMeasure == CorrelationGrid::BlockLengthMeasure::Weight
382 prefix.c_str(), forceCorrelationGrid().getNumBlocks());
386 void Bias::updateForceCorrelationGrid(gmx::ArrayRef<const double> probWeightNeighbor, double t)
388 if (forceCorrelationGrid_ == nullptr)
393 const std::vector<int>& neighbor = grid_.point(state_.coordState().gridpointIndex()).neighbor;
395 gmx::ArrayRef<double> forceFromNeighbor = tempForce_;
396 for (size_t n = 0; n < neighbor.size(); n++)
398 double weightNeighbor = probWeightNeighbor[n];
399 int indexNeighbor = neighbor[n];
401 /* Add the force data of this neighbor point. Note: the sum of these forces is the convolved force.
403 We actually add the force normalized by beta which has the units of 1/length. This means that the
404 resulting correlation time integral is directly in units of friction time/length^2 which is really what
405 we're interested in. */
406 state_.calcUmbrellaForceAndPotential(dimParams_, grid_, indexNeighbor, forceFromNeighbor);
408 /* Note: we might want to give a whole list of data to add instead and have this loop in the data adding function */
409 forceCorrelationGrid_->addData(indexNeighbor, weightNeighbor, forceFromNeighbor, t);
413 /* Return the number of data blocks that have been prepared for writing. */
414 int Bias::numEnergySubblocksToWrite() const
416 GMX_RELEASE_ASSERT(writer_ != nullptr, "Should only request data from an initialized writer");
418 return writer_->numBlocks();
421 /* Write bias data blocks to energy subblocks. */
422 int Bias::writeToEnergySubblocks(t_enxsubblock* subblock) const
424 GMX_RELEASE_ASSERT(writer_ != nullptr, "Should only request data from an initialized writer");
426 return writer_->writeToEnergySubblocks(*this, subblock);