Separate topology handling in selection tests.
authorTeemu Murtola <teemu.murtola@gmail.com>
Fri, 26 Apr 2013 12:16:18 +0000 (15:16 +0300)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Mon, 17 Jun 2013 15:15:37 +0000 (17:15 +0200)
Split common functionality for handling test topologies in selection
unit tests into a separate files (toputils.*).  This commit does not yet
add a lot of reusable code, but subsequent commits will (and they will
also introduce more files that use this functionality).

Part of #651.

Change-Id: Idda0099ba166014c9d6c36e0d8ba212e5fa8e940

src/gromacs/selection/tests/CMakeLists.txt
src/gromacs/selection/tests/selectioncollection.cpp
src/gromacs/selection/tests/selectionoption.cpp
src/gromacs/selection/tests/toputils.cpp [new file with mode: 0644]
src/gromacs/selection/tests/toputils.h [new file with mode: 0644]

index 3144b7157c7763f5e352655a7f51384f473eed64..231924ce197b5d1576dc93eba0aee078d280acd0 100644 (file)
@@ -1,7 +1,7 @@
 #
 # This file is part of the GROMACS molecular simulation package.
 #
-# Copyright (c) 2010,2011,2012, by the GROMACS development team, led by
+# Copyright (c) 2010,2011,2012,2013, by the GROMACS development team, led by
 # David van der Spoel, Berk Hess, 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,4 +34,5 @@
 
 gmx_add_unit_test(SelectionUnitTests selection-test
                   selectioncollection.cpp
-                  selectionoption.cpp)
+                  selectionoption.cpp
+                  toputils.cpp)
index 73b83a646e0bde484f3ed1f8ea2fa194be5f290a..32caba123cac5ae61d415f29540b5e48a23953c1 100644 (file)
  */
 #include <gtest/gtest.h>
 
-#include "gromacs/legacyheaders/smalloc.h"
-#include "gromacs/legacyheaders/statutil.h"
-#include "gromacs/legacyheaders/tpxio.h"
-#include "gromacs/legacyheaders/vec.h"
-
 #include "gromacs/options/basicoptions.h"
 #include "gromacs/options/options.h"
 #include "gromacs/selection/selectioncollection.h"
@@ -60,6 +55,8 @@
 #include "testutils/testfilemanager.h"
 #include "testutils/testoptions.h"
 
+#include "toputils.h"
+
 namespace
 {
 
@@ -75,7 +72,6 @@ class SelectionCollectionTest : public ::testing::Test
         static int               s_debugLevel;
 
         SelectionCollectionTest();
-        ~SelectionCollectionTest();
 
         void setAtomCount(int natoms)
         {
@@ -87,6 +83,9 @@ class SelectionCollectionTest : public ::testing::Test
         gmx::SelectionList       sel_;
         t_topology              *top_;
         t_trxframe              *frame_;
+
+    private:
+        gmx::test::TopologyManager  topManager_;
 };
 
 int SelectionCollectionTest::s_debugLevel = 0;
@@ -98,55 +97,21 @@ void SelectionCollectionTest::SetUpTestCase()
     gmx::test::parseTestOptions(&options);
 }
 
-
 SelectionCollectionTest::SelectionCollectionTest()
     : top_(NULL), frame_(NULL)
 {
+    topManager_.requestFrame();
     sc_.setDebugLevel(s_debugLevel);
     sc_.setReferencePosType("atom");
     sc_.setOutputPosType("atom");
 }
 
-
-SelectionCollectionTest::~SelectionCollectionTest()
-{
-    if (top_ != NULL)
-    {
-        free_t_atoms(&top_->atoms, TRUE);
-        done_top(top_);
-        sfree(top_);
-    }
-
-    if (frame_ != NULL)
-    {
-        sfree(frame_->x);
-        sfree(frame_);
-    }
-}
-
-
 void
 SelectionCollectionTest::loadTopology(const char *filename)
 {
-    char    title[STRLEN];
-    int     ePBC;
-    rvec   *xtop;
-    matrix  box;
-
-    snew(top_, 1);
-    read_tps_conf(gmx::test::TestFileManager::getInputFilePath(filename).c_str(),
-                  title, top_, &ePBC, &xtop, NULL, box, FALSE);
-
-    snew(frame_, 1);
-    frame_->flags  = TRX_NEED_X;
-    frame_->natoms = top_->atoms.nr;
-    frame_->bX     = TRUE;
-    snew(frame_->x, frame_->natoms);
-    memcpy(frame_->x, xtop, sizeof(*frame_->x) * frame_->natoms);
-    frame_->bBox   = TRUE;
-    copy_mat(box, frame_->box);
-
-    sfree(xtop);
+    topManager_.loadTopology(filename);
+    top_   = topManager_.topology();
+    frame_ = topManager_.frame();
 
     ASSERT_NO_THROW_GMX(sc_.setTopology(top_, -1));
 }
index 90f751170497a22d0392753a5f230025881359ec..00af92fac6e35d15785b4af3d364cf64fd1e753d 100644 (file)
@@ -41,9 +41,6 @@
  */
 #include <gtest/gtest.h>
 
-#include "gromacs/legacyheaders/smalloc.h"
-#include "gromacs/legacyheaders/tpxio.h"
-
 #include "gromacs/options/options.h"
 #include "gromacs/options/optionsassigner.h"
 #include "gromacs/selection/selection.h"
@@ -56,6 +53,8 @@
 #include "testutils/testasserts.h"
 #include "testutils/testfilemanager.h"
 
+#include "toputils.h"
+
 using gmx::test::TestFileManager;
 
 namespace
