Add tests for selection mass evaluation and fix it.
authorTeemu Murtola <teemu.murtola@gmail.com>
Fri, 8 Feb 2013 13:26:47 +0000 (15:26 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Wed, 13 Feb 2013 21:27:13 +0000 (22:27 +0100)
Removed unnecessary and at times incorrect precalculation of total
masses and charges for selected positions. Evaluate them again for each
frame if they may have changed.

Noticed when adding some tests for this behavior; the tests are also
added here.

Also store the masses and charges again in continuous arrays in
preparation for exposing those arrays.

Fixes #1145.

Change-Id: Ifa23d042706d9f2a0b8aa142338812f055f6f60a

src/gromacs/selection/evaluate.cpp
src/gromacs/selection/selection.cpp
src/gromacs/selection/selection.h
src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_ComputesMassesAndCharges.xml [new file with mode: 0644]
src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_ComputesMassesAndChargesWithoutTopology.xml [new file with mode: 0644]
src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesCharge.xml [new file with mode: 0644]
src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesMass.xml [new file with mode: 0644]
src/gromacs/selection/tests/refdata/referencedata.xsl
src/gromacs/selection/tests/selectioncollection.cpp

index 135790ba6d6a13f37de517621e33809abb47f75b..427d0d6f56d1fa949dc0dc78c85a45873d42c988 100644 (file)
@@ -426,7 +426,7 @@ SelectionEvaluator::evaluate(SelectionCollection *coll,
     for (isel = sc->sel.begin(); isel != sc->sel.end(); ++isel)
     {
         internal::SelectionData &sel = **isel;
-        sel.refreshMassesAndCharges();
+        sel.refreshMassesAndCharges(sc->top);
         sel.updateCoveredFractionForFrame();
     }
 }
@@ -444,7 +444,7 @@ SelectionEvaluator::evaluateFinal(SelectionCollection *coll, int nframes)
     for (isel = sc->sel.begin(); isel != sc->sel.end(); ++isel)
     {
         internal::SelectionData &sel = **isel;
-        sel.restoreOriginalPositions();
+        sel.restoreOriginalPositions(sc->top);
         sel.computeAverageCoveredFraction(nframes);
     }
 }
index 24c44db2f8913570ece8e7d8530f4265752876c9..ec75398083d5397067ca35597a9c3eeb2c432bfa 100644 (file)
@@ -124,47 +124,68 @@ SelectionData::initCoveredFraction(e_coverfrac_t type)
     return type == CFRAC_NONE || coveredFractionType_ != CFRAC_NONE;
 }
 
+namespace
+{
 
-void
-SelectionData::initializeMassesAndCharges(const t_topology *top)
+/*! \brief
+ * Helper function to compute total masses and charges for positions.
+ *
+ * \param[in]  top     Topology to take atom masses from.
+ * \param[in]  pos     Positions to compute masses and charges for.
+ * \param[out] masses  Output masses.
+ * \param[out] charges Output charges.
+ *
+ * Does not throw if enough space has been reserved for the output vectors.
+ */
+void computeMassesAndCharges(const t_topology *top, const gmx_ana_pos_t &pos,
+                             std::vector<real> *masses,
+                             std::vector<real> *charges)
 {
-    posInfo_.reserve(posCount());
-    for (int b = 0; b < posCount(); ++b)
+    GMX_ASSERT(top != NULL, "Should not have been called with NULL topology");
+    masses->clear();
+    charges->clear();
+    for (int b = 0; b < pos.nr; ++b)
     {
-        real mass   = 1.0;
+        real mass   = 0.0;
         real charge = 0.0;
-        if (top != NULL)
+        for (int i = pos.m.mapb.index[b]; i < pos.m.mapb.index[b+1]; ++i)
         {
-            mass = 0.0;
-            for (int i = rawPositions_.m.mapb.index[b];
-                 i < rawPositions_.m.mapb.index[b+1];
-                 ++i)
-            {
-                int index = rawPositions_.g->index[i];
-                mass   += top->atoms.atom[index].m;
-                charge += top->atoms.atom[index].q;
-            }
+            int index = pos.g->index[i];
+            mass   += top->atoms.atom[index].m;
+            charge += top->atoms.atom[index].q;
         }
-        posInfo_.push_back(PositionInfo(mass, charge));
+        masses->push_back(mass);
+        charges->push_back(charge);
     }
-    if (isDynamic() && !hasFlag(efSelection_DynamicMask))
+}
+
+} // namespace
+
+void
+SelectionData::initializeMassesAndCharges(const t_topology *top)
+{
+    GMX_ASSERT(posMass_.empty() && posCharge_.empty(),
+               "Should not be called more than once");
+    posMass_.reserve(posCount());
+    posCharge_.reserve(posCount());
+    if (top == NULL)
     {
-        originalPosInfo_ = posInfo_;
+        posMass_.resize(posCount(), 1.0);
+        posCharge_.resize(posCount(), 0.0);
+    }
+    else
+    {
+        computeMassesAndCharges(top, rawPositions_, &posMass_, &posCharge_);
     }
 }
 
 
 void
