Add unit test coverage for simple selection keywords.
authorTeemu Murtola <teemu.murtola@gmail.com>
Tue, 30 Apr 2013 03:51:34 +0000 (06:51 +0300)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Sun, 4 Aug 2013 19:45:34 +0000 (21:45 +0200)
"molindex" and "atomtype" keywords are now covered by separate tests.
Also updated one test to use the "z" keyword to have that covered.
Added a comment in sm_simple.cpp about unreachable code that is not
currently covered by tests.

Part of #651.

Change-Id: Id9ac9ea8e7d144214537265d9ed89335684934f1

src/gromacs/gmxlib/typedefs.c
src/gromacs/selection/sm_simple.cpp
src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesAtomtype.xml [new file with mode: 0644]
src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesDynamicAtomValuedParameters.xml
src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesMolIndex.xml [new file with mode: 0644]
src/gromacs/selection/tests/selectioncollection.cpp
src/gromacs/selection/tests/toputils.cpp
src/gromacs/selection/tests/toputils.h

index 5c3a1bef1c9c7025626b57dd27de54c0d5b4be0f..e18b161522ec6c41bf6b2e62f5ac2857d23d9eec 100644 (file)
@@ -870,17 +870,42 @@ void free_t_atoms(t_atoms *atoms, gmx_bool bFreeNames)
             }
         }
     }
+    if (bFreeNames && atoms->atomtype != NULL)
+    {
+        for (i = 0; i < atoms->nr; i++)
+        {
+            if (atoms->atomtype[i] != NULL)
+            {
+                sfree(*atoms->atomtype[i]);
+                *atoms->atomtype[i] = NULL;
+            }
+        }
+    }
+    if (bFreeNames && atoms->atomtypeB != NULL)
+    {
+        for (i = 0; i < atoms->nr; i++)
+        {
+            if (atoms->atomtypeB[i] != NULL)
+            {
+                sfree(*atoms->atomtypeB[i]);
+                *atoms->atomtypeB[i] = NULL;
+            }
+        }
+    }
     sfree(atoms->atomname);
-    /* Do we need to free atomtype and atomtypeB as well ? */
+    sfree(atoms->atomtype);
+    sfree(atoms->atomtypeB);
     sfree(atoms->resinfo);
     sfree(atoms->atom);
     sfree(atoms->pdbinfo);
-    atoms->nr       = 0;
-    atoms->nres     = 0;
-    atoms->atomname = NULL;
-    atoms->resinfo  = NULL;
-    atoms->atom     = NULL;
-    atoms->pdbinfo  = NULL;
+    atoms->nr        = 0;
+    atoms->nres      = 0;
+    atoms->atomname  = NULL;
+    atoms->atomtype  = NULL;
+    atoms->atomtypeB = NULL;
+    atoms->resinfo   = NULL;
+    atoms->atom      = NULL;
+    atoms->pdbinfo   = NULL;
 }
 
 real max_cutoff(real cutoff1, real cutoff2)
