Remove unnecessary includes of arrayref.h
[alexxy/gromacs.git] / src / gromacs / awh / tests / biasstate.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 "gromacs/awh/biasstate.h"
38
39 #include <cmath>
40
41 #include <memory>
42 #include <vector>
43
44 #include <gmock/gmock.h>
45 #include <gtest/gtest.h>
46
47 #include "gromacs/awh/biasgrid.h"
48 #include "gromacs/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"
53
54 #include "testutils/testasserts.h"
55 #include "testutils/testfilemanager.h"
56
57 namespace gmx
58 {
59
60 namespace test
61 {
62
63 /*! \internal \brief
64  * Struct that gathers all input for setting up and using a Bias
65  */
66 struct AwhTestParameters
67 {
68     double beta; //!< 1/(kB*T)
69
70     AwhDimParams  awhDimParams[2]; //!< Dimension parameters pointed to by \p awhBiasParams
71     AwhBiasParams awhBiasParams;   //!< Bias parameters pointed to by \[ awhParams
72     AwhParams     awhParams;       //!< AWH parameters, this is the struct to actually use
73 };
74
75 //! Helper function to set up the C-style AWH parameters for the test
76 static AwhTestParameters getAwhTestParameters()
77 {
78     AwhTestParameters params;
79
80     params.beta = 1.0;
81
82     AwhParams& awhParams = params.awhParams;
83     snew(params.awhParams.awhBiasParams, 1);
84     AwhBiasParams& awhBiasParams = params.awhParams.awhBiasParams[0];
85     snew(awhBiasParams.dimParams, 2);
86
87     AwhDimParams& awhDimParams0 = awhBiasParams.dimParams[0];
88
89     awhDimParams0.period         = 0;
90     awhDimParams0.diffusion      = 0.1;
91     awhDimParams0.origin         = 0.5;
92     awhDimParams0.end            = 1.5;
93     awhDimParams0.coordValueInit = awhDimParams0.origin;
94     awhDimParams0.coverDiameter  = 0;
95
96     AwhDimParams& awhDimParams1 = awhBiasParams.dimParams[1];
97
98     awhDimParams1.period         = 0;
99     awhDimParams1.diffusion      = 0.1;
100     awhDimParams1.origin         = 0.8;
101     awhDimParams1.end            = 1.3;
102     awhDimParams1.coordValueInit = awhDimParams1.origin;
103     awhDimParams1.coverDiameter  = 0;
104
105     awhBiasParams.ndim                 = 2;
106     awhBiasParams.eTarget              = eawhtargetCONSTANT;
107     awhBiasParams.targetBetaScaling    = 0;
108     awhBiasParams.targetCutoff         = 0;
109     awhBiasParams.eGrowth              = eawhgrowthLINEAR;
110     awhBiasParams.bUserData            = TRUE;
111     awhBiasParams.errorInitial         = 0.5;
112     awhBiasParams.shareGroup           = 0;
113     awhBiasParams.equilibrateHistogram = FALSE;
114
115     awhParams.numBias                    = 1;
116     awhParams.seed                       = 93471803;
117     awhParams.nstOut                     = 0;
118     awhParams.nstSampleCoord             = 1;
119     awhParams.numSamplesUpdateFreeEnergy = 10;
120     awhParams.ePotential                 = eawhpotentialCONVOLVED;
121     awhParams.shareBiasMultisim          = FALSE;
122
123     return params;
124 }
125
126 /*! \brief Test fixture for testing Bias updates
127  */
128 class BiasStateTest : public ::testing::TestWithParam<const char*>
129 {
130 public:
131     std::unique_ptr<BiasState> biasState_; //!< The bias state
132
133     BiasStateTest()
134     {
135         AwhTestParameters      params        = getAwhTestParameters();
136         const AwhParams&       awhParams     = params.awhParams;
137         const AwhBiasParams&   awhBiasParams = awhParams.awhBiasParams[0];
138         std::vector<DimParams> dimParams;
139         dimParams.emplace_back(1.0, 15.0, params.beta);
140         dimParams.emplace_back(1.0, 15.0, params.beta);
141         BiasGrid   grid(dimParams, awhBiasParams.dimParams);
142         BiasParams biasParams(awhParams, awhBiasParams, dimParams, 1.0, 1.0,
143                               BiasParams::DisableUpdateSkips::no, 1, grid.axis(), 0);
144         biasState_ = std::make_unique<BiasState>(awhBiasParams, 1.0, dimParams, grid);
145
146         // Here we initialize the grid point state using the input file
147         std::string filename = gmx::test::TestFileManager::getInputFilePath(GetParam());
148         biasState_->initGridPointState(awhBiasParams, dimParams, grid, biasParams, filename,
149                                        params.awhParams.numBias);
150
151         sfree(params.awhParams.awhBiasParams[0].dimParams);
152         sfree(params.awhParams.awhBiasParams);
153     }
154 };
155
156 TEST_P(BiasStateTest, InitializesFromFile)
157 {
158     gmx::ArrayRef<const PointState> points = biasState_->points();
159
160     /* Compute the mean square deviation from the expected values in the file.
161      * The PMF values are spaced by 0.5 per points and logPmfsum has opposite sign.
162      * The target is (index + 1)/120.
163      */
164     double msdPmf = 0;
165     for (index i = 0; i < points.ssize(); i++)
166     {
167         msdPmf += gmx::square(points[i].logPmfSum() - points[0].logPmfSum() + 0.5 * i) / points.size();
168         EXPECT_DOUBLE_EQ(points[i].target(), (i + 1) / 120.0);
169     }
170
171     EXPECT_NEAR(0.0, msdPmf, 1e-31);
172 }
173
174 // Test that Bias initialization open and reads the correct initialization
175 // files and the correct PMF and target distribution is set.
176 INSTANTIATE_TEST_CASE_P(WithParameters,
177                         BiasStateTest,
178                         ::testing::Values("pmf_target_format0.xvg", "pmf_target_format1.xvg"));
179
180 } // namespace test
181 } // namespace gmx