-SelectionData::refreshMassesAndCharges()
+SelectionData::refreshMassesAndCharges(const t_topology *top)
 {
-    if (!originalPosInfo_.empty())
+    if (top != NULL && isDynamic() && !hasFlag(efSelection_DynamicMask))
     {
-        posInfo_.clear();
-        for (int i = 0; i < posCount(); ++i)
-        {
-            int refid  = rawPositions_.m.refid[i];
-            posInfo_.push_back(originalPosInfo_[refid]);
-        }
+        computeMassesAndCharges(top, rawPositions_, &posMass_, &posCharge_);
     }
 }
 
@@ -192,7 +213,7 @@ SelectionData::computeAverageCoveredFraction(int nframes)
 
 
 void
-SelectionData::restoreOriginalPositions()
+SelectionData::restoreOriginalPositions(const t_topology *top)
 {
     if (isDynamic())
     {
@@ -201,7 +222,7 @@ SelectionData::restoreOriginalPositions()
         p.g->name = NULL;
         gmx_ana_indexmap_update(&p.m, p.g, hasFlag(gmx::efSelection_DynamicMask));
         p.nr = p.m.nr;
-        refreshMassesAndCharges();
+        refreshMassesAndCharges(top);
     }
 }
 
index 9c2d5ee894d34defec828504b2402df621357f08..f1a6801729df37ba3f0303cd078f6f1504fc01bc 100644 (file)
@@ -117,9 +117,9 @@ class SelectionData
          * \param[in] top   Topology information.
          * \throws    std::bad_alloc if out of memory.
          *
-         * Computed values are cached, and need to be updated for dynamic
-         * selections with refreshMassesAndCharges() after the selection has
-         * been evaluated.  This is done by SelectionEvaluator.
+         * For dynamic selections, the values need to be updated after each
+         * evaluation with refreshMassesAndCharges().
+         * This is done by SelectionEvaluator.
          *
          * This function is called by SelectionCompiler.
          *
@@ -130,9 +130,11 @@ class SelectionData
          * Updates masses and charges after dynamic selection has been
          * evaluated.
          *
+         * \param[in] top   Topology information.
+         *
          * Called by SelectionEvaluator.
          */
-        void refreshMassesAndCharges();
+        void refreshMassesAndCharges(const t_topology *top);
         /*! \brief
          * Updates the covered fraction after a selection has been evaluated.
          *
@@ -152,44 +154,25 @@ class SelectionData
         /*! \brief
          * Restores position information to state it was in after compilation.
          *
+         * \param[in] top   Topology information.
+         *
          * Depends on SelectionCompiler storing the original atoms in the
          * \a rootElement_ object.
          * Called by SelectionEvaluator::evaluateFinal().
          */
-        void restoreOriginalPositions();
+        void restoreOriginalPositions(const t_topology *top);
 
     private:
