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 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/compat/make_unique.h"
54 #include "gromacs/fileio/enxio.h"
55 #include "gromacs/trajectory/energyframe.h"
56 #include "gromacs/utility/exceptions.h"
57 #include "gromacs/utility/stringutil.h"
59 #include "testutils/testasserts.h"
67 openEnergyFileToReadFields(const std::string &filename,
68 const std::vector<std::string> &namesOfRequiredEnergyFields)
70 ener_file_ptr energyFile(open_enx(filename.c_str(), "r"));
74 GMX_THROW(FileIOError("Could not open energy file " + filename + " for reading"));
77 /* Read in the names of energy fields used in this file. The
78 * resulting data structure would leak if an exception was thrown,
79 * so transfer the contents that we actually need to a map we can
82 * TODO Technically, the insertions into the map could throw
83 * std::bad_alloc and we could leak memory allocated by
84 * do_enxnms(), but there's nothing we can do about this right
86 std::map<std::string, int> indicesOfEnergyFields;
89 gmx_enxnm_t *energyNames = nullptr;
90 do_enxnms(energyFile.get(), &numEnergyTerms, &energyNames);
91 for (int i = 0; i != numEnergyTerms; ++i)
93 const char *name = energyNames[i].name;
94 auto requiredEnergy = std::find_if(std::begin(namesOfRequiredEnergyFields),
95 std::end(namesOfRequiredEnergyFields),
96 [name](const std::string &n){
99 if (requiredEnergy != namesOfRequiredEnergyFields.end())
101 indicesOfEnergyFields[name] = i;
104 // Clean up old data structures
105 free_enxnms(numEnergyTerms, energyNames);
108 // Throw if we failed to find the fields we need
109 if (indicesOfEnergyFields.size() != namesOfRequiredEnergyFields.size())
111 std::string requiredEnergiesNotFound = "Did not find the following required energies in mdrun output:\n";
112 for (auto &name : namesOfRequiredEnergyFields)
114 auto possibleIndex = indicesOfEnergyFields.find(name);
115 if (possibleIndex == indicesOfEnergyFields.end())
117 requiredEnergiesNotFound += name + "\n";
120 GMX_THROW(APIError(requiredEnergiesNotFound));
123 return EnergyFrameReaderPtr(compat::make_unique<EnergyFrameReader>(indicesOfEnergyFields,
124 energyFile.release()));
127 //! Helper function to obtain resources
128 static t_enxframe *make_enxframe()
133 init_enxframe(frame);
138 //! Helper function to clean up resources
139 void done_enxframe(t_enxframe *fr)
141 // Free the contents, then the pointer itself
146 // === EnergyFrameReader ===
148 EnergyFrameReader::EnergyFrameReader(const std::map<std::string, int> &indicesOfEnergyFields,
149 ener_file *energyFile)
150 : indicesOfEnergyFields_(indicesOfEnergyFields),
151 energyFileGuard_(energyFile),
152 enxframeGuard_(make_enxframe()),
153 haveProbedForNextFrame_(false),
154 nextFrameExists_(false)
159 EnergyFrameReader::readNextFrame()
161 if (haveProbedForNextFrame_)
163 if (nextFrameExists_)
165 GMX_THROW(APIError("This frame has already been probed for, it should be used before probing again."));
169 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."));
172 haveProbedForNextFrame_ = true;
173 // If there's a next frame, read it into enxframe_, and report the result.
174 return nextFrameExists_ = do_enx(energyFileGuard_.get(), enxframeGuard_.get());
178 EnergyFrameReader::frame()
180 if (!haveProbedForNextFrame_)
184 if (!nextFrameExists_)
186 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()."));
189 // Prepare for reading future frames
190 haveProbedForNextFrame_ = false;
191 nextFrameExists_ = false;
193 // The probe filled enxframe_ with new data, so now we use that data to fill energyFrame
194 return EnergyFrame(*enxframeGuard_.get(), indicesOfEnergyFields_);