Improve interactive selection input
authorTeemu Murtola <teemu.murtola@gmail.com>
Wed, 12 Mar 2014 19:47:42 +0000 (21:47 +0200)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Tue, 18 Mar 2014 08:50:12 +0000 (09:50 +0100)
- Print the static index groups always at the start of the prompt.
- Add a hint for that an empty line provides information about the
  available index groups.
- Print the whole prompt again after an empty input line.
- Print only the selections from the current interactive prompt in the
  output.  This makes the listing clearer for commands like
    gmx sasa -surface -output
  when providing the latter set of selections interactively.
- Print the number of selections still required after the above list.
- Work better (on Unix) if stdin is from a pipe (using the same
  mechanism that was there in 4.5 and 4.6).

Change-Id: I38152fb4fccd142aeee318739ab3782fa6a775d8

src/gromacs/selection/parsetree.cpp
src/gromacs/selection/selectioncollection.cpp
src/gromacs/selection/selectioncollection.h
src/gromacs/selection/selectionoptionmanager.cpp
src/gromacs/trajectoryanalysis/cmdlinerunner.cpp
src/gromacs/utility/file.cpp
src/gromacs/utility/file.h

index 98627b2125c6fdbd554a64ec001c8f52b7719bd1..83b040cfdf246fdf0fa0b8f4666559e8909c5a40 100644 (file)
@@ -1164,38 +1164,12 @@ _gmx_sel_parser_should_finish(yyscan_t scanner)
     return (int)sc->sel.size() == _gmx_sel_lexer_exp_selcount(scanner);
 }
 
-/*!
- * \param[in] scanner Scanner data structure.
- */
 void
-_gmx_sel_handle_empty_cmd(yyscan_t scanner)
+_gmx_sel_handle_empty_cmd(yyscan_t /*scanner*/)
 {
-    gmx_ana_selcollection_t *sc   = _gmx_sel_lexer_selcollection(scanner);
-    gmx_ana_indexgrps_t     *grps = _gmx_sel_lexer_indexgrps(scanner);
-    int                      i;
-
-    if (!_gmx_sel_is_lexer_interactive(scanner))
-    {
-        return;
-    }
-
-    if (grps)
-    {
-        fprintf(stderr, "Available index groups:\n");
-        gmx_ana_indexgrps_print(stderr, _gmx_sel_lexer_indexgrps(scanner), 0);
-    }
-    if (sc->nvars > 0 || !sc->sel.empty())
-    {
-        fprintf(stderr, "Currently provided selections:\n");
-        for (i = 0; i < sc->nvars; ++i)
-        {
-            fprintf(stderr, "     %s\n", sc->varstrs[i]);
-        }
-        for (i = 0; i < (int)sc->sel.size(); ++i)
-        {
-            fprintf(stderr, " %2d. %s\n", i+1, sc->sel[i]->selectionText());
-        }
-    }
+    // This is now handled outside the actual parser, in
+    // selectioncollection.cpp.  This stub will be removed as part of later
+    // refactoring related to the selection parser.
 }
 
 /*!
index 0336378188690c00568446bcf732be02d4fbffa9..614dd01ae7babb1f44e4b86cf0a56df29bc544e7 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2010,2011,2012,2013, by the GROMACS development team, led by
+ * Copyright (c) 2010,2011,2012,2013,2014, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -206,6 +206,70 @@ int runParserLoop(yyscan_t scanner, _gmx_sel_yypstate *parserState,
     return status;
 }
 
+/*! \brief
+ * Print current status in response to empty line in interactive input.
+ *
+ * \param[in] sc             Selection collection data structure.
+ * \param[in] grps           Available index groups.
+ * \param[in] firstSelection Index of first selection from this interactive
+ *     session.
+ * \param[in] maxCount       Maximum number of selections.
+ * \param[in] context        Context to print for what the selections are for.
+ * \param[in] bFirst         Whether this is the header that is printed before
+ *     any user input.
+ *
+ * Prints the available index groups and currently provided selections.
+ */
+void printCurrentStatus(gmx_ana_selcollection_t *sc, gmx_ana_indexgrps_t *grps,
+                        size_t firstSelection, int maxCount,
+                        const std::string &context, bool bFirst)
+{
+    if (grps != NULL)
+    {
+        std::fprintf(stderr, "Available static index groups:\n");
+        gmx_ana_indexgrps_print(stderr, grps, 0);
+    }
+    std::fprintf(stderr, "Specify ");
+    if (maxCount < 0)
+    {
+        std::fprintf(stderr, "any number of selections");
+    }
+    else if (maxCount == 1)
+    {
+        std::fprintf(stderr, "a selection");
+    }
+    else
+    {
+        std::fprintf(stderr, "%d selections", maxCount);
+    }
+    std::fprintf(stderr, "%s%s:\n",
+                 context.empty() ? "" : " ", context.c_str());
+    std::fprintf(stderr,
+                 "(one per line, <enter> for status/groups, 'help' for help%s)\n",
+                 maxCount < 0 ? ", Ctrl-D to end" : "");
+    if (!bFirst && (sc->nvars > 0 || sc->sel.size() > firstSelection))
+    {
+        std::fprintf(stderr, "Currently provided selections:\n");
+        for (int i = 0; i < sc->nvars; ++i)
+        {
+            std::fprintf(stderr, "     %s\n", sc->varstrs[i]);
+        }
+        for (size_t i = firstSelection; i < sc->sel.size(); ++i)
+        {
+            std::fprintf(stderr, " %2d. %s\n",
+                         static_cast<int>(i - firstSelection + 1),
+                         sc->sel[i]->selectionText());
+        }
+        if (maxCount > 0)
+        {
+            const int remaining
+                = maxCount - static_cast<int>(sc->sel.size() - firstSelection);
+            std::fprintf(stderr, "(%d more selection%s required)\n",
+                         remaining, remaining > 1 ? "s" : "");
+        }
+    }
+}
+
 /*! \brief
  * Helper function that runs the parser once the tokenizer has been
  * initialized.
@@ -215,6 +279,7 @@ int runParserLoop(yyscan_t scanner, _gmx_sel_yypstate *parserState,
  *      algorithm designed for interactive input.
  * \param[in]     maxnr   Maximum number of selections to parse
  *      (if -1, parse as many as provided by the user).
+ * \param[in]     context Context to print for what the selections are for.
  * \returns       Vector of parsed selections.
  * \throws        std::bad_alloc if out of memory.
  * \throws        InvalidInputError if there is a parsing error.
@@ -222,27 +287,41 @@ int runParserLoop(yyscan_t scanner, _gmx_sel_yypstate *parserState,
  * Used internally to implement parseFromStdin(), parseFromFile() and
  * parseFromString().
  */