-        /*! \brief
-         * Additional information about positions.
-         *
-         * This structure contains information about positions in the
-         * selection that is not stored in ::gmx_ana_pos_t.
-         *
-         * Methods in this class do not throw.
-         */
-        struct PositionInfo
-        {
-            //! Construct position information with unit mass and no charge.
-            PositionInfo() : mass(1.0), charge(0.0) {}
-            //! Construct position information with the given information.
-            PositionInfo(real mass, real charge) : mass(mass), charge(charge) {}
-
-            //! Total mass of atoms that make up the position.
-            real                mass;
-            //! Total charge of atoms that make up the position.
-            real                charge;
-        };
-
         //! Name of the selection.
         std::string               name_;
         //! The actual selection string.
         std::string               selectionText_;
         //! Low-level representation of selected positions.
         gmx_ana_pos_t             rawPositions_;
-        //! Information associated with the current positions.
-        std::vector<PositionInfo> posInfo_;
-        //! Information for all possible positions.
-        std::vector<PositionInfo> originalPosInfo_;
+        //! Total masses for the current positions.
+        std::vector<real>         posMass_;
+        //! Total charges for the current positions.
+        std::vector<real>         posCharge_;
         SelectionFlags            flags_;
         //! Root of the selection evaluation tree.
         SelectionTreeElement     &rootElement_;
@@ -487,7 +470,7 @@ class SelectionPosition
          */
         real mass() const
         {
-            return sel_->posInfo_[i_].mass;
+            return sel_->posMass_[i_];
         }
         /*! \brief
          * Returns total charge for this position.
@@ -498,7 +481,7 @@ class SelectionPosition
          */
         real charge() const
         {
-            return sel_->posInfo_[i_].charge;
+            return sel_->posCharge_[i_];
         }
         //! Returns the number of atoms that make up this position.
         int atomCount() const
