Generalize command-line test output file checking
authorTeemu Murtola <teemu.murtola@gmail.com>
Mon, 24 Apr 2017 15:34:29 +0000 (18:34 +0300)
committerTeemu Murtola <teemu.murtola@gmail.com>
Thu, 11 May 2017 14:48:08 +0000 (16:48 +0200)
Add a separate interface that gets passed a path to the temporary file
created by the program under test, instead of only supporting checking a
TextInputStream.  This makes it possible 1) to use old I/O routines such
as read_tps_conf() for the checks, and 2) to support binary files.

For now, all existing matchers still use streams; any matcher using the
old interface can be wrapped in the new implementation, and overloads
for setOutputFile() exist that do just that.  At least conftest.h could
be easily generalized with the new approach quite a bit.

Change-Id: Ibe71e22e568a54044d5fe64f88c43bb60cc530cb

src/testutils/CMakeLists.txt
src/testutils/cmdlinetest.cpp
src/testutils/cmdlinetest.h
src/testutils/filematchers.cpp [new file with mode: 0644]
src/testutils/filematchers.h [new file with mode: 0644]

index e48e716e5560fc2966e51e2448ff5a377e3db845..8d48a3a3c3d0e1f1895081ffbf6543f76249b8fa 100644 (file)
@@ -41,6 +41,7 @@ include_directories(BEFORE SYSTEM ${GMOCK_INCLUDE_DIRS})
 set(TESTUTILS_SOURCES
     cmdlinetest.cpp
     conftest.cpp
+    filematchers.cpp
     interactivetest.cpp
     loggertest.cpp
     mpi-printer.cpp
index 0bc6376dc8b0c7f0b2469c4ff699a89f68ead0b2..5d2ec6e5f57f9aa715fc0b7a69ddcd4d4e263cb3 100644 (file)
@@ -55,7 +55,6 @@
 #include "gromacs/commandline/cmdlineoptionsmodule.h"
 #include "gromacs/commandline/cmdlineprogramcontext.h"
 #include "gromacs/utility/arrayref.h"
-#include "gromacs/utility/filestream.h"
 #include "gromacs/utility/gmxassert.h"
 #include "gromacs/utility/strconvert.h"
 #include "gromacs/utility/stringstream.h"
@@ -63,9 +62,9 @@
 #include "gromacs/utility/textreader.h"
 #include "gromacs/utility/textwriter.h"
 
+#include "testutils/filematchers.h"
 #include "testutils/refdata.h"
 #include "testutils/testfilemanager.h"
-#include "testutils/textblockmatchers.h"
 
 namespace gmx
 {
@@ -249,27 +248,14 @@ class CommandLineTestHelper::Impl
         struct OutputFileInfo
         {
             OutputFileInfo(const char *option, const std::string &path,
-                           TextBlockMatcherPointer matcher)
+                           FileMatcherPointer matcher)
                 : option(option), path(path), matcher(move(matcher))
             {
             }
-            OutputFileInfo(OutputFileInfo &&other)
-                : option(std::move(other.option)), path(std::move(other.path)),
-                  matcher(std::move(other.matcher))
-            {
-            }
-
-            OutputFileInfo &operator=(OutputFileInfo &&other)
-            {
-                option  = std::move(other.option);
-                path    = std::move(other.path);
-                matcher = std::move(other.matcher);
-                return *this;
-            }
 
             std::string              option;
             std::string              path;
-            TextBlockMatcherPointer  matcher;
+            FileMatcherPointer       matcher;
         };
 
         typedef std::vector<OutputFileInfo>        OutputFileList;
@@ -354,6 +340,13 @@ void CommandLineTestHelper::setInputFileContents(
 void CommandLineTestHelper::setOutputFile(
         CommandLine *args, const char *option, const char *filename,
         const ITextBlockMatcherSettings &matcher)
+{
+    setOutputFile(args, option, filename, TextFileMatch(matcher));
+}
+
+void CommandLineTestHelper::setOutputFile(
+        CommandLine *args, const char *option, const char *filename,
+        const IFileMatcherSettings &matcher)
 {
     std::string suffix(filename);
     if (startsWith(filename, "."))
@@ -362,7 +355,7 @@ void CommandLineTestHelper::setOutputFile(
     }
     std::string fullFilename = impl_->fileManager_.getTemporaryFilePath(suffix);
     args->addOption(option, fullFilename);
-    impl_->outputFiles_.emplace_back(option, fullFilename, matcher.createMatcher());
+    impl_->outputFiles_.emplace_back(option, fullFilename, matcher.createFileMatcher());
 }
 
 void CommandLineTestHelper::checkOutputFiles(TestReferenceChecker checker) const
@@ -372,15 +365,11 @@ void CommandLineTestHelper::checkOutputFiles(TestReferenceChecker checker) const
         TestReferenceChecker                 outputChecker(
                 checker.checkCompound("OutputFiles", "Files"));
         Impl::OutputFileList::const_iterator outfile;
-        for (outfile = impl_->outputFiles_.begin();
-             outfile != impl_->outputFiles_.end();
-             ++outfile)
+        for (const auto &outfile : impl_->outputFiles_)
         {
             TestReferenceChecker fileChecker(
-                    outputChecker.checkCompound("File", outfile->option.c_str()));
-            TextInputFile        stream(outfile->path);
-            outfile->matcher->checkStream(&stream, &fileChecker);
-            stream.close();
+                    outputChecker.checkCompound("File", outfile.option.c_str()));
+            outfile.matcher->checkFile(outfile.path, &fileChecker);
         }
     }
 }
