cda1460c8c83055b13f6f8898e5c55b72b1216dc
[alexxy/gromacs.git] / src / programs / mdrun / tests / trajectoryreader.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016, 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 Implementions of related classes for tests that want to
38  * inspect trajectories produced by mdrun.
39  *
40  * \author Mark Abraham <mark.j.abraham@gmail.com>
41  * \ingroup module_mdrun_integration_tests
42  */
43 #include "gmxpre.h"
44
45 #include "trajectoryreader.h"
46
47 #include <memory>
48 #include <string>
49
50 #include "gromacs/fileio/oenv.h"
51 #include "gromacs/fileio/trxio.h"
52 #include "gromacs/trajectory/trajectoryframe.h"
53 #include "gromacs/utility/exceptions.h"
54 #include "gromacs/utility/scoped_cptr.h"
55 #include "gromacs/utility/stringutil.h"
56
57 #include "testutils/testasserts.h"
58
59 namespace gmx
60 {
61 namespace test
62 {
63
64 //! Helper function to obtain resources
65 t_trxframe *make_trxframe()
66 {
67     t_trxframe *frame;
68
69     snew(frame, 1);
70     clear_trxframe(frame, true);
71
72     return frame;
73 }
74
75 //! Helper function to clean up resources
76 void done_trxframe(t_trxframe *fr)
77 {
78     // Free the contents, then the pointer itself
79     sfree(fr->x);
80     sfree(fr->v);
81     sfree(fr->f);
82     sfree(fr->index);
83     sfree(fr);
84 }
85
86 // === TrajectoryFrameReader ===
87
88 TrajectoryFrameReader::TrajectoryFrameReader(const std::string &filename)
89     : filename_(filename),
90       trajectoryFileGuard_(nullptr),
91       trxframeGuard_(make_trxframe()),
92       haveReadFirstFrame_(false),
93       haveProbedForNextFrame_(false),
94       nextFrameExists_(false)
95 {
96     gmx_output_env_t *oenv;
97     output_env_init_default(&oenv);
98     oenvGuard_.reset(oenv);
99 }
100
101 bool
102 TrajectoryFrameReader::readNextFrame()
103 {
104     if (haveProbedForNextFrame_)
105     {
106         if (nextFrameExists_)
107         {
108             GMX_THROW(APIError("This frame has already been probed for, it should be used before probing again."));
109         }
110         else
111         {
112             GMX_THROW(APIError("This frame has already been probed for, it doesn't exist, so there should not be subsequent attempts to probe for it."));
113         }
114     }
115     haveProbedForNextFrame_ = true;
116     // If there's a next frame, read it into trxframe_, and report the result.
117     if (!haveReadFirstFrame_)
118     {
119         t_trxstatus *trajectoryFile;
120         int          flags = TRX_READ_X | TRX_READ_V | TRX_READ_F;
121         nextFrameExists_ = read_first_frame(oenvGuard_.get(),
122                                             &trajectoryFile,
123                                             filename_.c_str(),
124                                             trxframeGuard_.get(),
125                                             flags);
126         if (!trajectoryFile)
127         {
128             GMX_THROW(FileIOError("Could not open trajectory file " + filename_ + " for reading"));
129         }
130         trajectoryFileGuard_.reset(trajectoryFile);
131         haveReadFirstFrame_ = true;
132     }
133     else
134     {
135         nextFrameExists_ = read_next_frame(oenvGuard_.get(),
136                                            trajectoryFileGuard_.get(),
137                                            trxframeGuard_.get());
138     }
139     return nextFrameExists_;
140 }
141
142 TrajectoryFrame
143 TrajectoryFrameReader::frame()
144 {
145     TrajectoryFrame frame;
146
147     if (!haveProbedForNextFrame_)
148     {
149         readNextFrame();
150     }
151     if (!nextFrameExists_)
152     {
153         GMX_THROW(APIError("There is no next frame, so there should have been no attempt to use the data, e.g. by reacting to a call to readNextFrame()."));
154     }
155
156     // Prepare for reading future frames
157     haveProbedForNextFrame_ = false;
158     nextFrameExists_        = false;
159
160     // The probe filled trxframeGuard_ with new data, so return it
161     frame.frame_ = trxframeGuard_.get();
162
163     if (!frame.frame_->bStep)
164     {
165         GMX_THROW(APIError("Cannot handle trajectory frame that lacks a step number"));
166     }
167
168     if (!frame.frame_->bTime)
169     {
170         GMX_THROW(APIError("Cannot handle trajectory frame that lacks a time"));
171     }
172
173     return frame;
174 }
175
176 // === TrajectoryFrame ===
177
178 TrajectoryFrame::TrajectoryFrame() : frame_(nullptr) {};
179
180 std::string TrajectoryFrame::getFrameName() const
181 {
182     GMX_RELEASE_ASSERT(frame_, "Cannot get name of invalid frame");
183     return formatString("Time %f Step %" GMX_PRId64, frame_->time, frame_->step);
184 }
185
186 void compareFrames(const std::pair<TrajectoryFrame, TrajectoryFrame> &frames,
187                    FloatingPointTolerance tolerance)
188 {
189     auto &reference = frames.first;
190     auto &test      = frames.second;
191
192     // NB We checked earlier for both frames that bStep and bTime are set
193
194     EXPECT_EQ(reference.frame_->step, test.frame_->step)
195     << "step didn't match between reference run " << reference.getFrameName() << " and test run " << test.getFrameName();
196
197     EXPECT_EQ(reference.frame_->time, test.frame_->time)
198     << "time didn't match between reference run " << reference.getFrameName() << " and test run " << test.getFrameName();
199
200     for (int i = 0; i < reference.frame_->natoms && i < test.frame_->natoms; ++i)
201     {
202         for (int d = 0; d < DIM; ++d)
203         {
204             if (reference.frame_->bX && test.frame_->bX)
205             {
206                 EXPECT_REAL_EQ_TOL(reference.frame_->x[i][d], test.frame_->x[i][d], tolerance)
207                 << " x[" << i << "][" << d <<"] didn't match between reference run " << reference.getFrameName() << " and test run " << test.getFrameName();
208             }
209             if (reference.frame_->bV && test.frame_->bV)
210             {
211                 EXPECT_REAL_EQ_TOL(reference.frame_->v[i][d], test.frame_->v[i][d], tolerance)
212                 << " v[" << i << "][" << d <<"] didn't match between reference run " << reference.getFrameName() << " and test run " << test.getFrameName();
213             }
214             if (reference.frame_->bF && test.frame_->bF)
215             {
216                 EXPECT_REAL_EQ_TOL(reference.frame_->f[i][d], test.frame_->f[i][d], tolerance)
217                 << " f[" << i << "][" << d <<"] didn't match between reference run " << reference.getFrameName() << " and test run " << test.getFrameName();
218             }
219         }
220     }
221 }
222
223 } // namespace
224 } // namespace