diff --git a/src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_ComputesMassesAndCharges.xml b/src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_ComputesMassesAndCharges.xml
new file mode 100644 (file)
index 0000000..63feebf
--- /dev/null
@@ -0,0 +1,478 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <ParsedSelections Name="Parsed">
+    <ParsedSelection Name="Selection1">
+      <String Name="Input">name CB</String>
+      <String Name="Name">name CB</String>
+      <String Name="Text">name CB</String>
+      <Bool Name="Dynamic">false</Bool>
+    </ParsedSelection>
+    <ParsedSelection Name="Selection2">
+      <String Name="Input">y &gt; 2</String>
+      <String Name="Name">y &gt; 2</String>
+      <String Name="Text">y &gt; 2</String>
+      <Bool Name="Dynamic">true</Bool>
+    </ParsedSelection>
+    <ParsedSelection Name="Selection3">
+      <String Name="Input">res_cog of y &gt; 2</String>
+      <String Name="Name">res_cog of y &gt; 2</String>
+      <String Name="Text">res_cog of y &gt; 2</String>
+      <Bool Name="Dynamic">true</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>
+      <Sequence Name="Positions">
+        <Int Name="Length">5</Int>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>0</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">-1.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>3</Int>
+          </Sequence>
+          <Real Name="Mass">1.030000</Real>
+          <Real Name="Charge">-1.030000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>6</Int>
+          </Sequence>
+          <Real Name="Mass">1.060000</Real>
+          <Real Name="Charge">-1.060000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>9</Int>
+          </Sequence>
+          <Real Name="Mass">1.090000</Real>
+          <Real Name="Charge">-1.090000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>12</Int>
+          </Sequence>
+          <Real Name="Mass">1.120000</Real>
+          <Real Name="Charge">-1.120000</Real>
+        </Position>
+      </Sequence>
+    </Selection>
+    <Selection Name="Selection2">
+      <Sequence Name="Atoms">
+        <Int Name="Length">15</Int>
+        <Int>0</Int>
+        <Int>1</Int>
+        <Int>2</Int>
+        <Int>3</Int>
+        <Int>4</Int>
+        <Int>5</Int>
+        <Int>6</Int>
+        <Int>7</Int>
+        <Int>8</Int>
+        <Int>9</Int>
+        <Int>10</Int>
+        <Int>11</Int>
+        <Int>12</Int>
+        <Int>13</Int>
+        <Int>14</Int>
+      </Sequence>
+      <Sequence Name="Positions">
+        <Int Name="Length">15</Int>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>0</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">-1.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>1</Int>
+          </Sequence>
+          <Real Name="Mass">1.010000</Real>
+          <Real Name="Charge">-1.010000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>2</Int>
+          </Sequence>
+          <Real Name="Mass">1.020000</Real>
+          <Real Name="Charge">-1.020000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>3</Int>
+          </Sequence>
+          <Real Name="Mass">1.030000</Real>
+          <Real Name="Charge">-1.030000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>4</Int>
+          </Sequence>
+          <Real Name="Mass">1.040000</Real>
+          <Real Name="Charge">-1.040000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>5</Int>
+          </Sequence>
+          <Real Name="Mass">1.050000</Real>
+          <Real Name="Charge">-1.050000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>6</Int>
+          </Sequence>
+          <Real Name="Mass">1.060000</Real>
+          <Real Name="Charge">-1.060000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>7</Int>
+          </Sequence>
+          <Real Name="Mass">1.070000</Real>
+          <Real Name="Charge">-1.070000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>8</Int>
+          </Sequence>
+          <Real Name="Mass">1.080000</Real>
+          <Real Name="Charge">-1.080000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>9</Int>
+          </Sequence>
+          <Real Name="Mass">1.090000</Real>
+          <Real Name="Charge">-1.090000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>10</Int>
+          </Sequence>
+          <Real Name="Mass">1.100000</Real>
+          <Real Name="Charge">-1.100000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>11</Int>
+          </Sequence>
+          <Real Name="Mass">1.110000</Real>
+          <Real Name="Charge">-1.110000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>12</Int>
+          </Sequence>
+          <Real Name="Mass">1.120000</Real>
+          <Real Name="Charge">-1.120000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>13</Int>
+          </Sequence>
+          <Real Name="Mass">1.130000</Real>
+          <Real Name="Charge">-1.130000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>14</Int>
+          </Sequence>
+          <Real Name="Mass">1.140000</Real>
+          <Real Name="Charge">-1.140000</Real>
+        </Position>
+      </Sequence>
+    </Selection>
+    <Selection Name="Selection3">
+      <Sequence Name="Atoms">
+        <Int Name="Length">15</Int>
+        <Int>0</Int>
+        <Int>1</Int>
+        <Int>2</Int>
+        <Int>3</Int>
+        <Int>4</Int>
+        <Int>5</Int>
+        <Int>6</Int>
+        <Int>7</Int>
+        <Int>8</Int>
+        <Int>9</Int>
+        <Int>10</Int>
+        <Int>11</Int>
+        <Int>12</Int>
+        <Int>13</Int>
+        <Int>14</Int>
+      </Sequence>
+      <Sequence Name="Positions">
+        <Int Name="Length">5</Int>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">3</Int>
+            <Int>0</Int>
+            <Int>1</Int>
+            <Int>2</Int>
+          </Sequence>
+          <Real Name="Mass">3.030000</Real>
+          <Real Name="Charge">-3.030000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">3</Int>
+            <Int>3</Int>
+            <Int>4</Int>
+            <Int>5</Int>
+          </Sequence>
+          <Real Name="Mass">3.120000</Real>
+          <Real Name="Charge">-3.120000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">3</Int>
+            <Int>6</Int>
+            <Int>7</Int>
+            <Int>8</Int>
+          </Sequence>
+          <Real Name="Mass">3.210000</Real>
+          <Real Name="Charge">-3.210000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">3</Int>
+            <Int>9</Int>
+            <Int>10</Int>
+            <Int>11</Int>
+          </Sequence>
+          <Real Name="Mass">3.300000</Real>
+          <Real Name="Charge">-3.300000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">3</Int>
+            <Int>12</Int>
+            <Int>13</Int>
+            <Int>14</Int>
+          </Sequence>
+          <Real Name="Mass">3.390000</Real>
+          <Real Name="Charge">-3.390000</Real>
+        </Position>
+      </Sequence>
+    </Selection>
+  </CompiledSelections>
+  <EvaluatedSelections Name="Frame1">
+    <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>
+      <Sequence Name="Positions">
+        <Int Name="Length">5</Int>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>0</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">-1.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>3</Int>
+          </Sequence>
+          <Real Name="Mass">1.030000</Real>
+          <Real Name="Charge">-1.030000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>6</Int>
+          </Sequence>
+          <Real Name="Mass">1.060000</Real>
+          <Real Name="Charge">-1.060000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>9</Int>
+          </Sequence>
+          <Real Name="Mass">1.090000</Real>
+          <Real Name="Charge">-1.090000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>12</Int>
+          </Sequence>
+          <Real Name="Mass">1.120000</Real>
+          <Real Name="Charge">-1.120000</Real>
+        </Position>
+      </Sequence>
+    </Selection>
+    <Selection Name="Selection2">
+      <Sequence Name="Atoms">
+        <Int Name="Length">7</Int>
+        <Int>2</Int>
+        <Int>3</Int>
+        <Int>6</Int>
+        <Int>7</Int>
+        <Int>10</Int>
+        <Int>11</Int>
+        <Int>14</Int>
+      </Sequence>
+      <Sequence Name="Positions">
+        <Int Name="Length">7</Int>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>2</Int>
+          </Sequence>
+          <Real Name="Mass">1.020000</Real>
+          <Real Name="Charge">-1.020000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>3</Int>
+          </Sequence>
+          <Real Name="Mass">1.030000</Real>
+          <Real Name="Charge">-1.030000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>6</Int>
+          </Sequence>
+          <Real Name="Mass">1.060000</Real>
+          <Real Name="Charge">-1.060000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>7</Int>
+          </Sequence>
+          <Real Name="Mass">1.070000</Real>
+          <Real Name="Charge">-1.070000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>10</Int>
+          </Sequence>
+          <Real Name="Mass">1.100000</Real>
+          <Real Name="Charge">-1.100000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>11</Int>
+          </Sequence>
+          <Real Name="Mass">1.110000</Real>
+          <Real Name="Charge">-1.110000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>14</Int>
+          </Sequence>
+          <Real Name="Mass">1.140000</Real>
+          <Real Name="Charge">-1.140000</Real>
+        </Position>
+      </Sequence>
+    </Selection>
+    <Selection Name="Selection3">
+      <Sequence Name="Atoms">
+        <Int Name="Length">7</Int>
+        <Int>2</Int>
+        <Int>3</Int>
+        <Int>6</Int>
+        <Int>7</Int>
+        <Int>10</Int>
+        <Int>11</Int>
+        <Int>14</Int>
+      </Sequence>
+      <Sequence Name="Positions">
+        <Int Name="Length">5</Int>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>2</Int>
+          </Sequence>
+          <Real Name="Mass">1.020000</Real>
+          <Real Name="Charge">-1.020000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>3</Int>
+          </Sequence>
+          <Real Name="Mass">1.030000</Real>
+          <Real Name="Charge">-1.030000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">2</Int>
+            <Int>6</Int>
+            <Int>7</Int>
+          </Sequence>
+          <Real Name="Mass">2.130000</Real>
+          <Real Name="Charge">-2.130000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">2</Int>
+            <Int>10</Int>
+            <Int>11</Int>
+          </Sequence>
+          <Real Name="Mass">2.210000</Real>
+          <Real Name="Charge">-2.210000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>14</Int>
+          </Sequence>
+          <Real Name="Mass">1.140000</Real>
+          <Real Name="Charge">-1.140000</Real>
+        </Position>
+      </Sequence>
+    </Selection>
+  </EvaluatedSelections>
+</ReferenceData>
diff --git a/src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_ComputesMassesAndChargesWithoutTopology.xml b/src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_ComputesMassesAndChargesWithoutTopology.xml
new file mode 100644 (file)
index 0000000..f8380a2
--- /dev/null
@@ -0,0 +1,212 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <ParsedSelections Name="Parsed">
+    <ParsedSelection Name="Selection1">
+      <String Name="Input">atomnr 1 to 3 8 to 9</String>
+      <String Name="Name">atomnr 1 to 3 8 to 9</String>
+      <String Name="Text">atomnr 1 to 3 8 to 9</String>
+      <Bool Name="Dynamic">false</Bool>
+    </ParsedSelection>
+    <ParsedSelection Name="Selection2">
+      <String Name="Input">y &gt; 2</String>
+      <String Name="Name">y &gt; 2</String>
+      <String Name="Text">y &gt; 2</String>
+      <Bool Name="Dynamic">true</Bool>
+    </ParsedSelection>
+    <ParsedSelection Name="Selection3">
+      <String Name="Input">cog of (y &gt; 2)</String>
+      <String Name="Name">cog of (y &gt; 2)</String>
+      <String Name="Text">cog of (y &gt; 2)</String>
+      <Bool Name="Dynamic">true</Bool>
+    </ParsedSelection>
+  </ParsedSelections>
+  <CompiledSelections Name="Compiled">
+    <Selection Name="Selection1">
+      <Sequence Name="Atoms">
+        <Int Name="Length">5</Int>
+        <Int>0</Int>
+        <Int>1</Int>
+        <Int>2</Int>
+        <Int>7</Int>
+        <Int>8</Int>
+      </Sequence>
+      <Sequence Name="Positions">
+        <Int Name="Length">5</Int>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>0</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>1</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>2</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>7</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>8</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+      </Sequence>
+    </Selection>
+    <Selection Name="Selection2">
+      <Sequence Name="Atoms">
+        <Int Name="Length">10</Int>
+        <Int>0</Int>
+        <Int>1</Int>
+        <Int>2</Int>
+        <Int>3</Int>
+        <Int>4</Int>
+        <Int>5</Int>
+        <Int>6</Int>
+        <Int>7</Int>
+        <Int>8</Int>
+        <Int>9</Int>
+      </Sequence>
+      <Sequence Name="Positions">
+        <Int Name="Length">10</Int>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>0</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>1</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>2</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>3</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>4</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>5</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>6</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>7</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>8</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">1</Int>
+            <Int>9</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+      </Sequence>
+    </Selection>
+    <Selection Name="Selection3">
+      <Sequence Name="Atoms">
+        <Int Name="Length">10</Int>
+        <Int>0</Int>
+        <Int>1</Int>
+        <Int>2</Int>
+        <Int>3</Int>
+        <Int>4</Int>
+        <Int>5</Int>
+        <Int>6</Int>
+        <Int>7</Int>
+        <Int>8</Int>
+        <Int>9</Int>
+      </Sequence>
+      <Sequence Name="Positions">
+        <Int Name="Length">1</Int>
+        <Position>
+          <Sequence Name="Atoms">
+            <Int Name="Length">10</Int>
+            <Int>0</Int>
+            <Int>1</Int>
+            <Int>2</Int>
+            <Int>3</Int>
+            <Int>4</Int>
+            <Int>5</Int>
+            <Int>6</Int>
+            <Int>7</Int>
+            <Int>8</Int>
+            <Int>9</Int>
+          </Sequence>
+          <Real Name="Mass">1.000000</Real>
+          <Real Name="Charge">0.000000</Real>
+        </Position>
+      </Sequence>
+    </Selection>
+  </CompiledSelections>
+</ReferenceData>
diff --git a/src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesCharge.xml b/src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesCharge.xml
new file mode 100644 (file)
index 0000000..62e2d06
--- /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">charge &lt; 0.5</String>
+      <String Name="Name">charge &lt; 0.5</String>
+      <String Name="Text">charge &lt; 0.5</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>1</Int>
+        <Int>2</Int>
+        <Int>3</Int>
+        <Int>4</Int>
+      </Sequence>
+    </Selection>
+  </CompiledSelections>
+</ReferenceData>
diff --git a/src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesMass.xml b/src/gromacs/selection/tests/refdata/SelectionCollectionDataTest_HandlesMass.xml
new file mode 100644 (file)
index 0000000..76b1e46
--- /dev/null
@@ -0,0 +1,29 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <ParsedSelections Name="Parsed">
+    <ParsedSelection Name="Selection1">
+      <String Name="Input">mass &gt; 5</String>
+      <String Name="Name">mass &gt; 5</String>
+      <String Name="Text">mass &gt; 5</String>
+      <Bool Name="Dynamic">false</Bool>
+    </ParsedSelection>
+  </ParsedSelections>
+  <CompiledSelections Name="Compiled">
+    <Selection Name="Selection1">
+      <Sequence Name="Atoms">
+        <Int Name="Length">10</Int>
+        <Int>5</Int>
+        <Int>6</Int>
+        <Int>7</Int>
+        <Int>8</Int>
+        <Int>9</Int>
+        <Int>10</Int>
+        <Int>11</Int>
+        <Int>12</Int>
+        <Int>13</Int>
+        <Int>14</Int>
+      </Sequence>
+    </Selection>
+  </CompiledSelections>
+</ReferenceData>
index ba7e14958b3a1bff7e17bd35d1c559014ca6e294..a624262800b294ad4c2f28986ebe8d99c336e664 100644 (file)
                 <xsl:if test="Position/Vector[@Name='Coordinates']">
                     <th>Coordinates</th>
                 </xsl:if>
