2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2018,2019, 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.
38 * Tests for the energy minimization functionality.
40 * \author Mark Abraham <mark.j.abraham@gmail.com>
41 * \ingroup module_mdrun_integration_tests
51 #include <gtest/gtest.h>
53 #include "gromacs/options/filenameoption.h"
54 #include "gromacs/topology/idef.h"
55 #include "gromacs/topology/ifunc.h"
56 #include "gromacs/trajectory/energyframe.h"
57 #include "gromacs/trajectory/trajectoryframe.h"
58 #include "gromacs/utility/basenetwork.h"
59 #include "gromacs/utility/stringutil.h"
61 #include "testutils/mpitest.h"
62 #include "testutils/refdata.h"
63 #include "testutils/simulationdatabase.h"
64 #include "testutils/testasserts.h"
66 #include "energycomparison.h"
67 #include "energyreader.h"
68 #include "moduletest.h"
77 /*! \brief Test fixture base for energy minimizaiton
79 * This test ensures mdrun can run an energy minimization, reaching
80 * a reproducible final energy.
82 * The choices for tolerance are arbitrary but sufficient. */
83 class EnergyMinimizationTest
:
84 public MdrunTestFixture
,
85 public ::testing::WithParamInterface
<std::tuple
<std::string
, std::string
>>
89 /*! \brief Database of empirical tolerances for EM integrators on the various systems. */
90 std::unordered_map
<std::string
, FloatingPointTolerance
> potentialEnergyToleranceForSystem_g
= {
91 { { "argon12", relativeToleranceAsPrecisionDependentUlp(-1, 10, 200) },
92 { "tip3p5", relativeToleranceAsPrecisionDependentUlp(-50, 200, 3800) },
93 { "glycine_vacuo", relativeToleranceAsPrecisionDependentUlp(1000, 100, 100) },
94 { "alanine_vsite_vacuo", relativeToleranceAsPrecisionDependentUlp(-160, 150, 400) },
95 { "glycine_no_constraints_vacuo", relativeToleranceAsPrecisionDependentUlp(2000, 100, 100) } }
98 TEST_P(EnergyMinimizationTest
, WithinTolerances
)
100 auto params
= GetParam();
101 auto simulationName
= std::get
<0>(params
);
102 auto minimizer
= std::get
<1>(params
);
103 SCOPED_TRACE(formatString("Comparing '%s' energy minimization for simulation '%s'",
104 minimizer
.c_str(), simulationName
.c_str()));
106 // TODO At some point we should also test PME-only ranks.
107 int numRanksAvailable
= getNumberOfTestMpiRanks();
108 if (!isNumberOfPpRanksSupported(simulationName
, numRanksAvailable
))
111 "Test system '%s' cannot run with %d ranks.\n"
112 "The supported numbers are: %s\n",
113 simulationName
.c_str(), numRanksAvailable
,
114 reportNumbersOfPpRanksSupported(simulationName
).c_str());
118 auto mdpFieldValues
=
119 prepareMdpFieldValues(simulationName
.c_str(), minimizer
.c_str(), "no", "no");
120 mdpFieldValues
["nsteps"] = "4";
122 int maxWarningsTolerated
= (minimizer
== "l-bfgs") ? 1 : 0;
123 // prepare the .tpr file
125 // TODO evolve grompp to report the number of warnings issued, so
126 // tests always expect the right number.
128 caller
.append("grompp");
129 caller
.addOption("-maxwarn", maxWarningsTolerated
);
130 runner_
.useTopGroAndNdxFromDatabase(simulationName
);
131 runner_
.useStringAsMdpFile(prepareMdpFileContents(mdpFieldValues
));
132 EXPECT_EQ(0, runner_
.callGrompp(caller
));
135 // do mdrun, preparing to check the energies later
136 runner_
.edrFileName_
= fileManager_
.getTemporaryFilePath("minimize.edr");
138 CommandLine mdrunCaller
;
139 mdrunCaller
.append("mdrun");
140 if (minimizer
== "l-bfgs" && getNumberOfTestMpiRanks() > 1)
142 // Ideally we would use this death test, but it is not
143 // stable enough in Jenkins, so we just skip it.
144 // EXPECT_DEATH_IF_SUPPORTED(runner_.callMdrun(mdrunCaller),
145 // "L-BFGS minimization only supports a single rank");
150 ASSERT_EQ(0, runner_
.callMdrun(mdrunCaller
));
154 EnergyTermsToCompare energyTermsToCompare
{ {
155 { interaction_function
[F_EPOT
].longname
, potentialEnergyToleranceForSystem_g
.at(simulationName
) },
158 TestReferenceData refData
;
159 auto checker
= refData
.rootChecker()
160 .checkCompound("Simulation", simulationName
)
161 .checkCompound("Minimizer", minimizer
);
162 checkEnergiesAgainstReferenceData(runner_
.edrFileName_
, energyTermsToCompare
, &checker
);
165 //! Containers of systems and integrators to test.
167 std::vector
<std::string
> unconstrainedSystemsToTest_g
= { "argon12",
168 "glycine_no_constraints_vacuo" };
169 std::vector
<std::string
> minimizersToTest_g
= { "steep", "cg", "l-bfgs" };
171 std::vector
<std::string
> constrainedSystemsToTest_g
= { "tip3p5", "glycine_vacuo",
172 "alanine_vsite_vacuo" };
173 std::vector
<std::string
> minimizersToTestWithConstraints_g
= { "steep", "cg" };
176 // The time for OpenCL kernel compilation means these tests might time
177 // out. If that proves to be a problem, these can be disabled for
178 // OpenCL builds. However, once that compilation is cached for the
179 // lifetime of the whole test binary process, these tests should run in
180 // such configurations.
181 INSTANTIATE_TEST_CASE_P(MinimizersWorkWithConstraints
,
182 EnergyMinimizationTest
,
183 ::testing::Combine(::testing::ValuesIn(constrainedSystemsToTest_g
),
184 ::testing::ValuesIn(minimizersToTestWithConstraints_g
)));
185 INSTANTIATE_TEST_CASE_P(MinimizersWork
,
186 EnergyMinimizationTest
,
187 ::testing::Combine(::testing::ValuesIn(unconstrainedSystemsToTest_g
),
188 ::testing::ValuesIn(minimizersToTest_g
)));