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