+                <xsl:if test="Position/Real[@Name='Mass']">
+                    <th>Mass</th>
+                </xsl:if>
+                <xsl:if test="Position/Real[@Name='Charge']">
+                    <th>Charge</th>
+                </xsl:if>
             </tr>
             <xsl:for-each select="Position">
             <tr>
                         <xsl:apply-templates select="Vector[@Name='Coordinates']"/>
                     </td>
                 </xsl:if>
+                <xsl:if test="Real[@Name='Mass']">
+                    <td><xsl:value-of select="Real[@Name='Mass']"/></td>
+                </xsl:if>
+                <xsl:if test="Real[@Name='Charge']">
+                    <td><xsl:value-of select="Real[@Name='Charge']"/></td>
+                </xsl:if>
             </tr>
             </xsl:for-each>
         </table>
index d03d9f4269998a9e19d4fdfdd89e368981d8c971..f0a6458d18ac38f7919b9cedebbfb59f0b519f0e 100644 (file)
@@ -160,6 +160,8 @@ class SelectionCollectionDataTest : public SelectionCollectionTest
             efTestPositionAtoms         = 1<<1,
             efTestPositionCoordinates   = 1<<2,
             efTestPositionMapping       = 1<<3,
+            efTestPositionMasses        = 1<<4,
+            efTestPositionCharges       = 1<<5,
             efDontTestCompiledAtoms     = 1<<8
         };
         typedef gmx::FlagsTemplate<TestFlag> TestFlags;
