Simplify selection parsing methods signature.
authorTeemu Murtola <teemu.murtola@gmail.com>
Sat, 19 May 2012 06:29:04 +0000 (09:29 +0300)
committerTeemu Murtola <teemu.murtola@gmail.com>
Sat, 19 May 2012 06:29:04 +0000 (09:29 +0300)
Selection parsing methods in SelectionCollection now return the list of
parsed selections instead of using an output variable.  This makes it
more obvious from the signature alone how the methods should be used.
Also, the ability to append to an existing vector was not used in any
code currently (except in tests), and the code to do that with the new
signature is still simple.

Change-Id: I344e4c15367e956b4294d18f56f4da34ccf42e27

src/gromacs/selection/selectioncollection-impl.h
src/gromacs/selection/selectioncollection.cpp
src/gromacs/selection/selectioncollection.h
src/gromacs/selection/selectionoptionmanager.cpp
src/gromacs/selection/tests/selectioncollection.cpp

index ddf51f7ec17cdda4e7f5f02177252eaacdb6d224..5b36531ed2f1d0836b0312df2fa2c23f71eb10e8 100644 (file)
@@ -129,19 +129,14 @@ class SelectionCollection::Impl
          * \param[in,out] scanner Scanner data structure.
          * \param[in]     maxnr   Maximum number of selections to parse
          *      (if -1, parse as many as provided by the user).
-         * \param[out]    output  Vector to which parsed selections are
-         *      appended.
+         * \returns       Vector of parsed selections.
          * \throws        std::bad_alloc if out of memory.
          * \throws        InvalidInputError if there is a parsing error.
          *
-         * Parsed selections are appended to \p output without clearing it
-         * first.  If parsing fails, \p output is not modified.
-         *
          * Used internally to implement parseFromStdin(), parseFromFile() and
          * parseFromString().
          */
-        void runParser(void *scanner, int maxnr,
-                       SelectionList *output);
+        SelectionList runParser(void *scanner, int maxnr);
         /*! \brief
          * Replace group references by group contents.
          *
index 8be14ae1045d7d565b44f9f07587c2e48558570e..9314e072a959ad8d7574926201644d740f8a3f8c 100644 (file)
@@ -123,9 +123,8 @@ SelectionCollection::Impl::clearSymbolTable()
 }
 
 
-void
-SelectionCollection::Impl::runParser(yyscan_t scanner, int maxnr,
-                                     SelectionList *output)
+SelectionList
+SelectionCollection::Impl::runParser(yyscan_t scanner, int maxnr)
 {
     gmx_ana_selcollection_t *sc = &_sc;
     GMX_ASSERT(sc == _gmx_sel_lexer_selcollection(scanner),
@@ -144,22 +143,21 @@ SelectionCollection::Impl::runParser(yyscan_t scanner, int maxnr,
         errors.append("Too few selections provided");
     }
 
-    if (bOk)
-    {
-        SelectionDataList::const_iterator i;
-        output->reserve(output->size() + nr);
-        for (i = _sc.sel.begin() + oldCount; i != _sc.sel.end(); ++i)
-        {
-            output->push_back(Selection(i->get()));
-        }
-    }
     // TODO: Remove added selections from the collection if parsing failed?
-
     if (!bOk || !errors.isEmpty())
     {
         GMX_ASSERT(!bOk && !errors.isEmpty(), "Inconsistent error reporting");
         GMX_THROW(InvalidInputError(errors.toString()));
     }
+
+    SelectionList result;
+    SelectionDataList::const_iterator i;
+    result.reserve(nr);
+    for (i = _sc.sel.begin() + oldCount; i != _sc.sel.end(); ++i)
+    {
+        result.push_back(Selection(i->get()));
+    }
+    return result;
 }
 
 
@@ -380,9 +378,8 @@ SelectionCollection::requiresTopology() const
 }
 
 
-void
-SelectionCollection::parseFromStdin(int nr, bool bInteractive,
-                                    SelectionList *output)
+SelectionList
+SelectionCollection::parseFromStdin(int nr, bool bInteractive)
 {
     yyscan_t scanner;
 
@@ -391,13 +388,12 @@ SelectionCollection::parseFromStdin(int nr, bool bInteractive,
                         _impl->_grps);
     /* We don't set the lexer input here, which causes it to use a special
      * internal implementation for reading from stdin. */
