2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,2018,2019,2020, 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.
37 * \brief Implemention of functions for comparing trajectories
40 * \author Mark Abraham <mark.j.abraham@gmail.com>
41 * \ingroup module_mdrun_integration_tests
45 #include "trajectorycomparison.h"
47 #include <gmock/gmock.h>
49 #include "gromacs/pbcutil/pbc.h"
50 #include "gromacs/trajectory/trajectoryframe.h"
52 #include "testutils/testasserts.h"
53 #include "testutils/testmatchers.h"
60 using ::testing::Pointwise;
62 /*! \brief Compares the box from \c reference and \c test
63 * according to the \c matchSettings and \c tolerance.
65 * \todo This could be streamlined when we have a proper 3D matrix
67 static void compareBox(const TrajectoryFrame& reference,
68 const TrajectoryFrame& test,
69 const TrajectoryFrameMatchSettings& matchSettings,
70 const FloatingPointTolerance tolerance)
72 if (!matchSettings.mustCompareBox)
76 bool canCompareBox = true;
77 if (!reference.hasBox())
79 ADD_FAILURE() << "Comparing the box was required, "
80 "but the reference frame did not have one";
81 canCompareBox = false;
85 ADD_FAILURE() << "Comparing the box was required, "
86 "but the test frame did not have one";
87 canCompareBox = false;
95 for (int d = 0; d < DIM; ++d)
97 for (int dd = 0; dd < DIM; ++dd)
99 EXPECT_REAL_EQ_TOL(reference.box()[d][dd], test.box()[d][dd], tolerance);
104 /*! \brief Help put all atom coordinates in \c frame into its box.
106 * This can perhaps go away when frame->x is a container. */
107 static std::vector<RVec> putAtomsInBox(const TrajectoryFrame& frame)
109 std::vector<RVec> x(frame.x().begin(), frame.x().end());
111 for (int d = 0; d < DIM; ++d)
113 for (int dd = 0; dd < DIM; ++dd)
115 box[d][dd] = frame.box()[d][dd];
118 // Note we don't need to compare bPBC because put_atoms_in_box
119 // implements a fallback if nothing specific was set in the
121 put_atoms_in_box(frame.pbc(), box, x);
125 /*! \brief Return whether the \c comparisonConditions and emptiness of
126 * reference and test frames means that a comparison should be
129 * This allows the framework to determine whether it is an error if a
130 * comparison cannot be made. */
131 static bool shouldDoComparison(const ComparisonConditions comparisonConditions,
132 const bool referenceIsEmpty,
133 const bool testIsEmpty)
135 if (comparisonConditions == ComparisonConditions::NoComparison)
140 bool doComparison = true;
143 if (comparisonConditions == ComparisonConditions::MustCompare
144 || comparisonConditions == ComparisonConditions::CompareIfBothFound
145 || (!referenceIsEmpty && comparisonConditions == ComparisonConditions::CompareIfReferenceFound))
147 ADD_FAILURE() << "Test frame lacked quantity for required comparison";
149 doComparison = false;
151 if (referenceIsEmpty)
153 if (comparisonConditions == ComparisonConditions::MustCompare
154 || comparisonConditions == ComparisonConditions::CompareIfBothFound
155 || (!testIsEmpty && comparisonConditions == ComparisonConditions::CompareIfTestFound))
157 ADD_FAILURE() << "Reference frame lacked quantity for required comparison";
159 doComparison = false;
164 /*! \brief Compares the position coordinates from \c reference and \c test
165 * according to the \c matchSettings and \c tolerance. */
166 static void compareCoordinates(const TrajectoryFrame& reference,
167 const TrajectoryFrame& test,
168 const TrajectoryFrameMatchSettings& matchSettings,
169 const FloatingPointTolerance tolerance)
171 SCOPED_TRACE("Comparing coordinates");
172 if (!shouldDoComparison(matchSettings.coordinatesComparison, reference.x().empty(), test.x().empty()))
177 bool canHandlePbc = true;
178 if (!reference.hasBox())
180 ADD_FAILURE() << "Comparing positions required PBC handling, "
181 "but the reference frame did not have a box";
182 canHandlePbc = false;
186 ADD_FAILURE() << "Comparing positions required PBC handling, "
187 "but the test frame did not have a box";
188 canHandlePbc = false;
191 if (matchSettings.requirePbcHandling && !canHandlePbc)
193 ADD_FAILURE() << "Cannot compare positions for the above reason(s)";
197 if ((matchSettings.handlePbcIfPossible || matchSettings.requirePbcHandling) && canHandlePbc)
199 EXPECT_THAT(putAtomsInBox(test), Pointwise(RVecEq(tolerance), putAtomsInBox(reference)));
203 EXPECT_THAT(test.x(), Pointwise(RVecEq(tolerance), reference.x()));
207 /*! \brief Compares the velocities from \c reference and \c test
208 * according to the \c matchSettings and \c tolerance. */
209 static void compareVelocities(const TrajectoryFrame& reference,
210 const TrajectoryFrame& test,
211 const TrajectoryFrameMatchSettings& matchSettings,
212 const FloatingPointTolerance tolerance)
214 SCOPED_TRACE("Comparing velocities");
215 if (!shouldDoComparison(matchSettings.velocitiesComparison, reference.v().empty(), test.v().empty()))
219 EXPECT_THAT(test.v(), Pointwise(RVecEq(tolerance), reference.v()));
222 /*! \brief Compares the forces from \c reference and \c test
223 * according to the \c matchSettings and \c tolerance. */
224 static void compareForces(const TrajectoryFrame& reference,
225 const TrajectoryFrame& test,
226 const TrajectoryFrameMatchSettings& matchSettings,
227 const FloatingPointTolerance tolerance)
229 SCOPED_TRACE("Comparing forces");
230 if (!shouldDoComparison(matchSettings.forcesComparison, reference.f().empty(), test.f().empty()))
234 EXPECT_THAT(test.f(), Pointwise(RVecEq(tolerance), reference.f()));
238 const TrajectoryTolerances TrajectoryComparison::s_defaultTrajectoryTolerances{
239 defaultRealTolerance(), // box
240 relativeToleranceAsFloatingPoint(1.0, 1.0e-3), // positions
241 defaultRealTolerance(), // velocities
242 relativeToleranceAsFloatingPoint(100.0, GMX_DOUBLE ? 1.0e-7 : 5.0e-5) // forces
246 TrajectoryComparison::TrajectoryComparison(const TrajectoryFrameMatchSettings& matchSettings,
247 const TrajectoryTolerances& tolerances) :
248 matchSettings_(matchSettings),
249 tolerances_(tolerances)
253 void TrajectoryComparison::operator()(const TrajectoryFrame& reference, const TrajectoryFrame& test) const
255 if (numComparedFrames_ >= matchSettings_.maxNumTrajectoryFrames)
257 // Nothing should be compared
260 SCOPED_TRACE("Comparing trajectory reference frame " + reference.frameName()
261 + " and test frame " + test.frameName());
262 EXPECT_EQ(reference.step(), test.step());
263 EXPECT_EQ(reference.time(), test.time());
264 compareBox(reference, test, matchSettings_, tolerances_.box);
265 compareCoordinates(reference, test, matchSettings_, tolerances_.coordinates);
266 compareVelocities(reference, test, matchSettings_, tolerances_.velocities);
267 compareForces(reference, test, matchSettings_, tolerances_.forces);
268 numComparedFrames_++;