Make PME OpenCL enabled only for AMD devices
[alexxy/gromacs.git] / src / programs / mdrun / tests / pmetest.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016,2017,2018, 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  * This implements basic PME sanity tests.
38  * It runs the input system with PME for several steps (on CPU and GPU, if available),
39  * and checks the reciprocal and conserved energies.
40  * As part of mdrun-test, this will always run single rank PME simulation.
41  * As part of mdrun-mpi-test, this will run same as above when a single rank is requested,
42  * or a simulation with a single separate PME rank ("-npme 1") when multiple ranks are requested.
43  * \todo Extend and generalize this for more multi-rank tests (-npme 0, -npme 2, etc).
44  * \todo Implement death tests (e.g. for PME GPU decomposition).
45  *
46  * \author Aleksei Iupinov <a.yupinov@gmail.com>
47  * \ingroup module_mdrun_integration_tests
48  */
49 #include "gmxpre.h"
50
51 #include <map>
52 #include <string>
53 #include <vector>
54
55 #include <gtest/gtest-spi.h>
56
57 #include "gromacs/ewald/pme.h"
58 #include "gromacs/gpu_utils/gpu_utils.h"
59 #include "gromacs/hardware/detecthardware.h"
60 #include "gromacs/hardware/gpu_hw_info.h"
61 #include "gromacs/trajectory/energyframe.h"
62 #include "gromacs/utility/cstringutil.h"
63 #include "gromacs/utility/gmxmpi.h"
64 #include "gromacs/utility/loggerbuilder.h"
65 #include "gromacs/utility/physicalnodecommunicator.h"
66 #include "gromacs/utility/stringutil.h"
67
68 #include "testutils/mpitest.h"
69 #include "testutils/refdata.h"
70
71 #include "energyreader.h"
72 #include "moduletest.h"
73
74 namespace gmx
75 {
76 namespace test
77 {
78 namespace
79 {
80
81 /*! \brief A basic PME runner
82  *
83  * \todo Consider also using GpuTest class. */
84 class PmeTest : public MdrunTestFixture
85 {
86     public:
87         //! Before any test is run, work out whether any compatible GPUs exist.
88         static void SetUpTestCase();
89         //! Store whether any compatible GPUs exist.
90         static bool s_hasCompatibleGpus;
91         //! Convenience typedef
92         using RunModesList = std::map < std::string, std::vector < const char *>>;
93         //! Runs the test with the given inputs
94         void runTest(const RunModesList &runModes);
95 };
96
97 bool PmeTest::s_hasCompatibleGpus = false;
98
99 void PmeTest::SetUpTestCase()
100 {
101     gmx_gpu_info_t gpuInfo {};
102     // It would be nicer to do this detection once and have mdrun
103     // re-use it, but this is OK. Note that this also caters for when
104     // there is no GPU support in the build.
105     //
106     // TODO report any error messages gracefully.
107     if (canDetectGpus(nullptr))
108     {
109         findGpus(&gpuInfo);
110         s_hasCompatibleGpus = (gpuInfo.n_dev_compatible > 0);
111     }
112     free_gpu_info(&gpuInfo);
113 }
114
115 void PmeTest::runTest(const RunModesList &runModes)
116 {
117     const std::string inputFile = "spc-and-methanol";
118     runner_.useTopGroAndNdxFromDatabase(inputFile);
119
120     // With single rank we can and will always test PP+PME as part of mdrun-test.
121     // With multiple ranks we can additionally test a single PME-only rank within mdrun-mpi-test.
122     const bool parallelRun    = (getNumberOfTestMpiRanks() > 1);
123     const bool useSeparatePme = parallelRun;
124
125     EXPECT_EQ(0, runner_.callGrompp());
126
127     TestReferenceData    refData;
128     TestReferenceChecker rootChecker(refData.rootChecker());
129     const bool           thisRankChecks = (gmx_node_rank() == 0);
130     if (!thisRankChecks)
131     {
132         EXPECT_NONFATAL_FAILURE(rootChecker.checkUnusedEntries(), ""); // skip checks on other ranks
133     }
134
135     auto hardwareInfo_ = gmx_detect_hardware(MDLogger {},
136                                              PhysicalNodeCommunicator(MPI_COMM_WORLD,
137                                                                       gmx_physicalnode_id_hash()));
138
139     for (const auto &mode : runModes)
140     {
141         auto modeTargetsGpus = (mode.first.find("Gpu") != std::string::npos);
142         if (modeTargetsGpus && !s_hasCompatibleGpus)
143         {
144             // This run mode will cause a fatal error from mdrun when
145             // it can't find GPUs, which is not something we're trying
146             // to test here.
147             continue;
148         }
149         auto modeTargetsPmeOnGpus = (mode.first.find("PmeOnGpu") != std::string::npos);
150         if (modeTargetsPmeOnGpus && !pme_gpu_supports_build(*hardwareInfo_, nullptr))
151         {
152             // This run mode will cause a fatal error from mdrun when
153             // it finds an unsuitable device, which is not something
154             // we're trying to test here.
155             continue;
156         }
157
158         runner_.edrFileName_ = fileManager_.getTemporaryFilePath(inputFile + "_" + mode.first + ".edr");
159
160         CommandLine commandLine(mode.second);
161
162         const bool  usePmeTuning = (mode.first.find("Tune") != std::string::npos);
163         if (usePmeTuning)
164         {
165             commandLine.append("-tunepme");
166             commandLine.addOption("-nstlist", 1); // a new grid every step
167         }
168         else
169         {
170             commandLine.append("-notunepme"); // for reciprocal energy reproducibility
171         }
172         if (useSeparatePme)
173         {
174             commandLine.addOption("-npme", 1);
175         }
176
177         ASSERT_EQ(0, runner_.callMdrun(commandLine));
178
179         if (thisRankChecks)
180         {
181             auto energyReader      = openEnergyFileToReadFields(runner_.edrFileName_, {"Coul. recip.", "Total Energy", "Kinetic En."});
182             auto conservedChecker  = rootChecker.checkCompound("Energy", "Conserved");
183             auto reciprocalChecker = rootChecker.checkCompound("Energy", "Reciprocal");
184             bool firstIteration    = true;
185             while (energyReader->readNextFrame())
186             {
187                 const EnergyFrame &frame            = energyReader->frame();
188                 const std::string  stepName         = frame.frameName();
189                 const real         conservedEnergy  = frame.at("Total Energy");
190                 const real         reciprocalEnergy = frame.at("Coul. recip.");
191                 if (firstIteration)
192                 {
193                     // use first step values as references for tolerance
194                     const real startingKineticEnergy = frame.at("Kinetic En.");
195                     const auto conservedTolerance    = relativeToleranceAsFloatingPoint(startingKineticEnergy, 2e-5);
196                     const auto reciprocalTolerance   = relativeToleranceAsFloatingPoint(reciprocalEnergy, 3e-5);
197                     reciprocalChecker.setDefaultTolerance(reciprocalTolerance);
198                     conservedChecker.setDefaultTolerance(conservedTolerance);
199                     firstIteration = false;
200                 }
201                 conservedChecker.checkReal(conservedEnergy, stepName.c_str());
202                 if (!usePmeTuning) // with PME tuning come differing grids and differing reciprocal energy
203                 {
204                     reciprocalChecker.checkReal(reciprocalEnergy, stepName.c_str());
205                 }
206             }
207         }
208     }
209 }
210
211 TEST_F(PmeTest, ReproducesEnergies)
212 {
213     const int         nsteps     = 20;
214     const std::string theMdpFile = formatString("coulombtype     = PME\n"
215                                                 "nstcalcenergy   = 1\n"
216                                                 "nstenergy       = 1\n"
217                                                 "pme-order       = 4\n"
218                                                 "nsteps          = %d\n",
219                                                 nsteps
220                                                 );
221
222     runner_.useStringAsMdpFile(theMdpFile);
223
224     //TODO test all proper/improper combinations in more thorough way?
225     RunModesList runModes;
226     runModes["PmeOnCpu"]         = {"-pme", "cpu"};
227     runModes["PmeAuto"]          = {"-pme", "auto"};
228     runModes["PmeOnGpuFftOnCpu"] = {"-pme", "gpu", "-pmefft", "cpu"};
229     runModes["PmeOnGpuFftOnGpu"] = {"-pme", "gpu", "-pmefft", "gpu"};
230     runModes["PmeOnGpuFftAuto"]  = {"-pme", "gpu", "-pmefft", "auto"};
231     // same manual modes but marked for PME tuning
232     runModes["PmeOnCpuTune"]         = {"-pme", "cpu"};
233     runModes["PmeOnGpuFftOnCpuTune"] = {"-pme", "gpu", "-pmefft", "cpu"};
234     runModes["PmeOnGpuFftOnGpuTune"] = {"-pme", "gpu", "-pmefft", "gpu"};
235
236     runTest(runModes);
237 }
238
239 TEST_F(PmeTest, ScalesTheBox)
240 {
241     const int         nsteps     = 0;
242     const std::string theMdpFile = formatString("coulombtype     = PME\n"
243                                                 "nstcalcenergy   = 1\n"
244                                                 "nstenergy       = 1\n"
245                                                 "pme-order       = 4\n"
246                                                 "pbc             = xy\n"
247                                                 "nwall           = 2\n"
248                                                 "ewald-geometry  = 3dc\n"
249                                                 "wall_atomtype   = OMet CMet\n"
250                                                 "wall_density    = 9 9.0\n"
251                                                 "wall-ewald-zfac = 5\n"
252                                                 "nsteps          = %d\n",
253                                                 nsteps
254                                                 );
255
256     runner_.useStringAsMdpFile(theMdpFile);
257
258     RunModesList runModes;
259     runModes["PmeOnCpu"]         = {"-pme", "cpu"};
260     runModes["PmeOnGpuFftOnCpu"] = {"-pme", "gpu", "-pmefft", "cpu"};
261     runModes["PmeOnGpuFftOnGpu"] = {"-pme", "gpu", "-pmefft", "gpu"};
262
263     runTest(runModes);
264 }
265
266 }  // namespace
267 }  // namespace test
268 }  // namespace gmx