Make some unit tests skip file system access
authorTeemu Murtola <teemu.murtola@gmail.com>
Wed, 24 Jun 2015 09:16:54 +0000 (12:16 +0300)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Mon, 29 Jun 2015 11:08:38 +0000 (13:08 +0200)
Make unit tests for FileNameOptionManager not use the file system.
Introduce a FileInputRedirectorInterface to support mocking file
existence checks, and use a mock implementation in the tests.

These particular tests did quite a bit of file system access, and the
speedup is only a few ms, although significant percentually (something
like 80%).  But there can be tests where this has more effect, and this
approach provides a starting point for more work on eliminating
unnecessary file system access from the tests.

The main benefit is clearer and more robust test code, as it is no
longer necessary to construct actual files and ensure that they do not
conflict with other tests or cause issues if the test crashes or such.

Change-Id: Ib9a171331e988fa7e74b16078164f477f8296c6e

src/gromacs/options/filenameoptionmanager.cpp
src/gromacs/options/filenameoptionmanager.h
src/gromacs/options/tests/filenameoptionmanager.cpp
src/gromacs/utility/fileredirector.cpp
src/gromacs/utility/fileredirector.h
src/testutils/CMakeLists.txt
src/testutils/testfileredirector.cpp [new file with mode: 0644]
src/testutils/testfileredirector.h [new file with mode: 0644]
src/testutils/testutils-doc.h

index fc630dfa25613435419e8dc505094bf2a29d2e1d..4f0767895be9a098fbc822e51808e2b3b130bed1 100644 (file)
@@ -53,7 +53,7 @@
 #include "gromacs/options/options.h"
 #include "gromacs/utility/arrayref.h"
 #include "gromacs/utility/exceptions.h"
-#include "gromacs/utility/file.h"
+#include "gromacs/utility/fileredirector.h"
 #include "gromacs/utility/path.h"
 #include "gromacs/utility/stringutil.h"
 
@@ -79,15 +79,16 @@ const char *const c_compressedExtensions[] =
  * The first match is returned.
  * Returns an empty string if no existing file is found.
  */
