2 * This file is part of the GROMACS molecular simulation package.
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.
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 Declares interface of a class for tests that want to inspect
38 * trajectories produced by mdrun.
40 * Intended usage is to create a TrajectoryFrameReader. Successive
41 * calls to its TrajectoryFrameReader::readNextFrame() and
42 * TrajectoryFrameReader::frame() methods will return a handle to a
43 * t_trxframe object for each frame.
45 * \author Mark Abraham <mark.j.abraham@gmail.com>
46 * \ingroup module_mdrun_integration_tests
48 #ifndef GMX_PROGRAMS_MDRUN_TESTS_TRAJECTORYREADER_H
49 #define GMX_PROGRAMS_MDRUN_TESTS_TRAJECTORYREADER_H
54 #include "gromacs/fileio/oenv.h"
55 #include "gromacs/fileio/trxio.h"
56 #include "gromacs/trajectory/trajectoryframe.h"
57 #include "gromacs/utility/unique_cptr.h"
59 #include "testutils/testasserts.h"
61 //! Forward declaration
62 struct gmx_output_env_t;
69 //! Forward declaration
70 class TrajectoryFrame;
72 //! Convenience smart pointer typedef
73 typedef unique_cptr<gmx_output_env_t, output_env_done> oenv_ptr;
74 //! Convenience smart pointer typedef
75 typedef unique_cptr<t_trxstatus, close_trx> trxstatus_file_ptr;
76 //! Helper function to free all resources
77 void done_trxframe(t_trxframe *fr);
78 //! Convenience smart pointer typedef
79 typedef unique_cptr<t_trxframe, done_trxframe> trxframe_ptr;
82 * \brief Manages returning a t_trxframe whose contents were read from
83 * successive frames of an trajectory file. */
84 class TrajectoryFrameReader
87 /*! \brief Attempt to read the next frame from the trajectory file.
89 * \return Whether a frame was available to read.
91 * This call wraps the read_first_frame()/read_next_frame()
92 * API, which does the file opening as a side effect of
93 * reading the first frame.
95 * If true is returned, then frame() should be called
96 * to get access to the data. If false is returned, then no
97 * further data exists and no further call to
98 * readNextFrame() or frame() should occur.
100 * \throws FileIOError upon reading the first frame, if the trajectory file cannot be opened
101 * \throws APIError if an earlier probe has not been properly handled
102 * (by calling frame(), or stopping trying to read
104 bool readNextFrame();
105 /*! \brief Return the next frame from the trajectory file.
107 * If the next frame has not been probed for, then probe for
108 * it. If no next frame exists, then throw APIError, because
109 * user code should have called readNextFrame() itself if this
110 * is possible. (This permits user code to avoid making calls
111 * to readNextFrame() in a case where it already knows that
114 * \throws APIError if no next frame exists, or if it lacks either time or step number. */
115 TrajectoryFrame frame();
116 /*! \brief Constructor
118 * \param[in] filename Name of trajectory file to open and read. */
119 explicit TrajectoryFrameReader(const std::string &filename);
121 //! Name of trajectory file to open and read
122 std::string filename_;
123 //! Owning handle of output environment object
125 //! Owning handle of an open trajectory file ready to read frames.
126 trxstatus_file_ptr trajectoryFileGuard_;
127 //! Owning handle of contents of trajectory file frame after reading.
128 const trxframe_ptr trxframeGuard_;
129 //! Whether the first frame has been read
130 bool haveReadFirstFrame_;
131 //! Whether the API has been used properly (ie. probe before reading).
132 bool haveProbedForNextFrame_;
133 //! Whether there has been a probe that found a next frame.
134 bool nextFrameExists_;
136 // Multiple owners of these resources isn't very sensible, so prevent it
137 GMX_DISALLOW_COPY_AND_ASSIGN(TrajectoryFrameReader);
140 //! Convenience smart pointer typedef
141 typedef std::unique_ptr<TrajectoryFrameReader> TrajectoryFrameReaderPtr;
143 /*! \brief Compare the fields of the two frames for equality within
146 * The two frames are required to have valid and matching values for
147 * time and step. Positions, velocities and/or forces will be compared
148 * when present in both frames, and expected to be equal within \c
150 void compareFrames(const std::pair<TrajectoryFrame, TrajectoryFrame> &frames,
151 FloatingPointTolerance tolerance);
154 * \brief Contains the content of a trajectory frame read by an TrajectoryFrameReader
156 * Objects of this type are intended to be constructed by
157 * TrajectoryFrameReader objects, and as such will always contain valid
158 * data from an trajectory file frame. */
159 class TrajectoryFrame
162 /*! \brief Return string that helps users identify this frame, containing time and step number.
164 * \throws std::bad_alloc when out of memory */
165 std::string getFrameName() const;
169 //! Handle to trajectory data