2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017,2018,2019,2020, 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>{},
144 &probWeightNeighbor);
146 if (isSampleCoordStep)
148 updateForceCorrelationGrid(probWeightNeighbor, neighborLambdaDhdl, t);
150 state_.sampleCoordAndPmf(dimParams_, grid_, probWeightNeighbor, convolvedBias);
152 /* Set the umbrella grid point (for the lambda axis) to the
153 * current grid point. */
154 if (params_.convolveForce && grid_.hasLambdaAxis())
156 state_.setUmbrellaGridpointToGridpoint();
160 /* Set the bias force and get the potential contribution from this bias.
161 * The potential jump occurs at different times depending on how
162 * the force is applied (and how the potential is normalized).
163 * For the convolved force it happens when the bias is updated,
164 * for the umbrella when the umbrella is moved.
168 if (params_.convolveForce)
170 state_.calcConvolvedForce(dimParams_, grid_, probWeightNeighbor,
171 moveUmbrella ? neighborLambdaDhdl : ArrayRef<const double>{},
172 tempForce_, biasForce_);
174 potential = -convolvedBias * params_.invBeta;
179 GMX_RELEASE_ASSERT(state_.points()[coordState.umbrellaGridpoint()].inTargetRegion(),
180 "AWH bias grid point for the umbrella reference value is outside of the "
182 potential = state_.calcUmbrellaForceAndPotential(
183 dimParams_, grid_, coordState.umbrellaGridpoint(),
184 moveUmbrella ? neighborLambdaDhdl : ArrayRef<const double>{}, biasForce_);
186 /* Moving the umbrella results in a force correction and
187 * a new potential. The umbrella center is sampled as often as
188 * the coordinate so we know the probability weights needed
189 * for moving the umbrella are up-to-date.
193 const bool onlySampleUmbrellaGridpoint = false;
194 double newPotential = state_.moveUmbrella(dimParams_, grid_, probWeightNeighbor,
195 neighborLambdaDhdl, biasForce_, step, seed,
196 params_.biasIndex, onlySampleUmbrellaGridpoint);
197 *potentialJump = newPotential - potential;
201 /* Update the free energy estimates and bias and other history dependent method parameters */
202 if (params_.isUpdateFreeEnergyStep(step))
204 state_.updateFreeEnergyAndAddSamplesToHistogram(dimParams_, grid_, params_, commRecord, ms,
205 t, step, fplog, &updateList_);
207 if (params_.convolveForce)
209 /* The update results in a potential jump, so we need the new convolved potential. */
210 double newPotential = -calcConvolvedBias(coordState.coordValue()) * params_.invBeta;
211 *potentialJump = newPotential - potential;
214 /* If there is a lambda axis it is still controlled using an umbrella even if the force
215 * is convolved in the other dimensions. */
216 if (moveUmbrella && params_.convolveForce && grid_.hasLambdaAxis())
218 const bool onlySampleUmbrellaGridpoint = true;
219 state_.moveUmbrella(dimParams_, grid_, probWeightNeighbor, neighborLambdaDhdl, biasForce_,
220 step, seed, params_.biasIndex, onlySampleUmbrellaGridpoint);
223 /* Return the potential. */
224 *awhPotential = potential;
226 /* Check the sampled histograms and potentially warn user if something is suspicious */
227 warnForHistogramAnomalies(t, step, fplog);
233 * Count the total number of samples / sample weight over all grid points.
235 * \param[in] pointState The state of the points in a bias.
236 * \returns the total sample count.
238 static int64_t countSamples(const std::vector<PointState>& pointState)
240 double numSamples = 0;
241 for (const PointState& point : pointState)
243 numSamples += point.weightSumTot();
246 return gmx::roundToInt64(numSamples);
250 * Check if the state (loaded from checkpoint) and the run are consistent.
252 * When the state and the run setup are inconsistent, an exception is thrown.
254 * \param[in] params The parameters of the bias.
255 * \param[in] state The state of the bias.
257 static void ensureStateAndRunConsistency(const BiasParams& params, const BiasState& state)
259 int64_t numSamples = countSamples(state.points());
260 int64_t numUpdatesFromSamples =
261 numSamples / (params.numSamplesUpdateFreeEnergy_ * params.numSharedUpdate);
262 int64_t numUpdatesExpected = state.histogramSize().numUpdates();
263 if (numUpdatesFromSamples != numUpdatesExpected)
265 std::string mesg = gmx::formatString(
266 "The number of AWH updates in the checkpoint file (%" PRId64
267 ") does not match the total number of AWH samples divided by the number of samples "
268 "per update for %d sharing AWH bias(es) (%" PRId64 "/%d=%" PRId64 ")",
269 numUpdatesExpected, params.numSharedUpdate, numSamples,
270 params.numSamplesUpdateFreeEnergy_ * params.numSharedUpdate, numUpdatesFromSamples);
271 mesg += " Maybe you changed AWH parameters.";
272 /* Unfortunately we currently do not store the number of simulations
273 * sharing the bias or the state to checkpoint. But we can hint at
274 * a mismatch in the number of sharing simulations.
276 if (numUpdatesFromSamples % state.histogramSize().numUpdates() == 0)
278 mesg += gmx::formatString(
279 " Or the run you continued from used %" PRId64
280 " sharing simulations, whereas you now specified %d sharing simulations.",
281 numUpdatesFromSamples / state.histogramSize().numUpdates(), params.numSharedUpdate);
283 GMX_THROW(InvalidInputError(mesg));
287 void Bias::restoreStateFromHistory(const AwhBiasHistory* biasHistory, const t_commrec* cr)
289 GMX_RELEASE_ASSERT(thisRankDoesIO_ == MASTER(cr),
290 "The master rank should do I/O, the other ranks should not");
294 GMX_RELEASE_ASSERT(biasHistory != nullptr,
295 "On the master rank we need a valid history object to restore from");
296 state_.restoreFromHistory(*biasHistory, grid_);
298 /* Ensure that the state is consistent with our current run setup,
299 * since the user can have changed AWH parameters or the number
300 * of simulations sharing the bias.
302 ensureStateAndRunConsistency(params_, state_);
304 if (forceCorrelationGrid_ != nullptr)
306 forceCorrelationGrid_->restoreStateFromHistory(biasHistory->forceCorrelationGrid);
312 state_.broadcast(cr);
316 void Bias::initHistoryFromState(AwhBiasHistory* biasHistory) const
318 GMX_RELEASE_ASSERT(biasHistory != nullptr, "Need a valid biasHistory");
320 state_.initHistoryFromState(biasHistory);
322 if (forceCorrelationGrid_ != nullptr)
324 biasHistory->forceCorrelationGrid = initCorrelationGridHistoryFromState(forceCorrelationGrid());
328 void Bias::updateHistory(AwhBiasHistory* biasHistory) const
330 GMX_RELEASE_ASSERT(biasHistory != nullptr, "Need a valid biasHistory");
332 state_.updateHistory(biasHistory, grid_);
334 if (forceCorrelationGrid_ != nullptr)
336 updateCorrelationGridHistory(&biasHistory->forceCorrelationGrid, forceCorrelationGrid());
340 Bias::Bias(int biasIndexInCollection,
341 const AwhParams& awhParams,
342 const AwhBiasParams& awhBiasParams,
343 const std::vector<DimParams>& dimParamsInit,
346 int numSharingSimulations,
347 const std::string& biasInitFilename,
348 ThisRankWillDoIO thisRankWillDoIO,
349 BiasParams::DisableUpdateSkips disableUpdateSkips) :
350 dimParams_(dimParamsInit),
351 grid_(dimParamsInit, awhBiasParams.dimParams),
358 numSharingSimulations,
360 biasIndexInCollection),
361 state_(awhBiasParams, params_.initialHistogramSize, dimParams_, grid_),
362 thisRankDoesIO_(thisRankWillDoIO == ThisRankWillDoIO::Yes),
366 numWarningsIssued_(0)
368 /* For a global update updateList covers all points, so reserve that */
369 updateList_.reserve(grid_.numPoints());
371 state_.initGridPointState(awhBiasParams, dimParams_, grid_, params_, biasInitFilename, awhParams.numBias);
375 /* Set up the force correlation object. */
377 /* We let the correlation init function set its parameters
378 * to something useful for now.
380 double blockLength = 0;
381 /* Construct the force correlation object. */
382 forceCorrelationGrid_ = std::make_unique<CorrelationGrid>(
383 state_.points().size(), ndim(), blockLength,
384 CorrelationGrid::BlockLengthMeasure::Time, awhParams.nstSampleCoord * mdTimeStep);
386 writer_ = std::make_unique<BiasWriter>(*this);
390 void Bias::printInitializationToLog(FILE* fplog) const
392 if (fplog != nullptr && forceCorrelationGrid_ != nullptr)
394 std::string prefix = gmx::formatString("\nawh%d:", params_.biasIndex + 1);
397 "%s initial force correlation block length = %g %s"
398 "%s force correlation number of blocks = %d",
399 prefix.c_str(), forceCorrelationGrid().getBlockLength(),
400 forceCorrelationGrid().blockLengthMeasure == CorrelationGrid::BlockLengthMeasure::Weight
403 prefix.c_str(), forceCorrelationGrid().getNumBlocks());
407 void Bias::updateForceCorrelationGrid(gmx::ArrayRef<const double> probWeightNeighbor,
408 ArrayRef<const double> neighborLambdaDhdl,
411 if (forceCorrelationGrid_ == nullptr)
416 const std::vector<int>& neighbor = grid_.point(state_.coordState().gridpointIndex()).neighbor;
418 gmx::ArrayRef<double> forceFromNeighbor = tempForce_;
419 for (size_t n = 0; n < neighbor.size(); n++)
421 double weightNeighbor = probWeightNeighbor[n];
422 int indexNeighbor = neighbor[n];
424 /* Add the force data of this neighbor point. Note: the sum of these forces is the convolved force.
426 We actually add the force normalized by beta which has the units of 1/length. This means that the
427 resulting correlation time integral is directly in units of friction time/length^2 which is really what
428 we're interested in. */
429 state_.calcUmbrellaForceAndPotential(dimParams_, grid_, indexNeighbor, neighborLambdaDhdl,
432 /* Note: we might want to give a whole list of data to add instead and have this loop in the data adding function */
433 forceCorrelationGrid_->addData(indexNeighbor, weightNeighbor, forceFromNeighbor, t);
437 /* Return the number of data blocks that have been prepared for writing. */
438 int Bias::numEnergySubblocksToWrite() const
440 GMX_RELEASE_ASSERT(writer_ != nullptr, "Should only request data from an initialized writer");
442 return writer_->numBlocks();
445 /* Write bias data blocks to energy subblocks. */
446 int Bias::writeToEnergySubblocks(t_enxsubblock* subblock) const
448 GMX_RELEASE_ASSERT(writer_ != nullptr, "Should only request data from an initialized writer");
450 return writer_->writeToEnergySubblocks(*this, subblock);
453 bool Bias::isFepLambdaDimension(int dim) const
455 GMX_ASSERT(dim < ndim(), "Requested dimension out of range.");
457 return dimParams_[dim].isFepLambdaDimension();
460 bool Bias::isSampleCoordStep(const int64_t step) const
462 return params_.isSampleCoordStep(step);