2 * This file is part of the GROMACS molecular simulation package.
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.
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.
37 #include "gromacs/applied_forces/awh/biasstate.h"
44 #include <gmock/gmock.h>
45 #include <gtest/gtest.h>
47 #include "gromacs/applied_forces/awh/biasgrid.h"
48 #include "gromacs/applied_forces/awh/pointstate.h"
49 #include "gromacs/math/functions.h"
50 #include "gromacs/mdtypes/awh_params.h"
51 #include "gromacs/utility/arrayref.h"
52 #include "gromacs/utility/smalloc.h"
54 #include "gromacs/applied_forces/awh/tests/awh_setup.h"
55 #include "testutils/testasserts.h"
56 #include "testutils/testfilemanager.h"
64 /*! \brief Test fixture for testing Bias updates
66 class BiasStateTest : public ::testing::TestWithParam<const char*>
69 std::unique_ptr<AwhTestParameters> params_;
72 std::unique_ptr<BiasState> biasState_; //!< The bias state
76 std::vector<std::vector<char>> awhDimParameters;
77 AwhCoordinateProviderType coordinateProvider = AwhCoordinateProviderType::Pull;
78 double diffusion = 0.1;
84 awhDimParameters.emplace_back(awhDimParamSerialized(
85 coordinateProvider, coordIndex, origin, end, period, diffusion));
92 awhDimParameters.emplace_back(awhDimParamSerialized(
93 coordinateProvider, coordIndex, origin, end, period, diffusion));
95 params_ = std::make_unique<AwhTestParameters>(getAwhTestParameters(
96 AwhHistogramGrowthType::Linear, AwhPotentialType::Convolved, awhDimParameters, 1, 1.0, false, 0.5, 0));
97 const AwhParams& awhParams = params_->awhParams;
98 const AwhBiasParams& awhBiasParams = awhParams.awhBiasParams()[0];
99 std::vector<DimParams> dimParams;
100 dimParams.push_back(DimParams::pullDimParams(1.0, 15.0, params_->beta));
101 dimParams.push_back(DimParams::pullDimParams(1.0, 15.0, params_->beta));
102 BiasGrid grid(dimParams, awhBiasParams.dimParams());
103 BiasParams biasParams(
104 awhParams, awhBiasParams, dimParams, 1.0, 1.0, BiasParams::DisableUpdateSkips::no, 1, grid.axis(), 0);
105 biasState_ = std::make_unique<BiasState>(awhBiasParams, 1.0, dimParams, grid);
107 // Here we initialize the grid point state using the input file
108 std::string filename = gmx::test::TestFileManager::getInputFilePath(GetParam());
109 biasState_->initGridPointState(
110 awhBiasParams, dimParams, grid, biasParams, filename, params_->awhParams.numBias());
114 TEST_P(BiasStateTest, InitializesFromFile)
116 gmx::ArrayRef<const PointState> points = biasState_->points();
118 /* Compute the mean square deviation from the expected values in the file.
119 * The PMF values are spaced by 0.5 per points and logPmfsum has opposite sign.
120 * The target is (index + 1)/120.
123 for (index i = 0; i < points.ssize(); i++)
125 msdPmf += gmx::square(points[i].logPmfSum() - points[0].logPmfSum() + 0.5 * i) / points.size();
126 EXPECT_DOUBLE_EQ(points[i].target(), (i + 1) / 120.0);
129 EXPECT_NEAR(0.0, msdPmf, 1e-31);
132 // Test that Bias initialization open and reads the correct initialization
133 // files and the correct PMF and target distribution is set.
134 INSTANTIATE_TEST_CASE_P(WithParameters,
136 ::testing::Values("pmf_target_format0.xvg", "pmf_target_format1.xvg"));