@@ -171,9 +173,6 @@ class SelectionCollectionDataTest : public SelectionCollectionTest
 
         void setFlags(TestFlags flags) { flags_ = flags; }
 
-        void runTest(int natoms, const char *const *selections, size_t count);
-        void runTest(const char *filename, const char *const *selections,
-                     size_t count);
         template <size_t count>
         void runTest(int natoms, const char *const (&selections)[count])
         {
@@ -185,15 +184,26 @@ class SelectionCollectionDataTest : public SelectionCollectionTest
             runTest(filename, selections, count);
         }
 
+        template <size_t count>
+        void runParser(const char *const (&selections)[count])
+        {
+            runParser(selections, count);
+        }
+
+        void runCompiler();
+        void runEvaluate();
+        void runEvaluateFinal();
+
     private:
         static void checkSelection(gmx::test::TestReferenceChecker *checker,
                                    const gmx::Selection &sel, TestFlags flags);
 
+        void runTest(int natoms, const char *const *selections, size_t count);
+        void runTest(const char *filename, const char *const *selections,
+                     size_t count);
         void runParser(const char *const *selections, size_t count);
-        void runCompiler();
+
         void checkCompiled();
-        void runEvaluate();
-        void runEvaluateFinal();
 
         gmx::test::TestReferenceData    data_;
         gmx::test::TestReferenceChecker checker_;
