Replace compat::make_unique with std::make_unique
[gromacs.git] / src / gromacs / awh / tests / biasstate.cpp
blob06c871608d41ec9c5283f8db5252250f88722364
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2017,2018,2019, 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.
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.
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.
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.
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.
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.
35 #include "gmxpre.h"
37 #include "gromacs/awh/biasstate.h"
39 #include <cmath>
41 #include <memory>
42 #include <vector>
44 #include <gmock/gmock.h>
45 #include <gtest/gtest.h>
47 #include "gromacs/awh/grid.h"
48 #include "gromacs/awh/pointstate.h"
49 #include "gromacs/math/functions.h"
50 #include "gromacs/mdtypes/awh-params.h"
51 #include "gromacs/utility/smalloc.h"
53 #include "testutils/testasserts.h"
54 #include "testutils/testfilemanager.h"
56 namespace gmx
59 namespace test
62 /*! \internal \brief
63 * Struct that gathers all input for setting up and using a Bias
65 struct AwhTestParameters
67 double beta; //!< 1/(kB*T)
69 AwhDimParams awhDimParams[2]; //!< Dimension parameters pointed to by \p awhBiasParams
70 AwhBiasParams awhBiasParams; //!< Bias parameters pointed to by \[ awhParams
71 AwhParams awhParams; //!< AWH parameters, this is the struct to actually use
74 //! Helper function to set up the C-style AWH parameters for the test
75 static AwhTestParameters getAwhTestParameters()
77 AwhTestParameters params;
79 params.beta = 1.0;
81 AwhParams &awhParams = params.awhParams;
82 snew(params.awhParams.awhBiasParams, 1);
83 AwhBiasParams &awhBiasParams = params.awhParams.awhBiasParams[0];
84 snew(awhBiasParams.dimParams, 2);
86 AwhDimParams &awhDimParams0 = awhBiasParams.dimParams[0];
88 awhDimParams0.period = 0;
89 awhDimParams0.diffusion = 0.1;
90 awhDimParams0.origin = 0.5;
91 awhDimParams0.end = 1.5;
92 awhDimParams0.coordValueInit = awhDimParams0.origin;
93 awhDimParams0.coverDiameter = 0;
95 AwhDimParams &awhDimParams1 = awhBiasParams.dimParams[1];
97 awhDimParams1.period = 0;
98 awhDimParams1.diffusion = 0.1;
99 awhDimParams1.origin = 0.8;
100 awhDimParams1.end = 1.3;
101 awhDimParams1.coordValueInit = awhDimParams1.origin;
102 awhDimParams1.coverDiameter = 0;
104 awhBiasParams.ndim = 2;
105 awhBiasParams.eTarget = eawhtargetCONSTANT;
106 awhBiasParams.targetBetaScaling = 0;
107 awhBiasParams.targetCutoff = 0;
108 awhBiasParams.eGrowth = eawhgrowthLINEAR;
109 awhBiasParams.bUserData = TRUE;
110 awhBiasParams.errorInitial = 0.5;
111 awhBiasParams.shareGroup = 0;
112 awhBiasParams.equilibrateHistogram = FALSE;
114 awhParams.numBias = 1;
115 awhParams.seed = 93471803;
116 awhParams.nstOut = 0;
117 awhParams.nstSampleCoord = 1;
118 awhParams.numSamplesUpdateFreeEnergy = 10;
119 awhParams.ePotential = eawhpotentialCONVOLVED;
120 awhParams.shareBiasMultisim = FALSE;
122 return params;
125 /*! \brief Test fixture for testing Bias updates
127 class BiasStateTest : public ::testing::TestWithParam<const char *>
129 public:
130 std::unique_ptr<BiasState> biasState_; //!< The bias state
132 BiasStateTest()
134 AwhTestParameters params = getAwhTestParameters();
135 const AwhParams &awhParams = params.awhParams;
136 const AwhBiasParams &awhBiasParams = awhParams.awhBiasParams[0];
137 std::vector<DimParams> dimParams;
138 dimParams.emplace_back(1.0, 15.0, params.beta);
139 dimParams.emplace_back(1.0, 15.0, params.beta);
140 Grid grid(dimParams, awhBiasParams.dimParams);
141 BiasParams biasParams(awhParams, awhBiasParams, dimParams, 1.0, 1.0, BiasParams::DisableUpdateSkips::no, 1, grid.axis(), 0);
142 biasState_ = std::make_unique<BiasState>(awhBiasParams, 1.0, dimParams, grid);
144 // Here we initialize the grid point state using the input file
145 std::string filename = gmx::test::TestFileManager::getInputFilePath(GetParam());
146 biasState_->initGridPointState(awhBiasParams, dimParams, grid, biasParams, filename, params.awhParams.numBias);
148 sfree(params.awhParams.awhBiasParams[0].dimParams);
149 sfree(params.awhParams.awhBiasParams);
153 TEST_P(BiasStateTest, InitializesFromFile)
155 gmx::ArrayRef<const PointState> points = biasState_->points();
157 /* Compute the mean square deviation from the expected values in the file.
158 * The PMF values are spaced by 0.5 per points and logPmfsum has opposite sign.
159 * The target is (index + 1)/120.
161 double msdPmf = 0;
162 for (index i = 0; i < points.size(); i++)
164 msdPmf += gmx::square(points[i].logPmfSum() - points[0].logPmfSum() + 0.5*i) / points.size();
165 EXPECT_DOUBLE_EQ(points[i].target(), (i + 1)/120.0);
168 EXPECT_NEAR(0.0, msdPmf, 1e-31);
171 // Test that Bias initialization open and reads the correct initialization
172 // files and the the correct PMF and target distribution is set.
173 INSTANTIATE_TEST_CASE_P(WithParameters, BiasStateTest,
174 ::testing::Values("pmf_target_format0.xvg",
175 "pmf_target_format1.xvg"));
177 } // namespace test
178 } // namespace gmx