-    _impl->runParser(scanner, nr, output);
+    return _impl->runParser(scanner, nr);
 }
 
 
-void
-SelectionCollection::parseFromFile(const std::string &filename,
-                                   SelectionList *output)
+SelectionList
+SelectionCollection::parseFromFile(const std::string &filename)
 {
     yyscan_t scanner;
 
@@ -407,14 +403,12 @@ SelectionCollection::parseFromFile(const std::string &filename,
                         _impl->_bExternalGroupsSet,
                         _impl->_grps);
     _gmx_sel_set_lex_input_file(scanner, file.handle());
-    _impl->runParser(scanner, -1, output);
-    file.close();
+    return _impl->runParser(scanner, -1);
 }
 
 
-void
-SelectionCollection::parseFromString(const std::string &str,
-                                     SelectionList *output)
+SelectionList
+SelectionCollection::parseFromString(const std::string &str)
 {
     yyscan_t scanner;
 
@@ -422,7 +416,7 @@ SelectionCollection::parseFromString(const std::string &str,
                         _impl->_bExternalGroupsSet,
                         _impl->_grps);
     _gmx_sel_set_lex_input_str(scanner, str.c_str());
-    _impl->runParser(scanner, -1, output);
+    return _impl->runParser(scanner, -1);
 }
 
 
index 6b74eeea00bca4a3a79e843a8d4b1e906811be82..5c19cf84ce12f5fe6674277109ed3dcec452ae76 100644 (file)
@@ -232,58 +232,46 @@ class SelectionCollection
          *      (if -1, parse as many as provided by the user).
          * \param[in]  bInteractive Whether the parser should behave
          *      interactively.
-         * \param[out] output   Vector to which parsed selections are appended.
+         * \returns    Vector of parsed selections.
          * \throws     std::bad_alloc if out of memory.
          * \throws     InvalidInputError if there is a parsing error
          *      (an interactive parser only throws this if too few selections
          *      are provided and the user forced the end of input).
          *
-         * Parsed selections are appended to \p output without clearing it
-         * first.  If parsing fails, \p output is not modified.
-         *
-         * The objects returned in \p output remain valid for the lifetime of
+         * The returned objects remain valid for the lifetime of
          * the selection collection.
          * Some information about the selections only becomes available once
          * compile() has been called; see \ref Selection.
          */
-        void parseFromStdin(int count, bool bInteractive,
-                            SelectionList *output);
+        SelectionList parseFromStdin(int count, bool bInteractive);
         /*! \brief
          * Parses selection(s) from a file.
          *
          * \param[in]  filename Name of the file to parse selections from.
-         * \param[out] output   Vector to which parsed selections are appended.
+         * \returns    Vector of parsed selections.
          * \throws     std::bad_alloc if out of memory.
          * \throws     InvalidInputError if there is a parsing error.
          *
-         * Parsed selections are appended to \p output without clearing it
-         * first.  If parsing fails, \p output is not modified.
-         *
-         * The objects returned in \p output remain valid for the lifetime of
+         * The returned objects remain valid for the lifetime of
          * the selection collection.
          * Some information about the selections only becomes available once
          * compile() has been called; see \ref Selection.
          */
-        void parseFromFile(const std::string &filename,
-                           SelectionList *output);
+        SelectionList parseFromFile(const std::string &filename);
         /*! \brief
          * Parses selection(s) from a string.
          *
          * \param[in]  str      String to parse selections from.
-         * \param[out] output   Vector to which parsed selections are appended.
+         * \returns    Vector of parsed selections.
          * \throws     std::bad_alloc if out of memory.
          * \throws     InvalidInputError if there is a parsing error.
          *
-         * Parsed selections are appended to \p output without clearing it
-         * first.  If parsing fails, \p output is not modified.
-         *
-         * The objects returned in \p output remain valid for the lifetime of
+         * The returned objects remain valid for the lifetime of
          * the selection collection.
          * Some information about the selections only becomes available once
          * compile() has been called; see \ref Selection.
          */
