Update instructions in containers.rst
[gromacs.git] / src / programs / mdrun / tests / energyreader.cpp
blob858051ffa53149874f8853a33aef4bebfa4295ce
1 /*
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.
36 /*! \internal \file
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
43 #include "gmxpre.h"
45 #include "energyreader.h"
47 #include <algorithm>
48 #include <map>
49 #include <memory>
50 #include <string>
51 #include <vector>
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"
60 namespace gmx
62 namespace test
65 EnergyFrameReaderPtr openEnergyFileToReadTerms(const std::string& filename,
66 const std::vector<std::string>& namesOfRequiredEnergyTerms)
68 ener_file_ptr energyFile(open_enx(filename.c_str(), "r"));
70 if (!energyFile)
72 GMX_THROW(FileIOError("Could not open energy file " + filename + " for reading"));
75 /* Read in the names of energy terms used in this file. The
76 * resulting data structure would leak if an exception was thrown,
77 * so transfer the contents that we actually need to a map we can
78 * keep.
80 * TODO Technically, the insertions into the map could throw
81 * std::bad_alloc and we could leak memory allocated by
82 * do_enxnms(), but there's nothing we can do about this right
83 * now. */
84 std::map<std::string, int> indicesOfEnergyTerms;
86 int numEnergyTerms;
87 gmx_enxnm_t* energyNames = nullptr;
88 do_enxnms(energyFile.get(), &numEnergyTerms, &energyNames);
89 for (int i = 0; i != numEnergyTerms; ++i)
91 const char* name = energyNames[i].name;
92 auto requiredEnergy = std::find_if(std::begin(namesOfRequiredEnergyTerms),
93 std::end(namesOfRequiredEnergyTerms),
94 [name](const std::string& n) { return name == n; });
95 if (requiredEnergy != namesOfRequiredEnergyTerms.end())
97 indicesOfEnergyTerms[name] = i;
100 // Clean up old data structures
101 free_enxnms(numEnergyTerms, energyNames);
104 // Throw if we failed to find the terms we need
105 if (indicesOfEnergyTerms.size() != namesOfRequiredEnergyTerms.size())
107 std::string requiredEnergiesNotFound =
108 "Did not find the following required energies in mdrun output:\n";
109 for (auto& name : namesOfRequiredEnergyTerms)
111 auto possibleIndex = indicesOfEnergyTerms.find(name);
112 if (possibleIndex == indicesOfEnergyTerms.end())
114 requiredEnergiesNotFound += name + "\n";
117 GMX_THROW(APIError(requiredEnergiesNotFound));
120 return EnergyFrameReaderPtr(
121 std::make_unique<EnergyFrameReader>(indicesOfEnergyTerms, energyFile.release()));
124 //! Helper function to obtain resources
125 static t_enxframe* make_enxframe()
127 t_enxframe* frame;
129 snew(frame, 1);
130 init_enxframe(frame);
132 return frame;
135 //! Helper function to clean up resources
136 void done_enxframe(t_enxframe* fr)
138 // Free the contents, then the pointer itself
139 free_enxframe(fr);
140 sfree(fr);
143 // === EnergyFrameReader ===
145 EnergyFrameReader::EnergyFrameReader(const std::map<std::string, int>& indicesOfEnergyTerms,
146 ener_file* energyFile) :
147 indicesOfEnergyTerms_(indicesOfEnergyTerms),
148 energyFileGuard_(energyFile),
149 enxframeGuard_(make_enxframe()),
150 haveProbedForNextFrame_(false),
151 nextFrameExists_(false)
155 bool EnergyFrameReader::readNextFrame()
157 if (haveProbedForNextFrame_)
159 if (nextFrameExists_)
161 GMX_THROW(
162 APIError("This frame has already been probed for, it should be used before "
163 "probing again."));
165 else
167 GMX_THROW(
168 APIError("This frame has already been probed for, it doesn't exist, so there "
169 "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());
177 EnergyFrame EnergyFrameReader::frame()
179 if (!haveProbedForNextFrame_)
181 readNextFrame();
183 if (!nextFrameExists_)
185 GMX_THROW(
186 APIError("There is no next frame, so there should have been no attempt to use the "
187 "data, e.g. by reacting to a call to readNextFrame()."));
190 // Prepare for reading future frames
191 haveProbedForNextFrame_ = false;
192 nextFrameExists_ = false;
194 // The probe filled enxframe_ with new data, so now we use that data to fill energyFrame
195 return EnergyFrame(*enxframeGuard_.get(), indicesOfEnergyTerms_);
198 } // namespace test
199 } // namespace gmx