Disable awh-potential=convolved when using AWH with pull and FEP
[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 "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 only eawhpotentialUMBRELLA (MC) for FEP lambda dimensions)
160          *   disableUpdateSkips (should not affect the results):
161          *     BiasParams::DisableUpdateSkips::yes: update the point state for every sample
162          *     BiasParams::DisableUpdateSkips::no:  update the point state at an interval > 1 sample
163          *
164          * Note: It would be nice to explicitly check that eawhpotential
165          *       and disableUpdateSkips do not affect the point state.
166          *       But the reference data will also ensure this.
167          */
168         int                            eawhgrowth;
169         int                            eawhpotential;
170         BiasParams::DisableUpdateSkips disableUpdateSkips;
171         std::tie(eawhgrowth, eawhpotential, disableUpdateSkips) = GetParam();
172
173         /* Set up a basic AWH setup with a single, 1D bias with parameters
174          * such that we can measure the effects of different parameters.
175          */
176         const AwhFepLambdaStateTestParameters params =
177                 getAwhFepLambdaTestParameters(eawhgrowth, eawhpotential);
178
179         seed_ = params.awhParams.seed;
180
181         double mdTimeStep = 0.1;
182
183         bias_ = std::make_unique<Bias>(-1, params.awhParams, params.awhBiasParams, params.dimParams,
184                                        params.beta, mdTimeStep, 1, "", Bias::ThisRankWillDoIO::No,
185                                        disableUpdateSkips);
186     }
187 };
188
189 TEST_P(BiasFepLambdaStateTest, ForcesBiasPmf)
190 {
191     gmx::test::TestReferenceData    data;
192     gmx::test::TestReferenceChecker checker(data.rootChecker());
193
194     Bias& bias = *bias_;
195
196     /* Make strings with the properties we expect to be different in the tests.
197      * These also helps to interpret the reference data.
198      */
199     std::vector<std::string> props;
200     props.push_back(formatString("stage:           %s", bias.state().inInitialStage() ? "initial" : "final"));
201     props.push_back(formatString("convolve forces: %s", bias.params().convolveForce ? "yes" : "no"));
202     props.push_back(formatString("skip updates:    %s", bias.params().skipUpdates() ? "yes" : "no"));
203
204     SCOPED_TRACE(gmx::formatString("%s, %s, %s", props[0].c_str(), props[1].c_str(), props[2].c_str()));
205
206     std::vector<double> force, pot;
207
208     double potentialJump = 0;
209     double mdTimeStep    = 0.1;
210     int    nSteps        = 501;
211
212     /* Some energies to use as base values (to which some noise is added later on). */
213     std::vector<double> neighborLambdaEnergies(numLambdaStates);
214     std::vector<double> neighborLambdaDhdl(numLambdaStates);
215     const double        magnitude = 12.0;
216     for (int i = 0; i < numLambdaStates; i++)
217     {
218         neighborLambdaEnergies[i] = magnitude * std::sin(i * 0.1);
219         neighborLambdaDhdl[i]     = magnitude * std::cos(i * 0.1);
220     }
221
222     for (int step = 0; step < nSteps; step++)
223     {
224         int      umbrellaGridpointIndex = bias.state().coordState().umbrellaGridpoint();
225         awh_dvec coordValue = { bias.getGridCoordValue(umbrellaGridpointIndex)[0], 0, 0, 0 };
226         double   potential  = 0;
227         gmx::ArrayRef<const double> biasForce = bias.calcForceAndUpdateBias(
228                 coordValue, neighborLambdaEnergies, neighborLambdaDhdl, &potential, &potentialJump,
229                 nullptr, nullptr, step * mdTimeStep, step, seed_, nullptr);
230
231         force.push_back(biasForce[0]);
232         pot.push_back(potential);
233     }
234
235     /* When skipping updates, ensure all skipped updates are performed here.
236      * This should result in the same bias state as at output in a normal run.
237      */
238     if (bias.params().skipUpdates())
239     {
240         bias.doSkippedUpdatesForAllPoints();
241     }
242
243     std::vector<double> pointBias, logPmfsum;
244     for (auto& point : bias.state().points())
245     {
246         pointBias.push_back(point.bias());
247         logPmfsum.push_back(point.logPmfSum());
248     }
249
250     constexpr int ulpTol = 10;
251
252     checker.checkSequence(props.begin(), props.end(), "Properties");
253     checker.setDefaultTolerance(absoluteTolerance(magnitude * GMX_DOUBLE_EPS * ulpTol));
254     checker.checkSequence(force.begin(), force.end(), "Force");
255     checker.checkSequence(pot.begin(), pot.end(), "Potential");
256     checker.setDefaultTolerance(relativeToleranceAsUlp(1.0, ulpTol));
257     checker.checkSequence(pointBias.begin(), pointBias.end(), "PointBias");
258     checker.checkSequence(logPmfsum.begin(), logPmfsum.end(), "PointLogPmfsum");
259 }
260
261 /* Scan initial/final phase, MC/convolved force and update skip (not) allowed
262  * Both the convolving and skipping should not affect the bias and PMF.
263  * It would be nice if the test would explicitly check for this.
264  * Currently this is tested through identical reference data.
265  */
266 INSTANTIATE_TEST_CASE_P(WithParameters,
267                         BiasFepLambdaStateTest,
268                         ::testing::Combine(::testing::Values(eawhgrowthLINEAR, eawhgrowthEXP_LINEAR),
269                                            ::testing::Values(eawhpotentialUMBRELLA),
270                                            ::testing::Values(BiasParams::DisableUpdateSkips::yes,
271                                                              BiasParams::DisableUpdateSkips::no)));
272
273 // Test that we detect coverings and exit the initial stage at the correct step
274 TEST(BiasFepLambdaStateTest, DetectsCovering)
275 {
276     const AwhFepLambdaStateTestParameters params =
277             getAwhFepLambdaTestParameters(eawhgrowthEXP_LINEAR, eawhpotentialUMBRELLA);
278
279     const double mdTimeStep = 0.1;
280
281     Bias bias(-1, params.awhParams, params.awhBiasParams, params.dimParams, params.beta, mdTimeStep,
282               1, "", Bias::ThisRankWillDoIO::No);
283
284     const int64_t exitStepRef = 320;
285
286     bool inInitialStage = bias.state().inInitialStage();
287
288     /* Some energies to use as base values (to which some noise is added later on). */
289     std::vector<double> neighborLambdaEnergies(numLambdaStates);
290     std::vector<double> neighborLambdaDhdl(numLambdaStates);
291     const double        magnitude = 12.0;
292     for (int i = 0; i < numLambdaStates; i++)
293     {
294         neighborLambdaEnergies[i] = magnitude * std::sin(i * 0.1);
295         neighborLambdaDhdl[i]     = magnitude * std::cos(i * 0.1);
296     }
297
298     int64_t step;
299     /* Normally this loop exits at exitStepRef, but we extend with failure */
300     for (step = 0; step <= 2 * exitStepRef; step++)
301     {
302         int      umbrellaGridpointIndex = bias.state().coordState().umbrellaGridpoint();
303         awh_dvec coordValue = { bias.getGridCoordValue(umbrellaGridpointIndex)[0], 0, 0, 0 };
304
305         double potential     = 0;
306         double potentialJump = 0;
307         bias.calcForceAndUpdateBias(coordValue, neighborLambdaEnergies, neighborLambdaDhdl,
308                                     &potential, &potentialJump, nullptr, nullptr, step, step,
309                                     params.awhParams.seed, nullptr);
310
311         inInitialStage = bias.state().inInitialStage();
312         if (!inInitialStage)
313         {
314             break;
315         }
316     }
317
318     EXPECT_EQ(false, inInitialStage);
319     if (!inInitialStage)
320     {
321         EXPECT_EQ(exitStepRef, step);
322     }
323 }
324
325 } // namespace test
326 } // namespace gmx