-        void parseFromString(const std::string &str,
-                             SelectionList *output);
+        SelectionList parseFromString(const std::string &str);
         /*! \brief
          * Prepares the selections for evaluation and performs optimizations.
          *
index fc905d6f01195e9dc2aab4878c1052e93761f652..5281dd4216df77936992c6e17f9f05b2c46dcf90 100644 (file)
@@ -248,8 +248,7 @@ void
 SelectionOptionManager::convertOptionValue(SelectionOptionStorage *storage,
                                            const std::string &value)
 {
-    SelectionList selections;
-    impl_->collection_.parseFromString(value, &selections);
+    SelectionList selections = impl_->collection_.parseFromString(value);
     storage->addSelections(selections, false);
 }
 
@@ -289,8 +288,8 @@ SelectionOptionManager::parseRequestedFromStdin(bool bInteractive)
             std::fprintf(stderr, "(one selection per line, 'help' for help%s)\n",
                          request.count() < 0 ? ", Ctrl-D to end" : "");
         }
-        SelectionList selections;
-        impl_->collection_.parseFromStdin(request.count(), bInteractive, &selections);
+        SelectionList selections
+            = impl_->collection_.parseFromStdin(request.count(), bInteractive);
         request.storage_->addSelections(selections, true);
     }
 }
@@ -298,16 +297,14 @@ SelectionOptionManager::parseRequestedFromStdin(bool bInteractive)
 void
 SelectionOptionManager::parseRequestedFromFile(const std::string &filename)
 {
-    SelectionList selections;
-    impl_->collection_.parseFromFile(filename, &selections);
+    SelectionList selections = impl_->collection_.parseFromFile(filename);
     impl_->placeSelectionsInRequests(selections);
 }
 
 void
 SelectionOptionManager::parseRequestedFromString(const std::string &str)
 {
-    SelectionList selections;
-    impl_->collection_.parseFromString(str, &selections);
+    SelectionList selections = impl_->collection_.parseFromString(str);
     impl_->placeSelectionsInRequests(selections);
 }
 
index 3eb3d597b17afe8d426c4abebdc3f5124e392cc6..5e7edb04f69f2f133a4967c44fd17e29b6307792 100644 (file)
@@ -232,7 +232,9 @@ SelectionCollectionDataTest::runParser(const char *const *selections)
     {
         SCOPED_TRACE(std::string("Parsing selection \"")
                      + selections[i] + "\"");
-        ASSERT_NO_THROW(_sc.parseFromString(selections[i], &_sel));
+        gmx::SelectionList result;
+        ASSERT_NO_THROW(result = _sc.parseFromString(selections[i]));
+        _sel.insert(_sel.end(), result.begin(), result.end());
         if (_sel.size() == _count)
         {
             std::string id = gmx::formatString("Variable%d", static_cast<int>(varcount + 1));
@@ -352,8 +354,7 @@ TEST_F(SelectionCollectionTest, HandlesNoSelections)
 
 TEST_F(SelectionCollectionTest, ParsesSelectionsFromFile)
 {
-    ASSERT_NO_THROW(_sc.parseFromFile(gmx::test::getTestFilePath("selfile.dat"),
-                                      &_sel));
+    ASSERT_NO_THROW(_sel = _sc.parseFromFile(gmx::test::getTestFilePath("selfile.dat")));
     // These should match the contents of selfile.dat
     ASSERT_EQ(2U, _sel.size());
     EXPECT_STREQ("resname RA RB", _sel[0].selectionText());
@@ -362,19 +363,19 @@ TEST_F(SelectionCollectionTest, ParsesSelectionsFromFile)
 
 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue)
 {
-    EXPECT_THROW(_sc.parseFromString("mindist from atomnr 1 cutoff", &_sel),
+    EXPECT_THROW(_sc.parseFromString("mindist from atomnr 1 cutoff"),
                  gmx::InvalidInputError);
 }
 
 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue2)
 {
-    EXPECT_THROW(_sc.parseFromString("within 1 of", &_sel),
+    EXPECT_THROW(_sc.parseFromString("within 1 of"),
                  gmx::InvalidInputError);
 }
 
 TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue3)
 {
-    EXPECT_THROW(_sc.parseFromString("within of atomnr 1", &_sel),
+    EXPECT_THROW(_sc.parseFromString("within of atomnr 1"),
                  gmx::InvalidInputError);
 }
 
@@ -382,7 +383,7 @@ TEST_F(SelectionCollectionTest, HandlesMissingMethodParamValue3)
 
 TEST_F(SelectionCollectionTest, RecoversFromUnknownGroupReference)
 {
-    ASSERT_NO_THROW(_sc.parseFromString("group \"foo\"", &_sel));
+    ASSERT_NO_THROW(_sc.parseFromString("group \"foo\""));
     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
     EXPECT_THROW(_sc.setIndexGroups(NULL), gmx::InvalidInputError);
     EXPECT_THROW(_sc.compile(), gmx::APIError);
@@ -390,42 +391,42 @@ TEST_F(SelectionCollectionTest, RecoversFromUnknownGroupReference)
 
 TEST_F(SelectionCollectionTest, RecoversFromMissingMoleculeInfo)
 {
-    ASSERT_NO_THROW(_sc.parseFromString("molindex 1 to 5", &_sel));
+    ASSERT_NO_THROW(_sc.parseFromString("molindex 1 to 5"));
     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
     EXPECT_THROW(_sc.compile(), gmx::InconsistentInputError);
 }
 
 TEST_F(SelectionCollectionTest, RecoversFromMissingAtomTypes)
 {
-    ASSERT_NO_THROW(_sc.parseFromString("type CA", &_sel));
+    ASSERT_NO_THROW(_sc.parseFromString("type CA"));
     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
     EXPECT_THROW(_sc.compile(), gmx::InconsistentInputError);
 }
 
 TEST_F(SelectionCollectionTest, RecoversFromMissingPDBInfo)
 {
-    ASSERT_NO_THROW(_sc.parseFromString("altloc A", &_sel));
+    ASSERT_NO_THROW(_sc.parseFromString("altloc A"));
     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
     EXPECT_THROW(_sc.compile(), gmx::InconsistentInputError);
 }
 
 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation)
 {
-    ASSERT_NO_THROW(_sc.parseFromString("all permute 1 1", &_sel));
+    ASSERT_NO_THROW(_sc.parseFromString("all permute 1 1"));
     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
     EXPECT_THROW(_sc.compile(), gmx::InvalidInputError);
 }
 
 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation2)
 {
-    ASSERT_NO_THROW(_sc.parseFromString("all permute 3 2 1", &_sel));
+    ASSERT_NO_THROW(_sc.parseFromString("all permute 3 2 1"));
     ASSERT_NO_FATAL_FAILURE(setAtomCount(10));
     EXPECT_THROW(_sc.compile(), gmx::InconsistentInputError);
 }
 
 TEST_F(SelectionCollectionTest, RecoversFromInvalidPermutation3)
 {
-    ASSERT_NO_THROW(_sc.parseFromString("x < 1.5 permute 3 2 1", &_sel));
+    ASSERT_NO_THROW(_sc.parseFromString("x < 1.5 permute 3 2 1"));
     ASSERT_NO_FATAL_FAILURE(loadTopology("simple.gro"));
     ASSERT_NO_THROW(_sc.compile());
     EXPECT_THROW(_sc.evaluate(_frame, NULL), gmx::InconsistentInputError);