2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2015,2016,2017,2018,2019, The GROMACS development team.
5 * Copyright (c) 2020,2021, by the GROMACS development team, led by
6 * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7 * and including many others, as listed in the AUTHORS file in the
8 * top-level source directory and at http://www.gromacs.org.
10 * GROMACS is free software; you can redistribute it and/or
11 * modify it under the terms of the GNU Lesser General Public License
12 * as published by the Free Software Foundation; either version 2.1
13 * of the License, or (at your option) any later version.
15 * GROMACS is distributed in the hope that it will be useful,
16 * but WITHOUT ANY WARRANTY; without even the implied warranty of
17 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
18 * Lesser General Public License for more details.
20 * You should have received a copy of the GNU Lesser General Public
21 * License along with GROMACS; if not, see
22 * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23 * Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301 USA.
25 * If you want to redistribute modifications to GROMACS, please
26 * consider that scientific software is very special. Version
27 * control is crucial - bugs must be traceable. We will be happy to
28 * consider code for inclusion in the official distribution, but
29 * derived work must not be called official GROMACS. Details are found
30 * in the README & COPYING files - if they are missing, get the
31 * official version at http://www.gromacs.org.
33 * To help us fund GROMACS development, we humbly ask that you cite
34 * the research papers on the package. Check out http://www.gromacs.org.
39 * Implements the BiasState class.
41 * \author Viveca Lindahl
42 * \author Berk Hess <hess@kth.se>
48 #include "biasstate.h"
59 #include "gromacs/fileio/gmxfio.h"
60 #include "gromacs/fileio/xvgr.h"
61 #include "gromacs/gmxlib/network.h"
62 #include "gromacs/math/utilities.h"
63 #include "gromacs/mdrunutility/multisim.h"
64 #include "gromacs/mdtypes/awh_history.h"
65 #include "gromacs/mdtypes/awh_params.h"
66 #include "gromacs/mdtypes/commrec.h"
67 #include "gromacs/simd/simd.h"
68 #include "gromacs/simd/simd_math.h"
69 #include "gromacs/utility/arrayref.h"
70 #include "gromacs/utility/exceptions.h"
71 #include "gromacs/utility/gmxassert.h"
72 #include "gromacs/utility/smalloc.h"
73 #include "gromacs/utility/stringutil.h"
76 #include "pointstate.h"
81 void BiasState::getPmf(gmx::ArrayRef<float> pmf) const
83 GMX_ASSERT(pmf.size() == points_.size(), "pmf should have the size of the bias grid");
85 /* The PMF is just the negative of the log of the sampled PMF histogram.
86 * Points with zero target weight are ignored, they will mostly contain noise.
88 for (size_t i = 0; i < points_.size(); i++)
90 pmf[i] = points_[i].inTargetRegion() ? -points_[i].logPmfSum() : GMX_FLOAT_MAX;
98 * Sum an array over all simulations on the master rank of each simulation.
100 * \param[in,out] arrayRef The data to sum.
101 * \param[in] multiSimComm Struct for multi-simulation communication.
103 void sumOverSimulations(gmx::ArrayRef<int> arrayRef, const gmx_multisim_t* multiSimComm)
105 gmx_sumi_sim(arrayRef.size(), arrayRef.data(), multiSimComm);
109 * Sum an array over all simulations on the master rank of each simulation.
111 * \param[in,out] arrayRef The data to sum.
112 * \param[in] multiSimComm Struct for multi-simulation communication.
114 void sumOverSimulations(gmx::ArrayRef<double> arrayRef, const gmx_multisim_t* multiSimComm)
116 gmx_sumd_sim(arrayRef.size(), arrayRef.data(), multiSimComm);
120 * Sum an array over all simulations on all ranks of each simulation.
122 * This assumes the data is identical on all ranks within each simulation.
124 * \param[in,out] arrayRef The data to sum.
125 * \param[in] commRecord Struct for intra-simulation communication.
126 * \param[in] multiSimComm Struct for multi-simulation communication.
129 void sumOverSimulations(gmx::ArrayRef<T> arrayRef, const t_commrec* commRecord, const gmx_multisim_t* multiSimComm)
131 if (MASTER(commRecord))
133 sumOverSimulations(arrayRef, multiSimComm);
135 if (commRecord->nnodes > 1)
137 gmx_bcast(arrayRef.size() * sizeof(T), arrayRef.data(), commRecord->mpi_comm_mygroup);
142 * Sum PMF over multiple simulations, when requested.
144 * \param[in,out] pointState The state of the points in the bias.
145 * \param[in] numSharedUpdate The number of biases sharing the histogram.
146 * \param[in] commRecord Struct for intra-simulation communication.
147 * \param[in] multiSimComm Struct for multi-simulation communication.
149 void sumPmf(gmx::ArrayRef<PointState> pointState,
151 const t_commrec* commRecord,
152 const gmx_multisim_t* multiSimComm)
154 if (numSharedUpdate == 1)
158 GMX_ASSERT(multiSimComm != nullptr && numSharedUpdate % multiSimComm->numSimulations_ == 0,
159 "numSharedUpdate should be a multiple of multiSimComm->numSimulations_");
160 GMX_ASSERT(numSharedUpdate == multiSimComm->numSimulations_,
161 "Sharing within a simulation is not implemented (yet)");
163 std::vector<double> buffer(pointState.size());
165 /* Need to temporarily exponentiate the log weights to sum over simulations */
166 for (size_t i = 0; i < buffer.size(); i++)
168 buffer[i] = pointState[i].inTargetRegion() ? std::exp(pointState[i].logPmfSum()) : 0;
171 sumOverSimulations(gmx::ArrayRef<double>(buffer), commRecord, multiSimComm);
173 /* Take log again to get (non-normalized) PMF */
174 double normFac = 1.0 / numSharedUpdate;
175 for (gmx::index i = 0; i < pointState.ssize(); i++)
177 if (pointState[i].inTargetRegion())
179 pointState[i].setLogPmfSum(std::log(buffer[i] * normFac));
185 * Find the minimum free energy value.
187 * \param[in] pointState The state of the points.
188 * \returns the minimum free energy value.
190 double freeEnergyMinimumValue(gmx::ArrayRef<const PointState> pointState)
192 double fMin = GMX_FLOAT_MAX;
194 for (auto const& ps : pointState)
196 if (ps.inTargetRegion() && ps.freeEnergy() < fMin)
198 fMin = ps.freeEnergy();
206 * Find and return the log of the probability weight of a point given a coordinate value.
208 * The unnormalized weight is given by
209 * w(point|value) = exp(bias(point) - U(value,point)),
210 * where U is a harmonic umbrella potential.
212 * \param[in] dimParams The bias dimensions parameters
213 * \param[in] points The point state.
214 * \param[in] grid The grid.
215 * \param[in] pointIndex Point to evaluate probability weight for.
216 * \param[in] pointBias Bias for the point (as a log weight).
217 * \param[in] value Coordinate value.
218 * \param[in] neighborLambdaEnergies The energy of the system in neighboring lambdas states. Can be
219 * empty when there are no free energy lambda state dimensions.
220 * \param[in] gridpointIndex The index of the current grid point.
221 * \returns the log of the biased probability weight.
223 double biasedLogWeightFromPoint(const std::vector<DimParams>& dimParams,
224 const std::vector<PointState>& points,
225 const BiasGrid& grid,
228 const awh_dvec value,
229 gmx::ArrayRef<const double> neighborLambdaEnergies,
232 double logWeight = detail::c_largeNegativeExponent;
234 /* Only points in the target region have non-zero weight */
235 if (points[pointIndex].inTargetRegion())
237 logWeight = pointBias;
239 /* Add potential for all parameter dimensions */
240 for (size_t d = 0; d < dimParams.size(); d++)
242 if (dimParams[d].isFepLambdaDimension())
244 /* If this is not a sampling step or if this function is called from
245 * calcConvolvedBias(), when writing energy subblocks, neighborLambdaEnergies will
246 * be empty. No convolution is required along the lambda dimension. */
247 if (!neighborLambdaEnergies.empty())
249 const int pointLambdaIndex = grid.point(pointIndex).coordValue[d];
250 const int gridpointLambdaIndex = grid.point(gridpointIndex).coordValue[d];
251 logWeight -= dimParams[d].fepDimParams().beta
252 * (neighborLambdaEnergies[pointLambdaIndex]
253 - neighborLambdaEnergies[gridpointLambdaIndex]);
258 double dev = getDeviationFromPointAlongGridAxis(grid, d, pointIndex, value[d]);
259 logWeight -= 0.5 * dimParams[d].pullDimParams().betak * dev * dev;
267 * Calculates the marginal distribution (marginal probability) for each value along
268 * a free energy lambda axis.
269 * The marginal distribution of one coordinate dimension value is the sum of the probability
270 * distribution of all values (herein all neighbor values) with the same value in the dimension
272 * \param[in] grid The bias grid.
273 * \param[in] neighbors The points to use for the calculation of the marginal distribution.
274 * \param[in] probWeightNeighbor Probability weights of the neighbors.
275 * \returns The calculated marginal distribution in a 1D array with
276 * as many elements as there are points along the axis of interest.
278 std::vector<double> calculateFELambdaMarginalDistribution(const BiasGrid& grid,
279 gmx::ArrayRef<const int> neighbors,
280 gmx::ArrayRef<const double> probWeightNeighbor)
282 const std::optional<int> lambdaAxisIndex = grid.lambdaAxisIndex();
283 GMX_RELEASE_ASSERT(lambdaAxisIndex,
284 "There must be a free energy lambda axis in order to calculate the free "
285 "energy lambda marginal distribution.");
286 const int numFepLambdaStates = grid.numFepLambdaStates();
287 std::vector<double> lambdaMarginalDistribution(numFepLambdaStates, 0);
289 for (size_t i = 0; i < neighbors.size(); i++)
291 const int neighbor = neighbors[i];
292 const int lambdaState = grid.point(neighbor).coordValue[lambdaAxisIndex.value()];
293 lambdaMarginalDistribution[lambdaState] += probWeightNeighbor[i];
295 return lambdaMarginalDistribution;
300 void BiasState::calcConvolvedPmf(const std::vector<DimParams>& dimParams,
301 const BiasGrid& grid,
302 std::vector<float>* convolvedPmf) const
304 size_t numPoints = grid.numPoints();
306 convolvedPmf->resize(numPoints);
308 /* Get the PMF to convolve. */
309 std::vector<float> pmf(numPoints);
312 for (size_t m = 0; m < numPoints; m++)
314 double freeEnergyWeights = 0;
315 const GridPoint& point = grid.point(m);
316 for (auto& neighbor : point.neighbor)
318 /* Do not convolve the bias along a lambda axis - only use the pmf from the current point */
319 if (!pointsHaveDifferentLambda(grid, m, neighbor))
321 /* The negative PMF is a positive bias. */
322 double biasNeighbor = -pmf[neighbor];
324 /* Add the convolved PMF weights for the neighbors of this point.
325 Note that this function only adds point within the target > 0 region.
326 Sum weights, take the logarithm last to get the free energy. */
327 double logWeight = biasedLogWeightFromPoint(dimParams, points_, grid, neighbor,
328 biasNeighbor, point.coordValue, {}, m);
329 freeEnergyWeights += std::exp(logWeight);
333 GMX_RELEASE_ASSERT(freeEnergyWeights > 0,
334 "Attempting to do log(<= 0) in AWH convolved PMF calculation.");
335 (*convolvedPmf)[m] = -std::log(static_cast<float>(freeEnergyWeights));
343 * Updates the target distribution for all points.
345 * The target distribution is always updated for all points
348 * \param[in,out] pointState The state of all points.
349 * \param[in] params The bias parameters.
351 void updateTargetDistribution(gmx::ArrayRef<PointState> pointState, const BiasParams& params)
353 double freeEnergyCutoff = 0;
354 if (params.eTarget == eawhtargetCUTOFF)
356 freeEnergyCutoff = freeEnergyMinimumValue(pointState) + params.freeEnergyCutoffInKT;
359 double sumTarget = 0;
360 for (PointState& ps : pointState)
362 sumTarget += ps.updateTargetWeight(params, freeEnergyCutoff);
364 GMX_RELEASE_ASSERT(sumTarget > 0, "We should have a non-zero distribution");
367 double invSum = 1.0 / sumTarget;
368 for (PointState& ps : pointState)
370 ps.scaleTarget(invSum);
375 * Puts together a string describing a grid point.
377 * \param[in] grid The grid.
378 * \param[in] point BiasGrid point index.
379 * \returns a string for the point.
381 std::string gridPointValueString(const BiasGrid& grid, int point)
383 std::string pointString;
387 for (int d = 0; d < grid.numDimensions(); d++)
389 pointString += gmx::formatString("%g", grid.point(point).coordValue[d]);
390 if (d < grid.numDimensions() - 1)
405 int BiasState::warnForHistogramAnomalies(const BiasGrid& grid, int biasIndex, double t, FILE* fplog, int maxNumWarnings) const
407 GMX_ASSERT(fplog != nullptr, "Warnings can only be issued if there is log file.");
408 const double maxHistogramRatio = 0.5; /* Tolerance for printing a warning about the histogram ratios */
410 /* Sum up the histograms and get their normalization */
411 double sumVisits = 0;
412 double sumWeights = 0;
413 for (auto& pointState : points_)
415 if (pointState.inTargetRegion())
417 sumVisits += pointState.numVisitsTot();
418 sumWeights += pointState.weightSumTot();
421 GMX_RELEASE_ASSERT(sumVisits > 0, "We should have visits");
422 GMX_RELEASE_ASSERT(sumWeights > 0, "We should have weight");
423 double invNormVisits = 1.0 / sumVisits;
424 double invNormWeight = 1.0 / sumWeights;
426 /* Check all points for warnings */
428 size_t numPoints = grid.numPoints();
429 for (size_t m = 0; m < numPoints; m++)
431 /* Skip points close to boundary or non-target region */
432 const GridPoint& gridPoint = grid.point(m);
433 bool skipPoint = false;
434 for (size_t n = 0; (n < gridPoint.neighbor.size()) && !skipPoint; n++)
436 int neighbor = gridPoint.neighbor[n];
437 skipPoint = !points_[neighbor].inTargetRegion();
438 for (int d = 0; (d < grid.numDimensions()) && !skipPoint; d++)
440 const GridPoint& neighborPoint = grid.point(neighbor);
441 skipPoint = neighborPoint.index[d] == 0
442 || neighborPoint.index[d] == grid.axis(d).numPoints() - 1;
446 /* Warn if the coordinate distribution is less than the target distribution with a certain fraction somewhere */
447 const double relativeWeight = points_[m].weightSumTot() * invNormWeight;
448 const double relativeVisits = points_[m].numVisitsTot() * invNormVisits;
449 if (!skipPoint && relativeVisits < relativeWeight * maxHistogramRatio)
451 std::string pointValueString = gridPointValueString(grid, m);
452 std::string warningMessage = gmx::formatString(
454 "at t = %g ps the obtained coordinate distribution at coordinate value %s "
455 "is less than a fraction %g of the reference distribution at that point. "
456 "If you are not certain about your settings you might want to increase your "
457 "pull force constant or "
458 "modify your sampling region.\n",
459 biasIndex + 1, t, pointValueString.c_str(), maxHistogramRatio);
460 gmx::TextLineWrapper wrapper;
461 wrapper.settings().setLineLength(c_linewidth);
462 fprintf(fplog, "%s", wrapper.wrapToString(warningMessage).c_str());
466 if (numWarnings >= maxNumWarnings)
475 double BiasState::calcUmbrellaForceAndPotential(const std::vector<DimParams>& dimParams,
476 const BiasGrid& grid,
478 ArrayRef<const double> neighborLambdaDhdl,
479 gmx::ArrayRef<double> force) const
481 double potential = 0;
482 for (size_t d = 0; d < dimParams.size(); d++)
484 if (dimParams[d].isFepLambdaDimension())
486 if (!neighborLambdaDhdl.empty())
488 const int coordpointLambdaIndex = grid.point(point).coordValue[d];
489 force[d] = neighborLambdaDhdl[coordpointLambdaIndex];
490 /* The potential should not be affected by the lambda dimension. */
496 getDeviationFromPointAlongGridAxis(grid, d, point, coordState_.coordValue()[d]);
497 double k = dimParams[d].pullDimParams().k;
499 /* Force from harmonic potential 0.5*k*dev^2 */
500 force[d] = -k * deviation;
501 potential += 0.5 * k * deviation * deviation;
508 void BiasState::calcConvolvedForce(const std::vector<DimParams>& dimParams,
509 const BiasGrid& grid,
510 gmx::ArrayRef<const double> probWeightNeighbor,
511 ArrayRef<const double> neighborLambdaDhdl,
512 gmx::ArrayRef<double> forceWorkBuffer,
513 gmx::ArrayRef<double> force) const
515 for (size_t d = 0; d < dimParams.size(); d++)
520 /* Only neighboring points have non-negligible contribution. */
521 const std::vector<int>& neighbor = grid.point(coordState_.gridpointIndex()).neighbor;
522 gmx::ArrayRef<double> forceFromNeighbor = forceWorkBuffer;
523 for (size_t n = 0; n < neighbor.size(); n++)
525 double weightNeighbor = probWeightNeighbor[n];
526 int indexNeighbor = neighbor[n];
528 /* Get the umbrella force from this point. The returned potential is ignored here. */
529 calcUmbrellaForceAndPotential(dimParams, grid, indexNeighbor, neighborLambdaDhdl, forceFromNeighbor);
531 /* Add the weighted umbrella force to the convolved force. */
532 for (size_t d = 0; d < dimParams.size(); d++)
534 force[d] += forceFromNeighbor[d] * weightNeighbor;
539 double BiasState::moveUmbrella(const std::vector<DimParams>& dimParams,
540 const BiasGrid& grid,
541 gmx::ArrayRef<const double> probWeightNeighbor,
542 ArrayRef<const double> neighborLambdaDhdl,
543 gmx::ArrayRef<double> biasForce,
548 /* Generate and set a new coordinate reference value */
549 coordState_.sampleUmbrellaGridpoint(grid, coordState_.gridpointIndex(), probWeightNeighbor,
550 step, seed, indexSeed);
552 std::vector<double> newForce(dimParams.size());
553 double newPotential = calcUmbrellaForceAndPotential(
554 dimParams, grid, coordState_.umbrellaGridpoint(), neighborLambdaDhdl, newForce);
556 /* A modification of the reference value at time t will lead to a different
557 force over t-dt/2 to t and over t to t+dt/2. For high switching rates
558 this means the force and velocity will change signs roughly as often.
559 To avoid any issues we take the average of the previous and new force
560 at steps when the reference value has been moved. E.g. if the ref. value
561 is set every step to (coord dvalue +/- delta) would give zero force.
563 for (gmx::index d = 0; d < biasForce.ssize(); d++)
565 /* Average of the current and new force */
566 biasForce[d] = 0.5 * (biasForce[d] + newForce[d]);
576 * Sets the histogram rescaling factors needed to control the histogram size.
578 * For sake of robustness, the reference weight histogram can grow at a rate
579 * different from the actual sampling rate. Typically this happens for a limited
580 * initial time, alternatively growth is scaled down by a constant factor for all
581 * times. Since the size of the reference histogram sets the size of the free
582 * energy update this should be reflected also in the PMF. Thus the PMF histogram
583 * needs to be rescaled too.
585 * This function should only be called by the bias update function or wrapped by a function that
586 * knows what scale factors should be applied when, e.g,
587 * getSkippedUpdateHistogramScaleFactors().
589 * \param[in] params The bias parameters.
590 * \param[in] newHistogramSize New reference weight histogram size.
591 * \param[in] oldHistogramSize Previous reference weight histogram size (before adding new samples).
592 * \param[out] weightHistScaling Scaling factor for the reference weight histogram.
593 * \param[out] logPmfSumScaling Log of the scaling factor for the PMF histogram.
595 void setHistogramUpdateScaleFactors(const BiasParams& params,
596 double newHistogramSize,
597 double oldHistogramSize,
598 double* weightHistScaling,
599 double* logPmfSumScaling)
602 /* The two scaling factors below are slightly different (ignoring the log factor) because the
603 reference and the PMF histogram apply weight scaling differently. The weight histogram
604 applies is locally, i.e. each sample is scaled down meaning all samples get equal weight.
605 It is done this way because that is what target type local Boltzmann (for which
606 target = weight histogram) needs. In contrast, the PMF histogram is rescaled globally
607 by repeatedly scaling down the whole histogram. The reasons for doing it this way are:
608 1) empirically this is necessary for converging the PMF; 2) since the extraction of
609 the PMF is theoretically only valid for a constant bias, new samples should get more
610 weight than old ones for which the bias is fluctuating more. */
612 newHistogramSize / (oldHistogramSize + params.updateWeight * params.localWeightScaling);
613 *logPmfSumScaling = std::log(newHistogramSize / (oldHistogramSize + params.updateWeight));
618 void BiasState::getSkippedUpdateHistogramScaleFactors(const BiasParams& params,
619 double* weightHistScaling,
620 double* logPmfSumScaling) const
622 GMX_ASSERT(params.skipUpdates(),
623 "Calling function for skipped updates when skipping updates is not allowed");
625 if (inInitialStage())
627 /* In between global updates the reference histogram size is kept constant so we trivially
628 know what the histogram size was at the time of the skipped update. */
629 double histogramSize = histogramSize_.histogramSize();
630 setHistogramUpdateScaleFactors(params, histogramSize, histogramSize, weightHistScaling,
635 /* In the final stage, the reference histogram grows at the sampling rate which gives trivial scale factors. */
636 *weightHistScaling = 1;
637 *logPmfSumScaling = 0;
641 void BiasState::doSkippedUpdatesForAllPoints(const BiasParams& params)
643 double weightHistScaling;
644 double logPmfsumScaling;
646 getSkippedUpdateHistogramScaleFactors(params, &weightHistScaling, &logPmfsumScaling);
648 for (auto& pointState : points_)
650 bool didUpdate = pointState.performPreviouslySkippedUpdates(
651 params, histogramSize_.numUpdates(), weightHistScaling, logPmfsumScaling);
653 /* Update the bias for this point only if there were skipped updates in the past to avoid calculating the log unneccessarily */
656 pointState.updateBias();
661 void BiasState::doSkippedUpdatesInNeighborhood(const BiasParams& params, const BiasGrid& grid)
663 double weightHistScaling;
664 double logPmfsumScaling;
666 getSkippedUpdateHistogramScaleFactors(params, &weightHistScaling, &logPmfsumScaling);
668 /* For each neighbor point of the center point, refresh its state by adding the results of all past, skipped updates. */
669 const std::vector<int>& neighbors = grid.point(coordState_.gridpointIndex()).neighbor;
670 for (auto& neighbor : neighbors)
672 bool didUpdate = points_[neighbor].performPreviouslySkippedUpdates(
673 params, histogramSize_.numUpdates(), weightHistScaling, logPmfsumScaling);
677 points_[neighbor].updateBias();
686 * Merge update lists from multiple sharing simulations.
688 * \param[in,out] updateList Update list for this simulation (assumed >= npoints long).
689 * \param[in] numPoints Total number of points.
690 * \param[in] commRecord Struct for intra-simulation communication.
691 * \param[in] multiSimComm Struct for multi-simulation communication.
693 void mergeSharedUpdateLists(std::vector<int>* updateList,
695 const t_commrec* commRecord,
696 const gmx_multisim_t* multiSimComm)
698 std::vector<int> numUpdatesOfPoint;
700 /* Flag the update points of this sim.
701 TODO: we can probably avoid allocating this array and just use the input array. */
702 numUpdatesOfPoint.resize(numPoints, 0);
703 for (auto& pointIndex : *updateList)
705 numUpdatesOfPoint[pointIndex] = 1;
708 /* Sum over the sims to get all the flagged points */
709 sumOverSimulations(arrayRefFromArray(numUpdatesOfPoint.data(), numPoints), commRecord, multiSimComm);
711 /* Collect the indices of the flagged points in place. The resulting array will be the merged update list.*/
713 for (int m = 0; m < numPoints; m++)
715 if (numUpdatesOfPoint[m] > 0)
717 updateList->push_back(m);
723 * Generate an update list of points sampled since the last update.
725 * \param[in] grid The AWH bias.
726 * \param[in] points The point state.
727 * \param[in] originUpdatelist The origin of the rectangular region that has been sampled since
728 * last update. \param[in] endUpdatelist The end of the rectangular that has been sampled since
729 * last update. \param[in,out] updateList Local update list to set (assumed >= npoints long).
731 void makeLocalUpdateList(const BiasGrid& grid,
732 const std::vector<PointState>& points,
733 const awh_ivec originUpdatelist,
734 const awh_ivec endUpdatelist,
735 std::vector<int>* updateList)
740 /* Define the update search grid */
741 for (int d = 0; d < grid.numDimensions(); d++)
743 origin[d] = originUpdatelist[d];
744 numPoints[d] = endUpdatelist[d] - originUpdatelist[d] + 1;
746 /* Because the end_updatelist is unwrapped it can be > (npoints - 1) so that numPoints can be > npoints in grid.
747 This helps for calculating the distance/number of points but should be removed and fixed when the way of
748 updating origin/end updatelist is changed (see sampleProbabilityWeights). */
749 numPoints[d] = std::min(grid.axis(d).numPoints(), numPoints[d]);
752 /* Make the update list */
755 bool pointExists = true;
758 pointExists = advancePointInSubgrid(grid, origin, numPoints, &pointIndex);
760 if (pointExists && points[pointIndex].inTargetRegion())
762 updateList->push_back(pointIndex);
769 void BiasState::resetLocalUpdateRange(const BiasGrid& grid)
771 const int gridpointIndex = coordState_.gridpointIndex();
772 for (int d = 0; d < grid.numDimensions(); d++)
774 /* This gives the minimum range consisting only of the current closest point. */
775 originUpdatelist_[d] = grid.point(gridpointIndex).index[d];
776 endUpdatelist_[d] = grid.point(gridpointIndex).index[d];
784 * Add partial histograms (accumulating between updates) to accumulating histograms.
786 * \param[in,out] pointState The state of the points in the bias.
787 * \param[in,out] weightSumCovering The weights for checking covering.
788 * \param[in] numSharedUpdate The number of biases sharing the histrogram.
789 * \param[in] commRecord Struct for intra-simulation communication.
790 * \param[in] multiSimComm Struct for multi-simulation communication.
791 * \param[in] localUpdateList List of points with data.
793 void sumHistograms(gmx::ArrayRef<PointState> pointState,
794 gmx::ArrayRef<double> weightSumCovering,
796 const t_commrec* commRecord,
797 const gmx_multisim_t* multiSimComm,
798 const std::vector<int>& localUpdateList)
800 /* The covering checking histograms are added before summing over simulations, so that the
801 weights from different simulations are kept distinguishable. */
802 for (int globalIndex : localUpdateList)
804 weightSumCovering[globalIndex] += pointState[globalIndex].weightSumIteration();
807 /* Sum histograms over multiple simulations if needed. */
808 if (numSharedUpdate > 1)
810 GMX_ASSERT(numSharedUpdate == multiSimComm->numSimulations_,
811 "Sharing within a simulation is not implemented (yet)");
813 /* Collect the weights and counts in linear arrays to be able to use gmx_sumd_sim. */
814 std::vector<double> weightSum;
815 std::vector<double> coordVisits;
817 weightSum.resize(localUpdateList.size());
818 coordVisits.resize(localUpdateList.size());
820 for (size_t localIndex = 0; localIndex < localUpdateList.size(); localIndex++)
822 const PointState& ps = pointState[localUpdateList[localIndex]];
824 weightSum[localIndex] = ps.weightSumIteration();
825 coordVisits[localIndex] = ps.numVisitsIteration();
828 sumOverSimulations(gmx::ArrayRef<double>(weightSum), commRecord, multiSimComm);
829 sumOverSimulations(gmx::ArrayRef<double>(coordVisits), commRecord, multiSimComm);
831 /* Transfer back the result */
832 for (size_t localIndex = 0; localIndex < localUpdateList.size(); localIndex++)
834 PointState& ps = pointState[localUpdateList[localIndex]];
836 ps.setPartialWeightAndCount(weightSum[localIndex], coordVisits[localIndex]);
840 /* Now add the partial counts and weights to the accumulating histograms.
841 Note: we still need to use the weights for the update so we wait
842 with resetting them until the end of the update. */
843 for (int globalIndex : localUpdateList)
845 pointState[globalIndex].addPartialWeightAndCount();
850 * Label points along an axis as covered or not.
852 * A point is covered if it is surrounded by visited points up to a radius = coverRadius.
854 * \param[in] visited Visited? For each point.
855 * \param[in] checkCovering Check for covering? For each point.
856 * \param[in] numPoints The number of grid points along this dimension.
857 * \param[in] period Period in number of points.
858 * \param[in] coverRadius Cover radius, in points, needed for defining a point as covered.
859 * \param[in,out] covered In this array elements are 1 for covered points and 0 for
860 * non-covered points, this routine assumes that \p covered has at least size \p numPoints.
862 void labelCoveredPoints(const std::vector<bool>& visited,
863 const std::vector<bool>& checkCovering,
867 gmx::ArrayRef<int> covered)
869 GMX_ASSERT(covered.ssize() >= numPoints, "covered should be at least as large as the grid");
871 bool haveFirstNotVisited = false;
872 int firstNotVisited = -1;
873 int notVisitedLow = -1;
874 int notVisitedHigh = -1;
876 for (int n = 0; n < numPoints; n++)
878 if (checkCovering[n] && !visited[n])
880 if (!haveFirstNotVisited)
884 haveFirstNotVisited = true;
890 /* Have now an interval I = [notVisitedLow,notVisitedHigh] of visited points bounded
891 by unvisited points. The unvisted end points affect the coveredness of the
892 visited with a reach equal to the cover radius. */
893 int notCoveredLow = notVisitedLow + coverRadius;
894 int notCoveredHigh = notVisitedHigh - coverRadius;
895 for (int i = notVisitedLow; i <= notVisitedHigh; i++)
897 covered[i] = static_cast<int>((i > notCoveredLow) && (i < notCoveredHigh));
900 /* Find a new interval to set covering for. Make the notVisitedHigh of this interval
901 the notVisitedLow of the next. */
902 notVisitedLow = notVisitedHigh;
907 /* Have labelled all the internal points. Now take care of the boundary regions. */
908 if (!haveFirstNotVisited)
910 /* No non-visited points <=> all points visited => all points covered. */
912 for (int n = 0; n < numPoints; n++)
919 int lastNotVisited = notVisitedLow;
921 /* For periodic boundaries, non-visited points can influence points
922 on the other side of the boundary so we need to wrap around. */
924 /* Lower end. For periodic boundaries the last upper end not visited point becomes the low-end not visited point.
925 For non-periodic boundaries there is no lower end point so a dummy value is used. */
926 int notVisitedHigh = firstNotVisited;
927 int notVisitedLow = period > 0 ? (lastNotVisited - period) : -(coverRadius + 1);
929 int notCoveredLow = notVisitedLow + coverRadius;
930 int notCoveredHigh = notVisitedHigh - coverRadius;
932 for (int i = 0; i <= notVisitedHigh; i++)
934 /* For non-periodic boundaries notCoveredLow = -1 will impose no restriction. */
935 covered[i] = static_cast<int>((i > notCoveredLow) && (i < notCoveredHigh));
938 /* Upper end. Same as for lower end but in the other direction. */
939 notVisitedHigh = period > 0 ? (firstNotVisited + period) : (numPoints + coverRadius);
940 notVisitedLow = lastNotVisited;
942 notCoveredLow = notVisitedLow + coverRadius;
943 notCoveredHigh = notVisitedHigh - coverRadius;
945 for (int i = notVisitedLow; i <= numPoints - 1; i++)
947 /* For non-periodic boundaries notCoveredHigh = numPoints will impose no restriction. */
948 covered[i] = static_cast<int>((i > notCoveredLow) && (i < notCoveredHigh));
955 bool BiasState::isSamplingRegionCovered(const BiasParams& params,
956 const std::vector<DimParams>& dimParams,
957 const BiasGrid& grid,
958 const t_commrec* commRecord,
959 const gmx_multisim_t* multiSimComm) const
961 /* Allocate and initialize arrays: one for checking visits along each dimension,
962 one for keeping track of which points to check and one for the covered points.
963 Possibly these could be kept as AWH variables to avoid these allocations. */
966 std::vector<bool> visited;
967 std::vector<bool> checkCovering;
968 // We use int for the covering array since we might use gmx_sumi_sim.
969 std::vector<int> covered;
972 std::vector<CheckDim> checkDim;
973 checkDim.resize(grid.numDimensions());
975 for (int d = 0; d < grid.numDimensions(); d++)
977 const size_t numPoints = grid.axis(d).numPoints();
978 checkDim[d].visited.resize(numPoints, false);
979 checkDim[d].checkCovering.resize(numPoints, false);
980 checkDim[d].covered.resize(numPoints, 0);
983 /* Set visited points along each dimension and which points should be checked for covering.
984 Specifically, points above the free energy cutoff (if there is one) or points outside
985 of the target region are ignored. */
987 /* Set the free energy cutoff */
988 double maxFreeEnergy = GMX_FLOAT_MAX;
990 if (params.eTarget == eawhtargetCUTOFF)
992 maxFreeEnergy = freeEnergyMinimumValue(points_) + params.freeEnergyCutoffInKT;
995 /* Set the threshold weight for a point to be considered visited. */
996 double weightThreshold = 1;
997 for (int d = 0; d < grid.numDimensions(); d++)
999 if (grid.axis(d).isFepLambdaAxis())
1001 /* Do not modify the weight threshold based on a FEP lambda axis. The spread
1002 * of the sampling weights is not depending on a Gaussian distribution (like
1004 weightThreshold *= 1.0;
1008 /* The spacing is proportional to 1/sqrt(betak). The weight threshold will be
1009 * approximately (given that the spacing can be modified if the dimension is periodic)
1010 * proportional to sqrt(1/(2*pi)). */
1011 weightThreshold *= grid.axis(d).spacing()
1012 * std::sqrt(dimParams[d].pullDimParams().betak * 0.5 * M_1_PI);
1016 /* Project the sampling weights onto each dimension */
1017 for (size_t m = 0; m < grid.numPoints(); m++)
1019 const PointState& pointState = points_[m];
1021 for (int d = 0; d < grid.numDimensions(); d++)
1023 int n = grid.point(m).index[d];
1025 /* Is visited if it was already visited or if there is enough weight at the current point */
1026 checkDim[d].visited[n] = checkDim[d].visited[n] || (weightSumCovering_[m] > weightThreshold);
1028 /* Check for covering if there is at least point in this slice that is in the target region and within the cutoff */
1029 checkDim[d].checkCovering[n] =
1030 checkDim[d].checkCovering[n]
1031 || (pointState.inTargetRegion() && pointState.freeEnergy() < maxFreeEnergy);
1035 /* Label each point along each dimension as covered or not. */
1036 for (int d = 0; d < grid.numDimensions(); d++)
1038 labelCoveredPoints(checkDim[d].visited, checkDim[d].checkCovering, grid.axis(d).numPoints(),
1039 grid.axis(d).numPointsInPeriod(), params.coverRadius()[d], checkDim[d].covered);
1042 /* Now check for global covering. Each dimension needs to be covered separately.
1043 A dimension is covered if each point is covered. Multiple simulations collectively
1044 cover the points, i.e. a point is covered if any of the simulations covered it.
1045 However, visited points are not shared, i.e. if a point is covered or not is
1046 determined by the visits of a single simulation. In general the covering criterion is
1047 all points covered => all points are surrounded by visited points up to a radius = coverRadius.
1048 For 1 simulation, all points covered <=> all points visited. For multiple simulations
1049 however, all points visited collectively !=> all points covered, except for coverRadius = 0.
1050 In the limit of large coverRadius, all points covered => all points visited by at least one
1051 simulation (since no point will be covered until all points have been visited by a
1052 single simulation). Basically coverRadius sets how much "connectedness" (or mixing) a point
1053 needs with surrounding points before sharing covering information with other simulations. */
1055 /* Communicate the covered points between sharing simulations if needed. */
1056 if (params.numSharedUpdate > 1)
1058 /* For multiple dimensions this may not be the best way to do it. */
1059 for (int d = 0; d < grid.numDimensions(); d++)
1062 gmx::arrayRefFromArray(checkDim[d].covered.data(), grid.axis(d).numPoints()),
1063 commRecord, multiSimComm);
1067 /* Now check if for each dimension all points are covered. Break if not true. */
1068 bool allPointsCovered = true;
1069 for (int d = 0; d < grid.numDimensions() && allPointsCovered; d++)
1071 for (int n = 0; n < grid.axis(d).numPoints() && allPointsCovered; n++)
1073 allPointsCovered = (checkDim[d].covered[n] != 0);
1077 return allPointsCovered;
1081 * Normalizes the free energy and PMF sum.
1083 * \param[in] pointState The state of the points.
1085 static void normalizeFreeEnergyAndPmfSum(std::vector<PointState>* pointState)
1087 double minF = freeEnergyMinimumValue(*pointState);
1089 for (PointState& ps : *pointState)
1091 ps.normalizeFreeEnergyAndPmfSum(minF);
1095 void BiasState::updateFreeEnergyAndAddSamplesToHistogram(const std::vector<DimParams>& dimParams,
1096 const BiasGrid& grid,
1097 const BiasParams& params,
1098 const t_commrec* commRecord,
1099 const gmx_multisim_t* multiSimComm,
1103 std::vector<int>* updateList)
1105 /* Note hat updateList is only used in this scope and is always
1106 * re-initialized. We do not use a local vector, because that would
1107 * cause reallocation every time this funtion is called and the vector
1108 * can be the size of the whole grid.
1111 /* Make a list of all local points, i.e. those that could have been touched since
1112 the last update. These are the points needed for summing histograms below
1113 (non-local points only add zeros). For local updates, this will also be the
1114 final update list. */
1115 makeLocalUpdateList(grid, points_, originUpdatelist_, endUpdatelist_, updateList);
1116 if (params.numSharedUpdate > 1)
1118 mergeSharedUpdateLists(updateList, points_.size(), commRecord, multiSimComm);
1121 /* Reset the range for the next update */
1122 resetLocalUpdateRange(grid);
1124 /* Add samples to histograms for all local points and sync simulations if needed */
1125 sumHistograms(points_, weightSumCovering_, params.numSharedUpdate, commRecord, multiSimComm, *updateList);
1127 sumPmf(points_, params.numSharedUpdate, commRecord, multiSimComm);
1129 /* Renormalize the free energy if values are too large. */
1130 bool needToNormalizeFreeEnergy = false;
1131 for (int& globalIndex : *updateList)
1133 /* We want to keep the absolute value of the free energies to be less
1134 c_largePositiveExponent to be able to safely pass these values to exp(). The check below
1135 ensures this as long as the free energy values grow less than 0.5*c_largePositiveExponent
1136 in a return time to this neighborhood. For reasonable update sizes it's unlikely that
1137 this requirement would be broken. */
1138 if (std::abs(points_[globalIndex].freeEnergy()) > 0.5 * detail::c_largePositiveExponent)
1140 needToNormalizeFreeEnergy = true;
1145 /* Update target distribution? */
1146 bool needToUpdateTargetDistribution =
1147 (params.eTarget != eawhtargetCONSTANT && params.isUpdateTargetStep(step));
1149 /* In the initial stage, the histogram grows dynamically as a function of the number of coverings. */
1150 bool detectedCovering = false;
1151 if (inInitialStage())
1154 (params.isCheckCoveringStep(step)
1155 && isSamplingRegionCovered(params, dimParams, grid, commRecord, multiSimComm));
1158 /* The weighthistogram size after this update. */
1159 double newHistogramSize = histogramSize_.newHistogramSize(params, t, detectedCovering, points_,
1160 weightSumCovering_, fplog);
1162 /* Make the update list. Usually we try to only update local points,
1163 * but if the update has non-trivial or non-deterministic effects
1164 * on non-local points a global update is needed. This is the case when:
1165 * 1) a covering occurred in the initial stage, leading to non-trivial
1166 * histogram rescaling factors; or
1167 * 2) the target distribution will be updated, since we don't make any
1168 * assumption on its form; or
1169 * 3) the AWH parameters are such that we never attempt to skip non-local
1171 * 4) the free energy values have grown so large that a renormalization
1174 if (needToUpdateTargetDistribution || detectedCovering || !params.skipUpdates() || needToNormalizeFreeEnergy)
1176 /* Global update, just add all points. */
1177 updateList->clear();
1178 for (size_t m = 0; m < points_.size(); m++)
1180 if (points_[m].inTargetRegion())
1182 updateList->push_back(m);
1187 /* Set histogram scale factors. */
1188 double weightHistScalingSkipped = 0;
1189 double logPmfsumScalingSkipped = 0;
1190 if (params.skipUpdates())
1192 getSkippedUpdateHistogramScaleFactors(params, &weightHistScalingSkipped, &logPmfsumScalingSkipped);
1194 double weightHistScalingNew;
1195 double logPmfsumScalingNew;
1196 setHistogramUpdateScaleFactors(params, newHistogramSize, histogramSize_.histogramSize(),
1197 &weightHistScalingNew, &logPmfsumScalingNew);
1199 /* Update free energy and reference weight histogram for points in the update list. */
1200 for (int pointIndex : *updateList)
1202 PointState* pointStateToUpdate = &points_[pointIndex];
1204 /* Do updates from previous update steps that were skipped because this point was at that time non-local. */
1205 if (params.skipUpdates())
1207 pointStateToUpdate->performPreviouslySkippedUpdates(params, histogramSize_.numUpdates(),
1208 weightHistScalingSkipped,
1209 logPmfsumScalingSkipped);
1212 /* Now do an update with new sampling data. */
1213 pointStateToUpdate->updateWithNewSampling(params, histogramSize_.numUpdates(),
1214 weightHistScalingNew, logPmfsumScalingNew);
1217 /* Only update the histogram size after we are done with the local point updates */
1218 histogramSize_.setHistogramSize(newHistogramSize, weightHistScalingNew);
1220 if (needToNormalizeFreeEnergy)
1222 normalizeFreeEnergyAndPmfSum(&points_);
1225 if (needToUpdateTargetDistribution)
1227 /* The target distribution is always updated for all points at once. */
1228 updateTargetDistribution(points_, params);
1231 /* Update the bias. The bias is updated separately and last since it simply a function of
1232 the free energy and the target distribution and we want to avoid doing extra work. */
1233 for (int pointIndex : *updateList)
1235 points_[pointIndex].updateBias();
1238 /* Increase the update counter. */
1239 histogramSize_.incrementNumUpdates();
1242 double BiasState::updateProbabilityWeightsAndConvolvedBias(const std::vector<DimParams>& dimParams,
1243 const BiasGrid& grid,
1244 gmx::ArrayRef<const double> neighborLambdaEnergies,
1245 std::vector<double, AlignedAllocator<double>>* weight) const
1247 /* Only neighbors of the current coordinate value will have a non-negligible chance of getting sampled */
1248 const std::vector<int>& neighbors = grid.point(coordState_.gridpointIndex()).neighbor;
1250 #if GMX_SIMD_HAVE_DOUBLE
1251 typedef SimdDouble PackType;
1252 constexpr int packSize = GMX_SIMD_DOUBLE_WIDTH;
1254 typedef double PackType;
1255 constexpr int packSize = 1;
1257 /* Round the size of the weight array up to packSize */
1258 const int weightSize = ((neighbors.size() + packSize - 1) / packSize) * packSize;
1259 weight->resize(weightSize);
1261 double* gmx_restrict weightData = weight->data();
1262 PackType weightSumPack(0.0);
1263 for (size_t i = 0; i < neighbors.size(); i += packSize)
1265 for (size_t n = i; n < i + packSize; n++)
1267 if (n < neighbors.size())
1269 const int neighbor = neighbors[n];
1270 (*weight)[n] = biasedLogWeightFromPoint(
1271 dimParams, points_, grid, neighbor, points_[neighbor].bias(),
1272 coordState_.coordValue(), neighborLambdaEnergies, coordState_.gridpointIndex());
1276 /* Pad with values that don't affect the result */
1277 (*weight)[n] = detail::c_largeNegativeExponent;
1280 PackType weightPack = load<PackType>(weightData + i);
1281 weightPack = gmx::exp(weightPack);
1282 weightSumPack = weightSumPack + weightPack;
1283 store(weightData + i, weightPack);
1285 /* Sum of probability weights */
1286 double weightSum = reduce(weightSumPack);
1287 GMX_RELEASE_ASSERT(weightSum > 0,
1288 "zero probability weight when updating AWH probability weights.");
1290 /* Normalize probabilities to sum to 1 */
1291 double invWeightSum = 1 / weightSum;
1293 /* When there is a free energy lambda state axis remove the convolved contributions along that
1294 * axis from the total bias. This must be done after calculating invWeightSum (since weightSum
1295 * will be modified), but before normalizing the weights (below). */
1296 if (grid.hasLambdaAxis())
1298 /* If there is only one axis the bias will not be convolved in any dimension. */
1299 if (grid.axis().size() == 1)
1301 weightSum = gmx::exp(points_[coordState_.gridpointIndex()].bias());
1305 for (size_t i = 0; i < neighbors.size(); i++)
1307 const int neighbor = neighbors[i];
1308 if (pointsHaveDifferentLambda(grid, coordState_.gridpointIndex(), neighbor))
1310 weightSum -= weightData[i];
1316 for (double& w : *weight)
1321 /* Return the convolved bias */
1322 return std::log(weightSum);
1325 double BiasState::calcConvolvedBias(const std::vector<DimParams>& dimParams,
1326 const BiasGrid& grid,
1327 const awh_dvec& coordValue) const
1329 int point = grid.nearestIndex(coordValue);
1330 const GridPoint& gridPoint = grid.point(point);
1332 /* Sum the probability weights from the neighborhood of the given point */
1333 double weightSum = 0;
1334 for (int neighbor : gridPoint.neighbor)
1336 /* No convolution is required along the lambda dimension. */
1337 if (pointsHaveDifferentLambda(grid, point, neighbor))
1341 double logWeight = biasedLogWeightFromPoint(dimParams, points_, grid, neighbor,
1342 points_[neighbor].bias(), coordValue, {}, point);
1343 weightSum += std::exp(logWeight);
1346 /* Returns -GMX_FLOAT_MAX if no neighboring points were in the target region. */
1347 return (weightSum > 0) ? std::log(weightSum) : -GMX_FLOAT_MAX;
1350 void BiasState::sampleProbabilityWeights(const BiasGrid& grid, gmx::ArrayRef<const double> probWeightNeighbor)
1352 const std::vector<int>& neighbor = grid.point(coordState_.gridpointIndex()).neighbor;
1354 /* Save weights for next update */
1355 for (size_t n = 0; n < neighbor.size(); n++)
1357 points_[neighbor[n]].increaseWeightSumIteration(probWeightNeighbor[n]);
1360 /* Update the local update range. Two corner points define this rectangular
1361 * domain. We need to choose two new corner points such that the new domain
1362 * contains both the old update range and the current neighborhood.
1363 * In the simplest case when an update is performed every sample,
1364 * the update range would simply equal the current neighborhood.
1366 int neighborStart = neighbor[0];
1367 int neighborLast = neighbor[neighbor.size() - 1];
1368 for (int d = 0; d < grid.numDimensions(); d++)
1370 int origin = grid.point(neighborStart).index[d];
1371 int last = grid.point(neighborLast).index[d];
1375 /* Unwrap if wrapped around the boundary (only happens for periodic
1376 * boundaries). This has been already for the stored index interval.
1378 /* TODO: what we want to do is to find the smallest the update
1379 * interval that contains all points that need to be updated.
1380 * This amounts to combining two intervals, the current
1381 * [origin, end] update interval and the new touched neighborhood
1382 * into a new interval that contains all points from both the old
1385 * For periodic boundaries it becomes slightly more complicated
1386 * than for closed boundaries because then it needs not be
1387 * true that origin < end (so one can't simply relate the origin/end
1388 * in the min()/max() below). The strategy here is to choose the
1389 * origin closest to a reference point (index 0) and then unwrap
1390 * the end index if needed and choose the largest end index.
1391 * This ensures that both intervals are in the new interval
1392 * but it's not necessarily the smallest.
1393 * Currently we solve this by going through each possibility
1394 * and checking them.
1396 last += grid.axis(d).numPointsInPeriod();
1399 originUpdatelist_[d] = std::min(originUpdatelist_[d], origin);
1400 endUpdatelist_[d] = std::max(endUpdatelist_[d], last);
1404 void BiasState::sampleCoordAndPmf(const std::vector<DimParams>& dimParams,
1405 const BiasGrid& grid,
1406 gmx::ArrayRef<const double> probWeightNeighbor,
1407 double convolvedBias)
1409 /* Sampling-based deconvolution extracting the PMF.
1410 * Update the PMF histogram with the current coordinate value.
1412 * Because of the finite width of the harmonic potential, the free energy
1413 * defined for each coordinate point does not exactly equal that of the
1414 * actual coordinate, the PMF. However, the PMF can be estimated by applying
1415 * the relation exp(-PMF) = exp(-bias_convolved)*P_biased/Z, i.e. by keeping a
1416 * reweighted histogram of the coordinate value. Strictly, this relies on
1417 * the unknown normalization constant Z being either constant or known. Here,
1418 * neither is true except in the long simulation time limit. Empirically however,
1419 * it works (mainly because how the PMF histogram is rescaled).
1422 const int gridPointIndex = coordState_.gridpointIndex();
1423 const std::optional<int> lambdaAxisIndex = grid.lambdaAxisIndex();
1425 /* Update the PMF of points along a lambda axis with their bias. */
1426 if (lambdaAxisIndex)
1428 const std::vector<int>& neighbors = grid.point(gridPointIndex).neighbor;
1430 std::vector<double> lambdaMarginalDistribution =
1431 calculateFELambdaMarginalDistribution(grid, neighbors, probWeightNeighbor);
1433 awh_dvec coordValueAlongLambda = { coordState_.coordValue()[0], coordState_.coordValue()[1],
1434 coordState_.coordValue()[2], coordState_.coordValue()[3] };
1435 for (size_t i = 0; i < neighbors.size(); i++)
1437 const int neighbor = neighbors[i];
1439 if (pointsAlongLambdaAxis(grid, gridPointIndex, neighbor))
1441 const double neighborLambda = grid.point(neighbor).coordValue[lambdaAxisIndex.value()];
1442 if (neighbor == gridPointIndex)
1444 bias = convolvedBias;
1448 coordValueAlongLambda[lambdaAxisIndex.value()] = neighborLambda;
1449 bias = calcConvolvedBias(dimParams, grid, coordValueAlongLambda);
1452 const double probWeight = lambdaMarginalDistribution[neighborLambda];
1453 const double weightedBias = bias - std::log(std::max(probWeight, GMX_DOUBLE_MIN)); // avoid log(0)
1455 if (neighbor == gridPointIndex && grid.covers(coordState_.coordValue()))
1457 points_[neighbor].samplePmf(weightedBias);
1461 points_[neighbor].updatePmfUnvisited(weightedBias);
1468 /* Only save coordinate data that is in range (the given index is always
1469 * in range even if the coordinate value is not).
1471 if (grid.covers(coordState_.coordValue()))
1473 /* Save PMF sum and keep a histogram of the sampled coordinate values */
1474 points_[gridPointIndex].samplePmf(convolvedBias);
1478 /* Save probability weights for the update */
1479 sampleProbabilityWeights(grid, probWeightNeighbor);
1482 void BiasState::initHistoryFromState(AwhBiasHistory* biasHistory) const
1484 biasHistory->pointState.resize(points_.size());
1487 void BiasState::updateHistory(AwhBiasHistory* biasHistory, const BiasGrid& grid) const
1489 GMX_RELEASE_ASSERT(biasHistory->pointState.size() == points_.size(),
1490 "The AWH history setup does not match the AWH state.");
1492 AwhBiasStateHistory* stateHistory = &biasHistory->state;
1493 stateHistory->umbrellaGridpoint = coordState_.umbrellaGridpoint();
1495 for (size_t m = 0; m < biasHistory->pointState.size(); m++)
1497 AwhPointStateHistory* psh = &biasHistory->pointState[m];
1499 points_[m].storeState(psh);
1501 psh->weightsum_covering = weightSumCovering_[m];
1504 histogramSize_.storeState(stateHistory);
1506 stateHistory->origin_index_updatelist = multiDimGridIndexToLinear(grid, originUpdatelist_);
1507 stateHistory->end_index_updatelist = multiDimGridIndexToLinear(grid, endUpdatelist_);
1510 void BiasState::restoreFromHistory(const AwhBiasHistory& biasHistory, const BiasGrid& grid)
1512 const AwhBiasStateHistory& stateHistory = biasHistory.state;
1514 coordState_.restoreFromHistory(stateHistory);
1516 if (biasHistory.pointState.size() != points_.size())
1519 InvalidInputError("Bias grid size in checkpoint and simulation do not match. "
1520 "Likely you provided a checkpoint from a different simulation."));
1522 for (size_t m = 0; m < points_.size(); m++)
1524 points_[m].setFromHistory(biasHistory.pointState[m]);
1527 for (size_t m = 0; m < weightSumCovering_.size(); m++)
1529 weightSumCovering_[m] = biasHistory.pointState[m].weightsum_covering;
1532 histogramSize_.restoreFromHistory(stateHistory);
1534 linearGridindexToMultiDim(grid, stateHistory.origin_index_updatelist, originUpdatelist_);
1535 linearGridindexToMultiDim(grid, stateHistory.end_index_updatelist, endUpdatelist_);
1538 void BiasState::broadcast(const t_commrec* commRecord)
1540 gmx_bcast(sizeof(coordState_), &coordState_, commRecord->mpi_comm_mygroup);
1542 gmx_bcast(points_.size() * sizeof(PointState), points_.data(), commRecord->mpi_comm_mygroup);
1544 gmx_bcast(weightSumCovering_.size() * sizeof(double), weightSumCovering_.data(),
1545 commRecord->mpi_comm_mygroup);
1547 gmx_bcast(sizeof(histogramSize_), &histogramSize_, commRecord->mpi_comm_mygroup);
1550 void BiasState::setFreeEnergyToConvolvedPmf(const std::vector<DimParams>& dimParams, const BiasGrid& grid)
1552 std::vector<float> convolvedPmf;
1554 calcConvolvedPmf(dimParams, grid, &convolvedPmf);
1556 for (size_t m = 0; m < points_.size(); m++)
1558 points_[m].setFreeEnergy(convolvedPmf[m]);
1563 * Count trailing data rows containing only zeros.
1565 * \param[in] data 2D data array.
1566 * \param[in] numRows Number of rows in array.
1567 * \param[in] numColumns Number of cols in array.
1568 * \returns the number of trailing zero rows.
1570 static int countTrailingZeroRows(const double* const* data, int numRows, int numColumns)
1572 int numZeroRows = 0;
1573 for (int m = numRows - 1; m >= 0; m--)
1575 bool rowIsZero = true;
1576 for (int d = 0; d < numColumns; d++)
1578 if (data[d][m] != 0)
1587 /* At a row with non-zero data */
1592 /* Still at a zero data row, keep checking rows higher up. */
1601 * Initializes the PMF and target with data read from an input table.
1603 * \param[in] dimParams The dimension parameters.
1604 * \param[in] grid The grid.
1605 * \param[in] filename The filename to read PMF and target from.
1606 * \param[in] numBias Number of biases.
1607 * \param[in] biasIndex The index of the bias.
1608 * \param[in,out] pointState The state of the points in this bias.
1610 static void readUserPmfAndTargetDistribution(const std::vector<DimParams>& dimParams,
1611 const BiasGrid& grid,
1612 const std::string& filename,
1615 std::vector<PointState>* pointState)
1617 /* Read the PMF and target distribution.
1618 From the PMF, the convolved PMF, or the reference value free energy, can be calculated
1619 base on the force constant. The free energy and target together determine the bias.
1621 std::string filenameModified(filename);
1624 size_t n = filenameModified.rfind('.');
1625 GMX_RELEASE_ASSERT(n != std::string::npos,
1626 "The filename should contain an extension starting with .");
1627 filenameModified.insert(n, formatString("%d", biasIndex));
1630 std::string correctFormatMessage = formatString(
1631 "%s is expected in the following format. "
1632 "The first ndim column(s) should contain the coordinate values for each point, "
1633 "each column containing values of one dimension (in ascending order). "
1634 "For a multidimensional coordinate, points should be listed "
1635 "in the order obtained by traversing lower dimensions first. "
1636 "E.g. for two-dimensional grid of size nxn: "
1637 "(1, 1), (1, 2),..., (1, n), (2, 1), (2, 2), ..., , (n, n - 1), (n, n). "
1638 "Column ndim + 1 should contain the PMF value for each coordinate value. "
1639 "The target distribution values should be in column ndim + 2 or column ndim + 5. "
1640 "Make sure the input file ends with a new line but has no trailing new lines.",
1642 gmx::TextLineWrapper wrapper;
1643 wrapper.settings().setLineLength(c_linewidth);
1644 correctFormatMessage = wrapper.wrapToString(correctFormatMessage);
1648 int numRows = read_xvg(filenameModified.c_str(), &data, &numColumns);
1650 /* Check basic data properties here. BiasGrid takes care of more complicated things. */
1654 std::string mesg = gmx::formatString("%s is empty!.\n\n%s", filename.c_str(),
1655 correctFormatMessage.c_str());
1656 GMX_THROW(InvalidInputError(mesg));
1659 /* Less than 2 points is not useful for PMF or target. */
1662 std::string mesg = gmx::formatString(
1663 "%s contains too few data points (%d)."
1664 "The minimum number of points is 2.",
1665 filename.c_str(), numRows);
1666 GMX_THROW(InvalidInputError(mesg));
1669 /* Make sure there are enough columns of data.
1671 Two formats are allowed. Either with columns {coords, PMF, target} or
1672 {coords, PMF, x, y, z, target, ...}. The latter format is allowed since that
1673 is how AWH output is written (x, y, z being other AWH variables). For this format,
1674 trailing columns are ignored.
1676 int columnIndexTarget;
1677 int numColumnsMin = dimParams.size() + 2;
1678 int columnIndexPmf = dimParams.size();
1679 if (numColumns == numColumnsMin)
1681 columnIndexTarget = columnIndexPmf + 1;
1685 columnIndexTarget = columnIndexPmf + 4;
1688 if (numColumns < numColumnsMin)
1690 std::string mesg = gmx::formatString(
1691 "The number of columns in %s should be at least %d."
1693 filename.c_str(), numColumnsMin, correctFormatMessage.c_str());
1694 GMX_THROW(InvalidInputError(mesg));
1697 /* read_xvg can give trailing zero data rows for trailing new lines in the input. We allow 1 zero row,
1698 since this could be real data. But multiple trailing zero rows cannot correspond to valid data. */
1699 int numZeroRows = countTrailingZeroRows(data, numRows, numColumns);
1700 if (numZeroRows > 1)
1702 std::string mesg = gmx::formatString(
1703 "Found %d trailing zero data rows in %s. Please remove trailing empty lines and "
1705 numZeroRows, filename.c_str());
1706 GMX_THROW(InvalidInputError(mesg));
1709 /* Convert from user units to internal units before sending the data of to grid. */
1710 for (size_t d = 0; d < dimParams.size(); d++)
1712 double scalingFactor = dimParams[d].scaleUserInputToInternal(1);
1713 if (scalingFactor == 1)
1717 for (size_t m = 0; m < pointState->size(); m++)
1719 data[d][m] *= scalingFactor;
1723 /* Get a data point for each AWH grid point so that they all get data. */
1724 std::vector<int> gridIndexToDataIndex(grid.numPoints());
1725 mapGridToDataGrid(&gridIndexToDataIndex, data, numRows, filename, grid, correctFormatMessage);
1727 /* Extract the data for each grid point.
1728 * We check if the target distribution is zero for all points.
1730 bool targetDistributionIsZero = true;
1731 for (size_t m = 0; m < pointState->size(); m++)
1733 (*pointState)[m].setLogPmfSum(-data[columnIndexPmf][gridIndexToDataIndex[m]]);
1734 double target = data[columnIndexTarget][gridIndexToDataIndex[m]];
1736 /* Check if the values are allowed. */
1739 std::string mesg = gmx::formatString(
1740 "Target distribution weight at point %zu (%g) in %s is negative.", m, target,
1742 GMX_THROW(InvalidInputError(mesg));
1746 targetDistributionIsZero = false;
1748 (*pointState)[m].setTargetConstantWeight(target);
1751 if (targetDistributionIsZero)
1754 gmx::formatString("The target weights given in column %d in %s are all 0",
1755 columnIndexTarget, filename.c_str());
1756 GMX_THROW(InvalidInputError(mesg));
1759 /* Free the arrays. */
1760 for (int m = 0; m < numColumns; m++)
1767 void BiasState::normalizePmf(int numSharingSims)
1769 /* The normalization of the PMF estimate matters because it determines how big effect the next sample has.
1770 Approximately (for large enough force constant) we should have:
1771 sum_x(exp(-pmf(x)) = nsamples*sum_xref(exp(-f(xref)).
1774 /* Calculate the normalization factor, i.e. divide by the pmf sum, multiply by the number of samples and the f sum */
1775 double expSumPmf = 0;
1777 for (const PointState& pointState : points_)
1779 if (pointState.inTargetRegion())
1781 expSumPmf += std::exp(pointState.logPmfSum());
1782 expSumF += std::exp(-pointState.freeEnergy());
1785 double numSamples = histogramSize_.histogramSize() / numSharingSims;
1788 double logRenorm = std::log(numSamples * expSumF / expSumPmf);
1789 for (PointState& pointState : points_)
1791 if (pointState.inTargetRegion())
1793 pointState.setLogPmfSum(pointState.logPmfSum() + logRenorm);
1798 void BiasState::initGridPointState(const AwhBiasParams& awhBiasParams,
1799 const std::vector<DimParams>& dimParams,
1800 const BiasGrid& grid,
1801 const BiasParams& params,
1802 const std::string& filename,
1805 /* Modify PMF, free energy and the constant target distribution factor
1806 * to user input values if there is data given.
1808 if (awhBiasParams.bUserData)
1810 readUserPmfAndTargetDistribution(dimParams, grid, filename, numBias, params.biasIndex, &points_);
1811 setFreeEnergyToConvolvedPmf(dimParams, grid);
1814 /* The local Boltzmann distribution is special because the target distribution is updated as a function of the reference weighthistogram. */
1815 GMX_RELEASE_ASSERT(params.eTarget != eawhtargetLOCALBOLTZMANN || points_[0].weightSumRef() != 0,
1816 "AWH reference weight histogram not initialized properly with local "
1817 "Boltzmann target distribution.");
1819 updateTargetDistribution(points_, params);
1821 for (PointState& pointState : points_)
1823 if (pointState.inTargetRegion())
1825 pointState.updateBias();
1829 /* Note that for zero target this is a value that represents -infinity but should not be used for biasing. */
1830 pointState.setTargetToZero();
1834 /* Set the initial reference weighthistogram. */
1835 const double histogramSize = histogramSize_.histogramSize();
1836 for (auto& pointState : points_)
1838 pointState.setInitialReferenceWeightHistogram(histogramSize);
1841 /* Make sure the pmf is normalized consistently with the histogram size.
1842 Note: the target distribution and free energy need to be set here. */
1843 normalizePmf(params.numSharedUpdate);
1846 BiasState::BiasState(const AwhBiasParams& awhBiasParams,
1847 double histogramSizeInitial,
1848 const std::vector<DimParams>& dimParams,
1849 const BiasGrid& grid) :
1850 coordState_(awhBiasParams, dimParams, grid),
1851 points_(grid.numPoints()),
1852 weightSumCovering_(grid.numPoints()),
1853 histogramSize_(awhBiasParams, histogramSizeInitial)
1855 /* The minimum and maximum multidimensional point indices that are affected by the next update */
1856 for (size_t d = 0; d < dimParams.size(); d++)
1858 int index = grid.point(coordState_.gridpointIndex()).index[d];
1859 originUpdatelist_[d] = index;
1860 endUpdatelist_[d] = index;