index b364904e48393693afb3b8534d6e3de08f8baab7..97c20e2d9b4c6ccf3c60c8d8d277532986603140 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2009,2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2009,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.
@@ -867,6 +867,9 @@ evaluate_coord(t_trxframe *fr, gmx_ana_index_t *g, real out[],
     }
     else
     {
+        // TODO: This loop is never reached in the current code.
+        // It would be useful to change the code such that it is, mostly for
+        // memory efficiency reasons.
         for (i = 0; i < g->isize; ++i)
         {
             out[i] = fr->x[g->index[i]][d];
diff --git a/src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesAtomtype.xml b/src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesAtomtype.xml
new file mode 100644 (file)
index 0000000..77c46a3
--- /dev/null
@@ -0,0 +1,24 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <ParsedSelections Name="Parsed">
+    <ParsedSelection Name="Selection1">
+      <String Name="Input">atomtype CA</String>
+      <String Name="Name">atomtype CA</String>
+      <String Name="Text">atomtype CA</String>
+      <Bool Name="Dynamic">false</Bool>
+    </ParsedSelection>
+  </ParsedSelections>
+  <CompiledSelections Name="Compiled">
+    <Selection Name="Selection1">
+      <Sequence Name="Atoms">
+        <Int Name="Length">5</Int>
+        <Int>0</Int>
+        <Int>3</Int>
+        <Int>6</Int>
+        <Int>9</Int>
+        <Int>12</Int>
+      </Sequence>
+    </Selection>
+  </CompiledSelections>
+</ReferenceData>
index bc1d89ded767b8dec59c6b6cbdcb084fc4e6658f..dfce299dd00322e765e5e61d5d56550be0efe29b 100644 (file)
@@ -9,9 +9,9 @@
       <Bool Name="Dynamic">true</Bool>
     </ParsedSelection>
     <ParsedSelection Name="Selection2">
-      <String Name="Input">(resnr 1 3 5 or x &gt; 10) and same residue as (atomnr 3 5 13 or y &gt; 5)</String>
-      <String Name="Name">(resnr 1 3 5 or x &gt; 10) and same residue as (atomnr 3 5 13 or y &gt; 5)</String>
-      <String Name="Text">(resnr 1 3 5 or x &gt; 10) and same residue as (atomnr 3 5 13 or y &gt; 5)</String>
+      <String Name="Input">(resnr 1 3 5 or x &gt; 10) and same residue as (atomnr 3 5 13 or z &gt; 5)</String>
+      <String Name="Name">(resnr 1 3 5 or x &gt; 10) and same residue as (atomnr 3 5 13 or z &gt; 5)</String>
+      <String Name="Text">(resnr 1 3 5 or x &gt; 10) and same residue as (atomnr 3 5 13 or z &gt; 5)</String>
       <Bool Name="Dynamic">true</Bool>
     </ParsedSelection>
   </ParsedSelections>
diff --git a/src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesMolIndex.xml b/src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesMolIndex.xml
new file mode 100644 (file)
index 0000000..bd56516
--- /dev/null
@@ -0,0 +1,45 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <ParsedSelections Name="Parsed">
+    <ParsedSelection Name="Selection1">
+      <String Name="Input">molindex 1 4</String>
+      <String Name="Name">molindex 1 4</String>
+      <String Name="Text">molindex 1 4</String>
+      <Bool Name="Dynamic">false</Bool>
+    </ParsedSelection>
+    <ParsedSelection Name="Selection2">
+      <String Name="Input">molecule 2 3 5</String>
+      <String Name="Name">molecule 2 3 5</String>
+      <String Name="Text">molecule 2 3 5</String>
+      <Bool Name="Dynamic">false</Bool>
+    </ParsedSelection>
+  </ParsedSelections>
+  <CompiledSelections Name="Compiled">
+    <Selection Name="Selection1">
+      <Sequence Name="Atoms">
+        <Int Name="Length">6</Int>
+        <Int>0</Int>
+        <Int>1</Int>
+        <Int>2</Int>
+        <Int>9</Int>
+        <Int>10</Int>
+        <Int>11</Int>
+      </Sequence>
+    </Selection>
+    <Selection Name="Selection2">
+      <Sequence Name="Atoms">
+        <Int Name="Length">9</Int>
+        <Int>3</Int>
+        <Int>4</Int>
+        <Int>5</Int>
+        <Int>6</Int>
+        <Int>7</Int>
+        <Int>8</Int>
+        <Int>12</Int>
+        <Int>13</Int>
+        <Int>14</Int>
+      </Sequence>
+    </Selection>
+  </CompiledSelections>
+</ReferenceData>
index 62b59ed6882e3628b5a1689ece4814ad5f5f3689..212cffc6303f02bf880a220b5ac36ed9f3168fb2 100644 (file)
@@ -76,14 +76,13 @@ class SelectionCollectionTest : public ::testing::Test
             ASSERT_NO_THROW_GMX(sc_.setTopology(NULL, natoms));
         }
         void loadTopology(const char *filename);
+        void setTopology();
 
-        gmx::SelectionCollection sc_;
-        gmx::SelectionList       sel_;
-        t_topology              *top_;
-        t_trxframe              *frame_;
-
-    private:
         gmx::test::TopologyManager  topManager_;
+        gmx::SelectionCollection    sc_;
+        gmx::SelectionList          sel_;
+        t_topology                 *top_;
+        t_trxframe                 *frame_;
 };
 
 int SelectionCollectionTest::s_debugLevel = 0;
@@ -108,6 +107,12 @@ void
 SelectionCollectionTest::loadTopology(const char *filename)
 {
     topManager_.loadTopology(filename);
+    setTopology();
+}
+
+void
+SelectionCollectionTest::setTopology()
+{
     top_   = topManager_.topology();
     frame_ = topManager_.frame();
 
@@ -535,7 +540,17 @@ TEST_F(SelectionCollectionDataTest, HandlesResIndex)
     runTest("simple.pdb", selections);
 }
 
