f69a45a650d741c3430a50b8383b3f737cfda3f5
[alexxy/gromacs.git] / src / gromacs / mdtypes / tests / multipletimestepping.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  * \brief
37  * Tests for the MultipleTimeStepping class and stand-alone functions.
38  *
39  * \author berk Hess <hess@kth.se>
40  * \ingroup module_mdtypes
41  */
42 #include "gmxpre.h"
43
44 #include "gromacs/mdtypes/multipletimestepping.h"
45
46 #include <gmock/gmock.h>
47 #include <gtest/gtest.h>
48
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/utility/gmxassert.h"
51 #include "gromacs/utility/smalloc.h"
52
53 #include "testutils/testasserts.h"
54
55 namespace gmx
56 {
57
58 namespace test
59 {
60
61 namespace
62 {
63
64 //! brief Sets up the MTS levels in \p ir and tests whether the number of errors matches \p numExpectedErrors
65 void setAndCheckMtsLevels(const GromppMtsOpts& mtsOpts, t_inputrec* ir, const int numExpectedErrors)
66 {
67     std::vector<std::string> errorMessages;
68     ir->useMts    = true;
69     ir->mtsLevels = setupMtsLevels(mtsOpts, &errorMessages);
70
71     if (haveValidMtsSetup(*ir))
72     {
73         std::vector<std::string> errorMessagesCheck = checkMtsRequirements(*ir);
74
75         // Concatenate the two lists with error messages
76         errorMessages.insert(errorMessages.end(), errorMessagesCheck.begin(), errorMessagesCheck.end());
77     }
78
79     EXPECT_EQ(errorMessages.size(), numExpectedErrors);
80 }
81
82 } // namespace
83
84 //! Checks that only numLevels = 2 does not produce an error
85 TEST(MultipleTimeStepping, ChecksNumLevels)
86 {
87     for (int numLevels = 0; numLevels <= 3; numLevels++)
88     {
89         GromppMtsOpts mtsOpts;
90         mtsOpts.numLevels    = numLevels;
91         mtsOpts.level2Factor = 2;
92
93         t_inputrec ir;
94
95         setAndCheckMtsLevels(mtsOpts, &ir, numLevels != 2 ? 1 : 0);
96     }
97 }
98
99 //! Test that each force group works
100 TEST(MultipleTimeStepping, SelectsForceGroups)
101 {
102     for (int forceGroupIndex = 0; forceGroupIndex < static_cast<int>(MtsForceGroups::Count);
103          forceGroupIndex++)
104     {
105         const MtsForceGroups forceGroup = static_cast<MtsForceGroups>(forceGroupIndex);
106         SCOPED_TRACE("Testing force group " + mtsForceGroupNames[forceGroup]);
107
108         GromppMtsOpts mtsOpts;
109         mtsOpts.numLevels    = 2;
110         mtsOpts.level2Forces = mtsForceGroupNames[forceGroup];
111         mtsOpts.level2Factor = 2;
112
113         t_inputrec ir;
114
115         setAndCheckMtsLevels(mtsOpts, &ir, 0);
116
117         EXPECT_EQ(ir.mtsLevels[1].forceGroups.count(), 1);
118         EXPECT_EQ(ir.mtsLevels[1].forceGroups[forceGroupIndex], true);
119     }
120 }
121
122 //! Checks that factor is checked
123 TEST(MultipleTimeStepping, ChecksStepFactor)
124 {
125     for (int stepFactor = 0; stepFactor <= 3; stepFactor++)
126     {
127         GromppMtsOpts mtsOpts;
128         mtsOpts.numLevels    = 2;
129         mtsOpts.level2Factor = stepFactor;
130
131         t_inputrec ir;
132
133         setAndCheckMtsLevels(mtsOpts, &ir, stepFactor < 2 ? 1 : 0);
134     }
135 }
136
137 namespace
138 {
139
140 GromppMtsOpts simpleMtsOpts()
141 {
142     GromppMtsOpts mtsOpts;
143     mtsOpts.numLevels    = 2;
144     mtsOpts.level2Forces = "nonbonded";
145     mtsOpts.level2Factor = 4;
146
147     return mtsOpts;
148 }
149
150 } // namespace
151
152 TEST(MultipleTimeStepping, ChecksPmeIsAtLastLevel)
153 {
154     const GromppMtsOpts mtsOpts = simpleMtsOpts();
155
156     t_inputrec ir;
157     ir.coulombtype = eelPME;
158
159     setAndCheckMtsLevels(mtsOpts, &ir, 1);
160 }
161
162 //! Test fixture base for parametrizing interval tests
163 using MtsIntervalTestParams = std::tuple<std::string, int>;
164 class MtsIntervalTest : public ::testing::Test, public ::testing::WithParamInterface<MtsIntervalTestParams>
165 {
166 public:
167     MtsIntervalTest()
168     {
169         const auto  params        = GetParam();
170         const auto& parameterName = std::get<0>(params);
171         const auto  interval      = std::get<1>(params);
172         numExpectedErrors_        = (interval == 4 ? 0 : 1);
173
174         if (parameterName == "nstcalcenergy")
175         {
176             ir_.nstcalcenergy = interval;
177         }
178         else if (parameterName == "nstenergy")
179         {
180             ir_.nstenergy = interval;
181         }
182         else if (parameterName == "nstfout")
183         {
184             ir_.nstfout = interval;
185         }
186         else if (parameterName == "nstlist")
187         {
188             ir_.nstlist = interval;
189         }
190         else if (parameterName == "nstdhdl")
191         {
192             ir_.efep             = efepYES;
193             ir_.fepvals->nstdhdl = interval;
194         }
195         else
196
197         {
198             GMX_RELEASE_ASSERT(false, "unknown parameter name");
199         }
200     }
201
202     t_inputrec ir_;
203     int        numExpectedErrors_;
204 };
205
206 TEST_P(MtsIntervalTest, Works)
207 {
208     const GromppMtsOpts mtsOpts = simpleMtsOpts();
209
210     setAndCheckMtsLevels(mtsOpts, &ir_, numExpectedErrors_);
211 }
212
213 INSTANTIATE_TEST_CASE_P(
214         ChecksStepInterval,
215         MtsIntervalTest,
216         ::testing::Combine(
217                 ::testing::Values("nstcalcenergy", "nstenergy", "nstfout", "nstlist", "nstdhdl"),
218                 ::testing::Values(3, 4, 5)));
219
220 // Check that correct input does not produce errors
221 TEST(MultipleTimeStepping, ChecksIntegrator)
222 {
223     const GromppMtsOpts mtsOpts = simpleMtsOpts();
224
225     t_inputrec ir;
226     ir.eI = eiBD;
227
228     setAndCheckMtsLevels(mtsOpts, &ir, 1);
229 }
230
231 } // namespace test
232 } // namespace gmx