2 * This file is part of the GROMACS molecular simulation package.
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.
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 ArrayRef<const double> neighborLambdaEnergies,
107 ArrayRef<const double> neighborLambdaDhdl,
108 double* awhPotential,
109 double* potentialJump,
110 const t_commrec* commRecord,
111 const gmx_multisim_t* ms,
119 GMX_THROW(InvalidInputError(
120 "The step number is negative which is not supported by the AWH code."));
123 GMX_RELEASE_ASSERT(!(params_.convolveForce && grid_.hasLambdaAxis()),
124 "When using AWH to sample an FEP lambda dimension the AWH potential cannot "
127 state_.setCoordValue(grid_, coordValue);
129 std::vector<double, AlignedAllocator<double>>& probWeightNeighbor = alignedTempWorkSpace_;
131 /* If the convolved force is needed or this is a sampling step,
132 * the bias in the current neighborhood needs to be up-to-date
133 * and the probablity weights need to be calculated.
135 const bool isSampleCoordStep = params_.isSampleCoordStep(step);
136 const bool moveUmbrella = (isSampleCoordStep || step == 0);
137 double convolvedBias = 0;
138 const CoordState& coordState = state_.coordState();
140 if (params_.convolveForce || moveUmbrella || isSampleCoordStep)
142 if (params_.skipUpdates())
144 state_.doSkippedUpdatesInNeighborhood(params_, grid_);
146 convolvedBias = state_.updateProbabilityWeightsAndConvolvedBias(
147 dimParams_, grid_, moveUmbrella ? neighborLambdaEnergies : ArrayRef<const double>{},
148 &probWeightNeighbor);
150 if (isSampleCoordStep)
152 updateForceCorrelationGrid(probWeightNeighbor, neighborLambdaDhdl, t);
154 state_.sampleCoordAndPmf(dimParams_, grid_, probWeightNeighbor, convolvedBias);
158 /* Set the bias force and get the potential contribution from this bias.
159 * The potential jump occurs at different times depending on how
160 * the force is applied (and how the potential is normalized).
161 * For the convolved force it happens when the bias is updated,
162 * for the umbrella when the umbrella is moved.
166 if (params_.convolveForce)
168 state_.calcConvolvedForce(dimParams_, grid_, probWeightNeighbor,
169 moveUmbrella ? neighborLambdaDhdl : ArrayRef<const double>{},
170 tempForce_, biasForce_);
172 potential = -convolvedBias * params_.invBeta;
177 GMX_RELEASE_ASSERT(state_.points()[coordState.umbrellaGridpoint()].inTargetRegion(),
178 "AWH bias grid point for the umbrella reference value is outside of the "
180 potential = state_.calcUmbrellaForceAndPotential(
181 dimParams_, grid_, coordState.umbrellaGridpoint(),
182 moveUmbrella ? neighborLambdaDhdl : ArrayRef<const double>{}, biasForce_);
184 /* Moving the umbrella results in a force correction and
185 * a new potential. The umbrella center is sampled as often as
186 * the coordinate so we know the probability weights needed
187 * for moving the umbrella are up-to-date.
191 double newPotential =
192 state_.moveUmbrella(dimParams_, grid_, probWeightNeighbor, neighborLambdaDhdl,
193 biasForce_, step, seed, params_.biasIndex);
194 *potentialJump = newPotential - potential;
198 /* Update the free energy estimates and bias and other history dependent method parameters */
199 if (params_.isUpdateFreeEnergyStep(step))
201 state_.updateFreeEnergyAndAddSamplesToHistogram(dimParams_, grid_, params_, commRecord, ms,
202 t, step, fplog, &updateList_);
204 if (params_.convolveForce)
206 /* The update results in a potential jump, so we need the new convolved potential. */
207 double newPotential = -calcConvolvedBias(coordState.coordValue()) * params_.invBeta;
208 *potentialJump = newPotential - potential;
212 /* Return the potential. */
213 *awhPotential = potential;
215 /* Check the sampled histograms and potentially warn user if something is suspicious */
216 warnForHistogramAnomalies(t, step, fplog);
222 * Count the total number of samples / sample weight over all grid points.
224 * \param[in] pointState The state of the points in a bias.
225 * \returns the total sample count.
227 static int64_t countSamples(const std::vector<PointState>& pointState)
229 double numSamples = 0;
230 for (const PointState& point : pointState)
232 numSamples += point.weightSumTot();
235 return gmx::roundToInt64(numSamples);
239 * Check if the state (loaded from checkpoint) and the run are consistent.
241 * When the state and the run setup are inconsistent, an exception is thrown.
243 * \param[in] params The parameters of the bias.
244 * \param[in] state The state of the bias.
246 static void ensureStateAndRunConsistency(const BiasParams& params, const BiasState& state)
248 int64_t numSamples = countSamples(state.points());
249 int64_t numUpdatesFromSamples =
250 numSamples / (params.numSamplesUpdateFreeEnergy_ * params.numSharedUpdate);
251 int64_t numUpdatesExpected = state.histogramSize().numUpdates();
252 if (numUpdatesFromSamples != numUpdatesExpected)
254 std::string mesg = gmx::formatString(
255 "The number of AWH updates in the checkpoint file (%" PRId64
256 ") does not match the total number of AWH samples divided by the number of samples "
257 "per update for %d sharing AWH bias(es) (%" PRId64 "/%d=%" PRId64 ")",
258 numUpdatesExpected, params.numSharedUpdate, numSamples,
259 params.numSamplesUpdateFreeEnergy_ * params.numSharedUpdate, numUpdatesFromSamples);
260 mesg += " Maybe you changed AWH parameters.";
261 /* Unfortunately we currently do not store the number of simulations
262 * sharing the bias or the state to checkpoint. But we can hint at
263 * a mismatch in the number of sharing simulations.
265 if (numUpdatesFromSamples % state.histogramSize().numUpdates() == 0)
267 mesg += gmx::formatString(
268 " Or the run you continued from used %" PRId64
269 " sharing simulations, whereas you now specified %d sharing simulations.",
270 numUpdatesFromSamples / state.histogramSize().numUpdates(), params.numSharedUpdate);
272 GMX_THROW(InvalidInputError(mesg));
276 void Bias::restoreStateFromHistory(const AwhBiasHistory* biasHistory, const t_commrec* cr)
278 GMX_RELEASE_ASSERT(thisRankDoesIO_ == MASTER(cr),
279 "The master rank should do I/O, the other ranks should not");
283 GMX_RELEASE_ASSERT(biasHistory != nullptr,
284 "On the master rank we need a valid history object to restore from");
285 state_.restoreFromHistory(*biasHistory, grid_);
287 /* Ensure that the state is consistent with our current run setup,
288 * since the user can have changed AWH parameters or the number
289 * of simulations sharing the bias.
291 ensureStateAndRunConsistency(params_, state_);
293 if (forceCorrelationGrid_ != nullptr)
295 forceCorrelationGrid_->restoreStateFromHistory(biasHistory->forceCorrelationGrid);
301 state_.broadcast(cr);
305 void Bias::initHistoryFromState(AwhBiasHistory* biasHistory) const
307 GMX_RELEASE_ASSERT(biasHistory != nullptr, "Need a valid biasHistory");
309 state_.initHistoryFromState(biasHistory);
311 if (forceCorrelationGrid_ != nullptr)
313 biasHistory->forceCorrelationGrid = initCorrelationGridHistoryFromState(forceCorrelationGrid());
317 void Bias::updateHistory(AwhBiasHistory* biasHistory) const
319 GMX_RELEASE_ASSERT(biasHistory != nullptr, "Need a valid biasHistory");
321 state_.updateHistory(biasHistory, grid_);
323 if (forceCorrelationGrid_ != nullptr)
325 updateCorrelationGridHistory(&biasHistory->forceCorrelationGrid, forceCorrelationGrid());
329 Bias::Bias(int biasIndexInCollection,
330 const AwhParams& awhParams,
331 const AwhBiasParams& awhBiasParams,
332 const std::vector<DimParams>& dimParamsInit,
335 int numSharingSimulations,
336 const std::string& biasInitFilename,
337 ThisRankWillDoIO thisRankWillDoIO,
338 BiasParams::DisableUpdateSkips disableUpdateSkips) :
339 dimParams_(dimParamsInit),
340 grid_(dimParamsInit, awhBiasParams.dimParams),
347 numSharingSimulations,
349 biasIndexInCollection),
350 state_(awhBiasParams, params_.initialHistogramSize, dimParams_, grid_),
351 thisRankDoesIO_(thisRankWillDoIO == ThisRankWillDoIO::Yes),
355 numWarningsIssued_(0)
357 /* For a global update updateList covers all points, so reserve that */
358 updateList_.reserve(grid_.numPoints());
360 state_.initGridPointState(awhBiasParams, dimParams_, grid_, params_, biasInitFilename, awhParams.numBias);
364 /* Set up the force correlation object. */
366 /* We let the correlation init function set its parameters
367 * to something useful for now.
369 double blockLength = 0;
370 /* Construct the force correlation object. */
371 forceCorrelationGrid_ = std::make_unique<CorrelationGrid>(
372 state_.points().size(), ndim(), blockLength,
373 CorrelationGrid::BlockLengthMeasure::Time, awhParams.nstSampleCoord * mdTimeStep);
375 writer_ = std::make_unique<BiasWriter>(*this);
379 void Bias::printInitializationToLog(FILE* fplog) const
381 if (fplog != nullptr && forceCorrelationGrid_ != nullptr)
383 std::string prefix = gmx::formatString("\nawh%d:", params_.biasIndex + 1);
386 "%s initial force correlation block length = %g %s"
387 "%s force correlation number of blocks = %d",
388 prefix.c_str(), forceCorrelationGrid().getBlockLength(),
389 forceCorrelationGrid().blockLengthMeasure == CorrelationGrid::BlockLengthMeasure::Weight
392 prefix.c_str(), forceCorrelationGrid().getNumBlocks());
396 void Bias::updateForceCorrelationGrid(gmx::ArrayRef<const double> probWeightNeighbor,
397 ArrayRef<const double> neighborLambdaDhdl,
400 if (forceCorrelationGrid_ == nullptr)
405 const std::vector<int>& neighbor = grid_.point(state_.coordState().gridpointIndex()).neighbor;
407 gmx::ArrayRef<double> forceFromNeighbor = tempForce_;
408 for (size_t n = 0; n < neighbor.size(); n++)
410 double weightNeighbor = probWeightNeighbor[n];
411 int indexNeighbor = neighbor[n];
413 /* Add the force data of this neighbor point. Note: the sum of these forces is the convolved force.
415 We actually add the force normalized by beta which has the units of 1/length. This means that the
416 resulting correlation time integral is directly in units of friction time/length^2 which is really what
417 we're interested in. */
418 state_.calcUmbrellaForceAndPotential(dimParams_, grid_, indexNeighbor, neighborLambdaDhdl,
421 /* Note: we might want to give a whole list of data to add instead and have this loop in the data adding function */
422 forceCorrelationGrid_->addData(indexNeighbor, weightNeighbor, forceFromNeighbor, t);
426 /* Return the number of data blocks that have been prepared for writing. */
427 int Bias::numEnergySubblocksToWrite() const
429 GMX_RELEASE_ASSERT(writer_ != nullptr, "Should only request data from an initialized writer");
431 return writer_->numBlocks();
434 /* Write bias data blocks to energy subblocks. */
435 int Bias::writeToEnergySubblocks(t_enxsubblock* subblock) const
437 GMX_RELEASE_ASSERT(writer_ != nullptr, "Should only request data from an initialized writer");
439 return writer_->writeToEnergySubblocks(*this, subblock);
442 bool Bias::isSampleCoordStep(const int64_t step) const
444 return params_.isSampleCoordStep(step);