Remove PCA_QUIET
authorTeemu Murtola <teemu.murtola@gmail.com>
Wed, 14 May 2014 18:37:35 +0000 (21:37 +0300)
committerGerrit Code Review <gerrit@gerrit.gromacs.org>
Tue, 8 Jul 2014 18:45:44 +0000 (20:45 +0200)
Instead of exposing this to all callers of parse_common_args(), make the
command-line runner only execute the module on the master node in the
case help is being printed.  It should be quite rare that anyone runs -h
with multiple nodes (which is the only case this nowadays affected), so
this doesn't need an extra flag that callers need to be aware of.

Change-Id: I0a0b76d358bd7db283e5244c49acbd4a40d7af6f

src/contrib/timefft.c
src/gromacs/commandline/cmdlinemodule.h
src/gromacs/commandline/cmdlinemodulemanager.cpp
src/gromacs/commandline/pargs.cpp
src/gromacs/commandline/pargs.h
src/gromacs/gmxana/gmx_pme_error.cpp
src/programs/mdrun/mdrun.cpp

index 55fd5fb3940023f82ad40630b7b5d51350eb80db..88811810af7d9208d6de84d642a0eebb8c0f16dc 100644 (file)
@@ -87,8 +87,7 @@ int main(int argc,char *argv[])
   cr = init_par(&argc,&argv);
   if (MASTER(cr))
     CopyRight(stdout,argv[0]);
-  parse_common_args(&argc,argv,
-                   PCA_CAN_SET_DEFFNM | (MASTER(cr) ? 0 : PCA_QUIET),
+  parse_common_args(&argc,argv, PCA_CAN_SET_DEFFNM,
                    NFILE,fnm,asize(pa),pa,0,NULL,0,NULL);
   gmx_log_open(ftp2fn(efLOG,NFILE,fnm),cr,1,0,&fplog);
 
index a482046f82b6c4c93e164ffabf19aece97276e6a..cfe908addb92125ac95c17b53b51b1e29a231730 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.
@@ -85,6 +85,9 @@ class CommandLineModuleInterface
          * \param[in] context  Context object for writing the help.
          * \throws    std::bad_alloc if out of memory.
          * \throws    FileIOError on any I/O error.
+         *
+         * Note that for MPI-enabled builds, this is called only on the master
+         * rank.
          */
         virtual void writeHelp(const CommandLineHelpContext &context) const = 0;
 };
index b88a3d3312dda96a92ef349c54b05e56a9a0018d..826d96fd9d9c24d86f7cc12764588dc3ccaef022 100644 (file)
@@ -515,7 +515,9 @@ void CommandLineModuleManager::addHelpTopic(HelpTopicPointer topic)
 int CommandLineModuleManager::run(int argc, char *argv[])
 {
     CommandLineModuleInterface    *module;
-    const bool                     bMaster = (gmx_node_rank() == 0);
+    // TODO: With thread-MPI, gmx_node_rank() returns random stuff here, so we
+    // need the first check.
+    const bool                     bMaster = (!gmx_mpi_initialized() || gmx_node_rank() == 0);
     bool                           bQuiet  = impl_->bQuiet_ || !bMaster;
     CommandLineCommonOptionsHolder optionsHolder;
     try
@@ -558,7 +560,11 @@ int CommandLineModuleManager::run(int argc, char *argv[])
         fprintf(stderr, "Will write debug log file: %s\n", filename.c_str());
         gmx_init_debug(optionsHolder.debugLevel(), filename.c_str());
     }