@@ -216,7 +226,9 @@ SelectionCollectionDataTest::checkSelection(
     }
     if (flags.test(efTestPositionAtoms)
         || flags.test(efTestPositionCoordinates)
-        || flags.test(efTestPositionMapping))
+        || flags.test(efTestPositionMapping)
+        || flags.test(efTestPositionMasses)
+        || flags.test(efTestPositionCharges))
     {
         TestReferenceChecker compound(
                 checker->checkSequenceCompound("Positions", sel.posCount()));
@@ -238,6 +250,14 @@ SelectionCollectionDataTest::checkSelection(
                 poscompound.checkInteger(p.refId(), "RefId");
                 poscompound.checkInteger(p.mappedId(), "MappedId");
             }
+            if (flags.test(efTestPositionMasses))
+            {
+                poscompound.checkReal(p.mass(), "Mass");
+            }
+            if (flags.test(efTestPositionCharges))
+            {
+                poscompound.checkReal(p.charge(), "Charge");
+            }
         }
     }
 }
@@ -567,8 +587,33 @@ TEST_F(SelectionCollectionDataTest, HandlesChain)
     runTest("simple.pdb", selections);
 }
 
-// TODO: Add test for mass
-// TODO: Add test for charge
+TEST_F(SelectionCollectionDataTest, HandlesMass)
+{
+    static const char * const selections[] = {
+        "mass > 5"
+    };
+    ASSERT_NO_FATAL_FAILURE(runParser(selections));
+    ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
+    for (int i = 0; i < top_->atoms.nr; ++i)
+    {
+        top_->atoms.atom[i].m = 1.0 + i;
+    }
+    ASSERT_NO_FATAL_FAILURE(runCompiler());
+}
+
+TEST_F(SelectionCollectionDataTest, HandlesCharge)
+{
+    static const char * const selections[] = {
+        "charge < 0.5"
+    };
+    ASSERT_NO_FATAL_FAILURE(runParser(selections));
+    ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
+    for (int i = 0; i < top_->atoms.nr; ++i)
+    {
+        top_->atoms.atom[i].q = i / 10.0;
+    }
+    ASSERT_NO_FATAL_FAILURE(runCompiler());
+}
 
 TEST_F(SelectionCollectionDataTest, HandlesAltLoc)
 {
@@ -741,6 +786,44 @@ TEST_F(SelectionCollectionDataTest, HandlesMergeModifier)
 }
 
 
