Reimplement constant acceleration groups
[alexxy/gromacs.git] / src / programs / mdrun / tests / constantacceleration.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2021, 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
36 /*! \internal \file
37  * \brief End-to-end tests checking sanity of results of simulations
38  *        containing constant acceleration groups
39  *
40  * \author Berk Hess <hess@kth.se>
41  * \ingroup module_mdrun_integration_tests
42  */
43 #include "gmxpre.h"
44
45 #include "config.h"
46
47 #include "gromacs/topology/ifunc.h"
48 #include "gromacs/utility/stringutil.h"
49
50 #include "testutils/mpitest.h"
51 #include "testutils/simulationdatabase.h"
52 #include "testutils/testmatchers.h"
53
54 #include "moduletest.h"
55 #include "simulatorcomparison.h"
56 #include "trajectoryreader.h"
57
58 namespace gmx::test
59 {
60 namespace
61 {
62 /*! \brief Test fixture checking the velocities of aomts
63  *
64  * This tests that velocities of non-accelerated atoms are zero
65  * and that velocities of accelarated atoms match acceleration*time.
66  * This code assumes the first half the atoms are non-accelerated
67  * and the second half are accelerated
68  */
69 using AccelerationGroupTestParams = std::tuple<std::string, std::string>;
70 class AccelerationGroupTest :
71     public MdrunTestFixture,
72     public ::testing::WithParamInterface<AccelerationGroupTestParams>
73 {
74 public:
75     //! Check that the velocities are zero or accelerated
76     static void checkVelocities(const std::string&            trajectoryName,
77                                 const RVec                    acceleration,
78                                 const FloatingPointTolerance& tolerance)
79     {
80         const size_t c_groupSize = 3;
81
82         const std::vector<RVec> zeroVelocities(c_groupSize, RVec{ 0, 0, 0 });
83
84         TrajectoryFrameReader trajectoryFrameReader(trajectoryName);
85         while (trajectoryFrameReader.readNextFrame())
86         {
87             const auto frame = trajectoryFrameReader.frame();
88             GMX_RELEASE_ASSERT(frame.v().size() == 2 * c_groupSize,
89                                "Expect velocities for both atom groups");
90
91             const RVec              referenceVelocity = real(frame.time()) * acceleration;
92             const std::vector<RVec> referenceVelocities(c_groupSize, referenceVelocity);
93
94             SCOPED_TRACE("Checking velocities of non-accelerated atoms");
95             ArrayRef<const RVec> nonAcceleratedVelocities = frame.v().subArray(0, c_groupSize);
96             EXPECT_THAT(nonAcceleratedVelocities, Pointwise(RVecEq(tolerance), zeroVelocities));
97
98             SCOPED_TRACE("Checking velocities of accelerated atoms");
99             ArrayRef<const RVec> acceleratedVelocities = frame.v().subArray(c_groupSize, c_groupSize);
100             EXPECT_THAT(acceleratedVelocities, Pointwise(RVecEq(tolerance), referenceVelocities));
101         }
102     }
103 };
104
105 TEST_P(AccelerationGroupTest, WithinTolerances)
106 {
107     const auto& params         = GetParam();
108     const auto& integrator     = std::get<0>(params);
109     const auto& tcoupling      = std::get<1>(params);
110     const auto& simulationName = "spc2";
111
112     // Prepare mdp input
113     auto mdpFieldValues       = prepareMdpFieldValues(simulationName, integrator, tcoupling, "no");
114     mdpFieldValues["nsteps"]  = "8";
115     mdpFieldValues["dt"]      = "0.002";
116     mdpFieldValues["nstxout"] = "0";
117     mdpFieldValues["nstvout"] = "8";
118     mdpFieldValues["nstfout"] = "0";
119     mdpFieldValues["comm-mode"] = "none";
120     // The two groups will not see each other when the cut-off is 0.9 nm
121     mdpFieldValues["coulombtype"]             = "reaction-field";
122     mdpFieldValues["rcoulomb"]                = "0.8";
123     mdpFieldValues["rvdw"]                    = "0.8";
124     mdpFieldValues["verlet-buffer-tolerance"] = "-1";
125     mdpFieldValues["rlist"]                   = "0.9";
126     // Couple the (non-)accelerated groups separately, so their velocties are independent
127     mdpFieldValues["tc-grps"] = "firstWaterMolecule secondWaterMolecule";
128     mdpFieldValues["ref-t"]   = "0.001 0.001";
129     // Use weak coupling so we can check vecolities of atoms with a tight tolerance
130     mdpFieldValues["tau-t"]      = "10.0 10.0";
131     const RVec c_acceleration    = { 2.0, 3.0, 1.5 };
132     mdpFieldValues["acc-grps"]   = "secondWaterMolecule";
133     mdpFieldValues["accelerate"] = "2.0 3.0 1.5";
134     // Set initial velocities to zero
135     mdpFieldValues["gen-vel"]  = "yes";
136     mdpFieldValues["gen-temp"] = "0";
137
138     // Run grompp
139     runner_.useTopGroAndNdxFromDatabase(simulationName);
140     runner_.useStringAsMdpFile(prepareMdpFileContents(mdpFieldValues));
141     runGrompp(&runner_);
142     // Run mdrun
143     runMdrun(&runner_);
144
145     // T-coupling causes changes in the velocities up to 1e-4
146     auto tolerance = absoluteTolerance((GMX_DOUBLE && tcoupling == "no") ? 1e-10 : 1e-4);
147
148     // Check the velocities of the non-accelerated and accelerated groups
149     checkVelocities(runner_.fullPrecisionTrajectoryFileName_, c_acceleration, tolerance);
150 }
151
152 // The v-rescale case check that Ekin computation and temperature coupling
153 // can be performed independently for atoms groups, so the accelerations
154 // are not affected. This can be useful in practice.
155 INSTANTIATE_TEST_SUITE_P(AccelerationWorks,
156                          AccelerationGroupTest,
157                          ::testing::Combine(::testing::Values("md", "md-vv"),
158                                             ::testing::Values("no", "v-rescale")));
159 } // namespace
160 } // namespace gmx::test