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