efd58c33977acce848e9d3b51a70303addba60bc
[alexxy/gromacs.git] / src / gromacs / awh / coordstate.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  * \brief
38  * Implements the CoordState class.
39  *
40  * \author Viveca Lindahl
41  * \author Berk Hess <hess@kth.se>
42  * \ingroup module_awh
43  */
44
45 #include "gmxpre.h"
46
47 #include "coordstate.h"
48
49 #include <algorithm>
50
51 #include "gromacs/math/utilities.h"
52 #include "gromacs/mdtypes/awh-history.h"
53 #include "gromacs/mdtypes/awh-params.h"
54 #include "gromacs/random/threefry.h"
55 #include "gromacs/random/uniformrealdistribution.h"
56 #include "gromacs/utility/arrayref.h"
57 #include "gromacs/utility/exceptions.h"
58 #include "gromacs/utility/gmxassert.h"
59 #include "gromacs/utility/stringutil.h"
60
61 #include "grid.h"
62
63 namespace gmx
64 {
65
66 CoordState::CoordState(const AwhBiasParams          &awhBiasParams,
67                        const std::vector<DimParams> &dimParams,
68                        const Grid                   &grid)
69 {
70     for (size_t d = 0; d < dimParams.size(); d++)
71     {
72         coordValue_[d] = dimParams[d].scaleUserInputToInternal(awhBiasParams.dimParams[d].coordValueInit);
73     }
74
75     /* The grid-point index is always the nearest point to the coordinate.
76      * We initialize the umbrella location to the same nearest point.
77      * More correctly one would sample from the biased distribution,
78      * but it doesn't really matter, as this will happen after a few steps.
79      */
80     gridpointIndex_    = grid.nearestIndex(coordValue_);
81     umbrellaGridpoint_ = gridpointIndex_;
82 }
83
84 namespace
85 {
86
87 /*! \brief Generate a sample from a discrete probability distribution defined on [0, distr.size() - 1].
88  *
89  * The pair (indexSeed0,indexSeed1) should be different for every invocation.
90  *
91  * \param[in] distr       Normalized probability distribution to generate a sample from.
92  * \param[in] seed        Random seed for initializing the random number generator.
93  * \param[in] indexSeed0  Random seed needed by the random number generator.
94  * \param[in] indexSeed1  Random seed needed by the random number generator.
95  * \returns a sample index in [0, distr.size() - 1]
96  */
97 int getSampleFromDistribution(ArrayRef<const double> distr,
98                               int64_t                seed,
99                               int64_t                indexSeed0,
100                               int64_t                indexSeed1)
101 {
102     gmx::ThreeFry2x64<0>               rng(seed, gmx::RandomDomain::AwhBiasing);
103     gmx::UniformRealDistribution<real> uniformRealDistr;
104
105     GMX_RELEASE_ASSERT(!distr.empty(), "We need a non-zero length distribution to sample from");
106
107     /* Generate the cumulative probability distribution function */
108     std::vector<double> cumulativeDistribution(distr.size());
109
110     cumulativeDistribution[0] = distr[0];
111
112     for (gmx::index i = 1; i < distr.size(); i++)
113     {
114         cumulativeDistribution[i] = cumulativeDistribution[i - 1] + distr[i];
115     }
116
117     GMX_RELEASE_ASSERT(gmx_within_tol(cumulativeDistribution.back(), 1.0, 0.01), "Attempt to get sample from non-normalized/zero distribution");
118
119     /* Use binary search to convert the real value to an integer in [0, ndistr - 1] distributed according to distr. */
120     rng.restart(indexSeed0, indexSeed1);
121
122     double value  = uniformRealDistr(rng);
123     int    sample = std::upper_bound(cumulativeDistribution.begin(), cumulativeDistribution.end() - 1, value) - cumulativeDistribution.begin();
124
125     return sample;
126 }
127
128 }   // namespace
129
130 void
131 CoordState::sampleUmbrellaGridpoint(const Grid                  &grid,
132                                     int                          gridpointIndex,
133                                     gmx::ArrayRef<const double>  probWeightNeighbor,
134                                     int64_t                      step,
135                                     int64_t                      seed,
136                                     int                          indexSeed)
137 {
138     /* Sample new umbrella reference value from the probability distribution
139      * which is defined for the neighboring points of the current coordinate.
140      */
141     const std::vector<int> &neighbor = grid.point(gridpointIndex).neighbor;
142
143     /* In order to use the same seed for all AWH biases and get independent
144        samples we use the index of the bias. */
145     int localIndex     = getSampleFromDistribution(probWeightNeighbor,
146                                                    seed, step, indexSeed);
147
148     umbrellaGridpoint_ = neighbor[localIndex];
149 }
150
151 void CoordState::setCoordValue(const Grid     &grid,
152                                const awh_dvec  coordValue)
153 {
154     /* We need to check for valid (probable) coordinate values, to give
155      * a clear error message instead of a low-level assertion failure.
156      * We allow values up to 10*sigma beyond the bounds. For points at
157      * the bounds this means a chance of less than 2*10-45 of a false positive
158      * in the usual case that the PMF is flat or goes up beyond the bounds.
159      * In the worst reasonable case of the PMF curving down with a curvature
160      * of half the harmonic force constant, the chance is 1.5*10-12.
161      */
162     constexpr int c_marginInSigma = 10;
163
164     for (int dim = 0; dim < grid.numDimensions(); dim++)
165     {
166         const GridAxis &axis = grid.axis(dim);
167         /* We do not check periodic coordinates, since that is more complicated
168          * and those cases are less likely to cause problems.
169          */
170         if (!axis.isPeriodic())
171         {
172             const double margin = axis.spacing()*c_marginInSigma/Grid::c_numPointsPerSigma;
173             if (coordValue[dim] < axis.origin() - margin ||
174                 coordValue[dim] > axis.origin() + axis.length() + margin)
175             {
176                 std::string mesg = gmx::formatString("Coordinate %d of an AWH bias has a value %f which is more than %d sigma out of the AWH range of [%f, %f]. You seem to have an unstable reaction coordinate setup or an unequilibrated system.",
177                                                      dim + 1,
178                                                      coordValue[dim],
179                                                      c_marginInSigma,
180                                                      axis.origin(),
181                                                      axis.origin() + axis.length());
182                 GMX_THROW(SimulationInstabilityError(mesg));
183             }
184         }
185
186         coordValue_[dim] = coordValue[dim];
187     }
188
189     /* The grid point closest to the coordinate value defines the current
190      * neighborhood of points. Besides at steps when global updates and/or
191      * checks are performed, only the neighborhood will be touched.
192      */
193     gridpointIndex_ = grid.nearestIndex(coordValue_);
194 }
195
196 void
197 CoordState::restoreFromHistory(const AwhBiasStateHistory &stateHistory)
198 {
199     umbrellaGridpoint_ = stateHistory.umbrellaGridpoint;
200 }
201
202 } // namespace gmx