a6ea79624dd7ebdc240935f1a250d29214c2ad1b
[alexxy/gromacs.git] / src / gromacs / mdlib / tests / freeenergyparameters.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 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 /*! \internal \file
36  *
37  * \brief Tests routines in freeenergyparameters.h .
38  *
39  * \author Christian Blau <blau@kth.se>
40  */
41
42 #include "gmxpre.h"
43
44 #include "gromacs/mdlib/freeenergyparameters.h"
45
46 #include <gtest/gtest.h>
47
48 #include "gromacs/mdtypes/inputrec.h"
49
50 #include "testutils/testasserts.h"
51 #include "testutils/testmatchers.h"
52
53
54 namespace gmx
55 {
56 namespace test
57 {
58 namespace
59 {
60
61
62 /*! \brief Parameters that will vary from test to test.
63  */
64 struct FreeEnergyParameterTestParameters
65 {
66     //! current state of lambda in the simulation, -1 if not set
67     int currentLambdaState = -1;
68     //! Fractional value of lambda to start from, -1 if not set
69     double initLambda = -1;
70     //! The initial number of the state, -1 if not set
71     int initFepState = -1;
72     //! Change of lambda per time step (fraction of (0.1)
73     double deltaLambda = 0;
74     //! number of lambda entries
75     int nLambda = 0;
76     //! the current simulation step
77     int64_t step = 0;
78     //! the expected lambda at the current simulation step
79     std::array<real, efptNR> expectedLambdas = { -1, -1, -1, -1, -1, -1, -1 };
80 };
81
82
83 /*! \brief Sets of parameters on which to run the tests.
84  */
85 const FreeEnergyParameterTestParameters freeEnergyParameterSets[] = {
86     // no parameters set at all
87     { -1, -1, -1, 0, 1, 0, { -1, -1, -1, -1, -1, -1, -1 } },
88     // setting current lambda state to 0, no other variables set, using eftpNR * [0.8] as lambda state matrix
89     { 0, -1, -1, 0, 1, 1, { 0.8, 0.8, 0.8, 0.8, 0.8, 0.8, 0.8 } },
90     // setting current lambda state to 1, no other variables set, using eftpNR * [0.2,0.8] as lambda state matrix
91     { 1, -1, -1, 0, 2, 1, { 0.8, 0.8, 0.8, 0.8, 0.8, 0.8, 0.8 } },
92     // test that non-zero deltaLambda trumps current lambda state
93     // setting current lambda state to 1, using deltaLambda 0.1 and setp 0 and eftpNR * [0.2,0.8] as lambda state matrix
94     { 1, -1, -1, 0.1, 2, 0, { 0.2, 0.2, 0.2, 0.2, 0.2, 0.2, 0.2 } },
95     // test that interpolating lambda values starting
96     // from lambda = 0, deltaLambda = 0.1, step = 10 results in values at end-state
97     { 1, 0, -1, 0.1, 2, 10, { 0.8, 0.8, 0.8, 0.8, 0.8, 0.8, 0.8 } },
98     // interpolating half the way
99     { 1, 0, -1, 0.1, 2, 5, { 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5 } },
100     // starting from end-state 1 and move lambda half-way with negative deltaLambda = -0.1
101     { -1, -1, 1, -0.1, 2, 5, { 0.5, 0.5, 0.5, 0.5, 0.5, 0.5, 0.5 } },
102     // starting from end-state 1 and move one step backwards with negative deltaLambda = -0.1
103     { -1, -1, 1, -0.1, 2, 1, { 0.74, 0.74, 0.74, 0.74, 0.74, 0.74, 0.74 } },
104     // three lambda states, the last two equal ([0.2,0.8,0.8]), move forward with deltaLambda = 0.1
105     { -1, -1, 0, 0.1, 3, 0, { 0.2, 0.2, 0.2, 0.2, 0.2, 0.2, 0.2 } },
106     { -1, -1, 0, 0.1, 3, 3, { 0.56, 0.56, 0.56, 0.56, 0.56, 0.56, 0.56 } },
107     { -1, -1, 0, 0.1, 3, 7, { 0.8, 0.8, 0.8, 0.8, 0.8, 0.8, 0.8 } },
108     { -1, -1, 0, 0.1, 3, 8, { 0.8, 0.8, 0.8, 0.8, 0.8, 0.8, 0.8 } },
109     { -1, -1, 0, 0.1, 3, 10, { 0.8, 0.8, 0.8, 0.8, 0.8, 0.8, 0.8 } },
110     // three lambda states, the last two equal ([0.2,0.8,0.8]), move backwards
111     { -1, -1, 2, -0.1, 3, 1, { 0.8, 0.8, 0.8, 0.8, 0.8, 0.8, 0.8 } },
112     { -1, -1, 2, -0.1, 3, 2, { 0.8, 0.8, 0.8, 0.8, 0.8, 0.8, 0.8 } },
113     { -1, -1, 2, -0.1, 3, 3, { 0.8, 0.8, 0.8, 0.8, 0.8, 0.8, 0.8 } },
114     { -1, -1, 2, -0.1, 3, 5, { 0.8, 0.8, 0.8, 0.8, 0.8, 0.8, 0.8 } },
115     { -1, -1, 2, -0.1, 3, 7, { 0.56, 0.56, 0.56, 0.56, 0.56, 0.56, 0.56 } },
116     { -1, -1, 2, -0.1, 3, 8, { 0.44, 0.44, 0.44, 0.44, 0.44, 0.44, 0.44 } },
117     { -1, -1, 2, -0.1, 3, 10, { 0.2, 0.2, 0.2, 0.2, 0.2, 0.2, 0.2 } },
118     // three lambda states, the last two equal ([0.2,0.8,0.8]), start in middle state, move backwards
119     { -1, -1, 1, -0.1, 3, 0, { 0.8, 0.8, 0.8, 0.8, 0.8, 0.8, 0.8 } },
120     { -1, -1, 1, -0.1, 3, 2, { 0.56, 0.56, 0.56, 0.56, 0.56, 0.56, 0.56 } },
121     { -1, -1, 1, -0.1, 3, 3, { 0.44, 0.44, 0.44, 0.44, 0.44, 0.44, 0.44 } },
122 };
123
124 /*! \brief Test for setting free energy parameters.
125  */
126 class FreeEnergyParameterTest : public ::testing::TestWithParam<FreeEnergyParameterTestParameters>
127 {
128 public:
129     //! Fills in the required FEP values for testing from the test parameters
130     t_lambda getFepVals()
131     {
132         t_lambda fepvals;
133         fepvals.init_fep_state = GetParam().initFepState;
134         fepvals.init_lambda    = GetParam().initLambda;
135         fepvals.delta_lambda   = GetParam().deltaLambda;
136         fepvals.all_lambda     = getLambdaMatrix(GetParam().nLambda);
137         fepvals.n_lambda       = GetParam().nLambda;
138         return fepvals;
139     }
140
141     /*! \brief construct a lambda matrix
142      *
143      * \param[in] nLambda
144      * \returns nLambda * eftpNR matrix with pre-defined values
145      */
146     double** getLambdaMatrix(int nLambda)
147     {
148         for (int i = 0; i < efptNR; ++i)
149         {
150             allLambda_[i] = defaultLambdaArrayForTest_[nLambda].data();
151         }
152         return allLambda_.data();
153     }
154
155 private:
156     //! Construction aide for double ** matrix without snew
157     std::array<double*, efptNR> allLambda_;
158     //! a set of default lambda arrays for different lengths
159     std::vector<std::vector<double>> defaultLambdaArrayForTest_ = { {}, { 0.8 }, { 0.2, 0.8 }, { 0.2, 0.8, 0.8 } };
160 };
161
162 TEST_P(FreeEnergyParameterTest, CorrectLambdas)
163 {
164     EXPECT_THAT(GetParam().expectedLambdas,
165                 Pointwise(RealEq(defaultRealTolerance()),
166                           currentLambdas(GetParam().step, getFepVals(), GetParam().currentLambdaState)));
167 }
168
169 INSTANTIATE_TEST_CASE_P(WithParameters, FreeEnergyParameterTest, ::testing::ValuesIn(freeEnergyParameterSets));
170
171 } // namespace
172
173 } // namespace test
174
175 } // namespace gmx