-// TODO: Add test for "molindex"
+TEST_F(SelectionCollectionDataTest, HandlesMolIndex)
+{
+    static const char * const selections[] = {
+        "molindex 1 4",
+        "molecule 2 3 5"
+    };
+    ASSERT_NO_FATAL_FAILURE(runParser(selections));
+    ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
+    topManager_.initUniformMolecules(3);
+    ASSERT_NO_FATAL_FAILURE(runCompiler());
+}
 
 TEST_F(SelectionCollectionDataTest, HandlesAtomname)
 {
@@ -557,7 +572,18 @@ TEST_F(SelectionCollectionDataTest, HandlesPdbAtomname)
     runTest("simple.pdb", selections);
 }
 
-// TODO: Add test for atomtype
+
+TEST_F(SelectionCollectionDataTest, HandlesAtomtype)
+{
+    static const char * const selections[] = {
+        "atomtype CA"
+    };
+    ASSERT_NO_FATAL_FAILURE(runParser(selections));
+    ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
+    const char *const types[] = { "CA", "SA", "SB" };
+    topManager_.initAtomTypes(types);
+    ASSERT_NO_FATAL_FAILURE(runCompiler());
+}
 
 TEST_F(SelectionCollectionDataTest, HandlesChain)
 {
@@ -879,7 +905,7 @@ TEST_F(SelectionCollectionDataTest, HandlesDynamicAtomValuedParameters)
 {
     static const char * const selections[] = {
         "same residue as (atomnr 3 5 13 or y > 5)",
-        "(resnr 1 3 5 or x > 10) and same residue as (atomnr 3 5 13 or y > 5)"
+        "(resnr 1 3 5 or x > 10) and same residue as (atomnr 3 5 13 or z > 5)"
     };
     setFlags(TestFlags() | efTestEvaluation);
     runTest("simple.gro", selections);
index c75e95ef2145415c7d90c791267b5f71645f67b1..36e8c97e9dc5c0f6952d9ebe9354491d95a2942b 100644 (file)
@@ -45,6 +45,7 @@
 
 #include "gromacs/legacyheaders/smalloc.h"
 #include "gromacs/legacyheaders/statutil.h"
+#include "gromacs/legacyheaders/string2.h"
 #include "gromacs/legacyheaders/tpxio.h"
 #include "gromacs/legacyheaders/typedefs.h"
 #include "gromacs/legacyheaders/vec.h"
@@ -166,6 +167,25 @@ void TopologyManager::initAtoms(int count)
     }
 }
 
+void TopologyManager::initAtomTypes(int count, const char *const types[])
+{
+    GMX_RELEASE_ASSERT(top_ != NULL, "Topology not initialized");
+    atomtypes_.reserve(count);
+    for (int i = 0; i < count; ++i)
+    {
+        atomtypes_.push_back(strdup(types[i]));
+    }
+    snew(top_->atoms.atomtype, top_->atoms.nr);
+    for (int i = 0, j = 0; i < top_->atoms.nr; ++i, ++j)
+    {
+        if (j == count)
+        {
+            j = 0;
+        }
+        top_->atoms.atomtype[i] = &atomtypes_[j];
+    }
+}
+
 void TopologyManager::initUniformResidues(int residueSize)
 {
     GMX_RELEASE_ASSERT(top_ != NULL, "Topology not initialized");
index 893cbc206fe6fe05505012a72c53e8b10baac256..a778bbd72c213eb8f46302c5ece7c1074e3f40b8 100644 (file)
@@ -42,6 +42,8 @@
 #ifndef GMX_SELECTION_TESTS_TOPUTILS_H
 #define GMX_SELECTION_TESTS_TOPUTILS_H
 
+#include <vector>
+
 #include "gromacs/legacyheaders/typedefs.h"
 
 namespace gmx
@@ -61,6 +63,12 @@ class TopologyManager
 
         void loadTopology(const char *filename);
         void initAtoms(int count);
+        void initAtomTypes(int count, const char *const types[]);
+        template <int count>
+        void initAtomTypes(const char *const (&types)[count])
+        {
+            initAtomTypes(count, types);
+        }
         void initUniformResidues(int residueSize);
         void initUniformMolecules(int moleculeSize);
 
@@ -68,8 +76,9 @@ class TopologyManager
         t_trxframe *frame() { return frame_; }
 
     private:
-        t_topology *top_;
-        t_trxframe *frame_;
+        t_topology             *top_;
+        t_trxframe             *frame_;
+        std::vector<char *>     atomtypes_;
 };
 
 } // namespace test