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 state_.setCoordValue(grid_, coordValue);
125 std::vector<double, AlignedAllocator<double>>& probWeightNeighbor = alignedTempWorkSpace_;
127 /* If the convolved force is needed or this is a sampling step,
128 * the bias in the current neighborhood needs to be up-to-date
129 * and the probablity weights need to be calculated.
131 const bool isSampleCoordStep = params_.isSampleCoordStep(step);
132 const bool moveUmbrella = (isSampleCoordStep || step == 0);
133 double convolvedBias = 0;
134 const CoordState& coordState = state_.coordState();
136 if (params_.convolveForce || moveUmbrella || isSampleCoordStep)
138 if (params_.skipUpdates())
140 state_.doSkippedUpdatesInNeighborhood(params_, grid_);
142 convolvedBias = state_.updateProbabilityWeightsAndConvolvedBias(
143 dimParams_, grid_, moveUmbrella ? neighborLambdaEnergies : ArrayRef<const double>{}, &probWeightNeighbor);
145 if (isSampleCoordStep)
147 updateForceCorrelationGrid(probWeightNeighbor, neighborLambdaDhdl, t);
149 state_.sampleCoordAndPmf(dimParams_, grid_, probWeightNeighbor, convolvedBias);
151 /* Set the umbrella grid point (for the lambda axis) to the
152 * current grid point. */
153 if (params_.convolveForce && grid_.hasLambdaAxis())
155 state_.setUmbrellaGridpointToGridpoint();
159 /* Set the bias force and get the potential contribution from this bias.
160 * The potential jump occurs at different times depending on how
161 * the force is applied (and how the potential is normalized).
162 * For the convolved force it happens when the bias is updated,
163 * for the umbrella when the umbrella is moved.
167 if (params_.convolveForce)
169 state_.calcConvolvedForce(dimParams_,
172 moveUmbrella ? neighborLambdaDhdl : ArrayRef<const double>{},
176 potential = -convolvedBias * params_.invBeta;
181 GMX_RELEASE_ASSERT(state_.points()[coordState.umbrellaGridpoint()].inTargetRegion(),
182 "AWH bias grid point for the umbrella reference value is outside of the "
184 potential = state_.calcUmbrellaForceAndPotential(
187 coordState.umbrellaGridpoint(),
188 moveUmbrella ? neighborLambdaDhdl : ArrayRef<const double>{},
191 /* Moving the umbrella results in a force correction and
192 * a new potential. The umbrella center is sampled as often as
193 * the coordinate so we know the probability weights needed
194 * for moving the umbrella are up-to-date.
198 const bool onlySampleUmbrellaGridpoint = false;
199 double newPotential = state_.moveUmbrella(dimParams_,
207 onlySampleUmbrellaGridpoint);
208 *potentialJump = newPotential - potential;
212 /* Update the free energy estimates and bias and other history dependent method parameters */
213 if (params_.isUpdateFreeEnergyStep(step))
215 state_.updateFreeEnergyAndAddSamplesToHistogram(
216 dimParams_, grid_, params_, commRecord, ms, t, step, fplog, &updateList_);
218 if (params_.convolveForce)
220 /* The update results in a potential jump, so we need the new convolved potential. */
221 double newPotential = -calcConvolvedBias(coordState.coordValue()) * params_.invBeta;
222 *potentialJump = newPotential - potential;
225 /* If there is a lambda axis it is still controlled using an umbrella even if the force
226 * is convolved in the other dimensions. */
227 if (moveUmbrella && params_.convolveForce && grid_.hasLambdaAxis())
229 const bool onlySampleUmbrellaGridpoint = true;
230 state_.moveUmbrella(dimParams_,
238 onlySampleUmbrellaGridpoint);
241 /* Return the potential. */
242 *awhPotential = potential;
244 /* Check the sampled histograms and potentially warn user if something is suspicious */
245 warnForHistogramAnomalies(t, step, fplog);
251 * Count the total number of samples / sample weight over all grid points.
253 * \param[in] pointState The state of the points in a bias.
254 * \returns the total sample count.
256 static int64_t countSamples(ArrayRef<const PointState> pointState)
258 double numSamples = 0;
259 for (const PointState& point : pointState)
261 numSamples += point.weightSumTot();
264 return gmx::roundToInt64(numSamples);
268 * Check if the state (loaded from checkpoint) and the run are consistent.
270 * When the state and the run setup are inconsistent, an exception is thrown.
272 * \param[in] params The parameters of the bias.
273 * \param[in] state The state of the bias.
275 static void ensureStateAndRunConsistency(const BiasParams& params, const BiasState& state)
277 int64_t numSamples = countSamples(state.points());
278 int64_t numUpdatesFromSamples =
279 numSamples / (params.numSamplesUpdateFreeEnergy_ * params.numSharedUpdate);
280 int64_t numUpdatesExpected = state.histogramSize().numUpdates();
281 if (numUpdatesFromSamples != numUpdatesExpected)
283 std::string mesg = gmx::formatString(
284 "The number of AWH updates in the checkpoint file (%" PRId64
285 ") does not match the total number of AWH samples divided by the number of samples "
286 "per update for %d sharing AWH bias(es) (%" PRId64 "/%d=%" PRId64 ")",
288 params.numSharedUpdate,
290 params.numSamplesUpdateFreeEnergy_ * params.numSharedUpdate,
291 numUpdatesFromSamples);
292 mesg += " Maybe you changed AWH parameters.";
293 /* Unfortunately we currently do not store the number of simulations
294 * sharing the bias or the state to checkpoint. But we can hint at
295 * a mismatch in the number of sharing simulations.
297 if (numUpdatesFromSamples % state.histogramSize().numUpdates() == 0)
299 mesg += gmx::formatString(
300 " Or the run you continued from used %" PRId64
301 " sharing simulations, whereas you now specified %d sharing simulations.",
302 numUpdatesFromSamples / state.histogramSize().numUpdates(),
303 params.numSharedUpdate);
305 GMX_THROW(InvalidInputError(mesg));
309 void Bias::restoreStateFromHistory(const AwhBiasHistory* biasHistory, const t_commrec* cr)
311 GMX_RELEASE_ASSERT(thisRankDoesIO_ == MASTER(cr),
312 "The master rank should do I/O, the other ranks should not");
316 GMX_RELEASE_ASSERT(biasHistory != nullptr,
317 "On the master rank we need a valid history object to restore from");
318 state_.restoreFromHistory(*biasHistory, grid_);
320 /* Ensure that the state is consistent with our current run setup,
321 * since the user can have changed AWH parameters or the number
322 * of simulations sharing the bias.
324 ensureStateAndRunConsistency(params_, state_);
326 if (forceCorrelationGrid_ != nullptr)
328 forceCorrelationGrid_->restoreStateFromHistory(biasHistory->forceCorrelationGrid);
334 state_.broadcast(cr);
338 void Bias::initHistoryFromState(AwhBiasHistory* biasHistory) const
340 GMX_RELEASE_ASSERT(biasHistory != nullptr, "Need a valid biasHistory");
342 state_.initHistoryFromState(biasHistory);
344 if (forceCorrelationGrid_ != nullptr)
346 biasHistory->forceCorrelationGrid = initCorrelationGridHistoryFromState(forceCorrelationGrid());
350 void Bias::updateHistory(AwhBiasHistory* biasHistory) const
352 GMX_RELEASE_ASSERT(biasHistory != nullptr, "Need a valid biasHistory");
354 state_.updateHistory(biasHistory, grid_);
356 if (forceCorrelationGrid_ != nullptr)
358 updateCorrelationGridHistory(&biasHistory->forceCorrelationGrid, forceCorrelationGrid());
362 Bias::Bias(int biasIndexInCollection,
363 const AwhParams& awhParams,
364 const AwhBiasParams& awhBiasParams,
365 ArrayRef<const DimParams> dimParamsInit,
368 int numSharingSimulations,
369 const std::string& biasInitFilename,
370 ThisRankWillDoIO thisRankWillDoIO,
371 BiasParams::DisableUpdateSkips disableUpdateSkips) :
372 dimParams_(dimParamsInit.begin(), dimParamsInit.end()),
373 grid_(dimParamsInit, awhBiasParams.dimParams),
380 numSharingSimulations,
382 biasIndexInCollection),
383 state_(awhBiasParams, params_.initialHistogramSize, dimParams_, grid_),
384 thisRankDoesIO_(thisRankWillDoIO == ThisRankWillDoIO::Yes),
388 numWarningsIssued_(0)
390 /* For a global update updateList covers all points, so reserve that */
391 updateList_.reserve(grid_.numPoints());
393 state_.initGridPointState(awhBiasParams, dimParams_, grid_, params_, biasInitFilename, awhParams.numBias);
397 /* Set up the force correlation object. */
399 /* We let the correlation init function set its parameters
400 * to something useful for now.
402 double blockLength = 0;
403 /* Construct the force correlation object. */
404 forceCorrelationGrid_ =
405 std::make_unique<CorrelationGrid>(state_.points().size(),
408 CorrelationGrid::BlockLengthMeasure::Time,
409 awhParams.nstSampleCoord * mdTimeStep);
411 writer_ = std::make_unique<BiasWriter>(*this);
415 void Bias::printInitializationToLog(FILE* fplog) const
417 if (fplog != nullptr && forceCorrelationGrid_ != nullptr)
419 std::string prefix = gmx::formatString("\nawh%d:", params_.biasIndex + 1);
422 "%s initial force correlation block length = %g %s"
423 "%s force correlation number of blocks = %d",
425 forceCorrelationGrid().getBlockLength(),
426 forceCorrelationGrid().blockLengthMeasure == CorrelationGrid::BlockLengthMeasure::Weight
430 forceCorrelationGrid().getNumBlocks());
434 void Bias::updateForceCorrelationGrid(gmx::ArrayRef<const double> probWeightNeighbor,
435 ArrayRef<const double> neighborLambdaDhdl,
438 if (forceCorrelationGrid_ == nullptr)
443 const std::vector<int>& neighbor = grid_.point(state_.coordState().gridpointIndex()).neighbor;
445 gmx::ArrayRef<double> forceFromNeighbor = tempForce_;
446 for (size_t n = 0; n < neighbor.size(); n++)
448 double weightNeighbor = probWeightNeighbor[n];
449 int indexNeighbor = neighbor[n];
451 /* Add the force data of this neighbor point. Note: the sum of these forces is the convolved force.
453 We actually add the force normalized by beta which has the units of 1/length. This means that the
454 resulting correlation time integral is directly in units of friction time/length^2 which is really what
455 we're interested in. */
456 state_.calcUmbrellaForceAndPotential(
457 dimParams_, grid_, indexNeighbor, neighborLambdaDhdl, forceFromNeighbor);
459 /* Note: we might want to give a whole list of data to add instead and have this loop in the data adding function */
460 forceCorrelationGrid_->addData(indexNeighbor, weightNeighbor, forceFromNeighbor, t);
464 /* Return the number of data blocks that have been prepared for writing. */
465 int Bias::numEnergySubblocksToWrite() const
467 GMX_RELEASE_ASSERT(writer_ != nullptr, "Should only request data from an initialized writer");
469 return writer_->numBlocks();
472 /* Write bias data blocks to energy subblocks. */
473 int Bias::writeToEnergySubblocks(t_enxsubblock* subblock) const
475 GMX_RELEASE_ASSERT(writer_ != nullptr, "Should only request data from an initialized writer");
477 return writer_->writeToEnergySubblocks(*this, subblock);
480 bool Bias::isSampleCoordStep(const int64_t step) const
482 return params_.isSampleCoordStep(step);