2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017,2018,2019, 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 energies produced by mdrun.
40 * \author Mark Abraham <mark.j.abraham@gmail.com>
41 * \ingroup module_mdrun_integration_tests
45 #include "energyreader.h"
53 #include "gromacs/fileio/enxio.h"
54 #include "gromacs/trajectory/energyframe.h"
55 #include "gromacs/utility/exceptions.h"
56 #include "gromacs/utility/stringutil.h"
58 #include "testutils/testasserts.h"
66 openEnergyFileToReadFields(const std::string &filename,
67 const std::vector<std::string> &namesOfRequiredEnergyFields)
69 ener_file_ptr energyFile(open_enx(filename.c_str(), "r"));
73 GMX_THROW(FileIOError("Could not open energy file " + filename + " for reading"));
76 /* Read in the names of energy fields used in this file. The
77 * resulting data structure would leak if an exception was thrown,
78 * so transfer the contents that we actually need to a map we can
81 * TODO Technically, the insertions into the map could throw
82 * std::bad_alloc and we could leak memory allocated by
83 * do_enxnms(), but there's nothing we can do about this right
85 std::map<std::string, int> indicesOfEnergyFields;
88 gmx_enxnm_t *energyNames = nullptr;
89 do_enxnms(energyFile.get(), &numEnergyTerms, &energyNames);
90 for (int i = 0; i != numEnergyTerms; ++i)
92 const char *name = energyNames[i].name;
93 auto requiredEnergy = std::find_if(std::begin(namesOfRequiredEnergyFields),
94 std::end(namesOfRequiredEnergyFields),
95 [name](const std::string &n){
98 if (requiredEnergy != namesOfRequiredEnergyFields.end())
100 indicesOfEnergyFields[name] = i;
103 // Clean up old data structures
104 free_enxnms(numEnergyTerms, energyNames);
107 // Throw if we failed to find the fields we need
108 if (indicesOfEnergyFields.size() != namesOfRequiredEnergyFields.size())
110 std::string requiredEnergiesNotFound = "Did not find the following required energies in mdrun output:\n";
111 for (auto &name : namesOfRequiredEnergyFields)
113 auto possibleIndex = indicesOfEnergyFields.find(name);
114 if (possibleIndex == indicesOfEnergyFields.end())
116 requiredEnergiesNotFound += name + "\n";
119 GMX_THROW(APIError(requiredEnergiesNotFound));
122 return EnergyFrameReaderPtr(std::make_unique<EnergyFrameReader>(indicesOfEnergyFields,
123 energyFile.release()));
126 //! Helper function to obtain resources
127 static t_enxframe *make_enxframe()
132 init_enxframe(frame);
137 //! Helper function to clean up resources
138 void done_enxframe(t_enxframe *fr)
140 // Free the contents, then the pointer itself
145 // === EnergyFrameReader ===
147 EnergyFrameReader::EnergyFrameReader(const std::map<std::string, int> &indicesOfEnergyFields,
148 ener_file *energyFile)
149 : indicesOfEnergyFields_(indicesOfEnergyFields),
150 energyFileGuard_(energyFile),
151 enxframeGuard_(make_enxframe()),
152 haveProbedForNextFrame_(false),
153 nextFrameExists_(false)
158 EnergyFrameReader::readNextFrame()
160 if (haveProbedForNextFrame_)
162 if (nextFrameExists_)
164 GMX_THROW(APIError("This frame has already been probed for, it should be used before probing again."));
168 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."));
171 haveProbedForNextFrame_ = true;
172 // If there's a next frame, read it into enxframe_, and report the result.
173 return nextFrameExists_ = do_enx(energyFileGuard_.get(), enxframeGuard_.get());
177 EnergyFrameReader::frame()
179 if (!haveProbedForNextFrame_)
183 if (!nextFrameExists_)
185 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()."));
188 // Prepare for reading future frames
189 haveProbedForNextFrame_ = false;
190 nextFrameExists_ = false;
192 // The probe filled enxframe_ with new data, so now we use that data to fill energyFrame
193 return EnergyFrame(*enxframeGuard_.get(), indicesOfEnergyFields_);