@@ -444,6 +433,13 @@ void CommandLineTestBase::setOutputFile(
     impl_->helper_.setOutputFile(&impl_->cmdline_, option, filename, matcher);
 }
 
+void CommandLineTestBase::setOutputFile(
+        const char *option, const char *filename,
+        const IFileMatcherSettings &matcher)
+{
+    impl_->helper_.setOutputFile(&impl_->cmdline_, option, filename, matcher);
+}
+
 CommandLine &CommandLineTestBase::commandLine()
 {
     return impl_->cmdline_;
index 589cdc22150325155f43b90bd786097bb75cc34c..5afb6c813e735d16e9708a2dc0f8cb63a0290c8c 100644 (file)
@@ -63,6 +63,7 @@ class ICommandLineOptionsModule;
 namespace test
 {
 
+class IFileMatcherSettings;
 class ITextBlockMatcherSettings;
 class TestFileManager;
 class TestReferenceChecker;
@@ -319,11 +320,16 @@ class CommandLineTestHelper
          *
          * If the output file is needed to trigger some computation, or is
          * unconditionally produced by the code under test, but the contents
-         * are not interesting for the test, use NoTextMatch as the matcher.
+         * are not interesting for the test, use NoContentsMatch as the matcher.
+         * Note that the existence of the output file is still verified.
          */
         void setOutputFile(CommandLine *args, const char *option,
                            const char *filename,
                            const ITextBlockMatcherSettings &matcher);
+        //! \copydoc setOutputFile(CommandLine *, const char *, const char *, const ITextBlockMatcherSettings &)
+        void setOutputFile(CommandLine *args, const char *option,
+                           const char *filename,
+                           const IFileMatcherSettings &matcher);
 
         /*! \brief
          * Checks output files added with setOutputFile() against reference
@@ -401,6 +407,13 @@ class CommandLineTestBase : public ::testing::Test
          */
         void setOutputFile(const char *option, const char *filename,
                            const ITextBlockMatcherSettings &matcher);
+        /*! \brief
+         * Sets an output file parameter and adds it to the set of tested files.
+         *
+         * \see CommandLineTestHelper::setOutputFile()
+         */
+        void setOutputFile(const char *option, const char *filename,
+                           const IFileMatcherSettings &matcher);
 
         /*! \brief
          * Returns the internal CommandLine object used to construct the
diff --git a/src/testutils/filematchers.cpp b/src/testutils/filematchers.cpp
new file mode 100644 (file)
index 0000000..cd68d38
--- /dev/null
@@ -0,0 +1,99 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2017, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \internal \file
+ * \brief
+ * Implements classes from filematchers.h.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \ingroup module_testutils
+ */
+#include "gmxpre.h"
+
+#include "filematchers.h"
+
+#include "gromacs/utility/filestream.h"
+
+#include "testutils/textblockmatchers.h"
+
+namespace gmx
+{
+namespace test
+{
+
+namespace
+{
+
+class TextFileMatcher : public IFileMatcher
+{
+    public:
+        explicit TextFileMatcher(TextBlockMatcherPointer matcher)
+            : matcher_(std::move(matcher))
+        {
+        }
+
+        virtual void checkFile(const std::string    &path,
+                               TestReferenceChecker *checker)
+        {
+            TextInputFile stream(path);
+            matcher_->checkStream(&stream, checker);
+            stream.close();
+        }
+
+    private:
+        TextBlockMatcherPointer matcher_;
+};
+
+}       // namespace
+
+IFileMatcher::~IFileMatcher()
+{
+}
+
+IFileMatcherSettings::~IFileMatcherSettings()
+{
+}
+
+FileMatcherPointer TextFileMatch::createFileMatcher() const
+{
+    return FileMatcherPointer(new TextFileMatcher(streamSettings_.createMatcher()));
+}
+
+FileMatcherPointer NoContentsMatch::createFileMatcher() const
+{
+    return TextFileMatch(NoTextMatch()).createFileMatcher();
+}
+
+} // namespace test
+} // namespace gmx
diff --git a/src/testutils/filematchers.h b/src/testutils/filematchers.h
new file mode 100644 (file)
index 0000000..df86c13
--- /dev/null
@@ -0,0 +1,156 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2017, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \libinternal \file
+ * \brief
+ * Declares utility classes for testing file contents.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \inlibraryapi
+ * \ingroup module_testutils
+ */
+#ifndef GMX_TESTUTILS_FILEMATCHERS_H
+#define GMX_TESTUTILS_FILEMATCHERS_H
+
+#include <memory>
+#include <string>
+
+namespace gmx
+{
+namespace test
+{
+
+class ITextBlockMatcherSettings;
+class TestReferenceChecker;
+
+/*! \libinternal \brief
+ * Represents a file matcher, matching file contents against reference (or
+ * other) data.
+ *
+ * Typical pattern of declaring such matchers is to
+ *  - Create a factory that implements IFileMatcherSettings,
+ *  - Make that factory provide any necessary parameters that the matcher needs,
+ *    using a "named parameter" idiom (see XvgMatch for an example), and
+ *  - Make the factory create and return an instance of an internal
+ *    implementation class that implements IFileMatcher and provides
+ *    the actual matching logic.
+ *
+ * Any method that then wants to accept a matcher can accept a
+ * IFileMatcherSettings.
+ *
+ * \inlibraryapi
+ * \ingroup module_testutils
+ */
+class IFileMatcher
+{
+    public:
+        virtual ~IFileMatcher();
+
+        /*! \brief
+         * Matches contents of a file.
+         *
+         * \param  path     Path to the file to match.
+         * \param  checker  Checker to use for matching.
+         *
+         * The method can change the state of the provided checker (e.g., by
+         * changing the default tolerance).
+         * The caller is responsible of providing a checker where such state
+         * changes do not matter.
+         */
+        virtual void checkFile(const std::string    &path,
+                               TestReferenceChecker *checker) = 0;
+};
+
+//! Smart pointer for managing a IFileMatcher.
+typedef std::unique_ptr<IFileMatcher> FileMatcherPointer;
+
+/*! \libinternal \brief
+ * Represents a factory for creating a file matcher.
+ *
+ * See derived classes for available matchers.  Each derived class represents
+ * one type of matcher (see IFileMatcher), and provides any methods
+ * necessary to pass parameters to such a matcher.  Methods that accept a
+ * matcher can then take in this interface, and call createFileMatcher() to use
+ * the matcher that the caller of the method specifies.
+ *
+ * \inlibraryapi
+ * \ingroup module_testutils
+ */
+class IFileMatcherSettings
+{
+    public:
+        //! Factory method that constructs the matcher after parameters are set.
+        virtual FileMatcherPointer createFileMatcher() const = 0;
+
+    protected:
+        virtual ~IFileMatcherSettings();
+};
+
+/*! \libinternal \brief
+ * Use a ITextStreamMatcher for matching the contents.
+ *
+ * \inlibraryapi
+ * \ingroup module_testutils
+ */
+class TextFileMatch : public IFileMatcherSettings
+{
+    public:
+        //! Creates a matcher to match contents with given text matcher.
+        explicit TextFileMatch(const ITextBlockMatcherSettings &streamSettings)
+            : streamSettings_(streamSettings)
+        {
+        }
+
+        virtual FileMatcherPointer createFileMatcher() const;
+
+    private:
+        const ITextBlockMatcherSettings &streamSettings_;
+};
+
+/*! \libinternal \brief
+ * Do not check the contents of the file.
+ *
+ * \inlibraryapi
+ * \ingroup module_testutils
+ */
+class NoContentsMatch : public IFileMatcherSettings
+{
+    public:
+        virtual FileMatcherPointer createFileMatcher() const;
+};
+
+} // namespace test
+} // namespace gmx
+
+#endif