2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,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.
37 * \brief Implementions of related classes for tests that want to
38 * inspect trajectories produced by mdrun.
40 * \author Mark Abraham <mark.j.abraham@gmail.com>
41 * \ingroup module_mdrun_integration_tests
45 #include "trajectoryreader.h"
50 #include "gromacs/fileio/oenv.h"
51 #include "gromacs/fileio/trxio.h"
52 #include "gromacs/trajectory/trajectoryframe.h"
53 #include "gromacs/utility/exceptions.h"
55 #include "testutils/testasserts.h"
56 #include "testutils/testmatchers.h"
63 //! Helper function to obtain resources
64 static t_trxframe
*make_trxframe()
69 clear_trxframe(frame
, true);
74 //! Helper function to clean up resources
75 void done_trxframe(t_trxframe
*fr
)
77 // Free the contents, then the pointer itself
85 TrajectoryFrameReader::TrajectoryFrameReader(const std::string
&filename
)
86 : filename_(filename
),
87 trajectoryFileGuard_(nullptr),
88 trxframeGuard_(make_trxframe()),
89 haveReadFirstFrame_(false),
90 haveProbedForNextFrame_(false),
91 nextFrameExists_(false)
93 gmx_output_env_t
*oenv
;
94 output_env_init_default(&oenv
);
95 oenvGuard_
.reset(oenv
);
99 TrajectoryFrameReader::readNextFrame()
101 if (haveProbedForNextFrame_
)
103 if (nextFrameExists_
)
105 GMX_THROW(APIError("This frame has already been probed for, it should be used before probing again."));
109 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."));
112 haveProbedForNextFrame_
= true;
113 // If there's a next frame, read it into trxframe_, and report the result.
114 if (!haveReadFirstFrame_
)
116 t_trxstatus
*trajectoryFile
;
117 int flags
= TRX_READ_X
| TRX_READ_V
| TRX_READ_F
;
118 nextFrameExists_
= read_first_frame(oenvGuard_
.get(),
121 trxframeGuard_
.get(),
125 GMX_THROW(FileIOError("Could not open trajectory file " + filename_
+ " for reading"));
127 trajectoryFileGuard_
.reset(trajectoryFile
);
128 haveReadFirstFrame_
= true;
132 nextFrameExists_
= read_next_frame(oenvGuard_
.get(),
133 trajectoryFileGuard_
.get(),
134 trxframeGuard_
.get());
136 return nextFrameExists_
;
140 TrajectoryFrameReader::frame()
142 if (!haveProbedForNextFrame_
)
146 if (!nextFrameExists_
)
148 GMX_THROW(APIError("There is no next frame, so there should have been no attempt to get it. Perhaps the return value of readNextFrame() was misused."));
151 // Prepare for reading future frames
152 haveProbedForNextFrame_
= false;
153 nextFrameExists_
= false;
155 // The probe filled trxframeGuard_ with new data, so return it
156 return TrajectoryFrame(*trxframeGuard_
.get());