Tests for md5 computation on open files
authorMark Abraham <mark.j.abraham@gmail.com>
Thu, 27 Sep 2018 09:24:13 +0000 (11:24 +0200)
committerMark Abraham <mark.j.abraham@gmail.com>
Mon, 1 Oct 2018 11:58:05 +0000 (13:58 +0200)
Needed to support planned refactoring.

Also some improvements to filename handling in futil to remove
arbitrary small limit, so that real-world filename lengths will work.

Change-Id: I8a353bfc76bde1bffa320ed398d4bc7ac6e3859b

src/gromacs/fileio/tests/CMakeLists.txt
src/gromacs/fileio/tests/filemd5.cpp [new file with mode: 0644]
src/gromacs/utility/futil.cpp

index d800b19cdf9efd81f57156cc8bc14a6a38d64cf2..e9883f837e46bd12f3a47408056b1cd9eb059d4c 100644 (file)
@@ -1,7 +1,7 @@
 #
 # This file is part of the GROMACS molecular simulation package.
 #
-# Copyright (c) 2013,2014,2015,2016, by the GROMACS development team, led by
+# Copyright (c) 2013,2014,2015,2016,2018, 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.
@@ -34,6 +34,7 @@
 
 set(test_sources
     confio.cpp
+    filemd5.cpp
     readinp.cpp
     )
 if (GMX_USE_TNG)
diff --git a/src/gromacs/fileio/tests/filemd5.cpp b/src/gromacs/fileio/tests/filemd5.cpp
new file mode 100644 (file)
index 0000000..7afad2e
--- /dev/null
@@ -0,0 +1,126 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2018, 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
+ * Tests for routines for computing MD5 sums on files.
+ *
+ * \author Mark Abraham <mark.j.abraham@gmail.com>
+ * \ingroup module_fileio
+ */
+#include "gmxpre.h"
+
+#include <cstdio>
+
+#include <algorithm>
+#include <numeric>
+#include <string>
+#include <vector>
+
+#include <gmock/gmock.h>
+#include <gtest/gtest.h>
+
+#include "gromacs/fileio/gmxfio.h"
+#include "gromacs/utility/arrayref.h"
+
+#include "testutils/refdata.h"
+#include "testutils/testfilemanager.h"
+
+namespace gmx
+{
+namespace test
+{
+namespace
+{
+
+class FileMD5Test : public ::testing::Test
+{
+    public:
+        void prepareFile(int lengthInBytes)
+        {
+            // Fill some memory with some arbitrary bits.
+            std::vector<char> data(lengthInBytes);
+            std::iota(data.begin(), data.end(), 1);
+            // Binary mode ensures it works the same on all OS
+            FILE *fp = fopen(filename_.c_str(), "wb");
+            fwrite(data.data(), sizeof(char), data.size(), fp);
+            fclose(fp);
+        }
+        ~FileMD5Test() override
+        {
+            if (file_)
+            {
+                gmx_fio_close(file_);
+            }
+        }
+        TestFileManager fileManager_;
+        // Make sure the file extension is one that gmx_fio_open will
+        // recognize to open as binary.
+        std::string     filename_ = fileManager_.getTemporaryFilePath("data.edr");
+        t_fileio       *file_     = nullptr;
+};
+
+TEST_F(FileMD5Test, CanComputeMD5)
+{
+    prepareFile(1000);
+    file_ = gmx_fio_open(filename_.c_str(), "r+");
+
+    unsigned char           digest[16] = {0};
+    ArrayRef<unsigned char> digestView(digest);
+    // Chosen to be less than the full file length
+    gmx_off_t               offset             = 64;
+    gmx_off_t               expectedLength     = 64;
+    gmx_off_t               lengthActuallyRead = gmx_fio_get_file_md5(file_, offset, digest);
+
+    EXPECT_EQ(expectedLength, lengthActuallyRead);
+    // Did we compute an actual reproducible checksum?
+    auto total = std::accumulate(digestView.begin(), digestView.end(), 0);
+    EXPECT_EQ(2111, total);
+}
+
+TEST_F(FileMD5Test, ReturnsErrorIfFileModeIsWrong)
+{
+    prepareFile(1000);
+    file_ = gmx_fio_open(filename_.c_str(), "r");
+
+    unsigned char           digest[16] = {0};
+    ArrayRef<unsigned char> digestView(digest);
+    gmx_off_t               offset             = 100;
+    gmx_off_t               lengthActuallyRead = gmx_fio_get_file_md5(file_, offset, digest);
+    EXPECT_EQ(-1, lengthActuallyRead);
+}
+
+} // namespace
+} // namespace test
+} // namespace gmx
index c2178ba974b6d5fc8c99e7ce129587ea7a92a49d..38533a71369774a9623b7b36f6495203a772f572 100644 (file)
@@ -441,8 +441,8 @@ FILE *gmx_ffopen(const char *file, const char *mode)
 #ifdef SKIP_FFOPS
     return fopen(file, mode);
 #else
-    FILE    *ff = nullptr;
-    char     buf[256], *bufsize = nullptr, *ptr;
+    FILE    *ff      = nullptr;
+    char    *bufsize = nullptr, *ptr;
     gmx_bool bRead;
     int      bs;
 
@@ -457,12 +457,11 @@ FILE *gmx_ffopen(const char *file, const char *mode)
     }
 
     bRead = (mode[0] == 'r' && mode[1] != '+');
-    strcpy(buf, file);
-    if (!bRead || gmx_fexist(buf))
+    if (!bRead || gmx_fexist(file))
     {
-        if ((ff = fopen(buf, mode)) == nullptr)
+        if ((ff = fopen(file, mode)) == nullptr)
         {
-            gmx_file(buf);
+            gmx_file(file);
         }
         /* Check whether we should be using buffering (default) or not
          * (for debugging)
@@ -494,17 +493,19 @@ FILE *gmx_ffopen(const char *file, const char *mode)
     }
     else
     {
-        sprintf(buf, "%s.Z", file);
-        if (gmx_fexist(buf))
+        std::string compressedFileName = file;
+        compressedFileName += ".Z";
+        if (gmx_fexist(compressedFileName.c_str()))
         {
-            ff = uncompress(buf, mode);
+            ff = uncompress(compressedFileName.c_str(), mode);
         }
         else
         {
-            sprintf(buf, "%s.gz", file);
-            if (gmx_fexist(buf))
+            compressedFileName  = file;
+            compressedFileName += ".gz";
+            if (gmx_fexist(compressedFileName.c_str()))
             {
-                ff = gunzip(buf, mode);
+                ff = gunzip(compressedFileName.c_str(), mode);
             }
             else
             {