57ae23e152f422767b151c97c2ca13806a624c37
[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, 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 "config.h"
52
53 #include <map>
54 #include <string>
55 #include <vector>
56
57 #include <gtest/gtest-spi.h>
58
59 #include "gromacs/gpu_utils/gpu_utils.h"
60 #include "gromacs/hardware/gpu_hw_info.h"
61 #include "gromacs/utility/cstringutil.h"
62 #include "gromacs/utility/gmxmpi.h"
63 #include "gromacs/utility/stringutil.h"
64
65 #include "testutils/mpitest.h"
66 #include "testutils/refdata.h"
67
68 #include "energyreader.h"
69 #include "moduletest.h"
70
71 namespace gmx
72 {
73 namespace test
74 {
75 namespace
76 {
77
78 //! A basic PME runner
79 class PmeTest : public MdrunTestFixture
80 {
81     public:
82         //! Before any test is run, work out whether any compatible GPUs exist.
83         static void SetUpTestCase();
84         //! Store whether any compatible GPUs exist.
85         static bool s_hasCompatibleCudaGpus;
86 };
87
88 bool PmeTest::s_hasCompatibleCudaGpus = false;
89
90 void PmeTest::SetUpTestCase()
91 {
92     gmx_gpu_info_t gpuInfo {};
93     char           detection_error[STRLEN];
94     GMX_UNUSED_VALUE(detection_error); //TODO
95     // It would be nicer to do this detection once and have mdrun
96     // re-use it, but this is OK. Note that this also caters for when
97     // there is no GPU support in the build.
98     if (GMX_GPU == GMX_GPU_CUDA &&
99         (detect_gpus(&gpuInfo, detection_error) >= 0) &&
100         gpuInfo.n_dev_compatible > 0)
101     {
102         s_hasCompatibleCudaGpus = true;
103     }
104     free_gpu_info(&gpuInfo);
105 }
106
107 TEST_F(PmeTest, ReproducesEnergies)
108 {
109     const int   nsteps     = 20;
110     std::string theMdpFile = formatString("coulombtype     = PME\n"
111                                           "nstcalcenergy   = 1\n"
112                                           "nstenergy       = 1\n"
113                                           "pme-order       = 4\n"
114                                           "nsteps          = %d\n",
115                                           nsteps);
116
117     runner_.useStringAsMdpFile(theMdpFile);
118
119     const std::string inputFile = "spc-and-methanol";
120     runner_.useTopGroAndNdxFromDatabase(inputFile.c_str());
121
122     // With single rank we can and will always test PP+PME as part of mdrun-test.
123     // With multiple ranks we can additionally test a single PME-only rank within mdrun-mpi-test.
124     const bool parallelRun    = (getNumberOfTestMpiRanks() > 1);
125     const bool useSeparatePme = parallelRun;
126
127     EXPECT_EQ(0, runner_.callGrompp());
128
129     //TODO test all proper/improper combinations in more thorough way?
130     std::map < std::string, std::vector < const char *>> runModes;
131     runModes["PmeOnCpu"]         = {"-pme", "cpu"};
132     runModes["PmeAuto"]          = {"-pme", "auto"};
133     runModes["PmeOnGpuFftOnCpu"] = {"-pme", "gpu", "-pmefft", "cpu"};
134     runModes["PmeOnGpuFftOnGpu"] = {"-pme", "gpu", "-pmefft", "gpu"};
135     runModes["PmeOnGpuFftAuto"]  = {"-pme", "gpu", "-pmefft", "auto"};
136     // same manual modes but marked for PME tuning
137     runModes["PmeOnCpuTune"]         = {"-pme", "cpu"};
138     runModes["PmeOnGpuFftOnCpuTune"] = {"-pme", "gpu", "-pmefft", "cpu"};
139     runModes["PmeOnGpuFftOnGpuTune"] = {"-pme", "gpu", "-pmefft", "gpu"};
140     TestReferenceData    refData;
141     TestReferenceChecker rootChecker(refData.rootChecker());
142     const bool           thisRankChecks = (gmx_node_rank() == 0);
143     if (!thisRankChecks)
144     {
145         EXPECT_NONFATAL_FAILURE(rootChecker.checkUnusedEntries(), ""); // skip checks on other ranks
146     }
147     for (const auto &mode : runModes)
148     {
149         auto modeTargetsGpus = (mode.first.find("Gpu") != std::string::npos);
150         if (modeTargetsGpus && !s_hasCompatibleCudaGpus)
151         {
152             // This run mode will cause a fatal error from mdrun when
153             // it can't find GPUs, which is not something we're trying
154             // 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         ASSERT_EQ(0, runner_.callMdrun(commandLine));
177
178         if (thisRankChecks)
179         {
180             auto energyReader      = openEnergyFileToReadFields(runner_.edrFileName_, {"Coul. recip.", "Total Energy", "Kinetic En."});
181             auto conservedChecker  = rootChecker.checkCompound("Energy", "Conserved");
182             auto reciprocalChecker = rootChecker.checkCompound("Energy", "Reciprocal");
183             for (int i = 0; i <= nsteps; i++)
184             {
185                 EnergyFrame frame            = energyReader->frame();
186                 std::string stepNum          = gmx::formatString("%d", i);
187                 const real  conservedEnergy  = frame.at("Total Energy");
188                 const real  reciprocalEnergy = frame.at("Coul. recip.");
189                 if (i == 0)
190                 {
191                     // use first step values as references for tolerance
192                     const real startingKineticEnergy = frame.at("Kinetic En.");
193                     const auto conservedTolerance    = relativeToleranceAsFloatingPoint(startingKineticEnergy, 2e-5);
194                     const auto reciprocalTolerance   = relativeToleranceAsFloatingPoint(reciprocalEnergy, 3e-5);
195                     reciprocalChecker.setDefaultTolerance(reciprocalTolerance);
196                     conservedChecker.setDefaultTolerance(conservedTolerance);
197                 }
198                 conservedChecker.checkReal(conservedEnergy, stepNum.c_str());
199                 if (!usePmeTuning) // with PME tuning come differing grids and differing reciprocal energy
200                 {
201                     reciprocalChecker.checkReal(reciprocalEnergy, stepNum.c_str());
202                 }
203             }
204         }
205         // FIXME: without this barrier, one of the mdruns was somehow having a non-PME inputrec (!)
206 #if GMX_LIB_MPI
207         if (parallelRun)
208         {
209             MPI_Barrier(MPI_COMM_WORLD);
210         }
211 #endif
212     }
213
214     // This is a workaround for the output files to not be deleted in a parallel run.
215     // TODO: consider moving the barrier into the file manager destructor.
216 #if GMX_LIB_MPI
217     if (parallelRun)
218     {
219         MPI_Barrier(MPI_COMM_WORLD);
220     }
221 #endif
222 }
223
224 }
225 }
226 }