Add basic tests for reading/writing structure files
authorTeemu Murtola <teemu.murtola@gmail.com>
Wed, 22 Jul 2015 04:14:14 +0000 (07:14 +0300)
committerTeemu Murtola <teemu.murtola@gmail.com>
Wed, 22 Jul 2015 13:52:16 +0000 (16:52 +0300)
For each of the supported format (except for tpr), make a simple
roundtrip test, starting from a dummy topology, and doing a
write-read-write cycle, checking that the output files remain the same.
Additional tests for preservation of various properties of the topology
could also be added later.

Test both read_tps_conf() and read_stx_conf() separately for now.

Added to support refactoring of the code.

Change-Id: Idbdfb4fc39b124d9236a54896ef4613ebfe80ceb

src/gromacs/fileio/confio.cpp
src/gromacs/fileio/filenm.h
src/gromacs/fileio/tests/CMakeLists.txt
src/gromacs/fileio/tests/confio.cpp [new file with mode: 0644]
src/testutils/stringtest.cpp
src/testutils/stringtest.h

index e7fbebf33306184f3082ea0796338ffb4197c688..6a022edff834039381419265552dcc40545af954 100644 (file)
@@ -599,7 +599,7 @@ const char *esp_prop[espNR] = {
     "molecule"
 };
 
