Update energy and trajectory frame types
[gromacs.git] / src / programs / mdrun / tests / initialconstraints.cpp
blob8441ef7a7344657725e592cbe0c385cc9c21f6a0
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 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.
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 * This implements basic initial constrains test (using single-rank mdrun)
39 * This test checks that the coordinates from file, which only satisfy
40 * the constraints up to gro precision, are constrained correctly and that
41 * the initial velocity of the center of mass does not contribute to the
42 * kinetic energy at step 0..
43 * It runs the input system for 1 step (no continuation), and compares the total energy.
45 * \author Aleksei Iupinov <a.yupinov@gmail.com>
46 * \ingroup module_mdrun_integration_tests
48 #include "gmxpre.h"
50 #include <string>
52 #include "gromacs/trajectory/energyframe.h"
53 #include "gromacs/utility/stringutil.h"
55 #include "energyreader.h"
56 #include "moduletest.h"
58 namespace gmx
60 namespace test
62 namespace
65 //! This type holds input integrators. Now it's holding names, but ei* enum values from md_enums.h could be used instead.
66 using EnergyIntegratorType = const char *;
68 //! Test fixture parametrized on integrators
69 class InitialConstraintsTest : public gmx::test::MdrunTestFixture,
70 public ::testing::WithParamInterface<EnergyIntegratorType>
74 TEST_P(InitialConstraintsTest, Works)
76 const int nsteps = 1;
77 const float timestep = 0.001;
78 auto integrator = GetParam();
79 const std::string integratorName(integrator);
80 SCOPED_TRACE("Integrating with " + integratorName);
81 const std::string theMdpFile = formatString("nstcalcenergy = 1\n"
82 "nstenergy = 1\n"
83 "comm-mode = linear\n"
84 "continuation = no\n"
85 "constraints = h-bonds\n"
86 "lincs_iter = 2\n"
87 "verlet-buffer-tolerance = 1e-4\n"
88 "nsttcouple = 1\n" // for md-vv-avek
89 "nstpcouple = 1\n" // for md-vv-avek
90 "integrator = %s\n"
91 "nsteps = %d\n"
92 "dt = %f\n",
93 integratorName.c_str(), nsteps, timestep);
95 runner_.useStringAsMdpFile(theMdpFile);
97 const std::string inputFile = "spc-and-methanol";
98 runner_.useTopGroAndNdxFromDatabase(inputFile.c_str());
99 EXPECT_EQ(0, runner_.callGrompp());
101 runner_.edrFileName_ = fileManager_.getTemporaryFilePath(inputFile + ".edr");
102 ASSERT_EQ(0, runner_.callMdrun());
104 auto energyReader = openEnergyFileToReadFields(runner_.edrFileName_, {"Total Energy", "Kinetic En."});
105 real totalEnergy = 0.0, prevTotalEnergy = 0.0;
106 auto tolerance = ulpTolerance(0); // The real value is set below from starting kinetic energy
107 for (int i = 0; i <= nsteps; i++)
109 EnergyFrame frame = energyReader->frame();
110 prevTotalEnergy = totalEnergy;
111 totalEnergy = frame.at("Total Energy");
112 if (i == 0)
114 // We set the tolerance for total energy based on magnitude of kinetic energy.
115 // The reason is that the other total energy component, the potential energy, can in theory have whatever magnitude.
116 const real startingKineticEnergy = frame.at("Kinetic En.");
117 tolerance = relativeToleranceAsFloatingPoint(startingKineticEnergy, 1e-5);
119 else
121 EXPECT_REAL_EQ_TOL(totalEnergy, prevTotalEnergy, tolerance);
126 //! Integrators with energy conservation to test
127 static const EnergyIntegratorType c_integratorsToTest [] = {"md", "md-vv", "md-vv-avek"};
129 INSTANTIATE_TEST_CASE_P(Checking, InitialConstraintsTest, ::testing::ValuesIn(c_integratorsToTest));