Merge branch release-2021 into merge-2021-into-master
[alexxy/gromacs.git] / src / gromacs / applied_forces / awh / tests / bias_fep_lambda_state.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 <cmath>
38
39 #include <memory>
40 #include <tuple>
41 #include <vector>
42
43 #include <gmock/gmock.h>
44 #include <gtest/gtest.h>
45
46 #include "gromacs/applied_forces/awh/bias.h"
47 #include "gromacs/applied_forces/awh/correlationgrid.h"
48 #include "gromacs/applied_forces/awh/pointstate.h"
49 #include "gromacs/mdtypes/awh_params.h"
50 #include "gromacs/utility/stringutil.h"
51
52 #include "gromacs/applied_forces/awh/tests/awh_setup.h"
53 #include "testutils/refdata.h"
54 #include "testutils/testasserts.h"
55
56 namespace gmx
57 {
58
59 namespace test
60 {
61
62 //! The number of lambda states to use in the tests.
63 static constexpr int c_numLambdaStates = 16;
64
65
66 //! Convenience typedef: growth type enum, potential type enum, disable update skips
67 typedef std::tuple<AwhHistogramGrowthType, AwhPotentialType, BiasParams::DisableUpdateSkips> BiasTestParameters;
68
69 /*! \brief Test fixture for testing Bias updates
70  */
71 class BiasFepLambdaStateTest : public ::testing::TestWithParam<BiasTestParameters>
72 {
73 private:
74     //! Storage for test parameters.
75     std::unique_ptr<AwhTestParameters> params_;
76
77 public:
78     //! Random seed for AWH MC sampling
79     int64_t seed_;
80
81     //! The awh Bias
82     std::unique_ptr<Bias> bias_;
83
84     BiasFepLambdaStateTest()
85     {
86         /* We test all combinations of:
87          *   eawhgrowth:
88          *     eawhgrowthLINEAR:     final, normal update phase
89          *     ewahgrowthEXP_LINEAR: intial phase, updated size is constant
90          *   eawhpotential (test only eawhpotentialUMBRELLA (MC) for FEP lambda dimensions)
91          *   disableUpdateSkips (should not affect the results):
92          *     BiasParams::DisableUpdateSkips::yes: update the point state for every sample
93          *     BiasParams::DisableUpdateSkips::no:  update the point state at an interval > 1 sample
94          *
95          * Note: It would be nice to explicitly check that eawhpotential
96          *       and disableUpdateSkips do not affect the point state.
97          *       But the reference data will also ensure this.
98          */
99         AwhHistogramGrowthType         eawhgrowth;
100         AwhPotentialType               eawhpotential;
101         BiasParams::DisableUpdateSkips disableUpdateSkips;
102         std::tie(eawhgrowth, eawhpotential, disableUpdateSkips) = GetParam();
103
104         /* Set up a basic AWH setup with a single, 1D bias with parameters
105          * such that we can measure the effects of different parameters.
106          */
107         constexpr AwhCoordinateProviderType coordinateProvider = AwhCoordinateProviderType::FreeEnergyLambda;
108         constexpr int                       coordIndex = 0;
109         constexpr double                    origin     = 0;
110         constexpr double                    end        = c_numLambdaStates - 1;
111         constexpr double                    period     = 0;
112         // Correction for removal of GaussianGeometryFactor/2 in histogram size
113         constexpr double diffusion = 1e-4 / (0.12927243028700 * 2);
114         const auto       awhDimBuffer =
115                 awhDimParamSerialized(coordinateProvider, coordIndex, origin, end, period, diffusion);
116         const auto awhDimArrayRef = gmx::arrayRefFromArray(&awhDimBuffer, 1);
117         params_                   = std::make_unique<AwhTestParameters>(getAwhTestParameters(
118                 eawhgrowth, eawhpotential, awhDimArrayRef, false, 0.4, true, 1.0, c_numLambdaStates));
119
120         seed_ = params_->awhParams.seed();
121
122         double mdTimeStep = 0.1;
123
124         bias_ = std::make_unique<Bias>(-1,
125                                        params_->awhParams,
126                                        params_->awhParams.awhBiasParams()[0],
127                                        params_->dimParams,
128                                        params_->beta,
129                                        mdTimeStep,
130                                        1,
131                                        "",
132                                        Bias::ThisRankWillDoIO::No,
133                                        disableUpdateSkips);
134     }
135 };
136
137 TEST_P(BiasFepLambdaStateTest, ForcesBiasPmf)
138 {
139     gmx::test::TestReferenceData    data;
140     gmx::test::TestReferenceChecker checker(data.rootChecker());
141
142     Bias& bias = *bias_;
143
144     /* Make strings with the properties we expect to be different in the tests.
145      * These also helps to interpret the reference data.
146      */
147     std::vector<std::string> props;
148     props.push_back(formatString("stage:           %s", bias.state().inInitialStage() ? "initial" : "final"));
149     props.push_back(formatString("convolve forces: %s", bias.params().convolveForce ? "yes" : "no"));
150     props.push_back(formatString("skip updates:    %s", bias.params().skipUpdates() ? "yes" : "no"));
151
152     SCOPED_TRACE(gmx::formatString("%s, %s, %s", props[0].c_str(), props[1].c_str(), props[2].c_str()));
153
154     std::vector<double> force, pot;
155
156     double potentialJump = 0;
157     double mdTimeStep    = 0.1;
158     int    nSteps        = 501;
159
160     /* Some energies to use as base values (to which some noise is added later on). */
161     std::vector<double> neighborLambdaEnergies(c_numLambdaStates);
162     std::vector<double> neighborLambdaDhdl(c_numLambdaStates);
163     const double        magnitude = 12.0;
164     for (int i = 0; i < c_numLambdaStates; i++)
165     {
166         neighborLambdaEnergies[i] = magnitude * std::sin(i * 0.1);
167         neighborLambdaDhdl[i]     = magnitude * std::cos(i * 0.1);
168     }
169
170     for (int step = 0; step < nSteps; step++)
171     {
172         int      umbrellaGridpointIndex = bias.state().coordState().umbrellaGridpoint();
173         awh_dvec coordValue = { bias.getGridCoordValue(umbrellaGridpointIndex)[0], 0, 0, 0 };
174         double   potential  = 0;
175         gmx::ArrayRef<const double> biasForce = bias.calcForceAndUpdateBias(coordValue,
176                                                                             neighborLambdaEnergies,
177                                                                             neighborLambdaDhdl,
178                                                                             &potential,
179                                                                             &potentialJump,
180                                                                             nullptr,
181                                                                             nullptr,
182                                                                             step * mdTimeStep,
183                                                                             step,
184                                                                             seed_,
185                                                                             nullptr);
186
187         force.push_back(biasForce[0]);
188         pot.push_back(potential);
189     }
190
191     /* When skipping updates, ensure all skipped updates are performed here.
192      * This should result in the same bias state as at output in a normal run.
193      */
194     if (bias.params().skipUpdates())
195     {
196         bias.doSkippedUpdatesForAllPoints();
197     }
198
199     std::vector<double> pointBias, logPmfsum;
200     for (const auto& point : bias.state().points())
201     {
202         pointBias.push_back(point.bias());
203         logPmfsum.push_back(point.logPmfSum());
204     }
205
206     constexpr int ulpTol = 10;
207
208     checker.checkSequence(props.begin(), props.end(), "Properties");
209     checker.setDefaultTolerance(absoluteTolerance(magnitude * GMX_DOUBLE_EPS * ulpTol));
210     checker.checkSequence(force.begin(), force.end(), "Force");
211     checker.checkSequence(pot.begin(), pot.end(), "Potential");
212     checker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpTol));
213     checker.checkSequence(pointBias.begin(), pointBias.end(), "PointBias");
214     checker.checkSequence(logPmfsum.begin(), logPmfsum.end(), "PointLogPmfsum");
215 }
216
217 /* Scan initial/final phase, MC/convolved force and update skip (not) allowed
218  * Both the convolving and skipping should not affect the bias and PMF.
219  * It would be nice if the test would explicitly check for this.
220  * Currently this is tested through identical reference data.
221  */
222 INSTANTIATE_TEST_CASE_P(WithParameters,
223                         BiasFepLambdaStateTest,
224                         ::testing::Combine(::testing::Values(AwhHistogramGrowthType::Linear,
225                                                              AwhHistogramGrowthType::ExponentialLinear),
226                                            ::testing::Values(AwhPotentialType::Umbrella),
227                                            ::testing::Values(BiasParams::DisableUpdateSkips::yes,
228                                                              BiasParams::DisableUpdateSkips::no)));
229
230 // Test that we detect coverings and exit the initial stage at the correct step
231 TEST(BiasFepLambdaStateTest, DetectsCovering)
232 {
233     constexpr AwhCoordinateProviderType coordinateProvider = AwhCoordinateProviderType::FreeEnergyLambda;
234     constexpr int                       coordIndex         = 0;
235     constexpr double                    origin             = 0;
236     constexpr double                    end                = c_numLambdaStates - 1;
237     constexpr double                    period             = 0;
238     constexpr double                    diffusion          = 1e-4 / (0.12927243028700 * 2);
239     auto                                awhDimBuffer =
240             awhDimParamSerialized(coordinateProvider, coordIndex, origin, end, period, diffusion);
241     auto                    awhDimArrayRef = gmx::arrayRefFromArray(&awhDimBuffer, 1);
242     const AwhTestParameters params(getAwhTestParameters(AwhHistogramGrowthType::ExponentialLinear,
243                                                         AwhPotentialType::Umbrella,
244                                                         awhDimArrayRef,
245                                                         false,
246                                                         0.4,
247                                                         true,
248                                                         1.0,
249                                                         c_numLambdaStates));
250
251     const double mdTimeStep = 0.1;
252
253     Bias bias(-1,
254               params.awhParams,
255               params.awhParams.awhBiasParams()[0],
256               params.dimParams,
257               params.beta,
258               mdTimeStep,
259               1,
260               "",
261               Bias::ThisRankWillDoIO::No);
262
263     const int64_t exitStepRef = 320;
264
265     bool inInitialStage = bias.state().inInitialStage();
266
267     /* Some energies to use as base values (to which some noise is added later on). */
268     std::vector<double> neighborLambdaEnergies(c_numLambdaStates);
269     std::vector<double> neighborLambdaDhdl(c_numLambdaStates);
270     const double        magnitude = 12.0;
271     for (int i = 0; i < c_numLambdaStates; i++)
272     {
273         neighborLambdaEnergies[i] = magnitude * std::sin(i * 0.1);
274         neighborLambdaDhdl[i]     = magnitude * std::cos(i * 0.1);
275     }
276
277     int64_t step;
278     /* Normally this loop exits at exitStepRef, but we extend with failure */
279     for (step = 0; step <= 2 * exitStepRef; step++)
280     {
281         int      umbrellaGridpointIndex = bias.state().coordState().umbrellaGridpoint();
282         awh_dvec coordValue = { bias.getGridCoordValue(umbrellaGridpointIndex)[0], 0, 0, 0 };
283
284         double potential     = 0;
285         double potentialJump = 0;
286         bias.calcForceAndUpdateBias(coordValue,
287                                     neighborLambdaEnergies,
288                                     neighborLambdaDhdl,
289                                     &potential,
290                                     &potentialJump,
291                                     nullptr,
292                                     nullptr,
293                                     step,
294                                     step,
295                                     params.awhParams.seed(),
296                                     nullptr);
297
298         inInitialStage = bias.state().inInitialStage();
299         if (!inInitialStage)
300         {
301             break;
302         }
303     }
304
305     EXPECT_EQ(false, inInitialStage);
306     if (!inInitialStage)
307     {
308         EXPECT_EQ(exitStepRef, step);
309     }
310 }
311
312 } // namespace test
313 } // namespace gmx