Update instructions in containers.rst
[gromacs.git] / src / gromacs / mdtypes / tests / multipletimestepping.cpp
blobf69a45a650d741c3430a50b8383b3f737cfda3f5
1 /*
2 * This file is part of the GROMACS molecular simulation package.
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.
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 /*! \internal \file
36 * \brief
37 * Tests for the MultipleTimeStepping class and stand-alone functions.
39 * \author berk Hess <hess@kth.se>
40 * \ingroup module_mdtypes
42 #include "gmxpre.h"
44 #include "gromacs/mdtypes/multipletimestepping.h"
46 #include <gmock/gmock.h>
47 #include <gtest/gtest.h>
49 #include "gromacs/mdtypes/inputrec.h"
50 #include "gromacs/utility/gmxassert.h"
51 #include "gromacs/utility/smalloc.h"
53 #include "testutils/testasserts.h"
55 namespace gmx
58 namespace test
61 namespace
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)
67 std::vector<std::string> errorMessages;
68 ir->useMts = true;
69 ir->mtsLevels = setupMtsLevels(mtsOpts, &errorMessages);
71 if (haveValidMtsSetup(*ir))
73 std::vector<std::string> errorMessagesCheck = checkMtsRequirements(*ir);
75 // Concatenate the two lists with error messages
76 errorMessages.insert(errorMessages.end(), errorMessagesCheck.begin(), errorMessagesCheck.end());
79 EXPECT_EQ(errorMessages.size(), numExpectedErrors);
82 } // namespace
84 //! Checks that only numLevels = 2 does not produce an error
85 TEST(MultipleTimeStepping, ChecksNumLevels)
87 for (int numLevels = 0; numLevels <= 3; numLevels++)
89 GromppMtsOpts mtsOpts;
90 mtsOpts.numLevels = numLevels;
91 mtsOpts.level2Factor = 2;
93 t_inputrec ir;
95 setAndCheckMtsLevels(mtsOpts, &ir, numLevels != 2 ? 1 : 0);
99 //! Test that each force group works
100 TEST(MultipleTimeStepping, SelectsForceGroups)
102 for (int forceGroupIndex = 0; forceGroupIndex < static_cast<int>(MtsForceGroups::Count);
103 forceGroupIndex++)
105 const MtsForceGroups forceGroup = static_cast<MtsForceGroups>(forceGroupIndex);
106 SCOPED_TRACE("Testing force group " + mtsForceGroupNames[forceGroup]);
108 GromppMtsOpts mtsOpts;
109 mtsOpts.numLevels = 2;
110 mtsOpts.level2Forces = mtsForceGroupNames[forceGroup];
111 mtsOpts.level2Factor = 2;
113 t_inputrec ir;
115 setAndCheckMtsLevels(mtsOpts, &ir, 0);
117 EXPECT_EQ(ir.mtsLevels[1].forceGroups.count(), 1);
118 EXPECT_EQ(ir.mtsLevels[1].forceGroups[forceGroupIndex], true);
122 //! Checks that factor is checked
123 TEST(MultipleTimeStepping, ChecksStepFactor)
125 for (int stepFactor = 0; stepFactor <= 3; stepFactor++)
127 GromppMtsOpts mtsOpts;
128 mtsOpts.numLevels = 2;
129 mtsOpts.level2Factor = stepFactor;
131 t_inputrec ir;
133 setAndCheckMtsLevels(mtsOpts, &ir, stepFactor < 2 ? 1 : 0);
137 namespace
140 GromppMtsOpts simpleMtsOpts()
142 GromppMtsOpts mtsOpts;
143 mtsOpts.numLevels = 2;
144 mtsOpts.level2Forces = "nonbonded";
145 mtsOpts.level2Factor = 4;
147 return mtsOpts;
150 } // namespace
152 TEST(MultipleTimeStepping, ChecksPmeIsAtLastLevel)
154 const GromppMtsOpts mtsOpts = simpleMtsOpts();
156 t_inputrec ir;
157 ir.coulombtype = eelPME;
159 setAndCheckMtsLevels(mtsOpts, &ir, 1);
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>
166 public:
167 MtsIntervalTest()
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);
174 if (parameterName == "nstcalcenergy")
176 ir_.nstcalcenergy = interval;
178 else if (parameterName == "nstenergy")
180 ir_.nstenergy = interval;
182 else if (parameterName == "nstfout")
184 ir_.nstfout = interval;
186 else if (parameterName == "nstlist")
188 ir_.nstlist = interval;
190 else if (parameterName == "nstdhdl")
192 ir_.efep = efepYES;
193 ir_.fepvals->nstdhdl = interval;
195 else
198 GMX_RELEASE_ASSERT(false, "unknown parameter name");
202 t_inputrec ir_;
203 int numExpectedErrors_;
206 TEST_P(MtsIntervalTest, Works)
208 const GromppMtsOpts mtsOpts = simpleMtsOpts();
210 setAndCheckMtsLevels(mtsOpts, &ir_, numExpectedErrors_);
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)));
220 // Check that correct input does not produce errors
221 TEST(MultipleTimeStepping, ChecksIntegrator)
223 const GromppMtsOpts mtsOpts = simpleMtsOpts();
225 t_inputrec ir;
226 ir.eI = eiBD;
228 setAndCheckMtsLevels(mtsOpts, &ir, 1);
231 } // namespace test
232 } // namespace gmx