3c72514481cab29dbfbffc189f2734696ae7376d
[alexxy/gromacs.git] / src / testutils / refdata.cpp
1 /*
2  * This file is part of the GROMACS molecular simulation package.
3  *
4  * Copyright (c) 2011-2018, The GROMACS development team.
5  * Copyright (c) 2019,2020,2021, by the GROMACS development team, led by
6  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
7  * and including many others, as listed in the AUTHORS file in the
8  * top-level source directory and at http://www.gromacs.org.
9  *
10  * GROMACS is free software; you can redistribute it and/or
11  * modify it under the terms of the GNU Lesser General Public License
12  * as published by the Free Software Foundation; either version 2.1
13  * of the License, or (at your option) any later version.
14  *
15  * GROMACS is distributed in the hope that it will be useful,
16  * but WITHOUT ANY WARRANTY; without even the implied warranty of
17  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
18  * Lesser General Public License for more details.
19  *
20  * You should have received a copy of the GNU Lesser General Public
21  * License along with GROMACS; if not, see
22  * http://www.gnu.org/licenses, or write to the Free Software Foundation,
23  * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
24  *
25  * If you want to redistribute modifications to GROMACS, please
26  * consider that scientific software is very special. Version
27  * control is crucial - bugs must be traceable. We will be happy to
28  * consider code for inclusion in the official distribution, but
29  * derived work must not be called official GROMACS. Details are found
30  * in the README & COPYING files - if they are missing, get the
31  * official version at http://www.gromacs.org.
32  *
33  * To help us fund GROMACS development, we humbly ask that you cite
34  * the research papers on the package. Check out http://www.gromacs.org.
35  */
36 /*! \internal \file
37  * \brief
38  * Implements classes and functions from refdata.h.
39  *
40  * \author Teemu Murtola <teemu.murtola@gmail.com>
41  * \author Mark Abraham <mark.j.abraham@gmail.com>
42  * \ingroup module_testutils
43  */
44 #include "gmxpre.h"
45
46 #include "testutils/refdata.h"
47
48 #include <cctype>
49 #include <cstdlib>
50
51 #include <algorithm>
52 #include <limits>
53 #include <optional>
54 #include <string>
55
56 #include <gtest/gtest.h>
57
58 #include "gromacs/math/vectypes.h"
59 #include "gromacs/options/basicoptions.h"
60 #include "gromacs/options/ioptionscontainer.h"
61 #include "gromacs/utility/any.h"
62 #include "gromacs/utility/exceptions.h"
63 #include "gromacs/utility/gmxassert.h"
64 #include "gromacs/utility/keyvaluetree.h"
65 #include "gromacs/utility/path.h"
66 #include "gromacs/utility/real.h"
67 #include "gromacs/utility/stringutil.h"
68
69 #include "testutils/testasserts.h"
70 #include "testutils/testexceptions.h"
71 #include "testutils/testfilemanager.h"
72
73 #include "refdata_checkers.h"
74 #include "refdata_impl.h"
75 #include "refdata_xml.h"
76
77 namespace gmx
78 {
79 namespace test
80 {
81
82 /********************************************************************
83  * TestReferenceData::Impl declaration
84  */
85
86 namespace internal
87 {
88
89 /*! \internal \brief
90  * Private implementation class for TestReferenceData.
91  *
92  * \ingroup module_testutils
93  */
94 class TestReferenceDataImpl
95 {
96 public:
97     //! Initializes a checker in the given mode.
98     TestReferenceDataImpl(ReferenceDataMode mode, bool bSelfTestMode, std::optional<std::string> testNameOverride);
99
100     //! Performs final reference data processing when test ends.
101     void onTestEnd(bool testPassed) const;
102
103     //! Full path of the reference data file.
104     std::string fullFilename_;
105     /*! \brief
106      * Root entry for comparing the reference data.
107      *
108      * Null after construction iff in compare mode and reference data was
109      * not loaded successfully.
110      * In all write modes, copies are present for nodes added to
111      * \a outputRootEntry_, and ReferenceDataEntry::correspondingOutputEntry()
112      * points to the copy in the output tree.
113      */
114     ReferenceDataEntry::EntryPointer compareRootEntry_;
115     /*! \brief
116      * Root entry for writing new reference data.
117      *
118      * Null if only comparing against existing data.  Otherwise, starts
119      * always as empty.
120      * When creating new reference data, this is maintained as a copy of
121      * \a compareRootEntry_.
122      * When updating existing data, entries are added either by copying
123      * from \a compareRootEntry_ (if they exist and comparison passes), or
124      * by creating new ones.
125      */
126     ReferenceDataEntry::EntryPointer outputRootEntry_;
127     /*! \brief
128      * Whether updating existing reference data.
129      */
130     bool updateMismatchingEntries_;
131     //! `true` if self-testing (enables extra failure messages).
132     bool bSelfTestMode_;
133     /*! \brief
134      * Whether any reference checkers have been created for this data.
135      */
136     bool bInUse_;
137 };
138
139 } // namespace internal
140
141 /********************************************************************
142  * Internal helpers
143  */
144
145 namespace
146 {
147
148 //! Convenience typedef for a smart pointer to TestReferenceDataImpl.
149 typedef std::shared_ptr<internal::TestReferenceDataImpl> TestReferenceDataImplPointer;
150
151 /*! \brief
152  * Global reference data instance.
153  *
154  * The object is created when the test creates a TestReferenceData, and the
155  * object is destructed (and other post-processing is done) at the end of each
156  * test by ReferenceDataTestEventListener (which is installed as a Google Test
157  * test listener).
158  */
159 TestReferenceDataImplPointer g_referenceData;
160 //! Global reference data mode set with setReferenceDataMode().
161 ReferenceDataMode g_referenceDataMode = ReferenceDataMode::Compare;
162
163 //! Returns the global reference data mode.
164 ReferenceDataMode getReferenceDataMode()
165 {
166     return g_referenceDataMode;
167 }
168
169 //! Returns a reference to the global reference data object.
170 TestReferenceDataImplPointer initReferenceDataInstance(std::optional<std::string> testNameOverride)
171 {
172     GMX_RELEASE_ASSERT(!g_referenceData, "Test cannot create multiple TestReferenceData instances");
173     g_referenceData.reset(new internal::TestReferenceDataImpl(
174             getReferenceDataMode(), false, std::move(testNameOverride)));
175     return g_referenceData;
176 }
177
178 //! Handles reference data creation for self-tests.
179 TestReferenceDataImplPointer initReferenceDataInstanceForSelfTest(ReferenceDataMode mode)
180 {
181     if (g_referenceData)
182     {
183         GMX_RELEASE_ASSERT(g_referenceData.unique(),
184                            "Test cannot create multiple TestReferenceData instances");
185         g_referenceData->onTestEnd(true);
186         g_referenceData.reset();
187     }
188     g_referenceData.reset(new internal::TestReferenceDataImpl(mode, true, std::nullopt));
189     return g_referenceData;
190 }
191
192 class ReferenceDataTestEventListener : public ::testing::EmptyTestEventListener
193 {
194 public:
195     void OnTestEnd(const ::testing::TestInfo& test_info) override
196     {
197         if (g_referenceData)
198         {
199             GMX_RELEASE_ASSERT(g_referenceData.unique(), "Test leaked TestRefeferenceData objects");
200             g_referenceData->onTestEnd(test_info.result()->Passed());
201             g_referenceData.reset();
202         }
203     }
204
205     void OnTestProgramEnd(const ::testing::UnitTest& /*unused*/) override
206     {
207         // Could be used e.g. to free internal buffers allocated by an XML parsing library
208     }
209 };
210
211 //! Formats a path to a reference data entry with a non-null id.
212 std::string formatEntryPath(const std::string& prefix, const std::string& id)
213 {
214     return prefix + "/" + id;
215 }
216
217 //! Formats a path to a reference data entry with a null id.
218 std::string formatSequenceEntryPath(const std::string& prefix, int seqIndex)
219 {
220     return formatString("%s/[%d]", prefix.c_str(), seqIndex + 1);
221 }
222
223 //! Finds all entries that have not been checked under a given root.
224 void gatherUnusedEntries(const ReferenceDataEntry& root,
225                          const std::string&        rootPath,
226                          std::vector<std::string>* unusedPaths)
227 {
228     if (!root.hasBeenChecked())
229     {
230         unusedPaths->push_back(rootPath);
231         return;
232     }
233     int seqIndex = 0;
234     for (const auto& child : root.children())
235     {
236         std::string path;
237         if (child->id().empty())
238         {
239             path = formatSequenceEntryPath(rootPath, seqIndex);
240             ++seqIndex;
241         }
242         else
243         {
244             path = formatEntryPath(rootPath, child->id());
245         }
246         gatherUnusedEntries(*child, path, unusedPaths);
247     }
248 }
249
250 //! Produces a GTest assertion of any entries under given root have not been checked.
251 void checkUnusedEntries(const ReferenceDataEntry& root, const std::string& rootPath)
252 {
253     std::vector<std::string> unusedPaths;
254     gatherUnusedEntries(root, rootPath, &unusedPaths);
255     if (!unusedPaths.empty())
256     {
257         std::string paths;
258         if (unusedPaths.size() > 5)
259         {
260             paths = joinStrings(unusedPaths.begin(), unusedPaths.begin() + 5, "\n  ");
261             paths = "  " + paths + "\n  ...";
262         }
263         else
264         {
265             paths = joinStrings(unusedPaths.begin(), unusedPaths.end(), "\n  ");
266             paths = "  " + paths;
267         }
268         ADD_FAILURE() << "Reference data items not used in test:" << std::endl << paths;
269     }
270 }
271
272 } // namespace
273
274 void initReferenceData(IOptionsContainer* options)
275 {
276     static const gmx::EnumerationArray<ReferenceDataMode, const char*> s_refDataNames = {
277         { "check", "create", "update-changed", "update-all" }
278     };
279     options->addOption(EnumOption<ReferenceDataMode>("ref-data")
280                                .enumValue(s_refDataNames)
281                                .store(&g_referenceDataMode)
282                                .description("Operation mode for tests that use reference data"));
283     ::testing::UnitTest::GetInstance()->listeners().Append(new ReferenceDataTestEventListener);
284 }
285
286 /********************************************************************
287  * TestReferenceDataImpl definition
288  */
289
290 namespace internal
291 {
292
293 TestReferenceDataImpl::TestReferenceDataImpl(ReferenceDataMode          mode,
294                                              bool                       bSelfTestMode,
295                                              std::optional<std::string> testNameOverride) :
296     updateMismatchingEntries_(false), bSelfTestMode_(bSelfTestMode), bInUse_(false)
297 {
298     const std::string dirname  = bSelfTestMode ? TestFileManager::getGlobalOutputTempDirectory()
299                                                : TestFileManager::getInputDataDirectory();
300     const std::string filename = testNameOverride.has_value()
301                                          ? testNameOverride.value()
302                                          : TestFileManager::getTestSpecificFileName(".xml");
303     fullFilename_              = Path::join(dirname, "refdata", filename);
304
305     switch (mode)
306     {
307         case ReferenceDataMode::Compare:
308             if (File::exists(fullFilename_, File::throwOnError))
309             {
310                 compareRootEntry_ = readReferenceDataFile(fullFilename_);
311             }
312             break;
313         case ReferenceDataMode::CreateMissing:
314             if (File::exists(fullFilename_, File::throwOnError))
315             {
316                 compareRootEntry_ = readReferenceDataFile(fullFilename_);
317             }
318             else
319             {
320                 compareRootEntry_ = ReferenceDataEntry::createRoot();
321                 outputRootEntry_  = ReferenceDataEntry::createRoot();
322             }
323             break;
324         case ReferenceDataMode::UpdateChanged:
325             if (File::exists(fullFilename_, File::throwOnError))
326             {
327                 compareRootEntry_ = readReferenceDataFile(fullFilename_);
328             }
329             else
330             {
331                 compareRootEntry_ = ReferenceDataEntry::createRoot();
332             }
333             outputRootEntry_          = ReferenceDataEntry::createRoot();
334             updateMismatchingEntries_ = true;
335             break;
336         case ReferenceDataMode::UpdateAll:
337             compareRootEntry_ = ReferenceDataEntry::createRoot();
338             outputRootEntry_  = ReferenceDataEntry::createRoot();
339             break;
340         case ReferenceDataMode::Count: GMX_THROW(InternalError("Invalid reference data mode"));
341     }
342 }
343
344 void TestReferenceDataImpl::onTestEnd(bool testPassed) const
345 {
346     if (!bInUse_)
347     {
348         return;
349     }
350     // TODO: Only write the file with update-changed if there were actual changes.
351     if (outputRootEntry_)
352     {
353         if (testPassed)
354         {
355             std::string dirname = Path::getParentPath(fullFilename_);
356             if (!Directory::exists(dirname))
357             {
358                 if (Directory::create(dirname) != 0)
359                 {
360                     GMX_THROW(TestException("Creation of reference data directory failed: " + dirname));
361                 }
362             }
363             writeReferenceDataFile(fullFilename_, *outputRootEntry_);
364         }
365     }
366     else if (compareRootEntry_)
367     {
368         checkUnusedEntries(*compareRootEntry_, "");
369     }
370 }
371
372 } // namespace internal
373
374
375 /********************************************************************
376  * TestReferenceChecker::Impl
377  */
378
379 /*! \internal \brief
380  * Private implementation class for TestReferenceChecker.
381  *
382  * \ingroup module_testutils
383  */
384 class TestReferenceChecker::Impl
385 {
386 public:
387     //! String constant for naming XML elements for boolean values.
388     static const char* const cBooleanNodeName;
389     //! String constant for naming XML elements for string values.
390     static const char* const cStringNodeName;
391     //! String constant for naming XML elements for unsigned char values.
392     static const char* const cUCharNodeName;
393     //! String constant for naming XML elements for integer values.
394     static const char* const cIntegerNodeName;
395     //! String constant for naming XML elements for int32 values.
396     static const char* const cInt32NodeName;
397     //! String constant for naming XML elements for unsigned int32 values.
398     static const char* const cUInt32NodeName;
399     //! String constant for naming XML elements for int32 values.
400     static const char* const cInt64NodeName;
401     //! String constant for naming XML elements for unsigned int64 values.
402     static const char* const cUInt64NodeName;
403     //! String constant for naming XML elements for floating-point values.
404     static const char* const cRealNodeName;
405     //! String constant for naming XML attribute for value identifiers.
406     static const char* const cIdAttrName;
407     //! String constant for naming compounds for vectors.
408     static const char* const cVectorType;
409     //! String constant for naming compounds for key-value tree objects.
410     static const char* const cObjectType;
411     //! String constant for naming compounds for sequences.
412     static const char* const cSequenceType;
413     //! String constant for value identifier for sequence length.
414     static const char* const cSequenceLengthName;
415
416     //! Creates a checker that does nothing.
417     explicit Impl(bool initialized);
418     //! Creates a checker with a given root entry.
419     Impl(const std::string&            path,
420          ReferenceDataEntry*           compareRootEntry,
421          ReferenceDataEntry*           outputRootEntry,
422          bool                          updateMismatchingEntries,
423          bool                          bSelfTestMode,
424          const FloatingPointTolerance& defaultTolerance);
425
426     //! Returns the path of this checker with \p id appended.
427     std::string appendPath(const char* id) const;
428
429     //! Creates an entry with given parameters and fills it with \p checker.
430     static ReferenceDataEntry::EntryPointer createEntry(const char*                       type,
431                                                         const char*                       id,
432                                                         const IReferenceDataEntryChecker& checker)
433     {
434         ReferenceDataEntry::EntryPointer entry(new ReferenceDataEntry(type, id));
435         checker.fillEntry(entry.get());
436         return entry;
437     }
438     //! Checks an entry for correct type and using \p checker.
439     static ::testing::AssertionResult checkEntry(const ReferenceDataEntry&         entry,
440                                                  const std::string&                fullId,
441                                                  const char*                       type,
442                                                  const IReferenceDataEntryChecker& checker)
443     {
444         if (entry.type() != type)
445         {
446             return ::testing::AssertionFailure() << "Mismatching reference data item type" << std::endl
447                                                  << "  In item: " << fullId << std::endl
448                                                  << "   Actual: " << type << std::endl
449                                                  << "Reference: " << entry.type();
450         }
451         return checker.checkEntry(entry, fullId);
452     }
453     //! Finds an entry by id and updates the last found entry pointer.
454     ReferenceDataEntry* findEntry(const char* id);
455     /*! \brief
456      * Finds/creates a reference data entry to match against.
457      *
458      * \param[in]  type   Type of entry to create.
459      * \param[in]  id     Unique identifier of the entry (can be NULL, in
460      *      which case the next entry without an id is matched).
461      * \param[out] checker  Checker to use for filling out created entries.
462      * \returns    Matching entry, or NULL if no matching entry found
463      *      (NULL is never returned in write mode; new entries are created
464      *      instead).
465      */
466     ReferenceDataEntry* findOrCreateEntry(const char*                       type,
467                                           const char*                       id,
468                                           const IReferenceDataEntryChecker& checker);
469     /*! \brief
470      * Helper method for checking a reference data value.
471      *
472      * \param[in]  name   Type of entry to find.
473      * \param[in]  id     Unique identifier of the entry (can be NULL, in
474      *     which case the next entry without an id is matched).
475      * \param[in]  checker  Checker that provides logic specific to the
476      *     type of the entry.
477      * \returns    Whether the reference data matched, including details
478      *     of the mismatch if the comparison failed.
479      * \throws     TestException if there is a problem parsing the
480      *     reference data.
481      *
482      * Performs common tasks in checking a reference value, such as
483      * finding or creating the correct entry.
484      * Caller needs to provide a checker object that provides the string
485      * value for a newly created entry and performs the actual comparison
486      * against a found entry.
487      */
488     ::testing::AssertionResult processItem(const char*                       name,
489                                            const char*                       id,
490                                            const IReferenceDataEntryChecker& checker);
491     /*! \brief
492      * Whether the checker is initialized.
493      */
494     bool initialized() const { return initialized_; }
495     /*! \brief
496      * Whether the checker should ignore all validation calls.
497      *
498      * This is used to ignore any calls within compounds for which
499      * reference data could not be found, such that only one error is
500      * issued for the missing compound, instead of every individual value.
501      */
502     bool shouldIgnore() const
503     {
504         GMX_RELEASE_ASSERT(initialized(), "Accessing uninitialized reference data checker.");
505         return compareRootEntry_ == nullptr;
506     }
507
508     //! Whether initialized with other means than the default constructor.
509     bool initialized_;
510     //! Default floating-point comparison tolerance.
511     FloatingPointTolerance defaultTolerance_;
512     /*! \brief
513      * Human-readable path to the root node of this checker.
514      *
515      * For the root checker, this will be "/", and for each compound, the
516      * id of the compound is added.  Used for reporting comparison
517      * mismatches.
518      */
519     std::string path_;
520     /*! \brief
521      * Current entry under which reference data is searched for comparison.
522      *
523      * Points to either the TestReferenceDataImpl::compareRootEntry_, or to
524      * a compound entry in the tree rooted at that entry.
525      *
526      * Can be NULL, in which case this checker does nothing (doesn't even
527      * report errors, see shouldIgnore()).
528      */
529     ReferenceDataEntry* compareRootEntry_;
530     /*! \brief
531      * Current entry under which entries for writing are created.
532      *
533      * Points to either the TestReferenceDataImpl::outputRootEntry_, or to
534      * a compound entry in the tree rooted at that entry.  NULL if only
535      * comparing, or if shouldIgnore() returns `false`.
536      */
537     ReferenceDataEntry* outputRootEntry_;
538     /*! \brief
539      * Iterator to a child of \a compareRootEntry_ that was last found.
540      *
541      * If `compareRootEntry_->isValidChild()` returns false, no entry has
542      * been found yet.
543      * After every check, is updated to point to the entry that was used
544      * for the check.
545      * Subsequent checks start the search for the matching node on this
546      * node.
547      */
548     ReferenceDataEntry::ChildIterator lastFoundEntry_;
549     /*! \brief
550      * Whether the reference data is being written (true) or compared
551      * (false).
552      */
553     bool updateMismatchingEntries_;
554     //! `true` if self-testing (enables extra failure messages).
555     bool bSelfTestMode_;
556     /*! \brief
557      * Current number of unnamed elements in a sequence.
558      *
559      * It is the index of the current unnamed element.
560      */
561     int seqIndex_;
562 };
563
564 const char* const TestReferenceChecker::Impl::cBooleanNodeName    = "Bool";
565 const char* const TestReferenceChecker::Impl::cStringNodeName     = "String";
566 const char* const TestReferenceChecker::Impl::cUCharNodeName      = "UChar";
567 const char* const TestReferenceChecker::Impl::cIntegerNodeName    = "Int";
568 const char* const TestReferenceChecker::Impl::cInt32NodeName      = "Int32";
569 const char* const TestReferenceChecker::Impl::cUInt32NodeName     = "UInt32";
570 const char* const TestReferenceChecker::Impl::cInt64NodeName      = "Int64";
571 const char* const TestReferenceChecker::Impl::cUInt64NodeName     = "UInt64";
572 const char* const TestReferenceChecker::Impl::cRealNodeName       = "Real";
573 const char* const TestReferenceChecker::Impl::cIdAttrName         = "Name";
574 const char* const TestReferenceChecker::Impl::cVectorType         = "Vector";
575 const char* const TestReferenceChecker::Impl::cObjectType         = "Object";
576 const char* const TestReferenceChecker::Impl::cSequenceType       = "Sequence";
577 const char* const TestReferenceChecker::Impl::cSequenceLengthName = "Length";
578
579
580 TestReferenceChecker::Impl::Impl(bool initialized) :
581     initialized_(initialized),
582     defaultTolerance_(defaultRealTolerance()),
583     compareRootEntry_(nullptr),
584     outputRootEntry_(nullptr),
585     updateMismatchingEntries_(false),
586     bSelfTestMode_(false),
587     seqIndex_(-1)
588 {
589 }
590
591
592 TestReferenceChecker::Impl::Impl(const std::string&            path,
593                                  ReferenceDataEntry*           compareRootEntry,
594                                  ReferenceDataEntry*           outputRootEntry,
595                                  bool                          updateMismatchingEntries,
596                                  bool                          bSelfTestMode,
597                                  const FloatingPointTolerance& defaultTolerance) :
598     initialized_(true),
599     defaultTolerance_(defaultTolerance),
600     path_(path),
601     compareRootEntry_(compareRootEntry),
602     outputRootEntry_(outputRootEntry),
603     lastFoundEntry_(compareRootEntry->children().end()),
604     updateMismatchingEntries_(updateMismatchingEntries),
605     bSelfTestMode_(bSelfTestMode),
606     seqIndex_(-1)
607 {
608 }
609
610
611 std::string TestReferenceChecker::Impl::appendPath(const char* id) const
612 {
613     return id != nullptr ? formatEntryPath(path_, id) : formatSequenceEntryPath(path_, seqIndex_);
614 }
615
616
617 ReferenceDataEntry* TestReferenceChecker::Impl::findEntry(const char* id)
618 {
619     ReferenceDataEntry::ChildIterator entry = compareRootEntry_->findChild(id, lastFoundEntry_);
620     seqIndex_                               = (id == nullptr) ? seqIndex_ + 1 : -1;
621     if (compareRootEntry_->isValidChild(entry))
622     {
623         lastFoundEntry_ = entry;
624         return entry->get();
625     }
626     return nullptr;
627 }
628
629 ReferenceDataEntry* TestReferenceChecker::Impl::findOrCreateEntry(const char* type,
630                                                                   const char* id,
631                                                                   const IReferenceDataEntryChecker& checker)
632 {
633     ReferenceDataEntry* entry = findEntry(id);
634     if (entry == nullptr && outputRootEntry_ != nullptr)
635     {
636         lastFoundEntry_ = compareRootEntry_->addChild(createEntry(type, id, checker));
637         entry           = lastFoundEntry_->get();
638     }
639     return entry;
640 }
641
642 ::testing::AssertionResult TestReferenceChecker::Impl::processItem(const char* type,
643                                                                    const char* id,
644                                                                    const IReferenceDataEntryChecker& checker)
645 {
646     if (shouldIgnore())
647     {
648         return ::testing::AssertionSuccess();
649     }
650     std::string         fullId = appendPath(id);
651     ReferenceDataEntry* entry  = findOrCreateEntry(type, id, checker);
652     if (entry == nullptr)
653     {
654         return ::testing::AssertionFailure() << "Reference data item " << fullId << " not found";
655     }
656     entry->setChecked();
657     ::testing::AssertionResult result(checkEntry(*entry, fullId, type, checker));
658     if (outputRootEntry_ != nullptr && entry->correspondingOutputEntry() == nullptr)
659     {
660         if (!updateMismatchingEntries_ || result)
661         {
662             outputRootEntry_->addChild(entry->cloneToOutputEntry());
663         }
664         else
665         {
666             ReferenceDataEntry::EntryPointer outputEntry(createEntry(type, id, checker));
667             entry->setCorrespondingOutputEntry(outputEntry.get());
668             outputRootEntry_->addChild(move(outputEntry));
669             return ::testing::AssertionSuccess();
670         }
671     }
672     if (bSelfTestMode_ && !result)
673     {
674         ReferenceDataEntry expected(type, id);
675         checker.fillEntry(&expected);
676         result << std::endl
677                << "String value: '" << expected.value() << "'" << std::endl
678                << " Ref. string: '" << entry->value() << "'";
679     }
680     return result;
681 }
682
683
684 /********************************************************************
685  * TestReferenceData
686  */
687
688 TestReferenceData::TestReferenceData() : impl_(initReferenceDataInstance(std::nullopt)) {}
689
690
691 TestReferenceData::TestReferenceData(std::string testNameOverride) :
692     impl_(initReferenceDataInstance(std::move(testNameOverride)))
693 {
694 }
695
696 TestReferenceData::TestReferenceData(ReferenceDataMode mode) :
697     impl_(initReferenceDataInstanceForSelfTest(mode))
698 {
699 }
700
701
702 TestReferenceData::~TestReferenceData() {}
703
704
705 TestReferenceChecker TestReferenceData::rootChecker()
706 {
707     if (!impl_->bInUse_ && !impl_->compareRootEntry_)
708     {
709         ADD_FAILURE() << "Reference data file not found: " << impl_->fullFilename_;
710     }
711     impl_->bInUse_ = true;
712     if (!impl_->compareRootEntry_)
713     {
714         return TestReferenceChecker(new TestReferenceChecker::Impl(true));
715     }
716     impl_->compareRootEntry_->setChecked();
717     return TestReferenceChecker(new TestReferenceChecker::Impl("",
718                                                                impl_->compareRootEntry_.get(),
719                                                                impl_->outputRootEntry_.get(),
720                                                                impl_->updateMismatchingEntries_,
721                                                                impl_->bSelfTestMode_,
722                                                                defaultRealTolerance()));
723 }
724
725
726 /********************************************************************
727  * TestReferenceChecker
728  */
729
730 TestReferenceChecker::TestReferenceChecker() : impl_(new Impl(false)) {}
731
732 TestReferenceChecker::TestReferenceChecker(Impl* impl) : impl_(impl) {}
733
734 TestReferenceChecker::TestReferenceChecker(const TestReferenceChecker& other) :
735     impl_(new Impl(*other.impl_))
736 {
737 }
738
739 TestReferenceChecker::TestReferenceChecker(TestReferenceChecker&& other) noexcept :
740     impl_(std::move(other.impl_))
741 {
742 }
743
744 TestReferenceChecker& TestReferenceChecker::operator=(TestReferenceChecker&& other) noexcept
745 {
746     impl_ = std::move(other.impl_);
747     return *this;
748 }
749
750 TestReferenceChecker::~TestReferenceChecker() {}
751
752 bool TestReferenceChecker::isValid() const
753 {
754     return impl_->initialized();
755 }
756
757
758 void TestReferenceChecker::setDefaultTolerance(const FloatingPointTolerance& tolerance)
759 {
760     impl_->defaultTolerance_ = tolerance;
761 }
762
763
764 void TestReferenceChecker::checkUnusedEntries()
765 {
766     if (impl_->compareRootEntry_)
767     {
768         gmx::test::checkUnusedEntries(*impl_->compareRootEntry_, impl_->path_);
769         // Mark them checked so that they are reported only once.
770         impl_->compareRootEntry_->setCheckedIncludingChildren();
771     }
772 }
773
774 void TestReferenceChecker::disableUnusedEntriesCheck()
775 {
776     if (impl_->compareRootEntry_)
777     {
778         impl_->compareRootEntry_->setCheckedIncludingChildren();
779     }
780 }
781
782
783 bool TestReferenceChecker::checkPresent(bool bPresent, const char* id)
784 {
785     if (impl_->shouldIgnore() || impl_->outputRootEntry_ != nullptr)
786     {
787         return bPresent;
788     }
789     ReferenceDataEntry::ChildIterator entry =
790             impl_->compareRootEntry_->findChild(id, impl_->lastFoundEntry_);
791     const bool bFound = impl_->compareRootEntry_->isValidChild(entry);
792     if (bFound != bPresent)
793     {
794         ADD_FAILURE() << "Mismatch while checking reference data item '" << impl_->appendPath(id) << "'\n"
795                       << "Expected: " << (bPresent ? "it is present.\n" : "it is absent.\n")
796                       << "  Actual: " << (bFound ? "it is present." : "it is absent.");
797     }
798     if (bFound && bPresent)
799     {
800         impl_->lastFoundEntry_ = entry;
801         return true;
802     }
803     return false;
804 }
805
806
807 TestReferenceChecker TestReferenceChecker::checkCompound(const char* type, const char* id)
808 {
809     if (impl_->shouldIgnore())
810     {
811         return TestReferenceChecker(new Impl(true));
812     }
813     std::string         fullId = impl_->appendPath(id);
814     NullChecker         checker;
815     ReferenceDataEntry* entry = impl_->findOrCreateEntry(type, id, checker);
816     if (entry == nullptr)
817     {
818         ADD_FAILURE() << "Reference data item " << fullId << " not found";
819         return TestReferenceChecker(new Impl(true));
820     }
821     entry->setChecked();
822     if (impl_->updateMismatchingEntries_)
823     {
824         entry->makeCompound(type);
825     }
826     else
827     {
828         ::testing::AssertionResult result(impl_->checkEntry(*entry, fullId, type, checker));
829         EXPECT_PLAIN(result);
830         if (!result)
831         {
832             return TestReferenceChecker(new Impl(true));
833         }
834     }
835     if (impl_->outputRootEntry_ != nullptr && entry->correspondingOutputEntry() == nullptr)
836     {
837         impl_->outputRootEntry_->addChild(entry->cloneToOutputEntry());
838     }
839     return TestReferenceChecker(new Impl(fullId,
840                                          entry,
841                                          entry->correspondingOutputEntry(),
842                                          impl_->updateMismatchingEntries_,
843                                          impl_->bSelfTestMode_,
844                                          impl_->defaultTolerance_));
845 }
846
847 TestReferenceChecker TestReferenceChecker::checkCompound(const char* type, const std::string& id)
848 {
849     return checkCompound(type, id.c_str());
850 }
851
852 /*! \brief Throw a TestException if the caller tries to write particular refdata that can't work.
853  *
854  * If the string to write is non-empty and has only whitespace,
855  * TinyXML2 can't read it correctly, so throw an exception for this
856  * case, so that we can't accidentally use it and run into mysterious
857  * problems.
858  *
859  * \todo Eliminate this limitation of TinyXML2. See
860  * e.g. https://github.com/leethomason/tinyxml2/issues/432
861  */
862 static void throwIfNonEmptyAndOnlyWhitespace(const std::string& s, const char* id)
863 {
864     if (!s.empty() && std::all_of(s.cbegin(), s.cend(), [](const char& c) { return std::isspace(c); }))
865     {
866         std::string message("String '" + s + "' with ");
867         message += (id != nullptr) ? "null " : "";
868         message += "ID ";
869         message += (id != nullptr) ? "" : id;
870         message +=
871                 " cannot be handled. We must refuse to write a refdata String"
872                 "field for a non-empty string that contains only whitespace, "
873                 "because it will not be read correctly by TinyXML2.";
874         GMX_THROW(TestException(message));
875     }
876 }
877
878 void TestReferenceChecker::checkBoolean(bool value, const char* id)
879 {
880     EXPECT_PLAIN(impl_->processItem(
881             Impl::cBooleanNodeName, id, ExactStringChecker(value ? "true" : "false")));
882 }
883
884
885 void TestReferenceChecker::checkString(const char* value, const char* id)
886 {
887     throwIfNonEmptyAndOnlyWhitespace(value, id);
888     EXPECT_PLAIN(impl_->processItem(Impl::cStringNodeName, id, ExactStringChecker(value)));
889 }
890
891
892 void TestReferenceChecker::checkString(const std::string& value, const char* id)
893 {
894     throwIfNonEmptyAndOnlyWhitespace(value, id);
895     EXPECT_PLAIN(impl_->processItem(Impl::cStringNodeName, id, ExactStringChecker(value)));
896 }
897
898
899 void TestReferenceChecker::checkTextBlock(const std::string& value, const char* id)
900 {
901     EXPECT_PLAIN(impl_->processItem(Impl::cStringNodeName, id, ExactStringBlockChecker(value)));
902 }
903
904
905 void TestReferenceChecker::checkUChar(unsigned char value, const char* id)
906 {
907     EXPECT_PLAIN(impl_->processItem(
908             Impl::cUCharNodeName, id, ExactStringChecker(formatString("%d", value))));
909 }
910
911 void TestReferenceChecker::checkInteger(int value, const char* id)
912 {
913     EXPECT_PLAIN(impl_->processItem(
914             Impl::cIntegerNodeName, id, ExactStringChecker(formatString("%d", value))));
915 }
916
917 void TestReferenceChecker::checkInt32(int32_t value, const char* id)
918 {
919     EXPECT_PLAIN(impl_->processItem(
920             Impl::cInt32NodeName, id, ExactStringChecker(formatString("%" PRId32, value))));
921 }
922
923 void TestReferenceChecker::checkUInt32(uint32_t value, const char* id)
924 {
925     EXPECT_PLAIN(impl_->processItem(
926             Impl::cUInt32NodeName, id, ExactStringChecker(formatString("%" PRIu32, value))));
927 }
928
929 void TestReferenceChecker::checkInt64(int64_t value, const char* id)
930 {
931     EXPECT_PLAIN(impl_->processItem(
932             Impl::cInt64NodeName, id, ExactStringChecker(formatString("%" PRId64, value))));
933 }
934
935 void TestReferenceChecker::checkUInt64(uint64_t value, const char* id)
936 {
937     EXPECT_PLAIN(impl_->processItem(
938             Impl::cUInt64NodeName, id, ExactStringChecker(formatString("%" PRIu64, value))));
939 }
940
941 void TestReferenceChecker::checkDouble(double value, const char* id)
942 {
943     FloatingPointChecker<double> checker(value, impl_->defaultTolerance_);
944     EXPECT_PLAIN(impl_->processItem(Impl::cRealNodeName, id, checker));
945 }
946
947
948 void TestReferenceChecker::checkFloat(float value, const char* id)
949 {
950     FloatingPointChecker<float> checker(value, impl_->defaultTolerance_);
951     EXPECT_PLAIN(impl_->processItem(Impl::cRealNodeName, id, checker));
952 }
953
954
955 void TestReferenceChecker::checkReal(float value, const char* id)
956 {
957     checkFloat(value, id);
958 }
959
960
961 void TestReferenceChecker::checkReal(double value, const char* id)
962 {
963     checkDouble(value, id);
964 }
965
966
967 void TestReferenceChecker::checkRealFromString(const std::string& value, const char* id)
968 {
969     FloatingPointFromStringChecker<real> checker(value, impl_->defaultTolerance_);
970     EXPECT_PLAIN(impl_->processItem(Impl::cRealNodeName, id, checker));
971 }
972
973
974 void TestReferenceChecker::checkVector(const BasicVector<int>& value, const char* id)
975 {
976     TestReferenceChecker compound(checkCompound(Impl::cVectorType, id));
977     compound.checkInteger(value[0], "X");
978     compound.checkInteger(value[1], "Y");
979     compound.checkInteger(value[2], "Z");
980 }
981
982
983 void TestReferenceChecker::checkVector(const BasicVector<float>& value, const char* id)
984 {
985     TestReferenceChecker compound(checkCompound(Impl::cVectorType, id));
986     compound.checkReal(value[0], "X");
987     compound.checkReal(value[1], "Y");
988     compound.checkReal(value[2], "Z");
989 }
990
991
992 void TestReferenceChecker::checkVector(const BasicVector<double>& value, const char* id)
993 {
994     TestReferenceChecker compound(checkCompound(Impl::cVectorType, id));
995     compound.checkReal(value[0], "X");
996     compound.checkReal(value[1], "Y");
997     compound.checkReal(value[2], "Z");
998 }
999
1000
1001 void TestReferenceChecker::checkVector(const int value[3], const char* id)
1002 {
1003     checkVector(BasicVector<int>(value), id);
1004 }
1005
1006
1007 void TestReferenceChecker::checkVector(const float value[3], const char* id)
1008 {
1009     checkVector(BasicVector<float>(value), id);
1010 }
1011
1012
1013 void TestReferenceChecker::checkVector(const double value[3], const char* id)
1014 {
1015     checkVector(BasicVector<double>(value), id);
1016 }
1017
1018
1019 void TestReferenceChecker::checkAny(const Any& any, const char* id)
1020 {
1021     if (any.isType<bool>())
1022     {
1023         checkBoolean(any.cast<bool>(), id);
1024     }
1025     else if (any.isType<int>())
1026     {
1027         checkInteger(any.cast<int>(), id);
1028     }
1029     else if (any.isType<int32_t>())
1030     {
1031         checkInt32(any.cast<int32_t>(), id);
1032     }
1033     else if (any.isType<uint32_t>())
1034     {
1035         checkInt32(any.cast<uint32_t>(), id);
1036     }
1037     else if (any.isType<int64_t>())
1038     {
1039         checkInt64(any.cast<int64_t>(), id);
1040     }
1041     else if (any.isType<uint64_t>())
1042     {
1043         checkInt64(any.cast<uint64_t>(), id);
1044     }
1045     else if (any.isType<float>())
1046     {
1047         checkFloat(any.cast<float>(), id);
1048     }
1049     else if (any.isType<double>())
1050     {
1051         checkDouble(any.cast<double>(), id);
1052     }
1053     else if (any.isType<std::string>())
1054     {
1055         checkString(any.cast<std::string>(), id);
1056     }
1057     else
1058     {
1059         GMX_THROW(TestException("Unsupported any type"));
1060     }
1061 }
1062
1063
1064 void TestReferenceChecker::checkKeyValueTreeObject(const KeyValueTreeObject& tree, const char* id)
1065 {
1066     TestReferenceChecker compound(checkCompound(Impl::cObjectType, id));
1067     for (const auto& prop : tree.properties())
1068     {
1069         compound.checkKeyValueTreeValue(prop.value(), prop.key().c_str());
1070     }
1071     compound.checkUnusedEntries();
1072 }
1073
1074
1075 void TestReferenceChecker::checkKeyValueTreeValue(const KeyValueTreeValue& value, const char* id)
1076 {
1077     if (value.isObject())
1078     {
1079         checkKeyValueTreeObject(value.asObject(), id);
1080     }
1081     else if (value.isArray())
1082     {
1083         const auto& values = value.asArray().values();
1084         checkSequence(values.begin(), values.end(), id);
1085     }
1086     else
1087     {
1088         checkAny(value.asAny(), id);
1089     }
1090 }
1091
1092
1093 TestReferenceChecker TestReferenceChecker::checkSequenceCompound(const char* id, size_t length)
1094 {
1095     TestReferenceChecker compound(checkCompound(Impl::cSequenceType, id));
1096     compound.checkInteger(static_cast<int>(length), Impl::cSequenceLengthName);
1097     return compound;
1098 }
1099
1100
1101 unsigned char TestReferenceChecker::readUChar(const char* id)
1102 {
1103     if (impl_->shouldIgnore())
1104     {
1105         GMX_THROW(TestException("Trying to read from non-existent reference data value"));
1106     }
1107     int value = 0;
1108     EXPECT_PLAIN(impl_->processItem(Impl::cUCharNodeName, id, ValueExtractor<int>(&value)));
1109     return value;
1110 }
1111
1112
1113 int TestReferenceChecker::readInteger(const char* id)
1114 {
1115     if (impl_->shouldIgnore())
1116     {
1117         GMX_THROW(TestException("Trying to read from non-existent reference data value"));
1118     }
1119     int value = 0;
1120     EXPECT_PLAIN(impl_->processItem(Impl::cIntegerNodeName, id, ValueExtractor<int>(&value)));
1121     return value;
1122 }
1123
1124
1125 int32_t TestReferenceChecker::readInt32(const char* id)
1126 {
1127     if (impl_->shouldIgnore())
1128     {
1129         GMX_THROW(TestException("Trying to read from non-existent reference data value"));
1130     }
1131     int32_t value = 0;
1132     EXPECT_PLAIN(impl_->processItem(Impl::cInt32NodeName, id, ValueExtractor<int32_t>(&value)));
1133     return value;
1134 }
1135
1136
1137 int64_t TestReferenceChecker::readInt64(const char* id)
1138 {
1139     if (impl_->shouldIgnore())
1140     {
1141         GMX_THROW(TestException("Trying to read from non-existent reference data value"));
1142     }
1143     int64_t value = 0;
1144     EXPECT_PLAIN(impl_->processItem(Impl::cInt64NodeName, id, ValueExtractor<int64_t>(&value)));
1145     return value;
1146 }
1147
1148
1149 float TestReferenceChecker::readFloat(const char* id)
1150 {
1151     if (impl_->shouldIgnore())
1152     {
1153         GMX_THROW(TestException("Trying to read from non-existent reference data value"));
1154     }
1155     float value = 0;
1156     EXPECT_PLAIN(impl_->processItem(Impl::cRealNodeName, id, ValueExtractor<float>(&value)));
1157     return value;
1158 }
1159
1160
1161 double TestReferenceChecker::readDouble(const char* id)
1162 {
1163     if (impl_->shouldIgnore())
1164     {
1165         GMX_THROW(TestException("Trying to read from non-existent reference data value"));
1166     }
1167     double value = 0;
1168     EXPECT_PLAIN(impl_->processItem(Impl::cRealNodeName, id, ValueExtractor<double>(&value)));
1169     return value;
1170 }
1171
1172
1173 std::string TestReferenceChecker::readString(const char* id)
1174 {
1175     if (impl_->shouldIgnore())
1176     {
1177         GMX_THROW(TestException("Trying to read from non-existent reference data value"));
1178     }
1179     std::string value;
1180     EXPECT_PLAIN(impl_->processItem(Impl::cStringNodeName, id, ValueExtractor<std::string>(&value)));
1181     return value;
1182 }
1183
1184 } // namespace test
1185 } // namespace gmx