-static void read_espresso_conf(const char *infile,
+static void read_espresso_conf(const char *infile, char *title,
                                t_atoms *atoms, rvec x[], rvec *v, matrix box)
 {
     t_symtab *symtab = NULL;
@@ -615,6 +615,8 @@ static void read_espresso_conf(const char *infile,
         snew(symtab, 1);
         open_symtab(symtab);
     }
+    // TODO: The code does not understand titles it writes...
+    title[0] = '\0';
 
     clear_mat(box);
 
@@ -1681,7 +1683,7 @@ void read_stx_conf(const char *infile, char *title, t_atoms *atoms,
             read_pdb_conf(infile, title, atoms, x, ePBC, box, TRUE, NULL);
             break;
         case efESP:
-            read_espresso_conf(infile, atoms, x, v, box);
+            read_espresso_conf(infile, title, atoms, x, v, box);
             break;
         case efTPR:
             snew(mtop, 1);
index 4d49c0e428fe9d73ae2d0d2cb1fddd0cf5b2b0c6..283f99f0436175c6e3e57fdfac1bf503d3395602 100644 (file)
@@ -45,8 +45,8 @@
 extern "C" {
 #endif
 
-/* this enum should correspond to the array deffile in gmxlib/filenm.c */
-enum {
+/* this enum should correspond to the array deffile in filenm.cpp */
+enum GromacsFileType {
     efMDP,
     efTRX, efTRO, efTRN, efTRR, efCOMPRESSED, efXTC, efTNG,
     efEDR,
index b0f62731b68418d10f8cf47e9d5f7b0bb073173a..8f83ed1cca72940202b559308d4b1eb2324ec14b 100644 (file)
@@ -1,7 +1,7 @@
 #
 # This file is part of the GROMACS molecular simulation package.
 #
-# Copyright (c) 2013,2014, by the GROMACS development team, led by
+# Copyright (c) 2013,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.
 # To help us fund GROMACS development, we humbly ask that you cite
 # the research papers on the package. Check out http://www.gromacs.org.
 
-if(GMX_USE_TNG)
-    gmx_add_unit_test(FileIOTests fileio-test
-        tngio.cpp)
+set(test_sources
+    confio.cpp
+    )
+if (GMX_USE_TNG)
+    list(APPEND test_sources tngio.cpp)
 endif()
+gmx_add_unit_test(FileIOTests fileio-test ${test_sources})
diff --git a/src/gromacs/fileio/tests/confio.cpp b/src/gromacs/fileio/tests/confio.cpp
new file mode 100644 (file)
index 0000000..00f971c
--- /dev/null
@@ -0,0 +1,238 @@
+/*
+ * 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
+ * Tests for reading/writing different structure file formats.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \ingroup module_fileio
+ */
+#include "gmxpre.h"
+
+#include "gromacs/fileio/confio.h"
+
+#include <iostream>
+#include <string>
+#include <vector>
+
+#include <gtest/gtest.h>
+
+#include "gromacs/fileio/filenm.h"
+#include "gromacs/fileio/tpxio.h"
+#include "gromacs/math/vec.h"
+#include "gromacs/math/vectypes.h"
+#include "gromacs/topology/atoms.h"
+#include "gromacs/topology/symtab.h"
+#include "gromacs/topology/topology.h"
+#include "gromacs/utility/cstringutil.h"
+#include "gromacs/utility/smalloc.h"
+
+#include "testutils/stringtest.h"
+#include "testutils/testfilemanager.h"
+
+// TODO: These should really appear somewhere centralized.
+/*! \brief
+ * Google Test formatter for GromacsFileType values.
+ */
+void PrintTo(const GromacsFileType &ftp, std::ostream *os)
+{
+    *os << "'" << ftp2ext(ftp) << "'";
+}
+
+namespace
+{
+
+class StructureIORoundtripTest : public gmx::test::StringTestBase,
+                                 public ::testing::WithParamInterface<GromacsFileType>
+{
+    public:
+        StructureIORoundtripTest()
+        {
+            generateReferenceTopology();
+            generateReferenceCoordinates();
+            testTop_ = NULL;
+            testX_   = NULL;
+            clear_mat(testBox_);
+            referenceFilename_ =
+                fileManager_.getTemporaryFilePath(getFileSuffix("ref"));
+            testFilename_ =
+                fileManager_.getTemporaryFilePath(getFileSuffix("test"));
+        }
+        ~StructureIORoundtripTest()
+        {
+            if (testTop_ != NULL)
+            {
+                free_t_atoms(&testTop_->atoms, TRUE);
+                done_top(testTop_);
+                sfree(testTop_);
+            }
+            sfree(testX_);
+            done_top(refTop_);
+            sfree(refTop_);
+        }
+
+        void writeReferenceFile()
+        {
+            write_sto_conf(referenceFilename_.c_str(), *refTop_->name,
+                           &refTop_->atoms, as_rvec_array(&refX_[0]), NULL, -1,
+                           refBox_);
+        }
+
+        void readReferenceFileStx()
+        {
+            int natoms = -1;
+            get_stx_coordnum(referenceFilename_.c_str(), &natoms);
+            ASSERT_EQ(refTop_->atoms.nr, natoms)
+            << "get_stx_coordnum() returned unexpected number of atoms";
+            char title[STRLEN];
+            snew(testTop_, 1);
+            init_t_atoms(&testTop_->atoms, natoms, GetParam() == efPDB);
+            snew(testX_, natoms);
+            read_stx_conf(referenceFilename_.c_str(), title, &testTop_->atoms,
+                          testX_, NULL, NULL, testBox_);
+            testTop_->name = put_symtab(&testTop_->symtab, title);
+        }
+
+        void readReferenceFileTps()
+        {
+            snew(testTop_, 1);
+            int  ePBC = -2;
+            char title[STRLEN];
+            read_tps_conf(referenceFilename_.c_str(), title, testTop_,
+                          &ePBC, &testX_, NULL, testBox_, FALSE);
+            testTop_->name = put_symtab(&testTop_->symtab, title);
+        }
+
+        void testTopologies()
+        {
+            // TODO: Compare the topologies.
+        }
+
+        void writeTestFileAndTest()
+        {
+            write_sto_conf(testFilename_.c_str(), *testTop_->name,
+                           &testTop_->atoms, testX_, NULL, -1, testBox_);
+            testFilesEqual(referenceFilename_, testFilename_);
+        }
+
+    private:
+        std::string getFileSuffix(const char *type)
+        {
+            return std::string(type) + "." + ftp2ext(GetParam());
+        }
+
+        void generateReferenceTopology()
+        {
+            snew(refTop_, 1);
+            open_symtab(&refTop_->symtab);
+            if (GetParam() == efESP)
+            {
+                // Titles cannot be read from an .esp file...
+                refTop_->name = put_symtab(&refTop_->symtab, "");
+            }
+            else
+            {
+                refTop_->name = put_symtab(&refTop_->symtab, "Test title");
+            }
+            const int atomCount = 10;
+            init_t_atoms(&refTop_->atoms, atomCount, FALSE);
+            for (int i = 0; i < atomCount; ++i)
+            {
+                char name[3];
+                name[0]                       = 'A';
+                name[1]                       = 'A' + i%3;
+                name[2]                       = '\0';
+                refTop_->atoms.atomname[i]    = put_symtab(&refTop_->symtab, name);
+                refTop_->atoms.atom[i].resind = i/3;
+                if (i%3 == 0)
+                {
+                    char resname[3];
+                    resname[0] = 'R';
+                    resname[1] = 'A' + i/3;
+                    resname[2] = '\0';
+                    t_atoms_set_resinfo(&refTop_->atoms, i, &refTop_->symtab,
+                                        resname, i/3 + 1, ' ', 0, ' ');
+                }
+            }
+            refTop_->atoms.nres = 4;
+            close_symtab(&refTop_->symtab);
+        }
+
+        void generateReferenceCoordinates()
+        {
+            clear_mat(refBox_);
+            refBox_[XX][XX] = 1;
+            refBox_[YY][YY] = 2;
+            refBox_[ZZ][ZZ] = 3;
+            const int atomCount = refTop_->atoms.nr;
+            refX_.reserve(atomCount);
+            for (int i = 0; i < atomCount; ++i)
+            {
+                refX_.push_back(gmx::RVec(i%4, i/4, (i/2)%3));
+            }
+        }
+
+        gmx::test::TestFileManager      fileManager_;
+        std::string                     referenceFilename_;
+        std::string                     testFilename_;
+        t_topology                     *refTop_;
+        std::vector<gmx::RVec>          refX_;
+        matrix                          refBox_;
+        t_topology                     *testTop_;
+        rvec                           *testX_;
+        matrix                          testBox_;
+};
+
+TEST_P(StructureIORoundtripTest, ReadWriteStxConf)
+{
+    writeReferenceFile();
+    readReferenceFileStx();
+    testTopologies();
+    writeTestFileAndTest();
+}
+
+TEST_P(StructureIORoundtripTest, ReadWriteTpsConf)
+{
+    writeReferenceFile();
+    readReferenceFileTps();
+    testTopologies();
+    writeTestFileAndTest();
+}
+
+INSTANTIATE_TEST_CASE_P(WithDifferentFormats,
+                        StructureIORoundtripTest,
+                            ::testing::Values(efGRO, efG96, efPDB, efESP));
+
+} // namespace
index ce37bdfee5f873631851795b37dd32500df8ea05..ad63a444f72a73bedd704209f4dcb0642478b510 100644 (file)
@@ -140,5 +140,18 @@ StringTestBase::checkFileContents(const std::string &filename, const char *id)
     checkText(text, id);
 }
 
+void
+StringTestBase::testFilesEqual(const std::string &refFilename,
+                               const std::string &testFilename)
+{
+    const std::string expectedContents = TextReader::readFileToString(refFilename);
+    const std::string contents         = TextReader::readFileToString(testFilename);
+    if (g_bWriteToStdOut)
+    {
+        printf("%s[END]\n", contents.c_str());
+    }
+    EXPECT_EQ(expectedContents, contents);
+}
+
 } // namespace test
 } // namespace gmx
index ce0f00af86084ab4b25caeef6283a45be62a1787..7f475b8e316c2b7b7de9b4644ed3c098863f7cc6 100644 (file)
@@ -110,6 +110,15 @@ class StringTestBase : public ::testing::Test
          */
         void checkFileContents(const std::string &filename, const char *id);
 
+        /*! \brief
+         * Tests that contents of two files are equal.
+         *
+         * \param[in] refFilename   File with the expected contents.
+         * \param[in] testFilename  File with the contents to be tested.
+         */
+        void testFilesEqual(const std::string &refFilename,
+                            const std::string &testFilename);
+
     private:
         class Impl;