@@ -69,7 +68,6 @@ class SelectionOptionTestBase : public ::testing::Test
 {
     public:
         SelectionOptionTestBase();
-        ~SelectionOptionTestBase();
 
         void setManager();
         void loadTopology(const char *filename);
@@ -77,43 +75,28 @@ class SelectionOptionTestBase : public ::testing::Test
         gmx::SelectionCollection    sc_;
         gmx::SelectionOptionManager manager_;
         gmx::Options                options_;
-        t_topology                 *top_;
+
+    private:
+        gmx::test::TopologyManager  topManager_;
 };
 
 SelectionOptionTestBase::SelectionOptionTestBase()
-    : manager_(&sc_), options_(NULL, NULL), top_(NULL)
+    : manager_(&sc_), options_(NULL, NULL)
 {
     sc_.setReferencePosType("atom");
     sc_.setOutputPosType("atom");
 }
 
-SelectionOptionTestBase::~SelectionOptionTestBase()
-{
-    if (top_ != NULL)
-    {
-        free_t_atoms(&top_->atoms, TRUE);
-        done_top(top_);
-        sfree(top_);
-    }
-}
-
 void SelectionOptionTestBase::setManager()
 {
     setManagerForSelectionOptions(&options_, &manager_);
 }
 
-
 void SelectionOptionTestBase::loadTopology(const char *filename)
 {
-    char    title[STRLEN];
-    int     ePBC;
-    matrix  box;
-
-    snew(top_, 1);
-    read_tps_conf(gmx::test::TestFileManager::getInputFilePath(filename).c_str(),
-                  title, top_, &ePBC, NULL, NULL, box, FALSE);
+    topManager_.loadTopology(filename);
 
-    ASSERT_NO_THROW_GMX(sc_.setTopology(top_, -1));
+    ASSERT_NO_THROW_GMX(sc_.setTopology(topManager_.topology(), -1));
 }
 
 
diff --git a/src/gromacs/selection/tests/toputils.cpp b/src/gromacs/selection/tests/toputils.cpp
new file mode 100644 (file)
index 0000000..0d5f285
--- /dev/null
@@ -0,0 +1,119 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2013, by the GROMACS development team, led by
+ * David van der Spoel, Berk Hess, 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 test helper routines from toputils.h.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \ingroup module_selection
+ */
+#include "toputils.h"
+
+#include <cstring>
+
+#include "gromacs/legacyheaders/smalloc.h"
+#include "gromacs/legacyheaders/statutil.h"
+#include "gromacs/legacyheaders/tpxio.h"
+#include "gromacs/legacyheaders/vec.h"
+
+#include "gromacs/utility/gmxassert.h"
+
+#include "testutils/testfilemanager.h"
+
+namespace gmx
+{
+namespace test
+{
+
+TopologyManager::TopologyManager()
+    : top_(NULL), frame_(NULL)
+{
+}
+
+TopologyManager::~TopologyManager()
+{
+    if (top_ != NULL)
+    {
+        free_t_atoms(&top_->atoms, TRUE);
+        done_top(top_);
+        sfree(top_);
+    }
+
+    if (frame_ != NULL)
+    {
+        sfree(frame_->x);
+        sfree(frame_);
+    }
+}
+
+void TopologyManager::requestFrame()
+{
+    GMX_RELEASE_ASSERT(top_ == NULL,
+                       "Frame must be requested before initializing topology");
+    if (frame_ == NULL)
+    {
+        snew(frame_, 1);
+    }
+}
+
+void TopologyManager::loadTopology(const char *filename)
+{
+    char    title[STRLEN];
+    int     ePBC;
+    rvec   *xtop = NULL;
+    matrix  box;
+
+    GMX_RELEASE_ASSERT(top_ == NULL, "Topology initialized more than once");
+    snew(top_, 1);
+    read_tps_conf(gmx::test::TestFileManager::getInputFilePath(filename).c_str(),
+                  title, top_, &ePBC, frame_ != NULL ? &xtop : NULL,
+                  NULL, box, FALSE);
+
+    if (frame_ != NULL)
+    {
+        frame_->flags  = TRX_NEED_X;
+        frame_->natoms = top_->atoms.nr;
+        frame_->bX     = TRUE;
+        snew(frame_->x, frame_->natoms);
+        std::memcpy(frame_->x, xtop, sizeof(*frame_->x) * frame_->natoms);
+        frame_->bBox   = TRUE;
+        copy_mat(box, frame_->box);
+    }
+
+    sfree(xtop);
+}
+
+} // namespace test
+} // namespace gmx
diff --git a/src/gromacs/selection/tests/toputils.h b/src/gromacs/selection/tests/toputils.h
new file mode 100644 (file)
index 0000000..8c74630
--- /dev/null
@@ -0,0 +1,73 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2013, by the GROMACS development team, led by
+ * David van der Spoel, Berk Hess, 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
+ * Helper routines for constructing topologies for tests.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \ingroup module_selection
+ */
+#ifndef GMX_SELECTION_TESTS_TOPUTILS_H
+#define GMX_SELECTION_TESTS_TOPUTILS_H
+
+#include "gromacs/legacyheaders/typedefs.h"
+
+namespace gmx
+{
+namespace test
+{
+
+class TopologyManager
+{
+    public:
+        TopologyManager();
+        ~TopologyManager();
+
+        void requestFrame();
+
+        void loadTopology(const char *filename);
+
+        t_topology *topology() { return top_; }
+        t_trxframe *frame() { return frame_; }
+
+    private:
+        t_topology *top_;
+        t_trxframe *frame_;
+};
+
+} // namespace test
+} // namespace gmx
+
+#endif