AWH: Remove random term from FEP tests
[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, 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 "testutils/refdata.h"
53 #include "testutils/testasserts.h"
54
55 namespace gmx
56 {
57
58 namespace test
59 {
60
61 //! The number of lambda states to use in the tests.
62 const int numLambdaStates = 16;
63
64 /*! \internal \brief
65  * Struct that gathers all input for setting up and using a Bias
66  */
67 struct AwhFepLambdaStateTestParameters
68 {
69     AwhFepLambdaStateTestParameters() = default;
70     //! Move constructor
71     AwhFepLambdaStateTestParameters(AwhFepLambdaStateTestParameters&& o) noexcept :
72         beta(o.beta),
73         awhDimParams(o.awhDimParams),
74         awhBiasParams(o.awhBiasParams),
75         awhParams(o.awhParams),
76         dimParams(std::move(o.dimParams))
77     {
78         awhBiasParams.dimParams = &awhDimParams;
79         awhParams.awhBiasParams = &awhBiasParams;
80     }
81     double beta; //!< 1/(kB*T)
82
83     AwhDimParams  awhDimParams;  //!< Dimension parameters pointed to by \p awhBiasParams
84     AwhBiasParams awhBiasParams; //!< Bias parameters pointed to by \[ awhParams
85     AwhParams     awhParams;     //!< AWH parameters, this is the struct to actually use
86
87     std::vector<DimParams> dimParams; //!< Dimension parameters for setting up Bias
88 };
89
90 //! Helper function to set up the C-style AWH parameters for the test
91 static AwhFepLambdaStateTestParameters getAwhFepLambdaTestParameters(int eawhgrowth, int eawhpotential)
92 {
93     AwhFepLambdaStateTestParameters params;
94
95     params.beta = 0.4;
96
97     AwhDimParams& awhDimParams = params.awhDimParams;
98
99     awhDimParams.period         = 0;
100     awhDimParams.diffusion      = 1e-4;
101     awhDimParams.origin         = 0;
102     awhDimParams.end            = numLambdaStates - 1;
103     awhDimParams.coordValueInit = awhDimParams.origin;
104     awhDimParams.coverDiameter  = 0;
105     awhDimParams.eCoordProvider = eawhcoordproviderFREE_ENERGY_LAMBDA;
106
107     AwhBiasParams& awhBiasParams = params.awhBiasParams;
108
109     awhBiasParams.ndim                 = 1;
110     awhBiasParams.dimParams            = &awhDimParams;
111     awhBiasParams.eTarget              = eawhtargetCONSTANT;
112     awhBiasParams.targetBetaScaling    = 0;
113     awhBiasParams.targetCutoff         = 0;
114     awhBiasParams.eGrowth              = eawhgrowth;
115     awhBiasParams.bUserData            = FALSE;
116     awhBiasParams.errorInitial         = 1.0 / params.beta;
117     awhBiasParams.shareGroup           = 0;
118     awhBiasParams.equilibrateHistogram = FALSE;
119
120     int64_t seed = 93471803;
121
122     params.dimParams.push_back(DimParams::fepLambdaDimParams(numLambdaStates, 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 //! Convenience typedef: growth type enum, potential type enum, disable update skips
139 typedef std::tuple<int, int, BiasParams::DisableUpdateSkips> BiasTestParameters;
140
141 /*! \brief Test fixture for testing Bias updates
142  */
143 class BiasFepLambdaStateTest : public ::testing::TestWithParam<BiasTestParameters>
144 {
145 public:
146     //! Random seed for AWH MC sampling
147     int64_t seed_;
148
149     //! The awh Bias
150     std::unique_ptr<Bias> bias_;
151
152     BiasFepLambdaStateTest()
153     {
154         /* We test all combinations of:
155          *   eawhgrowth:
156          *     eawhgrowthLINEAR:     final, normal update phase
157          *     ewahgrowthEXP_LINEAR: intial phase, updated size is constant
158          *   eawhpotential (test both, but for the FEP lambda state dimension MC will in practice be used,
159          *                  except that eawhpotentialCONVOLVED also gives a potential output):
160          *     eawhpotentialUMBRELLA:  MC on lambda state
161          *     eawhpotentialCONVOLVED: MD on a convolved potential landscape (falling back to MC on lambda state)
162          *   disableUpdateSkips (should not affect the results):
163          *     BiasParams::DisableUpdateSkips::yes: update the point state for every sample
164          *     BiasParams::DisableUpdateSkips::no:  update the point state at an interval > 1 sample
165          *
166          * Note: It would be nice to explicitly check that eawhpotential
167          *       and disableUpdateSkips do not affect the point state.
168          *       But the reference data will also ensure this.
169          */
170         int                            eawhgrowth;
171         int                            eawhpotential;
172         BiasParams::DisableUpdateSkips disableUpdateSkips;
173         std::tie(eawhgrowth, eawhpotential, disableUpdateSkips) = GetParam();
174
175         /* Set up a basic AWH setup with a single, 1D bias with parameters
176          * such that we can measure the effects of different parameters.
177          */
178         const AwhFepLambdaStateTestParameters params =
179                 getAwhFepLambdaTestParameters(eawhgrowth, eawhpotential);
180
181         seed_ = params.awhParams.seed;
182
183         double mdTimeStep = 0.1;
184
185         bias_ = std::make_unique<Bias>(-1, params.awhParams, params.awhBiasParams, params.dimParams,
186                                        params.beta, mdTimeStep, 1, "", Bias::ThisRankWillDoIO::No,
187                                        disableUpdateSkips);
188     }
189 };
190
191 TEST_P(BiasFepLambdaStateTest, ForcesBiasPmf)
192 {
193     gmx::test::TestReferenceData    data;
194     gmx::test::TestReferenceChecker checker(data.rootChecker());
195
196     Bias& bias = *bias_;
197
198     /* Make strings with the properties we expect to be different in the tests.
199      * These also helps to interpret the reference data.
200      */
201     std::vector<std::string> props;
202     props.push_back(formatString("stage:           %s", bias.state().inInitialStage() ? "initial" : "final"));
203     props.push_back(formatString("convolve forces: %s", bias.params().convolveForce ? "yes" : "no"));
204     props.push_back(formatString("skip updates:    %s", bias.params().skipUpdates() ? "yes" : "no"));
205
206     SCOPED_TRACE(gmx::formatString("%s, %s, %s", props[0].c_str(), props[1].c_str(), props[2].c_str()));
207
208     std::vector<double> force, pot;
209
210     double potentialJump = 0;
211     double mdTimeStep    = 0.1;
212     int    nSteps        = 501;
213
214     /* Some energies to use as base values (to which some noise is added later on). */
215     std::vector<double> neighborLambdaEnergies(numLambdaStates);
216     std::vector<double> neighborLambdaDhdl(numLambdaStates);
217     const double        magnitude = 12.0;
218     for (int i = 0; i < numLambdaStates; i++)
219     {
220         neighborLambdaEnergies[i] = magnitude * std::sin(i * 0.1);
221         neighborLambdaDhdl[i]     = magnitude * std::cos(i * 0.1);
222     }
223
224     for (int step = 0; step < nSteps; step++)
225     {
226         int      umbrellaGridpointIndex = bias.state().coordState().umbrellaGridpoint();
227         awh_dvec coordValue = { bias.getGridCoordValue(umbrellaGridpointIndex)[0], 0, 0, 0 };
228         double   potential  = 0;
229         gmx::ArrayRef<const double> biasForce = bias.calcForceAndUpdateBias(
230                 coordValue, neighborLambdaEnergies, neighborLambdaDhdl, &potential, &potentialJump,
231                 nullptr, nullptr, step * mdTimeStep, step, seed_, nullptr);
232
233         force.push_back(biasForce[0]);
234         pot.push_back(potential);
235     }
236
237     /* When skipping updates, ensure all skipped updates are performed here.
238      * This should result in the same bias state as at output in a normal run.
239      */
240     if (bias.params().skipUpdates())
241     {
242         bias.doSkippedUpdatesForAllPoints();
243     }
244
245     std::vector<double> pointBias, logPmfsum;
246     for (auto& point : bias.state().points())
247     {
248         pointBias.push_back(point.bias());
249         logPmfsum.push_back(point.logPmfSum());
250     }
251
252     constexpr int ulpTol = 10;
253
254     checker.checkSequence(props.begin(), props.end(), "Properties");
255     checker.setDefaultTolerance(absoluteTolerance(magnitude * GMX_DOUBLE_EPS * ulpTol));
256     checker.checkSequence(force.begin(), force.end(), "Force");
257     checker.checkSequence(pot.begin(), pot.end(), "Potential");
258     checker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpTol));
259     checker.checkSequence(pointBias.begin(), pointBias.end(), "PointBias");
260     checker.checkSequence(logPmfsum.begin(), logPmfsum.end(), "PointLogPmfsum");
261 }
262
263 /* Scan initial/final phase, MC/convolved force and update skip (not) allowed
264  * Both the convolving and skipping should not affect the bias and PMF.
265  * It would be nice if the test would explicitly check for this.
266  * Currently this is tested through identical reference data.
267  */
268 INSTANTIATE_TEST_CASE_P(WithParameters,
269                         BiasFepLambdaStateTest,
270                         ::testing::Combine(::testing::Values(eawhgrowthLINEAR, eawhgrowthEXP_LINEAR),
271                                            ::testing::Values(eawhpotentialUMBRELLA, eawhpotentialCONVOLVED),
272                                            ::testing::Values(BiasParams::DisableUpdateSkips::yes,
273                                                              BiasParams::DisableUpdateSkips::no)));
274
275 // Test that we detect coverings and exit the initial stage at the correct step
276 TEST(BiasFepLambdaStateTest, DetectsCovering)
277 {
278     const AwhFepLambdaStateTestParameters params =
279             getAwhFepLambdaTestParameters(eawhgrowthEXP_LINEAR, eawhpotentialCONVOLVED);
280
281     const double mdTimeStep = 0.1;
282
283     Bias bias(-1, params.awhParams, params.awhBiasParams, params.dimParams, params.beta, mdTimeStep,
284               1, "", Bias::ThisRankWillDoIO::No);
285
286     const int64_t exitStepRef = 320;
287
288     bool inInitialStage = bias.state().inInitialStage();
289
290     /* Some energies to use as base values (to which some noise is added later on). */
291     std::vector<double> neighborLambdaEnergies(numLambdaStates);
292     std::vector<double> neighborLambdaDhdl(numLambdaStates);
293     const double        magnitude = 12.0;
294     for (int i = 0; i < numLambdaStates; i++)
295     {
296         neighborLambdaEnergies[i] = magnitude * std::sin(i * 0.1);
297         neighborLambdaDhdl[i]     = magnitude * std::cos(i * 0.1);
298     }
299
300     int64_t step;
301     /* Normally this loop exits at exitStepRef, but we extend with failure */
302     for (step = 0; step <= 2 * exitStepRef; step++)
303     {
304         int      umbrellaGridpointIndex = bias.state().coordState().umbrellaGridpoint();
305         awh_dvec coordValue = { bias.getGridCoordValue(umbrellaGridpointIndex)[0], 0, 0, 0 };
306
307         double potential     = 0;
308         double potentialJump = 0;
309         bias.calcForceAndUpdateBias(coordValue, neighborLambdaEnergies, neighborLambdaDhdl,
310                                     &potential, &potentialJump, nullptr, nullptr, step, step,
311                                     params.awhParams.seed, nullptr);
312
313         inInitialStage = bias.state().inInitialStage();
314         if (!inInitialStage)
315         {
316             break;
317         }
318     }
319
320     EXPECT_EQ(false, inInitialStage);
321     if (!inInitialStage)
322     {
323         EXPECT_EQ(exitStepRef, step);
324     }
325 }
326
327 } // namespace test
328 } // namespace gmx