Update instructions in containers.rst
[gromacs.git] / src / testutils / refdata.cpp
blob6f44673a64e1e7c77eed24e133557b4ba39db988
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2011-2018, The GROMACS development team.
5 * Copyright (c) 2019,2020, 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.
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.
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.
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.
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.
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.
36 /*! \internal \file
37 * \brief
38 * Implements classes and functions from refdata.h.
40 * \author Teemu Murtola <teemu.murtola@gmail.com>
41 * \ingroup module_testutils
43 #include "gmxpre.h"
45 #include "testutils/refdata.h"
47 #include <cctype>
48 #include <cstdlib>
50 #include <algorithm>
51 #include <limits>
52 #include <string>
54 #include <gtest/gtest.h>
56 #include "gromacs/math/vectypes.h"
57 #include "gromacs/options/basicoptions.h"
58 #include "gromacs/options/ioptionscontainer.h"
59 #include "gromacs/utility/any.h"
60 #include "gromacs/utility/exceptions.h"
61 #include "gromacs/utility/gmxassert.h"
62 #include "gromacs/utility/keyvaluetree.h"
63 #include "gromacs/utility/path.h"
64 #include "gromacs/utility/real.h"
65 #include "gromacs/utility/stringutil.h"
67 #include "testutils/testasserts.h"
68 #include "testutils/testexceptions.h"
69 #include "testutils/testfilemanager.h"
71 #include "refdata_checkers.h"
72 #include "refdata_impl.h"
73 #include "refdata_xml.h"
75 namespace gmx
77 namespace test
80 /********************************************************************
81 * TestReferenceData::Impl declaration
84 namespace internal
87 /*! \internal \brief
88 * Private implementation class for TestReferenceData.
90 * \ingroup module_testutils
92 class TestReferenceDataImpl
94 public:
95 //! Initializes a checker in the given mode.
96 TestReferenceDataImpl(ReferenceDataMode mode, bool bSelfTestMode);
98 //! Performs final reference data processing when test ends.
99 void onTestEnd(bool testPassed);
101 //! Full path of the reference data file.
102 std::string fullFilename_;
103 /*! \brief
104 * Root entry for comparing the reference data.
106 * Null after construction iff in compare mode and reference data was
107 * not loaded successfully.
108 * In all write modes, copies are present for nodes added to
109 * \a outputRootEntry_, and ReferenceDataEntry::correspondingOutputEntry()
110 * points to the copy in the output tree.
112 ReferenceDataEntry::EntryPointer compareRootEntry_;
113 /*! \brief
114 * Root entry for writing new reference data.
116 * Null if only comparing against existing data. Otherwise, starts
117 * always as empty.
118 * When creating new reference data, this is maintained as a copy of
119 * \a compareRootEntry_.
120 * When updating existing data, entries are added either by copying
121 * from \a compareRootEntry_ (if they exist and comparison passes), or
122 * by creating new ones.
124 ReferenceDataEntry::EntryPointer outputRootEntry_;
125 /*! \brief
126 * Whether updating existing reference data.
128 bool updateMismatchingEntries_;
129 //! `true` if self-testing (enables extra failure messages).
130 bool bSelfTestMode_;
131 /*! \brief
132 * Whether any reference checkers have been created for this data.
134 bool bInUse_;
137 } // namespace internal
139 /********************************************************************
140 * Internal helpers
143 namespace
146 //! Convenience typedef for a smart pointer to TestReferenceDataImpl.
147 typedef std::shared_ptr<internal::TestReferenceDataImpl> TestReferenceDataImplPointer;
149 /*! \brief
150 * Global reference data instance.
152 * The object is created when the test creates a TestReferenceData, and the
153 * object is destructed (and other post-processing is done) at the end of each
154 * test by ReferenceDataTestEventListener (which is installed as a Google Test
155 * test listener).
157 TestReferenceDataImplPointer g_referenceData;
158 //! Global reference data mode set with setReferenceDataMode().
159 ReferenceDataMode g_referenceDataMode = ReferenceDataMode::Compare;
161 //! Returns the global reference data mode.
162 ReferenceDataMode getReferenceDataMode()
164 return g_referenceDataMode;
167 //! Returns a reference to the global reference data object.
168 TestReferenceDataImplPointer initReferenceDataInstance()
170 GMX_RELEASE_ASSERT(!g_referenceData, "Test cannot create multiple TestReferenceData instances");
171 g_referenceData.reset(new internal::TestReferenceDataImpl(getReferenceDataMode(), false));
172 return g_referenceData;
175 //! Handles reference data creation for self-tests.
176 TestReferenceDataImplPointer initReferenceDataInstanceForSelfTest(ReferenceDataMode mode)
178 if (g_referenceData)
180 GMX_RELEASE_ASSERT(g_referenceData.unique(),
181 "Test cannot create multiple TestReferenceData instances");
182 g_referenceData->onTestEnd(true);
183 g_referenceData.reset();
185 g_referenceData.reset(new internal::TestReferenceDataImpl(mode, true));
186 return g_referenceData;
189 class ReferenceDataTestEventListener : public ::testing::EmptyTestEventListener
191 public:
192 void OnTestEnd(const ::testing::TestInfo& test_info) override
194 if (g_referenceData)
196 GMX_RELEASE_ASSERT(g_referenceData.unique(), "Test leaked TestRefeferenceData objects");
197 g_referenceData->onTestEnd(test_info.result()->Passed());
198 g_referenceData.reset();
202 void OnTestProgramEnd(const ::testing::UnitTest& /*unused*/) override
204 // Could be used e.g. to free internal buffers allocated by an XML parsing library
208 //! Formats a path to a reference data entry with a non-null id.
209 std::string formatEntryPath(const std::string& prefix, const std::string& id)
211 return prefix + "/" + id;
214 //! Formats a path to a reference data entry with a null id.
215 std::string formatSequenceEntryPath(const std::string& prefix, int seqIndex)
217 return formatString("%s/[%d]", prefix.c_str(), seqIndex + 1);
220 //! Finds all entries that have not been checked under a given root.
221 void gatherUnusedEntries(const ReferenceDataEntry& root,
222 const std::string& rootPath,
223 std::vector<std::string>* unusedPaths)
225 if (!root.hasBeenChecked())
227 unusedPaths->push_back(rootPath);
228 return;
230 int seqIndex = 0;
231 for (const auto& child : root.children())
233 std::string path;
234 if (child->id().empty())
236 path = formatSequenceEntryPath(rootPath, seqIndex);
237 ++seqIndex;
239 else
241 path = formatEntryPath(rootPath, child->id());
243 gatherUnusedEntries(*child, path, unusedPaths);
247 //! Produces a GTest assertion of any entries under given root have not been checked.
248 void checkUnusedEntries(const ReferenceDataEntry& root, const std::string& rootPath)
250 std::vector<std::string> unusedPaths;
251 gatherUnusedEntries(root, rootPath, &unusedPaths);
252 if (!unusedPaths.empty())
254 std::string paths;
255 if (unusedPaths.size() > 5)
257 paths = joinStrings(unusedPaths.begin(), unusedPaths.begin() + 5, "\n ");
258 paths = " " + paths + "\n ...";
260 else
262 paths = joinStrings(unusedPaths.begin(), unusedPaths.end(), "\n ");
263 paths = " " + paths;
265 ADD_FAILURE() << "Reference data items not used in test:" << std::endl << paths;
269 } // namespace
271 void initReferenceData(IOptionsContainer* options)
273 static const gmx::EnumerationArray<ReferenceDataMode, const char*> s_refDataNames = {
274 { "check", "create", "update-changed", "update-all" }
276 options->addOption(EnumOption<ReferenceDataMode>("ref-data")
277 .enumValue(s_refDataNames)
278 .store(&g_referenceDataMode)
279 .description("Operation mode for tests that use reference data"));
280 ::testing::UnitTest::GetInstance()->listeners().Append(new ReferenceDataTestEventListener);
283 /********************************************************************
284 * TestReferenceDataImpl definition
287 namespace internal
290 TestReferenceDataImpl::TestReferenceDataImpl(ReferenceDataMode mode, bool bSelfTestMode) :
291 updateMismatchingEntries_(false),
292 bSelfTestMode_(bSelfTestMode),
293 bInUse_(false)
295 const std::string dirname = bSelfTestMode ? TestFileManager::getGlobalOutputTempDirectory()
296 : TestFileManager::getInputDataDirectory();
297 const std::string filename = TestFileManager::getTestSpecificFileName(".xml");
298 fullFilename_ = Path::join(dirname, "refdata", filename);
300 switch (mode)
302 case ReferenceDataMode::Compare:
303 if (File::exists(fullFilename_, File::throwOnError))
305 compareRootEntry_ = readReferenceDataFile(fullFilename_);
307 break;
308 case ReferenceDataMode::CreateMissing:
309 if (File::exists(fullFilename_, File::throwOnError))
311 compareRootEntry_ = readReferenceDataFile(fullFilename_);
313 else
315 compareRootEntry_ = ReferenceDataEntry::createRoot();
316 outputRootEntry_ = ReferenceDataEntry::createRoot();
318 break;
319 case ReferenceDataMode::UpdateChanged:
320 if (File::exists(fullFilename_, File::throwOnError))
322 compareRootEntry_ = readReferenceDataFile(fullFilename_);
324 else
326 compareRootEntry_ = ReferenceDataEntry::createRoot();
328 outputRootEntry_ = ReferenceDataEntry::createRoot();
329 updateMismatchingEntries_ = true;
330 break;
331 case ReferenceDataMode::UpdateAll:
332 compareRootEntry_ = ReferenceDataEntry::createRoot();
333 outputRootEntry_ = ReferenceDataEntry::createRoot();
334 break;
335 case ReferenceDataMode::Count: GMX_THROW(InternalError("Invalid reference data mode"));
339 void TestReferenceDataImpl::onTestEnd(bool testPassed)
341 if (!bInUse_)
343 return;
345 // TODO: Only write the file with update-changed if there were actual changes.
346 if (outputRootEntry_)
348 if (testPassed)
350 std::string dirname = Path::getParentPath(fullFilename_);
351 if (!Directory::exists(dirname))
353 if (Directory::create(dirname) != 0)
355 GMX_THROW(TestException("Creation of reference data directory failed: " + dirname));
358 writeReferenceDataFile(fullFilename_, *outputRootEntry_);
361 else if (compareRootEntry_)
363 checkUnusedEntries(*compareRootEntry_, "");
367 } // namespace internal
370 /********************************************************************
371 * TestReferenceChecker::Impl
374 /*! \internal \brief
375 * Private implementation class for TestReferenceChecker.
377 * \ingroup module_testutils
379 class TestReferenceChecker::Impl
381 public:
382 //! String constant for naming XML elements for boolean values.
383 static const char* const cBooleanNodeName;
384 //! String constant for naming XML elements for string values.
385 static const char* const cStringNodeName;
386 //! String constant for naming XML elements for unsigned char values.
387 static const char* const cUCharNodeName;
388 //! String constant for naming XML elements for integer values.
389 static const char* const cIntegerNodeName;
390 //! String constant for naming XML elements for int32 values.
391 static const char* const cInt32NodeName;
392 //! String constant for naming XML elements for unsigned int32 values.
393 static const char* const cUInt32NodeName;
394 //! String constant for naming XML elements for int32 values.
395 static const char* const cInt64NodeName;
396 //! String constant for naming XML elements for unsigned int64 values.
397 static const char* const cUInt64NodeName;
398 //! String constant for naming XML elements for floating-point values.
399 static const char* const cRealNodeName;
400 //! String constant for naming XML attribute for value identifiers.
401 static const char* const cIdAttrName;
402 //! String constant for naming compounds for vectors.
403 static const char* const cVectorType;
404 //! String constant for naming compounds for key-value tree objects.
405 static const char* const cObjectType;
406 //! String constant for naming compounds for sequences.
407 static const char* const cSequenceType;
408 //! String constant for value identifier for sequence length.
409 static const char* const cSequenceLengthName;
411 //! Creates a checker that does nothing.
412 explicit Impl(bool initialized);
413 //! Creates a checker with a given root entry.
414 Impl(const std::string& path,
415 ReferenceDataEntry* compareRootEntry,
416 ReferenceDataEntry* outputRootEntry,
417 bool updateMismatchingEntries,
418 bool bSelfTestMode,
419 const FloatingPointTolerance& defaultTolerance);
421 //! Returns the path of this checker with \p id appended.
422 std::string appendPath(const char* id) const;
424 //! Creates an entry with given parameters and fills it with \p checker.
425 static ReferenceDataEntry::EntryPointer createEntry(const char* type,
426 const char* id,
427 const IReferenceDataEntryChecker& checker)
429 ReferenceDataEntry::EntryPointer entry(new ReferenceDataEntry(type, id));
430 checker.fillEntry(entry.get());
431 return entry;
433 //! Checks an entry for correct type and using \p checker.
434 static ::testing::AssertionResult checkEntry(const ReferenceDataEntry& entry,
435 const std::string& fullId,
436 const char* type,
437 const IReferenceDataEntryChecker& checker)
439 if (entry.type() != type)
441 return ::testing::AssertionFailure() << "Mismatching reference data item type" << std::endl
442 << " In item: " << fullId << std::endl
443 << " Actual: " << type << std::endl
444 << "Reference: " << entry.type();
446 return checker.checkEntry(entry, fullId);
448 //! Finds an entry by id and updates the last found entry pointer.
449 ReferenceDataEntry* findEntry(const char* id);
450 /*! \brief
451 * Finds/creates a reference data entry to match against.
453 * \param[in] type Type of entry to create.
454 * \param[in] id Unique identifier of the entry (can be NULL, in
455 * which case the next entry without an id is matched).
456 * \param[out] checker Checker to use for filling out created entries.
457 * \returns Matching entry, or NULL if no matching entry found
458 * (NULL is never returned in write mode; new entries are created
459 * instead).
461 ReferenceDataEntry* findOrCreateEntry(const char* type,
462 const char* id,
463 const IReferenceDataEntryChecker& checker);
464 /*! \brief
465 * Helper method for checking a reference data value.
467 * \param[in] name Type of entry to find.
468 * \param[in] id Unique identifier of the entry (can be NULL, in
469 * which case the next entry without an id is matched).
470 * \param[in] checker Checker that provides logic specific to the
471 * type of the entry.
472 * \returns Whether the reference data matched, including details
473 * of the mismatch if the comparison failed.
474 * \throws TestException if there is a problem parsing the
475 * reference data.
477 * Performs common tasks in checking a reference value, such as
478 * finding or creating the correct entry.
479 * Caller needs to provide a checker object that provides the string
480 * value for a newly created entry and performs the actual comparison
481 * against a found entry.
483 ::testing::AssertionResult processItem(const char* name,
484 const char* id,
485 const IReferenceDataEntryChecker& checker);
486 /*! \brief
487 * Whether the checker is initialized.
489 bool initialized() const { return initialized_; }
490 /*! \brief
491 * Whether the checker should ignore all validation calls.
493 * This is used to ignore any calls within compounds for which
494 * reference data could not be found, such that only one error is
495 * issued for the missing compound, instead of every individual value.
497 bool shouldIgnore() const
499 GMX_RELEASE_ASSERT(initialized(), "Accessing uninitialized reference data checker.");
500 return compareRootEntry_ == nullptr;
503 //! Whether initialized with other means than the default constructor.
504 bool initialized_;
505 //! Default floating-point comparison tolerance.
506 FloatingPointTolerance defaultTolerance_;
507 /*! \brief
508 * Human-readable path to the root node of this checker.
510 * For the root checker, this will be "/", and for each compound, the
511 * id of the compound is added. Used for reporting comparison
512 * mismatches.
514 std::string path_;
515 /*! \brief
516 * Current entry under which reference data is searched for comparison.
518 * Points to either the TestReferenceDataImpl::compareRootEntry_, or to
519 * a compound entry in the tree rooted at that entry.
521 * Can be NULL, in which case this checker does nothing (doesn't even
522 * report errors, see shouldIgnore()).
524 ReferenceDataEntry* compareRootEntry_;
525 /*! \brief
526 * Current entry under which entries for writing are created.
528 * Points to either the TestReferenceDataImpl::outputRootEntry_, or to
529 * a compound entry in the tree rooted at that entry. NULL if only
530 * comparing, or if shouldIgnore() returns `false`.
532 ReferenceDataEntry* outputRootEntry_;
533 /*! \brief
534 * Iterator to a child of \a compareRootEntry_ that was last found.
536 * If `compareRootEntry_->isValidChild()` returns false, no entry has
537 * been found yet.
538 * After every check, is updated to point to the entry that was used
539 * for the check.
540 * Subsequent checks start the search for the matching node on this
541 * node.
543 ReferenceDataEntry::ChildIterator lastFoundEntry_;
544 /*! \brief
545 * Whether the reference data is being written (true) or compared
546 * (false).
548 bool updateMismatchingEntries_;
549 //! `true` if self-testing (enables extra failure messages).
550 bool bSelfTestMode_;
551 /*! \brief
552 * Current number of unnamed elements in a sequence.
554 * It is the index of the current unnamed element.
556 int seqIndex_;
559 const char* const TestReferenceChecker::Impl::cBooleanNodeName = "Bool";
560 const char* const TestReferenceChecker::Impl::cStringNodeName = "String";
561 const char* const TestReferenceChecker::Impl::cUCharNodeName = "UChar";
562 const char* const TestReferenceChecker::Impl::cIntegerNodeName = "Int";
563 const char* const TestReferenceChecker::Impl::cInt32NodeName = "Int32";
564 const char* const TestReferenceChecker::Impl::cUInt32NodeName = "UInt32";
565 const char* const TestReferenceChecker::Impl::cInt64NodeName = "Int64";
566 const char* const TestReferenceChecker::Impl::cUInt64NodeName = "UInt64";
567 const char* const TestReferenceChecker::Impl::cRealNodeName = "Real";
568 const char* const TestReferenceChecker::Impl::cIdAttrName = "Name";
569 const char* const TestReferenceChecker::Impl::cVectorType = "Vector";
570 const char* const TestReferenceChecker::Impl::cObjectType = "Object";
571 const char* const TestReferenceChecker::Impl::cSequenceType = "Sequence";
572 const char* const TestReferenceChecker::Impl::cSequenceLengthName = "Length";
575 TestReferenceChecker::Impl::Impl(bool initialized) :
576 initialized_(initialized),
577 defaultTolerance_(defaultRealTolerance()),
578 compareRootEntry_(nullptr),
579 outputRootEntry_(nullptr),
580 updateMismatchingEntries_(false),
581 bSelfTestMode_(false),
582 seqIndex_(-1)
587 TestReferenceChecker::Impl::Impl(const std::string& path,
588 ReferenceDataEntry* compareRootEntry,
589 ReferenceDataEntry* outputRootEntry,
590 bool updateMismatchingEntries,
591 bool bSelfTestMode,
592 const FloatingPointTolerance& defaultTolerance) :
593 initialized_(true),
594 defaultTolerance_(defaultTolerance),
595 path_(path),
596 compareRootEntry_(compareRootEntry),
597 outputRootEntry_(outputRootEntry),
598 lastFoundEntry_(compareRootEntry->children().end()),
599 updateMismatchingEntries_(updateMismatchingEntries),
600 bSelfTestMode_(bSelfTestMode),
601 seqIndex_(-1)
606 std::string TestReferenceChecker::Impl::appendPath(const char* id) const
608 return id != nullptr ? formatEntryPath(path_, id) : formatSequenceEntryPath(path_, seqIndex_);
612 ReferenceDataEntry* TestReferenceChecker::Impl::findEntry(const char* id)
614 ReferenceDataEntry::ChildIterator entry = compareRootEntry_->findChild(id, lastFoundEntry_);
615 seqIndex_ = (id == nullptr) ? seqIndex_ + 1 : -1;
616 if (compareRootEntry_->isValidChild(entry))
618 lastFoundEntry_ = entry;
619 return entry->get();
621 return nullptr;
624 ReferenceDataEntry* TestReferenceChecker::Impl::findOrCreateEntry(const char* type,
625 const char* id,
626 const IReferenceDataEntryChecker& checker)
628 ReferenceDataEntry* entry = findEntry(id);
629 if (entry == nullptr && outputRootEntry_ != nullptr)
631 lastFoundEntry_ = compareRootEntry_->addChild(createEntry(type, id, checker));
632 entry = lastFoundEntry_->get();
634 return entry;
637 ::testing::AssertionResult TestReferenceChecker::Impl::processItem(const char* type,
638 const char* id,
639 const IReferenceDataEntryChecker& checker)
641 if (shouldIgnore())
643 return ::testing::AssertionSuccess();
645 std::string fullId = appendPath(id);
646 ReferenceDataEntry* entry = findOrCreateEntry(type, id, checker);
647 if (entry == nullptr)
649 return ::testing::AssertionFailure() << "Reference data item " << fullId << " not found";
651 entry->setChecked();
652 ::testing::AssertionResult result(checkEntry(*entry, fullId, type, checker));
653 if (outputRootEntry_ != nullptr && entry->correspondingOutputEntry() == nullptr)
655 if (!updateMismatchingEntries_ || result)
657 outputRootEntry_->addChild(entry->cloneToOutputEntry());
659 else
661 ReferenceDataEntry::EntryPointer outputEntry(createEntry(type, id, checker));
662 entry->setCorrespondingOutputEntry(outputEntry.get());
663 outputRootEntry_->addChild(move(outputEntry));
664 return ::testing::AssertionSuccess();
667 if (bSelfTestMode_ && !result)
669 ReferenceDataEntry expected(type, id);
670 checker.fillEntry(&expected);
671 result << std::endl
672 << "String value: '" << expected.value() << "'" << std::endl
673 << " Ref. string: '" << entry->value() << "'";
675 return result;
679 /********************************************************************
680 * TestReferenceData
683 TestReferenceData::TestReferenceData() : impl_(initReferenceDataInstance()) {}
686 TestReferenceData::TestReferenceData(ReferenceDataMode mode) :
687 impl_(initReferenceDataInstanceForSelfTest(mode))
692 TestReferenceData::~TestReferenceData() {}
695 TestReferenceChecker TestReferenceData::rootChecker()
697 if (!impl_->bInUse_ && !impl_->compareRootEntry_)
699 ADD_FAILURE() << "Reference data file not found: " << impl_->fullFilename_;
701 impl_->bInUse_ = true;
702 if (!impl_->compareRootEntry_)
704 return TestReferenceChecker(new TestReferenceChecker::Impl(true));
706 impl_->compareRootEntry_->setChecked();
707 return TestReferenceChecker(new TestReferenceChecker::Impl(
708 "", impl_->compareRootEntry_.get(), impl_->outputRootEntry_.get(),
709 impl_->updateMismatchingEntries_, impl_->bSelfTestMode_, defaultRealTolerance()));
713 /********************************************************************
714 * TestReferenceChecker
717 TestReferenceChecker::TestReferenceChecker() : impl_(new Impl(false)) {}
719 TestReferenceChecker::TestReferenceChecker(Impl* impl) : impl_(impl) {}
721 TestReferenceChecker::TestReferenceChecker(const TestReferenceChecker& other) :
722 impl_(new Impl(*other.impl_))
726 TestReferenceChecker::TestReferenceChecker(TestReferenceChecker&& other) noexcept :
727 impl_(std::move(other.impl_))
731 TestReferenceChecker& TestReferenceChecker::operator=(TestReferenceChecker&& other) noexcept
733 impl_ = std::move(other.impl_);
734 return *this;
737 TestReferenceChecker::~TestReferenceChecker() {}
739 bool TestReferenceChecker::isValid() const
741 return impl_->initialized();
745 void TestReferenceChecker::setDefaultTolerance(const FloatingPointTolerance& tolerance)
747 impl_->defaultTolerance_ = tolerance;
751 void TestReferenceChecker::checkUnusedEntries()
753 if (impl_->compareRootEntry_)
755 gmx::test::checkUnusedEntries(*impl_->compareRootEntry_, impl_->path_);
756 // Mark them checked so that they are reported only once.
757 impl_->compareRootEntry_->setCheckedIncludingChildren();
761 void TestReferenceChecker::disableUnusedEntriesCheck()
763 if (impl_->compareRootEntry_)
765 impl_->compareRootEntry_->setCheckedIncludingChildren();
770 bool TestReferenceChecker::checkPresent(bool bPresent, const char* id)
772 if (impl_->shouldIgnore() || impl_->outputRootEntry_ != nullptr)
774 return bPresent;
776 ReferenceDataEntry::ChildIterator entry =
777 impl_->compareRootEntry_->findChild(id, impl_->lastFoundEntry_);
778 const bool bFound = impl_->compareRootEntry_->isValidChild(entry);
779 if (bFound != bPresent)
781 ADD_FAILURE() << "Mismatch while checking reference data item '" << impl_->appendPath(id) << "'\n"
782 << "Expected: " << (bPresent ? "it is present.\n" : "it is absent.\n")
783 << " Actual: " << (bFound ? "it is present." : "it is absent.");
785 if (bFound && bPresent)
787 impl_->lastFoundEntry_ = entry;
788 return true;
790 return false;
794 TestReferenceChecker TestReferenceChecker::checkCompound(const char* type, const char* id)
796 if (impl_->shouldIgnore())
798 return TestReferenceChecker(new Impl(true));
800 std::string fullId = impl_->appendPath(id);
801 NullChecker checker;
802 ReferenceDataEntry* entry = impl_->findOrCreateEntry(type, id, checker);
803 if (entry == nullptr)
805 ADD_FAILURE() << "Reference data item " << fullId << " not found";
806 return TestReferenceChecker(new Impl(true));
808 entry->setChecked();
809 if (impl_->updateMismatchingEntries_)
811 entry->makeCompound(type);
813 else
815 ::testing::AssertionResult result(impl_->checkEntry(*entry, fullId, type, checker));
816 EXPECT_PLAIN(result);
817 if (!result)
819 return TestReferenceChecker(new Impl(true));
822 if (impl_->outputRootEntry_ != nullptr && entry->correspondingOutputEntry() == nullptr)
824 impl_->outputRootEntry_->addChild(entry->cloneToOutputEntry());
826 return TestReferenceChecker(new Impl(fullId, entry, entry->correspondingOutputEntry(),
827 impl_->updateMismatchingEntries_, impl_->bSelfTestMode_,
828 impl_->defaultTolerance_));
831 TestReferenceChecker TestReferenceChecker::checkCompound(const char* type, const std::string& id)
833 return checkCompound(type, id.c_str());
836 /*! \brief Throw a TestException if the caller tries to write particular refdata that can't work.
838 * If the string to write is non-empty and has only whitespace,
839 * TinyXML2 can't read it correctly, so throw an exception for this
840 * case, so that we can't accidentally use it and run into mysterious
841 * problems.
843 * \todo Eliminate this limitation of TinyXML2. See
844 * e.g. https://github.com/leethomason/tinyxml2/issues/432
846 static void throwIfNonEmptyAndOnlyWhitespace(const std::string& s, const char* id)
848 if (!s.empty() && std::all_of(s.cbegin(), s.cend(), [](const char& c) { return std::isspace(c); }))
850 std::string message("String '" + s + "' with ");
851 message += (id != nullptr) ? "null " : "";
852 message += "ID ";
853 message += (id != nullptr) ? "" : id;
854 message +=
855 " cannot be handled. We must refuse to write a refdata String"
856 "field for a non-empty string that contains only whitespace, "
857 "because it will not be read correctly by TinyXML2.";
858 GMX_THROW(TestException(message));
862 void TestReferenceChecker::checkBoolean(bool value, const char* id)
864 EXPECT_PLAIN(impl_->processItem(Impl::cBooleanNodeName, id,
865 ExactStringChecker(value ? "true" : "false")));
869 void TestReferenceChecker::checkString(const char* value, const char* id)
871 throwIfNonEmptyAndOnlyWhitespace(value, id);
872 EXPECT_PLAIN(impl_->processItem(Impl::cStringNodeName, id, ExactStringChecker(value)));
876 void TestReferenceChecker::checkString(const std::string& value, const char* id)
878 throwIfNonEmptyAndOnlyWhitespace(value, id);
879 EXPECT_PLAIN(impl_->processItem(Impl::cStringNodeName, id, ExactStringChecker(value)));
883 void TestReferenceChecker::checkTextBlock(const std::string& value, const char* id)
885 EXPECT_PLAIN(impl_->processItem(Impl::cStringNodeName, id, ExactStringBlockChecker(value)));
889 void TestReferenceChecker::checkUChar(unsigned char value, const char* id)
891 EXPECT_PLAIN(impl_->processItem(Impl::cUCharNodeName, id,
892 ExactStringChecker(formatString("%d", value))));
895 void TestReferenceChecker::checkInteger(int value, const char* id)
897 EXPECT_PLAIN(impl_->processItem(Impl::cIntegerNodeName, id,
898 ExactStringChecker(formatString("%d", value))));
901 void TestReferenceChecker::checkInt32(int32_t value, const char* id)
903 EXPECT_PLAIN(impl_->processItem(Impl::cInt32NodeName, id,
904 ExactStringChecker(formatString("%" PRId32, value))));
907 void TestReferenceChecker::checkUInt32(uint32_t value, const char* id)
909 EXPECT_PLAIN(impl_->processItem(Impl::cUInt32NodeName, id,
910 ExactStringChecker(formatString("%" PRIu32, value))));
913 void TestReferenceChecker::checkInt64(int64_t value, const char* id)
915 EXPECT_PLAIN(impl_->processItem(Impl::cInt64NodeName, id,
916 ExactStringChecker(formatString("%" PRId64, value))));
919 void TestReferenceChecker::checkUInt64(uint64_t value, const char* id)
921 EXPECT_PLAIN(impl_->processItem(Impl::cUInt64NodeName, id,
922 ExactStringChecker(formatString("%" PRIu64, value))));
925 void TestReferenceChecker::checkDouble(double value, const char* id)
927 FloatingPointChecker<double> checker(value, impl_->defaultTolerance_);
928 EXPECT_PLAIN(impl_->processItem(Impl::cRealNodeName, id, checker));
932 void TestReferenceChecker::checkFloat(float value, const char* id)
934 FloatingPointChecker<float> checker(value, impl_->defaultTolerance_);
935 EXPECT_PLAIN(impl_->processItem(Impl::cRealNodeName, id, checker));
939 void TestReferenceChecker::checkReal(float value, const char* id)
941 checkFloat(value, id);
945 void TestReferenceChecker::checkReal(double value, const char* id)
947 checkDouble(value, id);
951 void TestReferenceChecker::checkRealFromString(const std::string& value, const char* id)
953 FloatingPointFromStringChecker<real> checker(value, impl_->defaultTolerance_);
954 EXPECT_PLAIN(impl_->processItem(Impl::cRealNodeName, id, checker));
958 void TestReferenceChecker::checkVector(const BasicVector<int>& value, const char* id)
960 TestReferenceChecker compound(checkCompound(Impl::cVectorType, id));
961 compound.checkInteger(value[0], "X");
962 compound.checkInteger(value[1], "Y");
963 compound.checkInteger(value[2], "Z");
967 void TestReferenceChecker::checkVector(const BasicVector<float>& value, const char* id)
969 TestReferenceChecker compound(checkCompound(Impl::cVectorType, id));
970 compound.checkReal(value[0], "X");
971 compound.checkReal(value[1], "Y");
972 compound.checkReal(value[2], "Z");
976 void TestReferenceChecker::checkVector(const BasicVector<double>& value, const char* id)
978 TestReferenceChecker compound(checkCompound(Impl::cVectorType, id));
979 compound.checkReal(value[0], "X");
980 compound.checkReal(value[1], "Y");
981 compound.checkReal(value[2], "Z");
985 void TestReferenceChecker::checkVector(const int value[3], const char* id)
987 checkVector(BasicVector<int>(value), id);
991 void TestReferenceChecker::checkVector(const float value[3], const char* id)
993 checkVector(BasicVector<float>(value), id);
997 void TestReferenceChecker::checkVector(const double value[3], const char* id)
999 checkVector(BasicVector<double>(value), id);
1003 void TestReferenceChecker::checkAny(const Any& any, const char* id)
1005 if (any.isType<bool>())
1007 checkBoolean(any.cast<bool>(), id);
1009 else if (any.isType<int>())
1011 checkInteger(any.cast<int>(), id);
1013 else if (any.isType<int32_t>())
1015 checkInt32(any.cast<int32_t>(), id);
1017 else if (any.isType<uint32_t>())
1019 checkInt32(any.cast<uint32_t>(), id);
1021 else if (any.isType<int64_t>())
1023 checkInt64(any.cast<int64_t>(), id);
1025 else if (any.isType<uint64_t>())
1027 checkInt64(any.cast<uint64_t>(), id);
1029 else if (any.isType<float>())
1031 checkFloat(any.cast<float>(), id);
1033 else if (any.isType<double>())
1035 checkDouble(any.cast<double>(), id);
1037 else if (any.isType<std::string>())
1039 checkString(any.cast<std::string>(), id);
1041 else
1043 GMX_THROW(TestException("Unsupported any type"));
1048 void TestReferenceChecker::checkKeyValueTreeObject(const KeyValueTreeObject& tree, const char* id)
1050 TestReferenceChecker compound(checkCompound(Impl::cObjectType, id));
1051 for (const auto& prop : tree.properties())
1053 compound.checkKeyValueTreeValue(prop.value(), prop.key().c_str());
1055 compound.checkUnusedEntries();
1059 void TestReferenceChecker::checkKeyValueTreeValue(const KeyValueTreeValue& value, const char* id)
1061 if (value.isObject())
1063 checkKeyValueTreeObject(value.asObject(), id);
1065 else if (value.isArray())
1067 const auto& values = value.asArray().values();
1068 checkSequence(values.begin(), values.end(), id);
1070 else
1072 checkAny(value.asAny(), id);
1077 TestReferenceChecker TestReferenceChecker::checkSequenceCompound(const char* id, size_t length)
1079 TestReferenceChecker compound(checkCompound(Impl::cSequenceType, id));
1080 compound.checkInteger(static_cast<int>(length), Impl::cSequenceLengthName);
1081 return compound;
1085 unsigned char TestReferenceChecker::readUChar(const char* id)
1087 if (impl_->shouldIgnore())
1089 GMX_THROW(TestException("Trying to read from non-existent reference data value"));
1091 int value = 0;
1092 EXPECT_PLAIN(impl_->processItem(Impl::cUCharNodeName, id, ValueExtractor<int>(&value)));
1093 return value;
1097 int TestReferenceChecker::readInteger(const char* id)
1099 if (impl_->shouldIgnore())
1101 GMX_THROW(TestException("Trying to read from non-existent reference data value"));
1103 int value = 0;
1104 EXPECT_PLAIN(impl_->processItem(Impl::cIntegerNodeName, id, ValueExtractor<int>(&value)));
1105 return value;
1109 int32_t TestReferenceChecker::readInt32(const char* id)
1111 if (impl_->shouldIgnore())
1113 GMX_THROW(TestException("Trying to read from non-existent reference data value"));
1115 int32_t value = 0;
1116 EXPECT_PLAIN(impl_->processItem(Impl::cInt32NodeName, id, ValueExtractor<int32_t>(&value)));
1117 return value;
1121 int64_t TestReferenceChecker::readInt64(const char* id)
1123 if (impl_->shouldIgnore())
1125 GMX_THROW(TestException("Trying to read from non-existent reference data value"));
1127 int64_t value = 0;
1128 EXPECT_PLAIN(impl_->processItem(Impl::cInt64NodeName, id, ValueExtractor<int64_t>(&value)));
1129 return value;
1133 float TestReferenceChecker::readFloat(const char* id)
1135 if (impl_->shouldIgnore())
1137 GMX_THROW(TestException("Trying to read from non-existent reference data value"));
1139 float value = 0;
1140 EXPECT_PLAIN(impl_->processItem(Impl::cRealNodeName, id, ValueExtractor<float>(&value)));
1141 return value;
1145 double TestReferenceChecker::readDouble(const char* id)
1147 if (impl_->shouldIgnore())
1149 GMX_THROW(TestException("Trying to read from non-existent reference data value"));
1151 double value = 0;
1152 EXPECT_PLAIN(impl_->processItem(Impl::cRealNodeName, id, ValueExtractor<double>(&value)));
1153 return value;
1157 std::string TestReferenceChecker::readString(const char* id)
1159 if (impl_->shouldIgnore())
1161 GMX_THROW(TestException("Trying to read from non-existent reference data value"));
1163 std::string value;
1164 EXPECT_PLAIN(impl_->processItem(Impl::cStringNodeName, id, ValueExtractor<std::string>(&value)));
1165 return value;
1168 } // namespace test
1169 } // namespace gmx