+/********************************************************************
+ * Tests for generic selection evaluation
+ */
+
+TEST_F(SelectionCollectionDataTest, ComputesMassesAndCharges)
+{
+    static const char * const selections[] = {
+        "name CB",
+        "y > 2",
+        "res_cog of y > 2"
+    };
+    setFlags(TestFlags() | efTestEvaluation | efTestPositionAtoms
+             | efTestPositionMasses | efTestPositionCharges);
+    ASSERT_NO_FATAL_FAILURE(runParser(selections));
+    ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
+    for (int i = 0; i < top_->atoms.nr; ++i)
+    {
+        top_->atoms.atom[i].m =   1.0 + i / 100.0;
+        top_->atoms.atom[i].q = -(1.0 + i / 100.0);
+    }
+    ASSERT_NO_FATAL_FAILURE(runCompiler());
+    ASSERT_NO_FATAL_FAILURE(runEvaluate());
+    ASSERT_NO_FATAL_FAILURE(runEvaluateFinal());
+}
+
+TEST_F(SelectionCollectionDataTest, ComputesMassesAndChargesWithoutTopology)
+{
+    static const char * const selections[] = {
+        "atomnr 1 to 3 8 to 9",
+        "y > 2",
+        "cog of (y > 2)"
+    };
+    setFlags(TestFlags() | efTestPositionAtoms
+             | efTestPositionMasses | efTestPositionCharges);
+    runTest(10, selections);
+}
+
+
 /********************************************************************
  * Tests for selection syntactic constructs
  */