-    int rc = module->run(argc, argv);
+    int rc = 0;
+    if (!(module == impl_->helpModule_ && !bMaster))
+    {
+        rc = module->run(argc, argv);
+    }
     if (!bQuiet)
     {
         gmx_thanx(stderr);
index 98377ca4b47a4faf4ff3dd89489841050e345db7..2d62788683cfad5339c64ddb801f387ab909ec53 100644 (file)
@@ -63,6 +63,7 @@
 #include "gromacs/options/options.h"
 #include "gromacs/options/timeunitmanager.h"
 #include "gromacs/utility/arrayref.h"
+#include "gromacs/utility/basenetwork.h"
 #include "gromacs/utility/common.h"
 #include "gromacs/utility/cstringutil.h"
 #include "gromacs/utility/exceptions.h"
@@ -568,14 +569,16 @@ gmx_bool parse_common_args(int *argc, char *argv[], unsigned long Flags,
             gmx::GlobalCommandLineHelpContext::get();
         if (context != NULL)
         {
-            if (!(FF(PCA_QUIET)))
-            {
-                gmx::CommandLineHelpWriter(options)
-                    .setShowDescriptions(true)
-                    .setTimeUnitString(timeUnitManager.timeUnitAsString())
-                    .setKnownIssues(gmx::ConstArrayRef<const char *>(bugs, nbugs))
-                    .writeHelp(*context);
-            }
+            // TODO: The first check should not be necessary, but with
+            // thread-MPI it is...
+            GMX_RELEASE_ASSERT(!gmx_mpi_initialized() || gmx_node_rank() == 0,
+                               "Help output should be handled higher up and "
+                               "only get called only on the master rank");
+            gmx::CommandLineHelpWriter(options)
+                .setShowDescriptions(true)
+                .setTimeUnitString(timeUnitManager.timeUnitAsString())
+                .setKnownIssues(gmx::ConstArrayRef<const char *>(bugs, nbugs))
+                .writeHelp(*context);
             return FALSE;
         }
 
index ac8fcf7e03df8de584c61140f71723b2ccd4a4e6..ee912a1399abdc65373bbe48c40fae13f32203de 100644 (file)
@@ -230,8 +230,6 @@ gmx_bool opt2parg_bSet(const char *option, int nparg, t_pargs pa[]);
 #define PCA_CAN_SET_DEFFNM (1<<10)
 /** Do not raise a fatal error when invalid options are encountered. */
 #define PCA_NOEXIT_ON_ARGS (1<<11)
-/** Do not print help output (used for non-master nodes). */
-#define PCA_QUIET          (1<<12)
 /** Default to low priority. */
 #define PCA_BE_NICE        (1<<13)
 /** Is this node not reading: for parallel all nodes but the master */
@@ -247,6 +245,17 @@ gmx_bool opt2parg_bSet(const char *option, int nparg, t_pargs pa[]);
  * by \p Flags.
  *
  * Recognized arguments are removed from the list.
+ *
+ * For full functionality, this function needs to be used within a function
+ * that is passed to gmx_run_cmain().  It should be called as the first thing in
+ * that function.  Initialization code can be executed before it, but you need
+ * to be aware that if the program is executed with -h and MPI, the code before
+ * parse_common_args() only executes on the master node.
+ *
+ * If the return value is `FALSE`, the program should return immediately (this
+ * is necessary for -h and a few other cases).
+ *
+ * \see gmx_run_cmain().
  */
 gmx_bool parse_common_args(int *argc, char *argv[], unsigned long Flags,
                            int nfile, t_filenm fnm[], int npargs, t_pargs *pa,
index 328c400006a8fe112464883b97827aa9e9171e8d..1d454bb84cbfdf26f8313fd87e144a561380a2d5 100644 (file)
@@ -1126,12 +1126,8 @@ int gmx_pme_error(int argc, char *argv[])
 #define NFILE asize(fnm)
 
     cr = init_commrec();
-#ifdef GMX_LIB_MPI
-    MPI_Barrier(MPI_COMM_WORLD);
-#endif
 
     PCA_Flags  = PCA_NOEXIT_ON_ARGS;
-    PCA_Flags |= (MASTER(cr) ? 0 : PCA_QUIET);
 
     if (!parse_common_args(&argc, argv, PCA_Flags,
                            NFILE, fnm, asize(pa), pa, asize(desc), desc,
index 9b0219e702c7e65d5bc9fad2073ebd711466efef..b01033e70c031dba112019dd276d056495dba71b 100644 (file)
@@ -609,7 +609,7 @@ int gmx_mdrun(int argc, char *argv[])
         { "-resethway", FALSE, etBOOL, {&bResetCountersHalfWay},
           "HIDDENReset the cycle counters after half the number of steps or halfway [TT]-maxh[tt]" }
     };
-    unsigned long   Flags, PCA_Flags;
+    unsigned long   Flags;
     ivec            ddxyz;
     int             dd_node_order;
     gmx_bool        bAddPart;
@@ -622,7 +622,7 @@ int gmx_mdrun(int argc, char *argv[])
 
     cr = init_commrec();
 
-    PCA_Flags = (PCA_CAN_SET_DEFFNM | (MASTER(cr) ? 0 : PCA_QUIET));
+    unsigned long PCA_Flags = PCA_CAN_SET_DEFFNM;
     // With -multi or -multidir, the file names are going to get processed
     // further (or the working directory changed), so we can't check for their
     // existence during parsing.  It isn't useful to do any completion based on