-SelectionList runParser(yyscan_t scanner, bool bStdIn, int maxnr)
+SelectionList runParser(yyscan_t scanner, bool bStdIn, int maxnr,
+                        const std::string &context)
 {
     boost::shared_ptr<void>  scannerGuard(scanner, &_gmx_sel_free_lexer);
-    gmx_ana_selcollection_t *sc = _gmx_sel_lexer_selcollection(scanner);
+    gmx_ana_selcollection_t *sc   = _gmx_sel_lexer_selcollection(scanner);
+    gmx_ana_indexgrps_t     *grps = _gmx_sel_lexer_indexgrps(scanner);
 
     MessageStringCollector   errors;
     _gmx_sel_set_lexer_error_reporter(scanner, &errors);
 
-    int  oldCount = sc->sel.size();
-    bool bOk      = false;
+    size_t oldCount = sc->sel.size();
+    bool   bOk      = false;
     {
         boost::shared_ptr<_gmx_sel_yypstate> parserState(
                 _gmx_sel_yypstate_new(), &_gmx_sel_yypstate_delete);
         if (bStdIn)
         {
             File       &stdinFile(File::standardInput());
-            bool        bInteractive = _gmx_sel_is_lexer_interactive(scanner);
+            const bool  bInteractive = _gmx_sel_is_lexer_interactive(scanner);
+            if (bInteractive)
+            {
+                printCurrentStatus(sc, grps, oldCount, maxnr, context, true);
+            }
             std::string line;
             int         status;
             while (promptLine(&stdinFile, bInteractive, &line))
             {
+                if (bInteractive)
+                {
+                    if (line.empty())
+                    {
+                        printCurrentStatus(sc, grps, oldCount, maxnr, context, false);
+                        continue;
+                    }
+                }
                 line.append("\n");
                 _gmx_sel_set_lex_input_str(scanner, line.c_str());
                 status = runParserLoop(scanner, parserState.get(), true);
@@ -500,14 +579,15 @@ SelectionCollection::requiresTopology() const
 
 
 SelectionList
-SelectionCollection::parseFromStdin(int nr, bool bInteractive)
+SelectionCollection::parseFromStdin(int nr, bool bInteractive,
+                                    const std::string &context)
 {
     yyscan_t scanner;
 
     _gmx_sel_init_lexer(&scanner, &impl_->sc_, bInteractive, nr,
                         impl_->bExternalGroupsSet_,
                         impl_->grps_);
-    return runParser(scanner, true, nr);
+    return runParser(scanner, true, nr, context);
 }
 
 
@@ -524,7 +604,7 @@ SelectionCollection::parseFromFile(const std::string &filename)
                             impl_->bExternalGroupsSet_,
                             impl_->grps_);
         _gmx_sel_set_lex_input_file(scanner, file.handle());
-        return runParser(scanner, false, -1);
+        return runParser(scanner, false, -1, std::string());
     }
     catch (GromacsException &ex)
     {
@@ -545,7 +625,7 @@ SelectionCollection::parseFromString(const std::string &str)
                         impl_->bExternalGroupsSet_,
                         impl_->grps_);
     _gmx_sel_set_lex_input_str(scanner, str.c_str());
-    return runParser(scanner, false, -1);
+    return runParser(scanner, false, -1, std::string());
 }
 
 
index 7da3d427d46c621efcf09b04e3d2aea69cc3cd13..cf79b660bd60cffe3dd87f6d2b16c2dc9ecc6ad2 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2010,2011,2012,2013, by the GROMACS development team, led by
+ * Copyright (c) 2010,2011,2012,2013,2014, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -245,6 +245,7 @@ class SelectionCollection
          *      (if -1, parse as many as provided by the user).
          * \param[in]  bInteractive Whether the parser should behave
          *      interactively.
+         * \param[in]  context  Context to print for interactive input.
          * \returns    Vector of parsed selections.
          * \throws     std::bad_alloc if out of memory.
          * \throws     InvalidInputError if there is a parsing error
@@ -255,8 +256,12 @@ class SelectionCollection
          * the selection collection.
          * Some information about the selections only becomes available once
          * compile() has been called; see \ref Selection.
+         *
+         * The string provided to \p context should be such that it can replace
+         * the three dots in "Specify selections ...:".  It can be empty.
          */
-        SelectionList parseFromStdin(int count, bool bInteractive);
+        SelectionList parseFromStdin(int count, bool bInteractive,
+                                     const std::string &context);
         /*! \brief
          * Parses selection(s) from a file.
          *
index 4094fc231d715f707104c89fc3a70b8bf83060b9..3600bd29233994a659941d064dcb21cd43be3b4b 100644 (file)
@@ -297,28 +297,12 @@ SelectionOptionManager::parseRequestedFromStdin(bool bInteractive)
     for (i = impl_->requests_.begin(); i != impl_->requests_.end(); ++i)
     {
         const Impl::SelectionRequest &request = *i;
-        if (bInteractive)
-        {
-            std::fprintf(stderr, "\nSpecify ");
-            if (request.count() < 0)
-            {
-                std::fprintf(stderr, "any number of selections");
-            }
-            else if (request.count() == 1)
-            {
-                std::fprintf(stderr, "a selection");
-            }
-            else
-            {
-                std::fprintf(stderr, "%d selections", request.count());
-            }
-            std::fprintf(stderr, " for option '%s' (%s):\n",
+        std::string                   context =
+            formatString("for option '%s'\n(%s)",
                          request.name().c_str(), request.description().c_str());
-            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);
+            = impl_->collection_.parseFromStdin(request.count(), bInteractive,
+                                                context);
         request.storage_->addSelections(selections, true);
     }
 }
index 9b6d7a79d755f08305e85cb4ca3d1506317e2962..2838db6d5ae15a7c2079555d4db91b976edafc8f 100644 (file)
@@ -138,8 +138,7 @@ TrajectoryAnalysisCommandLineRunner::Impl::parseOptions(
 
     common->initIndexGroups(selections, bUseDefaultGroups_);
 
-    // TODO: Check whether the input is a pipe.
-    const bool bInteractive = true;
+    const bool bInteractive = File::standardInput().isInteractive();
     seloptManager.parseRequestedFromStdin(bInteractive);
     common->doneIndexGroups(selections);
 
index b0adbce9d9050fde8c831c69a6710b2729ec5448..787cb72c5ab33702d54f68ba192bf17752712157 100644 (file)
 #include <string>
 #include <vector>
 
+#include "config.h"
+
 #include <sys/stat.h>
+#ifdef HAVE_UNISTD_H
+#include <unistd.h>
+#endif
 
 #include "gromacs/utility/exceptions.h"
 #include "gromacs/utility/gmxassert.h"
 #include "gromacs/utility/stringutil.h"
 
-#include "gmx_header_config.h"
-
 namespace gmx
 {
 
@@ -158,6 +161,17 @@ void File::close()
     }
 }
 
+bool File::isInteractive() const
+{
+    GMX_RELEASE_ASSERT(impl_->fp_ != NULL,
+                       "Attempted to access a file object that is not open");
+#ifdef HAVE_UNISTD_H
+    return isatty(fileno(impl_->fp_));
+#else
+    return true;
+#endif
+}
+
 FILE *File::handle()
 {
     GMX_RELEASE_ASSERT(impl_->fp_ != NULL,
index fbd6353c5f727b5167b371ff41ab81054a6f9409..15b48d16a559533fcb7507d625e3e008df17a6d2 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2012,2013, by the GROMACS development team, led by
+ * Copyright (c) 2012,2013,2014, by the GROMACS development team, led by
  * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
  * and including many others, as listed in the AUTHORS file in the
  * top-level source directory and at http://www.gromacs.org.
@@ -106,6 +106,17 @@ class File
          */
         void close();
 
+        /*! \brief
+         * Returns whether the file is an interactive terminal.
+         *
+         * Only works on Unix, otherwise always returns true.
+         * It only makes sense to call this for File::standardInput() and
+         * friends.
+         *
+         * Thie file must be open.
+         * Does not throw.
+         */
+        bool isInteractive() const;
         /*! \brief
          * Returns a file handle for interfacing with C functions.
          *