c59fb1327f92fbc44d08a4c13fc2898872fbaab0
[alexxy/gromacs.git] / src / gromacs / awh / tests / bias.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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.
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 #include "gmxpre.h"
36
37 #include "gromacs/awh/bias.h"
38
39 #include <cmath>
40
41 #include <memory>
42 #include <tuple>
43 #include <vector>
44
45 #include <gmock/gmock.h>
46 #include <gtest/gtest.h>
47
48 #include "gromacs/awh/correlationgrid.h"
49 #include "gromacs/awh/pointstate.h"
50 #include "gromacs/mdtypes/awh_params.h"
51 #include "gromacs/utility/stringutil.h"
52
53 #include "testutils/refdata.h"
54 #include "testutils/testasserts.h"
55
56 namespace gmx
57 {
58
59 namespace test
60 {
61
62 /*! \internal \brief
63  * Struct that gathers all input for setting up and using a Bias
64  */
65 struct AwhTestParameters
66 {
67     AwhTestParameters() = default;
68     //! Move constructor
69     AwhTestParameters(AwhTestParameters&& o) noexcept :
70         beta(o.beta),
71         awhDimParams(o.awhDimParams),
72         awhBiasParams(o.awhBiasParams),
73         awhParams(o.awhParams),
74         dimParams(std::move(o.dimParams))
75     {
76         awhBiasParams.dimParams = &awhDimParams;
77         awhParams.awhBiasParams = &awhBiasParams;
78     }
79     double beta; //!< 1/(kB*T)
80
81     AwhDimParams  awhDimParams;  //!< Dimension parameters pointed to by \p awhBiasParams
82     AwhBiasParams awhBiasParams; //!< Bias parameters pointed to by \[ awhParams
83     AwhParams     awhParams;     //!< AWH parameters, this is the struct to actually use
84
85     std::vector<DimParams> dimParams; //!< Dimension parameters for setting up Bias
86 };
87
88 //! Helper function to set up the C-style AWH parameters for the test
89 static AwhTestParameters getAwhTestParameters(int eawhgrowth, int eawhpotential)
90 {
91     AwhTestParameters params;
92
93     params.beta = 0.4;
94
95     AwhDimParams& awhDimParams = params.awhDimParams;
96
97     awhDimParams.period         = 0;
98     awhDimParams.diffusion      = 0.1;
99     awhDimParams.origin         = 0.5;
100     awhDimParams.end            = 1.5;
101     awhDimParams.coordValueInit = awhDimParams.origin;
102     awhDimParams.coverDiameter  = 0;
103     awhDimParams.eCoordProvider = eawhcoordproviderPULL;
104
105     AwhBiasParams& awhBiasParams = params.awhBiasParams;
106
107     awhBiasParams.ndim                 = 1;
108     awhBiasParams.dimParams            = &awhDimParams;
109     awhBiasParams.eTarget              = eawhtargetCONSTANT;
110     awhBiasParams.targetBetaScaling    = 0;
111     awhBiasParams.targetCutoff         = 0;
112     awhBiasParams.eGrowth              = eawhgrowth;
113     awhBiasParams.bUserData            = FALSE;
114     awhBiasParams.errorInitial         = 0.5 / params.beta;
115     awhBiasParams.shareGroup           = 0;
116     awhBiasParams.equilibrateHistogram = FALSE;
117
118     double  convFactor = 1;
119     double  k          = 1000;
120     int64_t seed       = 93471803;
121
122     params.dimParams.emplace_back(convFactor, k, params.beta);
123
124     AwhParams& awhParams = params.awhParams;
125
126     awhParams.numBias                    = 1;
127     awhParams.awhBiasParams              = &awhBiasParams;
128     awhParams.seed                       = seed;
129     awhParams.nstOut                     = 0;
130     awhParams.nstSampleCoord             = 1;
131     awhParams.numSamplesUpdateFreeEnergy = 10;
132     awhParams.ePotential                 = eawhpotential;
133     awhParams.shareBiasMultisim          = FALSE;
134
135     return params;
136 }
137
138 //! Database of 21 test coordinates that represent a trajectory */
139 const double g_coords[] = { 0.62, 0.70, 0.68, 0.80, 0.93, 0.87, 1.16, 1.14, 0.95, 0.89, 0.91,
140                             0.86, 0.88, 0.79, 0.75, 0.82, 0.74, 0.70, 0.68, 0.71, 0.73 };
141
142 //! Convenience typedef: growth type enum, potential type enum, disable update skips
143 typedef std::tuple<int, int, BiasParams::DisableUpdateSkips> BiasTestParameters;
144
145 /*! \brief Test fixture for testing Bias updates
146  */
147 class BiasTest : public ::testing::TestWithParam<BiasTestParameters>
148 {
149 public:
150     //! Random seed for AWH MC sampling
151     int64_t seed_;
152
153     //! Coordinates representing a trajectory in time
154     std::vector<double> coordinates_;
155     //! The awh Bias
156     std::unique_ptr<Bias> bias_;
157
158     BiasTest() : coordinates_(std::begin(g_coords), std::end(g_coords))
159     {
160         /* We test all combinations of:
161          *   eawhgrowth:
162          *     eawhgrowthLINEAR:     final, normal update phase
163          *     ewahgrowthEXP_LINEAR: intial phase, updated size is constant
164          *   eawhpotential (should only affect the force output):
165          *     eawhpotentialUMBRELLA:  MC on lambda (umbrella potential location)
166          *     eawhpotentialCONVOLVED: MD on a convolved potential landscape
167          *   disableUpdateSkips (should not affect the results):
168          *     BiasParams::DisableUpdateSkips::yes: update the point state for every sample
169          *     BiasParams::DisableUpdateSkips::no:  update the point state at an interval > 1 sample
170          *
171          * Note: It would be nice to explicitly check that eawhpotential
172          *       and disableUpdateSkips do not affect the point state.
173          *       But the reference data will also ensure this.
174          */
175         int                            eawhgrowth;
176         int                            eawhpotential;
177         BiasParams::DisableUpdateSkips disableUpdateSkips;
178         std::tie(eawhgrowth, eawhpotential, disableUpdateSkips) = GetParam();
179
180         /* Set up a basic AWH setup with a single, 1D bias with parameters
181          * such that we can measure the effects of different parameters.
182          * The idea is to, among other things, have part of the interval
183          * not covered by samples.
184          */
185         const AwhTestParameters params = getAwhTestParameters(eawhgrowth, eawhpotential);
186
187         seed_ = params.awhParams.seed;
188
189         double mdTimeStep = 0.1;
190
191         int numSamples = coordinates_.size() - 1; // No sample taken at step 0
192         GMX_RELEASE_ASSERT(numSamples % params.awhParams.numSamplesUpdateFreeEnergy == 0,
193                            "This test is intended to reproduce the situation when the might need "
194                            "to write output during a normal AWH run, therefore the number of "
195                            "samples should be a multiple of the free-energy update interval (but "
196                            "the test should also runs fine without this condition).");
197
198         bias_ = std::make_unique<Bias>(-1, params.awhParams, params.awhBiasParams, params.dimParams,
199                                        params.beta, mdTimeStep, 1, "", Bias::ThisRankWillDoIO::No,
200                                        disableUpdateSkips);
201     }
202 };
203
204 TEST_P(BiasTest, ForcesBiasPmf)
205 {
206     gmx::test::TestReferenceData    data;
207     gmx::test::TestReferenceChecker checker(data.rootChecker());
208
209     Bias& bias = *bias_;
210
211     /* Make strings with the properties we expect to be different in the tests.
212      * These also helps to interpret the reference data.
213      */
214     std::vector<std::string> props;
215     props.push_back(formatString("stage:           %s", bias.state().inInitialStage() ? "initial" : "final"));
216     props.push_back(formatString("convolve forces: %s", bias.params().convolveForce ? "yes" : "no"));
217     props.push_back(formatString("skip updates:    %s", bias.params().skipUpdates() ? "yes" : "no"));
218
219     SCOPED_TRACE(gmx::formatString("%s, %s, %s", props[0].c_str(), props[1].c_str(), props[2].c_str()));
220
221     std::vector<double> force, pot, potJump;
222
223     double  coordMaxValue = 0;
224     double  potentialJump = 0;
225     int64_t step          = 0;
226     for (auto& coord : coordinates_)
227     {
228         coordMaxValue = std::max(coordMaxValue, std::abs(coord));
229
230         awh_dvec                    coordValue = { coord, 0, 0, 0 };
231         double                      potential  = 0;
232         gmx::ArrayRef<const double> biasForce =
233                 bias.calcForceAndUpdateBias(coordValue, {}, {}, &potential, &potentialJump, nullptr,
234                                             nullptr, step, step, seed_, nullptr);
235
236         force.push_back(biasForce[0]);
237         pot.push_back(potential);
238         potJump.push_back(potentialJump);
239
240         step++;
241     }
242
243     /* When skipping updates, ensure all skipped updates are performed here.
244      * This should result in the same bias state as at output in a normal run.
245      */
246     if (bias.params().skipUpdates())
247     {
248         bias.doSkippedUpdatesForAllPoints();
249     }
250
251     std::vector<double> pointBias, logPmfsum;
252     for (auto& point : bias.state().points())
253     {
254         pointBias.push_back(point.bias());
255         logPmfsum.push_back(point.logPmfSum());
256     }
257
258     /* The umbrella force is computed from the coordinate deviation.
259      * In taking this deviation we lose a lot of precision, so we should
260      * compare against k*max(coord) instead of the instantaneous force.
261      */
262     const double kCoordMax = bias.dimParams()[0].k * coordMaxValue;
263
264     constexpr int ulpTol = 10;
265
266     checker.checkSequence(props.begin(), props.end(), "Properties");
267     checker.setDefaultTolerance(absoluteTolerance(kCoordMax * GMX_DOUBLE_EPS * ulpTol));
268     checker.checkSequence(force.begin(), force.end(), "Force");
269     checker.checkSequence(pot.begin(), pot.end(), "Potential");
270     checker.checkSequence(potJump.begin(), potJump.end(), "PotentialJump");
271     checker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpTol));
272     checker.checkSequence(pointBias.begin(), pointBias.end(), "PointBias");
273     checker.checkSequence(logPmfsum.begin(), logPmfsum.end(), "PointLogPmfsum");
274 }
275
276 /* Scan initial/final phase, MC/convolved force and update skip (not) allowed
277  * Both the convolving and skipping should not affect the bias and PMF.
278  * It would be nice if the test would explicitly check for this.
279  * Currently this is tested through identical reference data.
280  */
281 INSTANTIATE_TEST_CASE_P(WithParameters,
282                         BiasTest,
283                         ::testing::Combine(::testing::Values(eawhgrowthLINEAR, eawhgrowthEXP_LINEAR),
284                                            ::testing::Values(eawhpotentialUMBRELLA, eawhpotentialCONVOLVED),
285                                            ::testing::Values(BiasParams::DisableUpdateSkips::yes,
286                                                              BiasParams::DisableUpdateSkips::no)));
287
288 // Test that we detect coverings and exit the initial stage at the correct step
289 TEST(BiasTest, DetectsCovering)
290 {
291     const AwhTestParameters params = getAwhTestParameters(eawhgrowthEXP_LINEAR, eawhpotentialCONVOLVED);
292     const AwhDimParams&     awhDimParams = params.awhParams.awhBiasParams[0].dimParams[0];
293
294     const double mdTimeStep = 0.1;
295
296     Bias bias(-1, params.awhParams, params.awhBiasParams, params.dimParams, params.beta, mdTimeStep,
297               1, "", Bias::ThisRankWillDoIO::No);
298
299     /* We use a trajectory of the sum of two sines to cover the reaction
300      * coordinate range in a semi-realistic way. The period is 4*pi=12.57.
301      * We get out of the initial stage after 4 coverings at step 300.
302      */
303     const int64_t exitStepRef = 300;
304     const double  midPoint    = 0.5 * (awhDimParams.end + awhDimParams.origin);
305     const double  halfWidth   = 0.5 * (awhDimParams.end - awhDimParams.origin);
306
307     bool inInitialStage = bias.state().inInitialStage();
308     /* Normally this loop exits at exitStepRef, but we extend with failure */
309     int64_t step;
310     for (step = 0; step <= 2 * exitStepRef; step++)
311     {
312         double t     = step * mdTimeStep;
313         double coord = midPoint + halfWidth * (0.5 * std::sin(t) + 0.55 * std::sin(1.5 * t));
314
315         awh_dvec coordValue    = { coord, 0, 0, 0 };
316         double   potential     = 0;
317         double   potentialJump = 0;
318         bias.calcForceAndUpdateBias(coordValue, {}, {}, &potential, &potentialJump, nullptr,
319                                     nullptr, step, step, params.awhParams.seed, nullptr);
320
321         inInitialStage = bias.state().inInitialStage();
322         if (!inInitialStage)
323         {
324             break;
325         }
326     }
327
328     EXPECT_EQ(false, inInitialStage);
329     if (!inInitialStage)
330     {
331         EXPECT_EQ(exitStepRef, step);
332     }
333 }
334
335 } // namespace test
336 } // namespace gmx