47e0414959d67602e129a5d613cb24643a0a4494
[alexxy/gromacs.git] / src / programs / mdrun / tests / energyreader.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2016, 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.
8  *
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.
13  *
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.
18  *
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.
23  *
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.
31  *
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.
34  */
35
36 /*! \internal \file
37  * \brief Implementions of related classes for tests that want to
38  * inspect energies produced by mdrun.
39  *
40  * \author Mark Abraham <mark.j.abraham@gmail.com>
41  * \ingroup module_mdrun_integration_tests
42  */
43 #include "gmxpre.h"
44
45 #include "energyreader.h"
46
47 #include <algorithm>
48 #include <map>
49 #include <memory>
50 #include <string>
51 #include <vector>
52
53 #include "gromacs/fileio/enxio.h"
54 #include "gromacs/utility/exceptions.h"
55 #include "gromacs/utility/scoped_cptr.h"
56 #include "gromacs/utility/stringutil.h"
57
58 #include "testutils/testasserts.h"
59
60 namespace gmx
61 {
62 namespace test
63 {
64
65 EnergyFrameReaderPtr
66 openEnergyFileToReadFields(const std::string              &filename,
67                            const std::vector<std::string> &namesOfRequiredEnergyFields)
68 {
69     ener_file_ptr energyFile(open_enx(filename.c_str(), "r"));
70
71     if (!energyFile)
72     {
73         GMX_THROW(FileIOError("Could not open energy file " + filename + " for reading"));
74     }
75
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
79      * keep.
80      *
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
84      * now. */
85     std::map<std::string, int> indicesOfEnergyFields;
86     {
87         int          numEnergyTerms;
88         gmx_enxnm_t *energyNames = nullptr;
89         do_enxnms(energyFile.get(), &numEnergyTerms, &energyNames);
90         for (int i = 0; i != numEnergyTerms; ++i)
91         {
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){
96                                                           return 0 == n.compare(name);
97                                                       });
98             if (requiredEnergy != namesOfRequiredEnergyFields.end())
99             {
100                 indicesOfEnergyFields[name] = i;
101             }
102         }
103         // Clean up old data structures
104         free_enxnms(numEnergyTerms, energyNames);
105     }
106
107     // Throw if we failed to find the fields we need
108     if (indicesOfEnergyFields.size() != namesOfRequiredEnergyFields.size())
109     {
110         std::string requiredEnergiesNotFound = "Did not find the following required energies in mdrun output:\n";
111         for (auto &name : namesOfRequiredEnergyFields)
112         {
113             auto possibleIndex = indicesOfEnergyFields.find(name);
114             if (possibleIndex == indicesOfEnergyFields.end())
115             {
116                 requiredEnergiesNotFound += name + "\n";
117             }
118         }
119         GMX_THROW(APIError(requiredEnergiesNotFound));
120     }
121
122     return EnergyFrameReaderPtr(new EnergyFrameReader(indicesOfEnergyFields, energyFile.release()));
123 }
124
125 //! Helper function to obtain resources
126 t_enxframe *make_enxframe()
127 {
128     t_enxframe *frame;
129
130     snew(frame, 1);
131     init_enxframe(frame);
132
133     return frame;
134 }
135
136 //! Helper function to clean up resources
137 void done_enxframe(t_enxframe *fr)
138 {
139     // Free the contents, then the pointer itself
140     free_enxframe(fr);
141     sfree(fr);
142 }
143
144 // === EnergyFrameReader ===
145
146 EnergyFrameReader::EnergyFrameReader(const std::map<std::string, int> &indicesOfEnergyFields,
147                                      ener_file *energyFile)
148     : indicesOfEnergyFields_(indicesOfEnergyFields),
149       energyFileGuard_(energyFile),
150       enxframeGuard_(make_enxframe()),
151       haveProbedForNextFrame_(false),
152       nextFrameExists_(false)
153 {
154 }
155
156 bool
157 EnergyFrameReader::readNextFrame()
158 {
159     if (haveProbedForNextFrame_)
160     {
161         if (nextFrameExists_)
162         {
163             GMX_THROW(APIError("This frame has already been probed for, it should be used before probing again."));
164         }
165         else
166         {
167             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."));
168         }
169     }
170     haveProbedForNextFrame_ = true;
171     // If there's a next frame, read it into enxframe_, and report the result.
172     return nextFrameExists_ = do_enx(energyFileGuard_.get(), enxframeGuard_.get());
173 }
174
175 EnergyFrame
176 EnergyFrameReader::frame()
177 {
178     EnergyFrame energyFrame;
179
180     if (!haveProbedForNextFrame_)
181     {
182         readNextFrame();
183     }
184     if (!nextFrameExists_)
185     {
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()."));
187     }
188
189     // The probe filled enxframe_ with new data, so now we use that data to fill energyFrame
190     t_enxframe *enxframe = enxframeGuard_.get();
191     energyFrame.time_ = enxframe->t;
192     energyFrame.step_ = enxframe->step;
193     for (auto &index : indicesOfEnergyFields_)
194     {
195         if (index.second >= enxframe->nre)
196         {
197             GMX_THROW(InternalError(formatString("Index %d for energy %s not present in energy frame with %d energies",
198                                                  index.second, index.first.c_str(), enxframe->nre)));
199         }
200         energyFrame.values_[index.first] = enxframe->ener[index.second].e;
201     }
202
203     // Prepare for reading future frames
204     haveProbedForNextFrame_ = false;
205     nextFrameExists_        = false;
206
207     return energyFrame;
208 }
209
210 // === EnergyFrame ===
211
212 EnergyFrame::EnergyFrame() : values_(), step_(), time_() {};
213
214 std::string EnergyFrame::getFrameName() const
215 {
216     return formatString("Time %f Step %" GMX_PRId64, time_, step_);
217 }
218
219 const real &EnergyFrame::at(const std::string &name) const
220 {
221     auto valueIterator = values_.find(name);
222     if (valueIterator == values_.end())
223     {
224         GMX_THROW(APIError("Cannot get energy value " + name + " unless previously registered when constructing EnergyFrameReader"));
225     }
226     return valueIterator->second;
227 }
228
229 void compareFrames(const std::pair<EnergyFrame, EnergyFrame> &frames,
230                    FloatingPointTolerance tolerance)
231 {
232     auto &reference = frames.first;
233     auto &test      = frames.second;
234
235     for (auto referenceIt = reference.values_.begin(); referenceIt != reference.values_.end(); ++referenceIt)
236     {
237         auto testIt = test.values_.find(referenceIt->first);
238         if (testIt != test.values_.end())
239         {
240             auto energyFieldInReference = referenceIt->second;
241             auto energyFieldInTest      = testIt->second;
242             EXPECT_REAL_EQ_TOL(energyFieldInReference, energyFieldInTest, tolerance)
243             << referenceIt->first << " didn't match between reference run " << reference.getFrameName() << " and test run " << test.getFrameName();
244         }
245     }
246 }
247
248 } // namespace
249 } // namespace