Fix velocity and force evaluation for empty selections.
authorTeemu Murtola <teemu.murtola@gmail.com>
Sat, 9 Feb 2013 15:36:07 +0000 (17:36 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Wed, 13 Feb 2013 21:27:26 +0000 (22:27 +0100)
Requesting velocity or force evaluation and then passing an empty
selection resulted in an assert (in multiple places in the code).
Added a test for the behavior in this case.

Change-Id: I70cfb43c652f709fd555b4f23467f71fdb4bd3e7

src/gromacs/selection/position.cpp
src/gromacs/selection/selection.h
src/gromacs/selection/selectionoption.h
src/gromacs/selection/tests/selectioncollection.cpp

index b2c61669096e962d115e5a1768d893e4dd6360c9..93f3f545b5e3f33850bc7221dcec6adda8462aa2 100644 (file)
@@ -73,6 +73,14 @@ gmx_ana_pos_clear(gmx_ana_pos_t *pos)
 void
 gmx_ana_pos_reserve(gmx_ana_pos_t *pos, int n, int isize)
 {
+    GMX_RELEASE_ASSERT(n >= 0, "Invalid position allocation count");
+    // Always reserve at least one entry to make NULL checks against pos->x
+    // and gmx_ana_pos_reserve_velocities/forces() work as expected in the case
+    // that there are actually no positions.
+    if (n == 0)
+    {
+        n = 1;
+    }
     if (pos->nalloc_x < n)
     {
         pos->nalloc_x = n;
@@ -96,7 +104,7 @@ gmx_ana_pos_reserve(gmx_ana_pos_t *pos, int n, int isize)
  * \param[in,out] pos   Position data structure.
  *
  * Currently, this function can only be called after gmx_ana_pos_reserve()
- * has been called at least once with a \p n > 0.
+ * has been called at least once with a \p n >= 0.
  */
 void
 gmx_ana_pos_reserve_velocities(gmx_ana_pos_t *pos)
@@ -113,7 +121,7 @@ gmx_ana_pos_reserve_velocities(gmx_ana_pos_t *pos)
  * \param[in,out] pos   Position data structure.
  *
  * Currently, this function can only be called after gmx_ana_pos_reserve()
- * has been called at least once with a \p n > 0.
+ * has been called at least once with a \p n >= 0.
  */
 void
 gmx_ana_pos_reserve_forces(gmx_ana_pos_t *pos)
index c87f5b40585e774a06e2346e049703dfbd383575..00d5b63c48170e2d77fbe462abccb70e37bc85f4 100644 (file)
@@ -210,6 +210,11 @@ class SelectionData
  * any time after the selection has been parsed, and type() can be accessed
  * after the selection has been compiled.
  *
+ * There are a few methods that can be used to change the behavior of the
+ * selection.  setEvaluateVelocities() and setEvaluateForces() can be called
+ * before the selection is compiled to request evaluation of velocities and/or
+ * forces in addition to coordinates.
+ *
  * Each selection is made of a set of positions.  Each position has associated
  * coordinates, and possibly velocities and forces if they have been requested
  * and are available.  It also has a set of atoms associated with it; typically
@@ -370,22 +375,6 @@ class Selection
         {
             return ConstArrayRef<int>(data().rawPositions_.m.mapid, posCount());
         }
-        /*! \brief
-         * Sets the ID for the \p i'th position for use with
-         * SelectionPosition::mappedId().
-         *
-         * \param[in] i  Zero-based index
-         * \param[in] id Identifier to set.
-         *
-         * This method is not part of SelectionPosition because that interface
-         * only provides access to const data by design.
-         *
-         * This method can only be called after compilation, before the
-         * selection has been evaluated for any frame.
-         *
-         * \see SelectionPosition::mappedId()
-         */
-        void setOriginalId(int i, int id) { data().rawPositions_.m.orgid[i] = id; }
 
         //! Deprecated method for direct access to position data.
         const gmx_ana_pos_t *positions() const { return &data().rawPositions_; }
@@ -394,6 +383,7 @@ class Selection
         bool isCoveredFractionDynamic() const { return data().isCoveredFractionDynamic(); }
         //! Returns the covered fraction for the current frame.
         real coveredFraction() const { return data().coveredFraction_; }
+
         /*! \brief
          * Initializes information about covered fractions.
          *
@@ -405,6 +395,57 @@ class Selection
         {
             return data().initCoveredFraction(type);
         }
+        /*! \brief
+         * Sets whether this selection evaluates velocities for positions.
+         *
+         * \param[in] bEnabled  If true, velocities are evaluated.
+         *
+         * If you request the evaluation, but then evaluate the selection for
+         * a frame that does not contain velocity information, results are
+         * undefined.
+         *
+         * \todo
+         * Implement it such that in the above case, hasVelocities() will
+         * return false for such frames.
+         *
+         * Does not throw.
+         */
+        void setEvaluateVelocities(bool bEnabled)
+        {
+            data().flags_.set(efSelection_EvaluateVelocities, bEnabled);
+        }
+        /*! \brief
+         * Sets whether this selection evaluates forces for positions.
+         *
+         * \param[in] bEnabled  If true, forces are evaluated.
+         *
+         * If you request the evaluation, but then evaluate the selection for
+         * a frame that does not contain force information, results are
+         * undefined.
+         *
+         * Does not throw.
+         */
+        void setEvaluateForces(bool bEnabled)
+        {
+            data().flags_.set(efSelection_EvaluateForces, bEnabled);
+        }
+
+        /*! \brief
+         * Sets the ID for the \p i'th position for use with
+         * SelectionPosition::mappedId().
+         *
+         * \param[in] i  Zero-based index
+         * \param[in] id Identifier to set.
+         *
+         * This method is not part of SelectionPosition because that interface
+         * only provides access to const data by design.
+         *
+         * This method can only be called after compilation, before the
+         * selection has been evaluated for any frame.
+         *
+         * \see SelectionPosition::mappedId()
+         */
+        void setOriginalId(int i, int id) { data().rawPositions_.m.orgid[i] = id; }
 
         /*! \brief
          * Prints out one-line description of the selection.
index e5aa27c35ef92f596cfb93a38031822af36e67b9..c457b2773c1207dec1ea93a6155e0aca51def609 100644 (file)
@@ -79,16 +79,14 @@ class SelectionOption : public OptionTemplate<Selection, SelectionOption>
         /*! \brief
          * Request velocity evaluation for output positions.
          *
-         * Note that even with this flag set, velocities may not be available,
-         * in which case Selection::hasVelocities() returns false.
+         * \see Selection::setEvaluateVelocities()
          */
         MyClass &evaluateVelocities()
         { selectionFlags_.set(efSelection_EvaluateVelocities); return me(); }
         /*! \brief
          * Request force evaluation for output positions.
          *
-         * Note that even with this flag set, forces may not be available,
-         * in which case Selection::hasForces() returns false.
+         * \see Selection::setEvaluateForces()
          */
         MyClass &evaluateForces()
         { selectionFlags_.set(efSelection_EvaluateForces); return me(); }
@@ -227,7 +225,7 @@ class SelectionOptionInfo : public OptionInfo
          *
          * Does not throw.
          *
-         * \see SelectionOption::evaluateVelocities()
+         * \see Selection::setEvaluateVelocities()
          */
         void setEvaluateVelocities(bool bEnabled);
         /*! \brief
@@ -237,7 +235,7 @@ class SelectionOptionInfo : public OptionInfo
          *
          * Does not throw.
          *
-         * \see SelectionOption::evaluateForces()
+         * \see Selection::setEvaluateForces()
          */
         void setEvaluateForces(bool bEnabled);
         /*! \brief
index f0a6458d18ac38f7919b9cedebbfb59f0b519f0e..cf5be54c2f46749a4232d0cfa346aff863eae43f 100644 (file)
@@ -402,6 +402,22 @@ TEST_F(SelectionCollectionTest, HandlesNoSelections)
     EXPECT_NO_THROW(sc_.compile());
 }
 
+TEST_F(SelectionCollectionTest, HandlesVelocityAndForceRequests)
+{
+    ASSERT_NO_THROW(sel_ = sc_.parseFromString("atomnr 1 to 10; none"));
+    ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
+    ASSERT_EQ(2U, sel_.size());
+    ASSERT_NO_THROW(sel_[0].setEvaluateVelocities(true));
+    ASSERT_NO_THROW(sel_[1].setEvaluateVelocities(true));
+    ASSERT_NO_THROW(sel_[0].setEvaluateForces(true));
+    ASSERT_NO_THROW(sel_[1].setEvaluateForces(true));
+    ASSERT_NO_THROW(sc_.compile());
+    EXPECT_TRUE(sel_[0].hasVelocities());
+    EXPECT_TRUE(sel_[1].hasVelocities());
+    EXPECT_TRUE(sel_[0].hasForces());
+    EXPECT_TRUE(sel_[1].hasForces());
+}
+
 TEST_F(SelectionCollectionTest, ParsesSelectionsFromFile)
 {
     ASSERT_NO_THROW(sel_ = sc_.parseFromFile(