-std::string findExistingExtension(const std::string        &prefix,
-                                  const FileNameOptionInfo &option)
+std::string findExistingExtension(const std::string                  &prefix,
+                                  const FileNameOptionInfo           &option,
+                                  const FileInputRedirectorInterface *redirector)
 {
     ConstArrayRef<int>                 types = option.fileTypes();
     ConstArrayRef<int>::const_iterator i;
     for (i = types.begin(); i != types.end(); ++i)
     {
         std::string testFilename(prefix + ftp2ext_with_dot(*i));
-        if (File::exists(testFilename))
+        if (redirector->fileExists(testFilename))
         {
             return testFilename;
         }
@@ -109,12 +110,18 @@ std::string findExistingExtension(const std::string        &prefix,
 class FileNameOptionManager::Impl
 {
     public:
-        Impl() : bInputCheckingDisabled_(false) {}
+        Impl()
+            : redirector_(&defaultFileInputRedirector()),
+              bInputCheckingDisabled_(false)
+        {
+        }
 
+        //! Redirector for file existence checks.
+        const FileInputRedirectorInterface *redirector_;
         //! Global default file name, if set.
-        std::string     defaultFileName_;
+        std::string                         defaultFileName_;
         //! Whether input option processing has been disabled.
-        bool            bInputCheckingDisabled_;
+        bool                                bInputCheckingDisabled_;
 };
 
 /********************************************************************
@@ -130,6 +137,12 @@ FileNameOptionManager::~FileNameOptionManager()
 {
 }
 
+void FileNameOptionManager::setInputRedirector(
+        const FileInputRedirectorInterface *redirector)
+{
+    impl_->redirector_ = redirector;
+}
+
 void FileNameOptionManager::disableInputOptionChecking(bool bDisable)
 {
     impl_->bInputCheckingDisabled_ = bDisable;
@@ -169,7 +182,7 @@ std::string FileNameOptionManager::completeFileName(
     const int fileType = fn2ftp(value.c_str());
     if (bInput && !impl_->bInputCheckingDisabled_)
     {
-        if (fileType == efNR && File::exists(value))
+        if (fileType == efNR && impl_->redirector_->fileExists(value))
         {
             ConstArrayRef<const char *>                 compressedExtensions(c_compressedExtensions);
             ConstArrayRef<const char *>::const_iterator ext;
@@ -196,7 +209,8 @@ std::string FileNameOptionManager::completeFileName(
         }
         else if (fileType == efNR)
         {
-            std::string processedValue = findExistingExtension(value, option);
+            const std::string processedValue
+                = findExistingExtension(value, option, impl_->redirector_);
             if (!processedValue.empty())
             {
                 return processedValue;
@@ -225,7 +239,7 @@ std::string FileNameOptionManager::completeFileName(
             {
                 // TODO: Treat also library files.
             }
-            else if (!bAllowMissing && !File::exists(value))
+            else if (!bAllowMissing && !impl_->redirector_->fileExists(value))
             {
                 std::string message
                     = formatString("File '%s' does not exist or is not accessible.",
@@ -264,7 +278,8 @@ std::string FileNameOptionManager::completeDefaultFileName(
     const bool        bAllowMissing = option.allowMissing();
     if (bInput)
     {
-        std::string completedName = findExistingExtension(realPrefix, option);
+        const std::string completedName
+            = findExistingExtension(realPrefix, option, impl_->redirector_);
         if (!completedName.empty())
         {
             return completedName;
index ac588a00606dde9f0f530f12880fe684d510e7c0..37f9f7c3bb880673bc18ee49e617deacc4da3ab8 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2014, by the GROMACS development team, led by
+ * Copyright (c) 2014,2015, 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.
@@ -51,6 +51,7 @@
 namespace gmx
 {
 
+class FileInputRedirectorInterface;
 class FileNameOptionInfo;
 class Options;
 
@@ -84,6 +85,24 @@ class FileNameOptionManager : public OptionManagerInterface
         FileNameOptionManager();
         virtual ~FileNameOptionManager();
 
+        /*! \brief
+         * Redirects file existence checks.
+         *
+         * \param[in] redirector  File redirector to use for existence checks.
+         *
+         * The manager checks for existence of various files on the file system
+         * to complete file extensions.  This method can be used to redirect
+         * those checks to an alternative implementation.
+         *
+         * This is used for unit tests to more easily control the result of the
+         * checks and to keep the tests as fast as possible by avoiding real
+         * file system access.  To keep implementation options open, behavior
+         * with `redirector == NULL` is undefined and should not be relied on.
+         * For tests, there should only be need to call this a single time,
+         * right after creating the manager.
+         */
+        void setInputRedirector(const FileInputRedirectorInterface *redirector);
+
         /*! \brief
          * Disables special input file option handling.
          *
index bb01302e8629ef62e00e298f70615066759efeaa..30b6697e2f5ac260b7d410ec622782718a040294 100644 (file)
 #include "gromacs/options/options.h"
 #include "gromacs/options/optionsassigner.h"
 #include "gromacs/utility/exceptions.h"
-#include "gromacs/utility/file.h"
-#include "gromacs/utility/path.h"
 
 #include "testutils/testasserts.h"
-#include "testutils/testfilemanager.h"
+#include "testutils/testfileredirector.h"
 
 namespace
 {
@@ -67,19 +65,18 @@ class FileNameOptionManagerTest : public ::testing::Test
         FileNameOptionManagerTest()
             : options_(NULL, NULL)
         {
+            manager_.setInputRedirector(&redirector_);
             options_.addManager(&manager_);
         }
 
-        std::string createDummyFile(const char *suffix)
+        void addExistingFile(const char *filename)
         {
-            std::string filename(tempFiles_.getTemporaryFilePath(suffix));
-            gmx::File::writeFileFromString(filename, "Dummy file");
-            return filename;
+            redirector_.addExistingFile(filename);
         }
 
-        gmx::FileNameOptionManager manager_;
-        gmx::Options               options_;
-        gmx::test::TestFileManager tempFiles_;
+        gmx::test::TestFileInputRedirector redirector_;
+        gmx::FileNameOptionManager         manager_;
+        gmx::Options                       options_;
 };
 
 /********************************************************************
@@ -251,8 +248,7 @@ TEST_F(FileNameOptionManagerTest, AcceptsMissingRequiredInputFileIfSpecified)
 
 TEST_F(FileNameOptionManagerTest, AddsMissingExtensionBasedOnExistingFile)
 {
-    std::string filename(createDummyFile(".trr"));
-    std::string inputValue(gmx::Path::stripExtension(filename));
+    addExistingFile("testfile.trr");
 
     std::string value;
     ASSERT_NO_THROW_GMX(options_.addOption(
@@ -262,26 +258,25 @@ TEST_F(FileNameOptionManagerTest, AddsMissingExtensionBasedOnExistingFile)
     gmx::OptionsAssigner assigner(&options_);
     EXPECT_NO_THROW_GMX(assigner.start());
     EXPECT_NO_THROW_GMX(assigner.startOption("f"));
-    EXPECT_NO_THROW_GMX(assigner.appendValue(inputValue));
+    EXPECT_NO_THROW_GMX(assigner.appendValue("testfile"));
     EXPECT_NO_THROW_GMX(assigner.finishOption());
     EXPECT_NO_THROW_GMX(assigner.finish());
     EXPECT_NO_THROW_GMX(options_.finish());
 
-    EXPECT_EQ(filename, value);
+    EXPECT_EQ("testfile.trr", value);
 }
 
 TEST_F(FileNameOptionManagerTest,
        AddsMissingExtensionForRequiredDefaultNameBasedOnExistingFile)
 {
-    std::string filename(createDummyFile(".trr"));
-    std::string inputValue(gmx::Path::stripExtension(filename));
+    addExistingFile("testfile.trr");
 
     std::string value;
     ASSERT_NO_THROW_GMX(options_.addOption(
                                 FileNameOption("f").store(&value).required()
                                     .filetype(gmx::eftTrajectory).inputFile()
-                                    .defaultBasename(inputValue.c_str())));
-    EXPECT_EQ(inputValue + ".xtc", value);
+                                    .defaultBasename("testfile")));
+    EXPECT_EQ("testfile.xtc", value);
 
     gmx::OptionsAssigner assigner(&options_);
     EXPECT_NO_THROW_GMX(assigner.start());
@@ -290,20 +285,19 @@ TEST_F(FileNameOptionManagerTest,
     EXPECT_NO_THROW_GMX(assigner.finish());
     EXPECT_NO_THROW_GMX(options_.finish());
 
-    EXPECT_EQ(filename, value);
+    EXPECT_EQ("testfile.trr", value);
 }
 
 TEST_F(FileNameOptionManagerTest,
        AddsMissingExtensionForOptionalDefaultNameBasedOnExistingFile)
 {
-    std::string filename(createDummyFile(".trr"));
-    std::string inputValue(gmx::Path::stripExtension(filename));
+    addExistingFile("testfile.trr");
 
     std::string value;
     ASSERT_NO_THROW_GMX(options_.addOption(
                                 FileNameOption("f").store(&value)
                                     .filetype(gmx::eftTrajectory).inputFile()
-                                    .defaultBasename(inputValue.c_str())));
+                                    .defaultBasename("testfile")));
 
     gmx::OptionsAssigner assigner(&options_);
     EXPECT_NO_THROW_GMX(assigner.start());
@@ -312,14 +306,13 @@ TEST_F(FileNameOptionManagerTest,
     EXPECT_NO_THROW_GMX(assigner.finish());
     EXPECT_NO_THROW_GMX(options_.finish());
 
-    EXPECT_EQ(filename, value);
+    EXPECT_EQ("testfile.trr", value);
 }
 
 TEST_F(FileNameOptionManagerTest,
        AddsMissingExtensionForRequiredFromDefaultNameOptionBasedOnExistingFile)
 {
-    std::string filename(createDummyFile(".trr"));
-    std::string inputValue(gmx::Path::stripExtension(filename));
+    addExistingFile("testfile.trr");
 
     std::string value;
     ASSERT_NO_THROW_GMX(options_.addOption(
@@ -332,19 +325,18 @@ TEST_F(FileNameOptionManagerTest,
     gmx::OptionsAssigner assigner(&options_);
     EXPECT_NO_THROW_GMX(assigner.start());
     EXPECT_NO_THROW_GMX(assigner.startOption("deffnm"));
-    EXPECT_NO_THROW_GMX(assigner.appendValue(inputValue));
+    EXPECT_NO_THROW_GMX(assigner.appendValue("testfile"));
     EXPECT_NO_THROW_GMX(assigner.finishOption());
     EXPECT_NO_THROW_GMX(assigner.finish());
     EXPECT_NO_THROW_GMX(options_.finish());
 
-    EXPECT_EQ(filename, value);
+    EXPECT_EQ("testfile.trr", value);
 }
 
 TEST_F(FileNameOptionManagerTest,
        AddsMissingExtensionForOptionalFromDefaultNameOptionBasedOnExistingFile)
 {
-    std::string filename(createDummyFile(".trr"));
-    std::string inputValue(gmx::Path::stripExtension(filename));
+    addExistingFile("testfile.trr");
 
     std::string value;
     ASSERT_NO_THROW_GMX(options_.addOption(
@@ -356,14 +348,14 @@ TEST_F(FileNameOptionManagerTest,
     gmx::OptionsAssigner assigner(&options_);
     EXPECT_NO_THROW_GMX(assigner.start());
     EXPECT_NO_THROW_GMX(assigner.startOption("deffnm"));
-    EXPECT_NO_THROW_GMX(assigner.appendValue(inputValue));
+    EXPECT_NO_THROW_GMX(assigner.appendValue("testfile"));
     EXPECT_NO_THROW_GMX(assigner.finishOption());
     EXPECT_NO_THROW_GMX(assigner.startOption("f"));
     EXPECT_NO_THROW_GMX(assigner.finishOption());
     EXPECT_NO_THROW_GMX(assigner.finish());
     EXPECT_NO_THROW_GMX(options_.finish());
 
-    EXPECT_EQ(filename, value);
+    EXPECT_EQ("testfile.trr", value);
 }
 
 } // namespace
index 9a1406122d1d51a44d3b6b38067cfa3f5fa2efe4..405702bd7280e180f17002a0abd6713d32adab96 100644 (file)
 namespace gmx
 {
 
+FileInputRedirectorInterface::~FileInputRedirectorInterface()
+{
+}
+
 FileOutputRedirectorInterface::~FileOutputRedirectorInterface()
 {
 }
@@ -55,6 +59,23 @@ FileOutputRedirectorInterface::~FileOutputRedirectorInterface()
 namespace
 {
 
+/*! \internal
+ * \brief
+ * Implements the redirector returned by defaultFileInputRedirector().
+ *
+ * Does not redirect anything, but uses the file system as requested.
+ *
+ * \ingroup module_utility
+ */
+class DefaultInputRedirector : public FileInputRedirectorInterface
+{
+    public:
+        virtual bool fileExists(const char *filename) const
+        {
+            return File::exists(filename);
+        }
+};
+
 /*! \internal
  * \brief
  * Implements the redirector returned by defaultFileOutputRedirector().
@@ -80,6 +101,12 @@ class DefaultOutputRedirector : public FileOutputRedirectorInterface
 }   // namespace
 
 //! \cond libapi
+FileInputRedirectorInterface &defaultFileInputRedirector()
+{
+    static DefaultInputRedirector instance;
+    return instance;
+}
+
 FileOutputRedirectorInterface &defaultFileOutputRedirector()
 {
     static DefaultOutputRedirector instance;
index f65ce72c78aa9d6c53308ccd9e770848a0e733cf..d0a7c6e47edf55ffaac4a34be22d3f100ce53acb 100644 (file)
 namespace gmx
 {
 
+/*! \libinternal \brief
+ * Allows overriding file existence checks from code that supports it.
+ *
+ * The calling code should take in this interface and use the methods in it
+ * all file system operations that need to support this redirection.
+ * By default, the code can then use defaultFileInputRedirector() in case no
+ * redirection is needed.
+ *
+ * This allows tests to override the file existence checks without actually
+ * using the file system.
+ *
+ * With some further refactoring of the File class, this could also support
+ * redirecting input files from in-memory buffers as well, but for now the
+ * current capabilities are sufficient.
+ *
+ * \inlibraryapi
+ * \ingroup module_utility
+ */
+class FileInputRedirectorInterface
+{
+    public:
+        virtual ~FileInputRedirectorInterface();
+
+        /*! \brief
+         * Checks whether the provided path exists (and is a file).
+         */
+        virtual bool fileExists(const char *filename) const = 0;
+
+        //! Convenience method to check file existence using an std::string path.
+        bool fileExists(const std::string &filename) const
+        {
+            return fileExists(filename.c_str());
+        }
+};
+
 /*! \libinternal \brief
  * Allows capturing `stdout` and file output from code that supports it.
  *
@@ -91,12 +126,25 @@ class FileOutputRedirectorInterface
 };
 
 //! \cond libapi
+/*! \brief
+ * Returns default implementation for FileInputRedirectorInterface.
+ *
+ * The returned implementation does not redirect anything, but just uses the
+ * file system normally.
+ *
+ * Does not throw.
+ *
+ * \ingroup module_utility
+ */
+FileInputRedirectorInterface &defaultFileInputRedirector();
 /*! \brief
  * Returns default implementation for FileOutputRedirectorInterface.
  *
  * The returned implementation does not redirect anything, but just opens the
  * files at requested locations.
  *
+ * Does not throw.
+ *
  * \ingroup module_utility
  */
 FileOutputRedirectorInterface &defaultFileOutputRedirector();
index 381039d8c1bf0e51a6a451b6e6b969dd96194b21..131d9351f5de5ed471f404d874e6d7fa806807b4 100644 (file)
 
 include_directories(BEFORE SYSTEM ${GMOCK_INCLUDE_DIRS})
 include_directories(${LIBXML2_INCLUDE_DIR})
-file(GLOB TESTUTILS_SOURCES *.cpp)
+set(TESTUTILS_SOURCES
+    cmdlinetest.cpp
+    integrationtests.cpp
+    mpi-printer.cpp
+    refdata.cpp
+    stringtest.cpp
+    testasserts.cpp
+    testfilemanager.cpp
+    testfileredirector.cpp
+    testinit.cpp
+    testoptions.cpp)
 
 add_library(testutils STATIC ${UNITTEST_TARGET_OPTIONS} ${TESTUTILS_SOURCES})
 set(TESTUTILS_LIBS testutils)
diff --git a/src/testutils/testfileredirector.cpp b/src/testutils/testfileredirector.cpp
new file mode 100644 (file)
index 0000000..62ea648
--- /dev/null
@@ -0,0 +1,73 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2015, 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 testfileredirector.h.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \ingroup module_testutils
+ */
+#include "gmxpre.h"
+
+#include "testfileredirector.h"
+
+#include <set>
+#include <string>
+
+namespace gmx
+{
+namespace test
+{
+
+TestFileInputRedirector::TestFileInputRedirector()
+{
+}
+
+TestFileInputRedirector::~TestFileInputRedirector()
+{
+}
+
+void TestFileInputRedirector::addExistingFile(const char *filename)
+{
+    existingFiles_.insert(filename);
+}
+
+bool TestFileInputRedirector::fileExists(const char *filename) const
+{
+    return existingFiles_.count(filename) > 0;
+}
+
+} // namespace test
+} // namespace gmx
diff --git a/src/testutils/testfileredirector.h b/src/testutils/testfileredirector.h
new file mode 100644 (file)
index 0000000..711c6c8
--- /dev/null
@@ -0,0 +1,93 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2015, 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 generic mock implementations for interfaces in fileredirector.h.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \inlibraryapi
+ * \ingroup module_testutils
+ */
+#ifndef GMX_TESTUTILS_TESTFILEREDIRECTOR_H
+#define GMX_TESTUTILS_TESTFILEREDIRECTOR_H
+
+#include <set>
+#include <string>
+
+#include "gromacs/utility/classhelpers.h"
+#include "gromacs/utility/fileredirector.h"
+
+namespace gmx
+{
+namespace test
+{
+
+/*! \libinternal \brief
+ * In-memory implementation for FileInputRedirectorInterface for tests.
+ *
+ * By default, this implementation will return `false` for all file existence
+ * checks.  To return `true` for a specific path, use addExistingFile().
+ *
+ * \inlibraryapi
+ * \ingroup module_testutils
+ */
+class TestFileInputRedirector : public FileInputRedirectorInterface
+{
+    public:
+        TestFileInputRedirector();
+        virtual ~TestFileInputRedirector();
+
+        /*! \brief
+         * Marks the provided path as an existing file.
+         *
+         * \throws std::bad_alloc if out of memory.
+         *
+         * Further checks for existence of the given path will return `true`.
+         */
+        void addExistingFile(const char *filename);
+
+        // From FileInputRedirectorInterface
+        virtual bool fileExists(const char *filename) const;
+
+    private:
+        std::set<std::string> existingFiles_;
+
+        GMX_DISALLOW_COPY_AND_ASSIGN(TestFileInputRedirector);
+};
+
+} // namespace test
+} // namespace gmx
+
+#endif
index 50c622ab53d49986682dc0ccb270f922c66fe4ce..a99fccd9dac05911488d62b1c6d4d85fa27b274a 100644 (file)
@@ -51,6 +51,9 @@
  *  - gmx::test::TestFileManager (in testfilemanager.h) provides functionality
  *    for locating test input files from the source directory and managing
  *    temporary files that need to be created during the test.
+ *  - gmx::test::TestFileInputRedirector (in testfileredirector.h) provides
+ *    functionality for capturing file existence checks in code that uses
+ *    gmx::FileInputRedirectorInterface.
  *  - #GMX_TEST_OPTIONS macro provides facilities for adding custom command
  *    line options for the test binary.
  *  - testasserts.h provides several custom test assertions for better