2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017,2018, 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>
58 #include "gromacs/compat/make_unique.h"
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,
90 maxNumWarningsInCheck);
92 if (numWarningsIssued_ >= maxNumWarningsInRun)
94 fprintf(fplog, "\nawh%d: suppressing future AWH warnings.\n", biasIndex() + 1);
98 void Bias::doSkippedUpdatesForAllPoints()
100 if (params_.skipUpdates())
102 state_.doSkippedUpdatesForAllPoints(params_);
106 gmx::ArrayRef<const double>
107 Bias::calcForceAndUpdateBias(const awh_dvec coordValue,
108 double *awhPotential,
109 double *potentialJump,
110 const t_commrec *commRecord,
111 const gmx_multisim_t *ms,
119 GMX_THROW(InvalidInputError("The step number is negative which is not supported by the AWH code."));
122 state_.setCoordValue(grid_, coordValue);
124 std::vector < double, AlignedAllocator < double>> &probWeightNeighbor = alignedTempWorkSpace_;
126 /* If the convolved force is needed or this is a sampling step,
127 * the bias in the current neighborhood needs to be up-to-date
128 * and the probablity weights need to be calculated.
130 const bool sampleCoord = params_.isSampleCoordStep(step);
131 const bool moveUmbrella = (sampleCoord || step == 0);
132 double convolvedBias = 0;
133 if (params_.convolveForce || moveUmbrella || sampleCoord)
135 if (params_.skipUpdates())
137 state_.doSkippedUpdatesInNeighborhood(params_, grid_);
140 convolvedBias = 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,
163 tempForce_, biasForce_);
165 potential = -convolvedBias*params_.invBeta;
170 GMX_RELEASE_ASSERT(state_.points()[coordState.umbrellaGridpoint()].inTargetRegion(),
171 "AWH bias grid point for the umbrella reference value is outside of the target region.");
173 state_.calcUmbrellaForceAndPotential(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, biasForce_, step, seed, params_.biasIndex);
183 *potentialJump = newPotential - potential;
187 /* Update the free energy estimates and bias and other history dependent method parameters */
188 if (params_.isUpdateFreeEnergyStep(step))
190 state_.updateFreeEnergyAndAddSamplesToHistogram(dimParams_, grid_,
196 if (params_.convolveForce)
198 /* The update results in a potential jump, so we need the new convolved potential. */
199 double newPotential = -calcConvolvedBias(coordState.coordValue())*params_.invBeta;
200 *potentialJump = newPotential - potential;
204 /* Return the potential. */
205 *awhPotential = potential;
207 /* Check the sampled histograms and potentially warn user if something is suspicious */
208 warnForHistogramAnomalies(t, step, fplog);
214 * Count the total number of samples / sample weight over all grid points.
216 * \param[in] pointState The state of the points in a bias.
217 * \returns the total sample count.
219 static int64_t countSamples(const std::vector<PointState> &pointState)
221 double numSamples = 0;
222 for (const PointState &point : pointState)
224 numSamples += point.weightSumTot();
227 return gmx::roundToInt64(numSamples);
231 * Check if the state (loaded from checkpoint) and the run are consistent.
233 * When the state and the run setup are inconsistent, an exception is thrown.
235 * \param[in] params The parameters of the bias.
236 * \param[in] state The state of the bias.
238 static void ensureStateAndRunConsistency(const BiasParams ¶ms,
239 const BiasState &state)
241 int64_t numSamples = countSamples(state.points());
242 int64_t numUpdatesFromSamples = numSamples/(params.numSamplesUpdateFreeEnergy_*params.numSharedUpdate);
243 int64_t numUpdatesExpected = state.histogramSize().numUpdates();
244 if (numUpdatesFromSamples != numUpdatesExpected)
246 std::string mesg = gmx::formatString("The number of AWH updates in the checkpoint file (%" PRId64 ") does not match the total number of AWH samples divided by the number of samples per update for %d sharing AWH bias(es) (%" PRId64 "/%d=%" PRId64 ")",
248 params.numSharedUpdate,
250 params.numSamplesUpdateFreeEnergy_*params.numSharedUpdate,
251 numUpdatesFromSamples);
252 mesg += " Maybe you changed AWH parameters.";
253 /* Unfortunately we currently do not store the number of simulations
254 * sharing the bias or the state to checkpoint. But we can hint at
255 * a mismatch in the number of sharing simulations.
257 if (numUpdatesFromSamples % state.histogramSize().numUpdates() == 0)
259 mesg += gmx::formatString(" Or the run you continued from used %" PRId64 " sharing simulations, whereas you now specified %d sharing simulations.",
260 numUpdatesFromSamples/state.histogramSize().numUpdates(),
261 params.numSharedUpdate);
263 GMX_THROW(InvalidInputError(mesg));
267 void Bias::restoreStateFromHistory(const AwhBiasHistory *biasHistory,
270 GMX_RELEASE_ASSERT(thisRankDoesIO_ == MASTER(cr), "The master rank should do I/O, the other ranks should not");
274 GMX_RELEASE_ASSERT(biasHistory != nullptr, "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),
331 params_(awhParams, awhBiasParams, dimParams_, beta, mdTimeStep, disableUpdateSkips, numSharingSimulations, grid_.axis(), biasIndexInCollection),
332 state_(awhBiasParams, params_.initialHistogramSize, dimParams_, grid_),
333 thisRankDoesIO_(thisRankWillDoIO == ThisRankWillDoIO::Yes),
337 numWarningsIssued_(0)
339 /* For a global update updateList covers all points, so reserve that */
340 updateList_.reserve(grid_.numPoints());
342 state_.initGridPointState(awhBiasParams, dimParams_, grid_, params_, biasInitFilename, awhParams.numBias);
346 /* Set up the force correlation object. */
348 /* We let the correlation init function set its parameters
349 * to something useful for now.
351 double blockLength = 0;
352 /* Construct the force correlation object. */
353 forceCorrelationGrid_ =
354 compat::make_unique<CorrelationGrid>(state_.points().size(), ndim(),
355 blockLength, CorrelationGrid::BlockLengthMeasure::Time,
356 awhParams.nstSampleCoord*mdTimeStep);
358 writer_ = compat::make_unique<BiasWriter>(*this);
362 void Bias::printInitializationToLog(FILE *fplog) const
364 if (fplog != nullptr && forceCorrelationGrid_ != nullptr)
367 gmx::formatString("\nawh%d:", params_.biasIndex + 1);
370 "%s initial force correlation block length = %g %s"
371 "%s force correlation number of blocks = %d",
372 prefix.c_str(), forceCorrelationGrid().getBlockLength(),
373 forceCorrelationGrid().blockLengthMeasure == CorrelationGrid::BlockLengthMeasure::Weight ? "" : "ps",
374 prefix.c_str(), forceCorrelationGrid().getNumBlocks());
378 void Bias::updateForceCorrelationGrid(gmx::ArrayRef<const double> probWeightNeighbor,
381 if (forceCorrelationGrid_ == nullptr)
386 const std::vector<int> &neighbor = grid_.point(state_.coordState().gridpointIndex()).neighbor;
388 gmx::ArrayRef<double> forceFromNeighbor = tempForce_;
389 for (size_t n = 0; n < neighbor.size(); n++)
391 double weightNeighbor = probWeightNeighbor[n];
392 int indexNeighbor = neighbor[n];
394 /* Add the force data of this neighbor point. Note: the sum of these forces is the convolved force.
396 We actually add the force normalized by beta which has the units of 1/length. This means that the
397 resulting correlation time integral is directly in units of friction time/length^2 which is really what
398 we're interested in. */
399 state_.calcUmbrellaForceAndPotential(dimParams_, grid_, indexNeighbor, forceFromNeighbor);
401 /* Note: we might want to give a whole list of data to add instead and have this loop in the data adding function */
402 forceCorrelationGrid_->addData(indexNeighbor, weightNeighbor, forceFromNeighbor, t);
406 /* Return the number of data blocks that have been prepared for writing. */
407 int Bias::numEnergySubblocksToWrite() const
409 GMX_RELEASE_ASSERT(writer_ != nullptr, "Should only request data from an initialized writer");
411 return writer_->numBlocks();
414 /* Write bias data blocks to energy subblocks. */
415 int Bias::writeToEnergySubblocks(t_enxsubblock *subblock) const
417 GMX_RELEASE_ASSERT(writer_ != nullptr, "Should only request data from an initialized writer");
419 return writer_->writeToEnergySubblocks(*this, subblock);