Add support for int64 in refdata
[gromacs.git] / src / testutils / refdata.cpp
blob61e1786517e496107b120575ce32e05dd8d0c19b
1 /*
2 * This file is part of the GROMACS molecular simulation package.
4 * Copyright (c) 2011,2012,2013,2014, 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.
35 /*! \internal \file
36 * \brief
37 * Implements classes and functions from refdata.h.
39 * \author Teemu Murtola <teemu.murtola@gmail.com>
40 * \ingroup module_testutils
42 #include "refdata.h"
44 #include <cstdio>
45 #include <cstdlib>
47 #include <limits>
48 #include <string>
50 #include <gtest/gtest.h>
51 #include <libxml/parser.h>
52 #include <libxml/xmlmemory.h>
54 #include "gromacs/fileio/path.h"
55 #include "gromacs/options/basicoptions.h"
56 #include "gromacs/options/options.h"
57 #include "gromacs/utility/exceptions.h"
58 #include "gromacs/utility/gmxassert.h"
59 #include "gromacs/utility/stringutil.h"
61 #include "testutils/testasserts.h"
62 #include "testutils/testexceptions.h"
63 #include "testutils/testfilemanager.h"
65 namespace
68 /*! \internal \brief
69 * Global test environment for freeing up libxml2 internal buffers.
71 class TestReferenceDataEnvironment : public ::testing::Environment
73 public:
74 //! Frees internal buffers allocated by libxml2.
75 virtual void TearDown()
77 xmlCleanupParser();
81 //! Global reference data mode set with gmx::test::setReferenceDataMode().
82 // TODO: Make this a real enum; requires solving a TODO in StringOption.
83 int g_referenceDataMode = gmx::test::erefdataCompare;
85 } // namespace
87 namespace gmx
89 namespace test
92 ReferenceDataMode getReferenceDataMode()
94 return static_cast<ReferenceDataMode>(g_referenceDataMode);
97 void setReferenceDataMode(ReferenceDataMode mode)
99 g_referenceDataMode = mode;
102 std::string getReferenceDataPath()
104 return TestFileManager::getInputFilePath("refdata");
107 void initReferenceData(Options *options)
109 // Needs to correspond to the enum order in refdata.h.
110 const char *const refDataEnum[] = { "check", "create", "update" };
111 options->addOption(
112 StringOption("ref-data").enumValue(refDataEnum)
113 .defaultEnumIndex(0)
114 .storeEnumIndex(&g_referenceDataMode)
115 .description("Operation mode for tests that use reference data"));
116 ::testing::AddGlobalTestEnvironment(new TestReferenceDataEnvironment);
119 /********************************************************************
120 * TestReferenceData::Impl
123 /*! \internal \brief
124 * Private implementation class for TestReferenceData.
126 * \ingroup module_testutils
128 class TestReferenceData::Impl
130 public:
131 //! String constant for output XML version string.
132 static const xmlChar * const cXmlVersion;
133 //! String constant for XML stylesheet processing instruction name.
134 static const xmlChar * const cXmlStyleSheetNodeName;
135 //! String constant for XML stylesheet reference.
136 static const xmlChar * const cXmlStyleSheetContent;
137 //! String constant for naming the root XML element.
138 static const xmlChar * const cRootNodeName;
140 //! Initializes a checker in the given mode.
141 Impl(ReferenceDataMode mode, bool bSelfTestMode);
142 ~Impl();
144 //! Full path of the reference data file.
145 std::string fullFilename_;
146 /*! \brief
147 * XML document for the reference data.
149 * May be NULL if there was an I/O error in initialization.
151 xmlDocPtr refDoc_;
152 /*! \brief
153 * Whether the reference data is being written (true) or compared
154 * (false).
156 bool bWrite_;
157 //! `true` if self-testing (enables extra failure messages).
158 bool bSelfTestMode_;
159 /*! \brief
160 * Whether any reference checkers have been created for this data.
162 bool bInUse_;
165 const xmlChar * const TestReferenceData::Impl::cXmlVersion =
166 (const xmlChar *)"1.0";
167 const xmlChar * const TestReferenceData::Impl::cXmlStyleSheetNodeName =
168 (const xmlChar *)"xml-stylesheet";
169 const xmlChar * const TestReferenceData::Impl::cXmlStyleSheetContent =
170 (const xmlChar *)"type=\"text/xsl\" href=\"referencedata.xsl\"";
171 const xmlChar * const TestReferenceData::Impl::cRootNodeName =
172 (const xmlChar *)"ReferenceData";
175 TestReferenceData::Impl::Impl(ReferenceDataMode mode, bool bSelfTestMode)
176 : refDoc_(NULL), bWrite_(false), bSelfTestMode_(bSelfTestMode), bInUse_(false)
178 std::string dirname = getReferenceDataPath();
179 std::string filename = TestFileManager::getTestSpecificFileName(".xml");
180 fullFilename_ = Path::join(dirname, filename);
182 bWrite_ = true;
183 if (mode != erefdataUpdateAll)
185 FILE *fp = std::fopen(fullFilename_.c_str(), "r");
186 if (fp != NULL)
188 bWrite_ = false;
189 fclose(fp);
191 else if (mode == erefdataCompare)
193 bWrite_ = false;
194 return;
197 if (bWrite_)
199 // TODO: Error checking
200 refDoc_ = xmlNewDoc(cXmlVersion);
201 xmlNodePtr rootNode = xmlNewDocNode(refDoc_, NULL, cRootNodeName, NULL);
202 xmlDocSetRootElement(refDoc_, rootNode);
203 xmlNodePtr xslNode = xmlNewDocPI(refDoc_, cXmlStyleSheetNodeName,
204 cXmlStyleSheetContent);
205 xmlAddPrevSibling(rootNode, xslNode);
207 else
209 refDoc_ = xmlParseFile(fullFilename_.c_str());
210 if (refDoc_ == NULL)
212 GMX_THROW(TestException("Reference data not parsed successfully: " + fullFilename_));
214 xmlNodePtr rootNode = xmlDocGetRootElement(refDoc_);
215 if (rootNode == NULL)
217 xmlFreeDoc(refDoc_);
218 GMX_THROW(TestException("Reference data is empty: " + fullFilename_));
220 if (xmlStrcmp(rootNode->name, cRootNodeName) != 0)
222 xmlFreeDoc(refDoc_);
223 GMX_THROW(TestException("Invalid root node type in " + fullFilename_));
229 TestReferenceData::Impl::~Impl()
231 if (bWrite_ && bInUse_ && refDoc_ != NULL)
233 std::string dirname = getReferenceDataPath();
234 if (!Directory::exists(dirname))
236 if (Directory::create(dirname) != 0)
238 ADD_FAILURE() << "Creation of reference data directory failed for " << dirname;
241 if (xmlSaveFormatFile(fullFilename_.c_str(), refDoc_, 1) == -1)
243 ADD_FAILURE() << "Saving reference data failed for " + fullFilename_;
246 if (refDoc_ != NULL)
248 xmlFreeDoc(refDoc_);
253 /********************************************************************
254 * TestReferenceChecker::Impl
257 /*! \internal \brief
258 * Private implementation class for TestReferenceChecker.
260 * \ingroup module_testutils
262 class TestReferenceChecker::Impl
264 public:
265 //! String constant for naming XML elements for boolean values.
266 static const xmlChar * const cBooleanNodeName;
267 //! String constant for naming XML elements for string values.
268 static const xmlChar * const cStringNodeName;
269 //! String constant for naming XML elements for integer values.
270 static const xmlChar * const cIntegerNodeName;
271 //! String constant for naming XML elements for int64 values.
272 static const xmlChar * const cInt64NodeName;
273 //! String constant for naming XML elements for unsigned int64 values.
274 static const xmlChar * const cUInt64NodeName;
275 //! String constant for naming XML elements for floating-point values.
276 static const xmlChar * const cRealNodeName;
277 //! String constant for naming XML attribute for value identifiers.
278 static const xmlChar * const cIdAttrName;
279 //! String constant for naming compounds for vectors.
280 static const char * const cVectorType;
281 //! String constant for naming compounds for sequences.
282 static const char * const cSequenceType;
283 //! String constant for value identifier for sequence length.
284 static const char * const cSequenceLengthName;
286 //! Creates a checker that does nothing.
287 explicit Impl(bool bWrite);
288 //! Creates a checker with a given root node.
289 Impl(const std::string &path, xmlNodePtr rootNode, bool bWrite,
290 bool bSelfTestMode, const FloatingPointTolerance &defaultTolerance);
292 //! Returns a string for SCOPED_TRACE() for checking element \p id.
293 std::string traceString(const char *id) const;
294 //! Returns the path of this checker with \p id appended.
295 std::string appendPath(const char *id) const;
297 /*! \brief
298 * Finds a reference data node.
300 * \param[in] name Type of node to find (can be NULL, in which case
301 * any type is matched).
302 * \param[in] id Unique identifier of the node (can be NULL, in
303 * which case the next node without an id is matched).
304 * \returns Matching node, or NULL if no matching node found.
306 * Searches for a node in the reference data that matches the given
307 * \p name and \p id. Searching starts from the node that follows the
308 * previously matched node (relevant for performance, and if there are
309 * duplicate ids or nodes without ids). Note that the match pointer is
310 * not updated by this method.
312 xmlNodePtr findNode(const xmlChar *name, const char *id) const;
313 /*! \brief
314 * Finds/creates a reference data node to match against.
316 * \param[in] name Type of node to find.
317 * \param[in] id Unique identifier of the node (can be NULL, in
318 * which case the next node without an id is matched).
319 * \param[out] bFound Whether the node was found (false if the node was
320 * created in write mode).
321 * \returns Matching node, or NULL if no matching node found
322 * (NULL is never returned in write mode).
323 * \throws TestException if node creation fails in write mode.
325 * Finds a node using findNode() and updates the match pointer is a
326 * match is found. If a match is not found, the method returns NULL in
327 * read mode and creates a new node in write mode. If the creation
328 * fails in write mode, throws.
330 xmlNodePtr findOrCreateNode(const xmlChar *name, const char *id,
331 bool *bFound);
332 /*! \brief
333 * Helper method for checking a reference data value.
335 * \param[in] name Type of node to find.
336 * \param[in] id Unique identifier of the node (can be NULL, in
337 * which case the next node without an id is matched).
338 * \param[in] value String value of the value to be compared.
339 * \param[out] bFound true if a matchin value was found.
340 * \returns String value for the reference value.
341 * \throws TestException if node creation fails in write mode.
343 * Performs common tasks in checking a reference value:
344 * finding/creating the correct XML node and reading/writing its string
345 * value. Caller is responsible for converting the value to and from
346 * string where necessary and performing the actual comparison.
348 * In read mode, if a value is not found, adds a Google Test failure
349 * and returns an empty string. If the reference value is found,
350 * returns it (\p value is not used in this case).
352 * In write mode, creates the node if it is not found, sets its value
353 * as \p value and returns \p value.
355 std::string processItem(const xmlChar *name, const char *id,
356 const char *value, bool *bFound);
357 //! Convenience wrapper that takes a std::string.
358 std::string processItem(const xmlChar *name, const char *id,
359 const std::string &value, bool *bFound);
360 /*! \brief
361 * Whether the checker should ignore all validation calls.
363 * This is used to ignore any calls within compounds for which
364 * reference data could not be found, such that only one error is
365 * issued for the missing compound, instead of every individual value.
367 bool shouldIgnore() const;
369 //! Default floating-point comparison tolerance.
370 FloatingPointTolerance defaultTolerance_;
371 /*! \brief
372 * Human-readable path to the root node of this checker.
374 * For the root checker, this will be "/", and for each compound, the
375 * id of the compound is added. Used for reporting comparison
376 * mismatches.
378 std::string path_;
379 /*! \brief
380 * Current node under which reference data is searched.
382 * Points to either the root of TestReferenceData::Impl::refDoc_, or to
383 * a compound node.
385 * Can be NULL, in which case this checker does nothing (doesn't even
386 * report errors, see shouldIgnore()).
388 xmlNodePtr currNode_;
389 /*! \brief
390 * Points to a child of \a currNode_ that was last found.
392 * On initialization, is initialized to NULL. After every check, is
393 * updated to point to the node that was used for the check.
394 * Subsequent checks start the search for the matching node on this
395 * node.
397 * Is NULL if \a currNode_ contains no children or if no checks have
398 * yet been made.
399 * Otherwise, always points to a direct child of \a currNode_.
401 xmlNodePtr prevFoundNode_;
402 /*! \brief
403 * Whether the reference data is being written (true) or compared
404 * (false).
406 bool bWrite_;
407 //! `true` if self-testing (enables extra failure messages).
408 bool bSelfTestMode_;
409 /*! \brief
410 * Current number of unnamed elements in a sequence.
412 * It is the index of the next added unnamed element.
414 int seqIndex_;
417 const xmlChar * const TestReferenceChecker::Impl::cBooleanNodeName =
418 (const xmlChar *)"Bool";
419 const xmlChar * const TestReferenceChecker::Impl::cStringNodeName =
420 (const xmlChar *)"String";
421 const xmlChar * const TestReferenceChecker::Impl::cIntegerNodeName =
422 (const xmlChar *)"Int";
423 const xmlChar * const TestReferenceChecker::Impl::cInt64NodeName =
424 (const xmlChar *)"Int64";
425 const xmlChar * const TestReferenceChecker::Impl::cUInt64NodeName =
426 (const xmlChar *)"UInt64";
427 const xmlChar * const TestReferenceChecker::Impl::cRealNodeName =
428 (const xmlChar *)"Real";
429 const xmlChar * const TestReferenceChecker::Impl::cIdAttrName =
430 (const xmlChar *)"Name";
431 const char * const TestReferenceChecker::Impl::cVectorType =
432 "Vector";
433 const char * const TestReferenceChecker::Impl::cSequenceType =
434 "Sequence";
435 const char * const TestReferenceChecker::Impl::cSequenceLengthName =
436 "Length";
439 TestReferenceChecker::Impl::Impl(bool bWrite)
440 : defaultTolerance_(defaultRealTolerance()),
441 currNode_(NULL), prevFoundNode_(NULL), bWrite_(bWrite),
442 bSelfTestMode_(false), seqIndex_(0)
447 TestReferenceChecker::Impl::Impl(const std::string &path, xmlNodePtr rootNode,
448 bool bWrite, bool bSelfTestMode,
449 const FloatingPointTolerance &defaultTolerance)
450 : defaultTolerance_(defaultTolerance), path_(path + "/"),
451 currNode_(rootNode), prevFoundNode_(NULL), bWrite_(bWrite),
452 bSelfTestMode_(bSelfTestMode), seqIndex_(0)
457 std::string
458 TestReferenceChecker::Impl::traceString(const char *id) const
460 return "Checking '" + appendPath(id) + "'";
464 std::string
465 TestReferenceChecker::Impl::appendPath(const char *id) const
467 std::string printId = (id != NULL) ? id : formatString("[%d]", seqIndex_);
468 return path_ + printId;
472 xmlNodePtr
473 TestReferenceChecker::Impl::findNode(const xmlChar *name, const char *id) const
475 if (currNode_ == NULL || currNode_->children == NULL)
477 return NULL;
479 const xmlChar *xmlId = reinterpret_cast<const xmlChar *>(id);
480 xmlNodePtr node = prevFoundNode_;
481 bool bWrap = true;
482 if (node != NULL)
484 if (id == NULL)
486 xmlChar *refId = xmlGetProp(node, cIdAttrName);
487 if (refId == NULL)
489 if (name == NULL || xmlStrcmp(node->name, name) == 0)
491 bWrap = false;
492 node = node->next;
493 if (node == NULL)
495 return NULL;
499 else
501 xmlFree(refId);
505 else
507 node = currNode_->children;
508 bWrap = false;
512 if (name == NULL || xmlStrcmp(node->name, name) == 0)
514 xmlChar *refId = xmlGetProp(node, cIdAttrName);
515 if (xmlId == NULL && refId == NULL)
517 return node;
519 if (refId != NULL)
521 if (xmlId != NULL && xmlStrcmp(refId, xmlId) == 0)
523 xmlFree(refId);
524 return node;
526 xmlFree(refId);
529 node = node->next;
530 if (bWrap && node == NULL)
532 node = currNode_->children;
535 while (node != NULL && node != prevFoundNode_);
536 return NULL;
540 xmlNodePtr
541 TestReferenceChecker::Impl::findOrCreateNode(const xmlChar *name,
542 const char *id,
543 bool *bFound)
545 *bFound = false;
546 xmlNodePtr node = findNode(name, id);
547 if (node != NULL)
549 *bFound = true;
550 prevFoundNode_ = node;
552 if (node == NULL)
554 if (bWrite_)
556 node = xmlNewTextChild(currNode_, NULL, name, NULL);
557 if (node != NULL && id != NULL)
559 const xmlChar *xmlId = reinterpret_cast<const xmlChar *>(id);
560 xmlAttrPtr prop = xmlNewProp(node, cIdAttrName, xmlId);
561 if (prop == NULL)
563 xmlFreeNode(node);
564 node = NULL;
567 if (node == NULL)
569 GMX_THROW(TestException("XML node creation failed"));
571 prevFoundNode_ = node;
573 else
575 ADD_FAILURE() << "Reference data item not found";
578 seqIndex_ = (id == NULL) ? seqIndex_+1 : 0;
580 return node;
584 std::string
585 TestReferenceChecker::Impl::processItem(const xmlChar *name, const char *id,
586 const char *value, bool *bFound)
588 xmlNodePtr node = findOrCreateNode(name, id, bFound);
589 if (node == NULL)
591 return std::string();
593 if (bWrite_ && !*bFound)
595 xmlNodeAddContent(node, reinterpret_cast<const xmlChar *>(value));
596 *bFound = true;
597 return std::string(value);
599 else
601 xmlChar *refXmlValue = xmlNodeGetContent(node);
602 std::string refValue(reinterpret_cast<const char *>(refXmlValue));
603 xmlFree(refXmlValue);
604 return refValue;
609 std::string
610 TestReferenceChecker::Impl::processItem(const xmlChar *name, const char *id,
611 const std::string &value, bool *bFound)
613 return processItem(name, id, value.c_str(), bFound);
617 bool
618 TestReferenceChecker::Impl::shouldIgnore() const
620 return currNode_ == NULL;
624 /********************************************************************
625 * TestReferenceData
628 TestReferenceData::TestReferenceData()
629 : impl_(new Impl(getReferenceDataMode(), false))
634 TestReferenceData::TestReferenceData(ReferenceDataMode mode)
635 : impl_(new Impl(mode, true))
640 TestReferenceData::~TestReferenceData()
645 bool TestReferenceData::isWriteMode() const
647 return impl_->bWrite_;
651 TestReferenceChecker TestReferenceData::rootChecker()
653 if (!isWriteMode() && !impl_->bInUse_ && impl_->refDoc_ == NULL)
655 ADD_FAILURE() << "Reference data file not found: "
656 << impl_->fullFilename_;
658 impl_->bInUse_ = true;
659 if (impl_->refDoc_ == NULL)
661 return TestReferenceChecker(new TestReferenceChecker::Impl(isWriteMode()));
663 xmlNodePtr rootNode = xmlDocGetRootElement(impl_->refDoc_);
664 // TODO: The default tolerance for double-precision builds that explicitly
665 // call checkFloat() may not be ideal.
666 return TestReferenceChecker(
667 new TestReferenceChecker::Impl("", rootNode, isWriteMode(),
668 impl_->bSelfTestMode_,
669 defaultRealTolerance()));
673 /********************************************************************
674 * TestReferenceChecker
677 TestReferenceChecker::TestReferenceChecker(Impl *impl)
678 : impl_(impl)
683 TestReferenceChecker::TestReferenceChecker(const TestReferenceChecker &other)
684 : impl_(new Impl(*other.impl_))
689 TestReferenceChecker &
690 TestReferenceChecker::operator=(const TestReferenceChecker &other)
692 impl_.reset(new Impl(*other.impl_));
693 return *this;
697 TestReferenceChecker::~TestReferenceChecker()
702 bool TestReferenceChecker::isWriteMode() const
704 return impl_->bWrite_;
708 void TestReferenceChecker::setDefaultTolerance(
709 const FloatingPointTolerance &tolerance)
711 impl_->defaultTolerance_ = tolerance;
715 bool TestReferenceChecker::checkPresent(bool bPresent, const char *id)
717 if (isWriteMode() || impl_->shouldIgnore())
719 return bPresent;
721 xmlNodePtr node = impl_->findNode(NULL, id);
722 bool bFound = (node != NULL);
723 if (bFound != bPresent)
725 ADD_FAILURE() << "Mismatch while checking reference data item '"
726 << impl_->appendPath(id) << "'\n"
727 << "Expected: " << (bPresent ? "it is present.\n" : "it is absent.\n")
728 << " Actual: " << (bFound ? "it is present." : "it is absent.");
730 if (bFound && bPresent)
732 impl_->prevFoundNode_ = node;
733 return true;
735 return false;
739 TestReferenceChecker TestReferenceChecker::checkCompound(const char *type, const char *id)
741 SCOPED_TRACE(impl_->traceString(id));
742 if (impl_->shouldIgnore())
744 return TestReferenceChecker(new Impl(isWriteMode()));
746 const xmlChar *xmlNodeName = reinterpret_cast<const xmlChar *>(type);
747 bool bFound;
748 xmlNodePtr newNode = impl_->findOrCreateNode(xmlNodeName, id, &bFound);
749 if (newNode == NULL)
751 return TestReferenceChecker(new Impl(isWriteMode()));
753 return TestReferenceChecker(
754 new Impl(impl_->appendPath(id), newNode, isWriteMode(),
755 impl_->bSelfTestMode_, impl_->defaultTolerance_));
759 void TestReferenceChecker::checkBoolean(bool value, const char *id)
761 if (impl_->shouldIgnore())
763 return;
765 SCOPED_TRACE(impl_->traceString(id));
766 bool bFound = false;
767 const char *strValue = value ? "true" : "false";
768 std::string refStrValue =
769 impl_->processItem(Impl::cBooleanNodeName, id, strValue, &bFound);
770 if (bFound)
772 EXPECT_EQ(refStrValue, strValue);
777 void TestReferenceChecker::checkString(const char *value, const char *id)
779 if (impl_->shouldIgnore())
781 return;
783 SCOPED_TRACE(impl_->traceString(id));
784 bool bFound = false;
785 std::string refStrValue =
786 impl_->processItem(Impl::cStringNodeName, id, value, &bFound);
787 if (bFound)
789 EXPECT_EQ(refStrValue, value);
794 void TestReferenceChecker::checkString(const std::string &value, const char *id)
796 checkString(value.c_str(), id);
800 void TestReferenceChecker::checkStringBlock(const std::string &value,
801 const char *id)
803 if (impl_->shouldIgnore())
805 return;
807 SCOPED_TRACE(impl_->traceString(id));
808 bool bFound;
809 xmlNodePtr node = impl_->findOrCreateNode(Impl::cStringNodeName, id, &bFound);
810 if (node == NULL)
812 return;
814 // An extra newline is written in the beginning to make lines align
815 // in the output xml (otherwise, the first line would be off by the length
816 // of the starting CDATA tag).
817 if (isWriteMode() && !bFound)
819 std::string adjustedValue = "\n" + value;
820 const xmlChar *xmlValue
821 = reinterpret_cast<const xmlChar *>(adjustedValue.c_str());
822 // TODO: Figure out if \r and \r\n can be handled without them changing
823 // to \n in the roundtrip
824 xmlNodePtr cdata
825 = xmlNewCDataBlock(node->doc, xmlValue,
826 static_cast<int>(adjustedValue.length()));
827 xmlAddChild(node, cdata);
829 else
831 xmlNodePtr cdata = node->children;
832 while (cdata != NULL && cdata->type != XML_CDATA_SECTION_NODE)
834 cdata = cdata->next;
836 if (cdata == NULL)
838 ADD_FAILURE() << "Invalid string block element";
839 return;
841 xmlChar *refXmlValue = xmlNodeGetContent(cdata);
842 if (refXmlValue[0] != '\n')
844 ADD_FAILURE() << "Invalid string block element";
845 xmlFree(refXmlValue);
846 return;
848 std::string refValue(reinterpret_cast<const char *>(refXmlValue + 1));
849 xmlFree(refXmlValue);
850 EXPECT_EQ(refValue, value);
855 void TestReferenceChecker::checkInteger(int value, const char *id)
857 if (impl_->shouldIgnore())
859 return;
861 SCOPED_TRACE(impl_->traceString(id));
862 bool bFound = false;
863 std::string strValue = formatString("%d", value);
864 std::string refStrValue =
865 impl_->processItem(Impl::cIntegerNodeName, id, strValue, &bFound);
866 if (bFound)
868 EXPECT_EQ(refStrValue, strValue);
872 void TestReferenceChecker::checkInt64(gmx_int64_t value, const char *id)
874 if (impl_->shouldIgnore())
876 return;
878 SCOPED_TRACE(impl_->traceString(id));
879 bool bFound = false;
880 std::string strValue = formatString("%" GMX_PRId64, value);
881 std::string refStrValue =
882 impl_->processItem(Impl::cInt64NodeName, id, strValue, &bFound);
883 if (bFound)
885 EXPECT_EQ(refStrValue, strValue);
889 void TestReferenceChecker::checkUInt64(gmx_uint64_t value, const char *id)
891 if (impl_->shouldIgnore())
893 return;
895 SCOPED_TRACE(impl_->traceString(id));
896 bool bFound = false;
897 std::string strValue = formatString("%" GMX_PRIu64, value);
898 std::string refStrValue =
899 impl_->processItem(Impl::cUInt64NodeName, id, strValue, &bFound);
900 if (bFound)
902 EXPECT_EQ(refStrValue, strValue);
906 void TestReferenceChecker::checkDouble(double value, const char *id)
908 if (impl_->shouldIgnore())
910 return;
912 SCOPED_TRACE(impl_->traceString(id));
913 bool bFound = false;
914 const int prec = std::numeric_limits<double>::digits10 + 2;
915 std::string strValue = formatString("%.*g", prec, value);
916 std::string refStrValue =
917 impl_->processItem(Impl::cRealNodeName, id, strValue, &bFound);
918 if (bFound)
920 char *endptr;
921 double refValue = std::strtod(refStrValue.c_str(), &endptr);
922 EXPECT_EQ('\0', *endptr);
923 if (impl_->bSelfTestMode_)
925 EXPECT_DOUBLE_EQ_TOL(refValue, value, impl_->defaultTolerance_)
926 << "String value: " << strValue << std::endl
927 << " Ref. string: " << refStrValue;
929 else
931 EXPECT_DOUBLE_EQ_TOL(refValue, value, impl_->defaultTolerance_);
937 void TestReferenceChecker::checkFloat(float value, const char *id)
939 if (impl_->shouldIgnore())
941 return;
943 SCOPED_TRACE(impl_->traceString(id));
944 bool bFound = false;
945 const int prec = std::numeric_limits<float>::digits10 + 2;
946 std::string strValue = formatString("%.*g", prec, value);
947 std::string refStrValue =
948 impl_->processItem(Impl::cRealNodeName, id, strValue, &bFound);
949 if (bFound)
951 char *endptr;
952 float refValue = static_cast<float>(std::strtod(refStrValue.c_str(), &endptr));
953 EXPECT_EQ('\0', *endptr);
954 if (impl_->bSelfTestMode_)
956 EXPECT_FLOAT_EQ_TOL(refValue, value, impl_->defaultTolerance_)
957 << "String value: " << strValue << std::endl
958 << " Ref. string: " << refStrValue;
960 else
962 EXPECT_FLOAT_EQ_TOL(refValue, value, impl_->defaultTolerance_);
968 void TestReferenceChecker::checkReal(float value, const char *id)
970 checkFloat(value, id);
974 void TestReferenceChecker::checkReal(double value, const char *id)
976 checkDouble(value, id);
980 void TestReferenceChecker::checkVector(const int value[3], const char *id)
982 TestReferenceChecker compound(checkCompound(Impl::cVectorType, id));
983 compound.checkInteger(value[0], "X");
984 compound.checkInteger(value[1], "Y");
985 compound.checkInteger(value[2], "Z");
989 void TestReferenceChecker::checkVector(const float value[3], const char *id)
991 TestReferenceChecker compound(checkCompound(Impl::cVectorType, id));
992 compound.checkReal(value[0], "X");
993 compound.checkReal(value[1], "Y");
994 compound.checkReal(value[2], "Z");
998 void TestReferenceChecker::checkVector(const double value[3], const char *id)
1000 TestReferenceChecker compound(checkCompound(Impl::cVectorType, id));
1001 compound.checkReal(value[0], "X");
1002 compound.checkReal(value[1], "Y");
1003 compound.checkReal(value[2], "Z");
1007 TestReferenceChecker
1008 TestReferenceChecker::checkSequenceCompound(const char *id, size_t length)
1010 TestReferenceChecker compound(checkCompound(Impl::cSequenceType, id));
1011 compound.checkInteger(static_cast<int>(length), Impl::cSequenceLengthName);
1012 return compound;
1015 } // namespace test
1016 } // namespace gmx