Only issue FFT warning messages on changes
[gromacs/AngularHB.git] / src / programs / mdrun / tests / energyreader.cpp
blob0fe91c7e1e78b9ea6ed13df16313b9bc1fe79b9b
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2016,2017, 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/utility/exceptions.h"
55 #include "gromacs/utility/stringutil.h"
57 #include "testutils/testasserts.h"
59 namespace gmx
61 namespace test
64 EnergyFrameReaderPtr
65 openEnergyFileToReadFields(const std::string &filename,
66 const std::vector<std::string> &namesOfRequiredEnergyFields)
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 fields 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> indicesOfEnergyFields;
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(namesOfRequiredEnergyFields),
93 std::end(namesOfRequiredEnergyFields),
94 [name](const std::string &n){
95 return 0 == n.compare(name);
96 });
97 if (requiredEnergy != namesOfRequiredEnergyFields.end())
99 indicesOfEnergyFields[name] = i;
102 // Clean up old data structures
103 free_enxnms(numEnergyTerms, energyNames);
106 // Throw if we failed to find the fields we need
107 if (indicesOfEnergyFields.size() != namesOfRequiredEnergyFields.size())
109 std::string requiredEnergiesNotFound = "Did not find the following required energies in mdrun output:\n";
110 for (auto &name : namesOfRequiredEnergyFields)
112 auto possibleIndex = indicesOfEnergyFields.find(name);
113 if (possibleIndex == indicesOfEnergyFields.end())
115 requiredEnergiesNotFound += name + "\n";
118 GMX_THROW(APIError(requiredEnergiesNotFound));
121 return EnergyFrameReaderPtr(new EnergyFrameReader(indicesOfEnergyFields, 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> &indicesOfEnergyFields,
146 ener_file *energyFile)
147 : indicesOfEnergyFields_(indicesOfEnergyFields),
148 energyFileGuard_(energyFile),
149 enxframeGuard_(make_enxframe()),
150 haveProbedForNextFrame_(false),
151 nextFrameExists_(false)
155 bool
156 EnergyFrameReader::readNextFrame()
158 if (haveProbedForNextFrame_)
160 if (nextFrameExists_)
162 GMX_THROW(APIError("This frame has already been probed for, it should be used before probing again."));
164 else
166 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."));
169 haveProbedForNextFrame_ = true;
170 // If there's a next frame, read it into enxframe_, and report the result.
171 return nextFrameExists_ = do_enx(energyFileGuard_.get(), enxframeGuard_.get());
174 EnergyFrame
175 EnergyFrameReader::frame()
177 EnergyFrame energyFrame;
179 if (!haveProbedForNextFrame_)
181 readNextFrame();
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 // The probe filled enxframe_ with new data, so now we use that data to fill energyFrame
189 t_enxframe *enxframe = enxframeGuard_.get();
190 energyFrame.time_ = enxframe->t;
191 energyFrame.step_ = enxframe->step;
192 for (auto &index : indicesOfEnergyFields_)
194 if (index.second >= enxframe->nre)
196 GMX_THROW(InternalError(formatString("Index %d for energy %s not present in energy frame with %d energies",
197 index.second, index.first.c_str(), enxframe->nre)));
199 energyFrame.values_[index.first] = enxframe->ener[index.second].e;
202 // Prepare for reading future frames
203 haveProbedForNextFrame_ = false;
204 nextFrameExists_ = false;
206 return energyFrame;
209 // === EnergyFrame ===
211 EnergyFrame::EnergyFrame() : values_(), step_(), time_() {};
213 std::string EnergyFrame::getFrameName() const
215 return formatString("Time %f Step %" GMX_PRId64, time_, step_);
218 const real &EnergyFrame::at(const std::string &name) const
220 auto valueIterator = values_.find(name);
221 if (valueIterator == values_.end())
223 GMX_THROW(APIError("Cannot get energy value " + name + " unless previously registered when constructing EnergyFrameReader"));
225 return valueIterator->second;
228 void compareFrames(const std::pair<EnergyFrame, EnergyFrame> &frames,
229 FloatingPointTolerance tolerance)
231 auto &reference = frames.first;
232 auto &test = frames.second;
234 for (auto referenceIt = reference.values_.begin(); referenceIt != reference.values_.end(); ++referenceIt)
236 auto testIt = test.values_.find(referenceIt->first);
237 if (testIt != test.values_.end())
239 auto energyFieldInReference = referenceIt->second;
240 auto energyFieldInTest = testIt->second;
241 EXPECT_REAL_EQ_TOL(energyFieldInReference, energyFieldInTest, tolerance)
242 << referenceIt->first << " didn't match between reference run " << reference.getFrameName() << " and test run " << test.getFrameName();
247 } // namespace
248 } // namespace