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