Quiet stderr output, particularly for multi-simulations
authorMark Abraham <mark.j.abraham@gmail.com>
Fri, 11 Jan 2013 17:06:33 +0000 (18:06 +0100)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 16 Jan 2013 01:16:01 +0000 (02:16 +0100)
* removed printing of DD info to stderr
* printed multi-simulation information only once in places where
  repetition is clearly redundant (Some repetition remains,
  but only from each simulation master.)
* add option to not print result of multi-simulation check
  if it passed, so that we don't have to print things to
  stderr/stdout just because the .log file is not yet open
* printing of diagnostics about the number of MPI
  processes present when mdrun starts only goes (once) to
  each debug file, and not to stderr
* reduced printing of diagnostics about the number of OpenMP
  threads; now goes to stderr only on SIMMASTER, or
  once to each debug file
* clarified errors and informational messages about selecting
  the number of OpenMP threads

Fixes #1078, refs #1083

Change-Id: If782259bcd62ddd9be325393930080b70c5cfb4e

include/gmx_omp_nthreads.h
include/main.h
src/gmxlib/disre.c
src/gmxlib/gmx_omp_nthreads.c
src/gmxlib/main.c
src/gmxlib/network.c
src/gmxlib/orires.c
src/kernel/mdrun.c
src/kernel/repl_ex.c
src/kernel/runner.c
src/mdlib/domdec.c

index 4f126f45f94a53835fd691b871101df6a49d2e18..9cc36db39096a9099ff22ed231356fef6ea60e8a 100644 (file)
@@ -69,6 +69,7 @@ int gmx_omp_nthreads_get(int mod);
 
 /*! Read the OMP_NUM_THREADS env. var. and check against the value set on the command line. */
 GMX_LIBGMX_EXPORT
-void gmx_omp_nthreads_read_env(int *nthreads_omp);
+void gmx_omp_nthreads_read_env(int *nthreads_omp,
+                               gmx_bool bIsSimMaster);
 
 #endif /* GMX_OMP_NTHREADS */
index 64d64466aaa0239de32221304ca942ec2e8fa641..cb51a88a4d13b5b4257564e3f23e4ef813961e2e 100644 (file)
@@ -67,14 +67,18 @@ void gmx_log_close(FILE *fp);
 
 GMX_LIBGMX_EXPORT
 void check_multi_int(FILE *log,const gmx_multisim_t *ms,
-                           int val,const char *name);
+                     int val,const char *name,
+                     gmx_bool bQuiet);
 GMX_LIBGMX_EXPORT
 void check_multi_large_int(FILE *log,const gmx_multisim_t *ms,
-                           gmx_large_int_t val,const char *name);
+                           gmx_large_int_t val,const char *name,
+                           gmx_bool bQuiet);
 /* Check if val is the same on all processors for a mdrun -multi run
  * The string name is used to print to the log file and in a fatal error
- * if the val's don't match.
+ * if the val's don't match. If bQuiet is true and the check passes,
+ * no output is written.
  */
+
 GMX_LIBGMX_EXPORT
 void init_multisystem(t_commrec *cr, int nsim, char **multidirs,
                       int nfile, const t_filenm fnm[], gmx_bool bParFn);
index 21732cdd3b2d34e068e83d4ab83059549dcdf9fe..2124384c90448b3f2846e70c5df8a228801fa870 100644 (file)
@@ -193,7 +193,8 @@ void init_disres(FILE *fplog,const gmx_mtop_t *mtop,
             fprintf(fplog,"Found GMX_DISRE_ENSEMBLE_SIZE set to %d systems per ensemble\n",dd->nsystems);
         }
         check_multi_int(fplog,cr->ms,dd->nsystems,
-                        "the number of systems per ensemble");
+                        "the number of systems per ensemble",
+                        FALSE);
         /* We use to allow any value of nsystems which was a divisor
          * of ms->nsim. But this required an extra communicator which
          * was stored in t_fcdata. This pulled in mpi.h in nearly all C files.
@@ -230,7 +231,8 @@ void init_disres(FILE *fplog,const gmx_mtop_t *mtop,
         if (cr && cr->ms)
         {
             check_multi_int(fplog,cr->ms,fcd->disres.nres,
-                            "the number of distance restraints");
+                            "the number of distance restraints",
+                            FALSE);
         }
         please_cite(fplog,"Tropp80a");
         please_cite(fplog,"Torda89a");
index f05ce9e895361c79e6936c209e3c77092cb80479..f4f320917ba259922a5fa4b23943c1dfcc209378 100644 (file)
@@ -189,9 +189,12 @@ static int pick_module_nthreads(FILE *fplog, int m,
     return modth.nth[m] = nth;
 }
 
-void gmx_omp_nthreads_read_env(int *nthreads_omp)
+void gmx_omp_nthreads_read_env(int *nthreads_omp,
+                               gmx_bool bIsSimMaster)
 {
     char *env;
+    gmx_bool bCommandLineSetNthreadsOMP = *nthreads_omp > 0;
+    char buffer[STRLEN];
 
     assert(nthreads_omp);
 
@@ -205,17 +208,33 @@ void gmx_omp_nthreads_read_env(int *nthreads_omp)
             gmx_fatal(FARGS,"OMP_NUM_THREADS is invalid: '%s'",env);
         }
 
-        if (*nthreads_omp > 0 && nt_omp != *nthreads_omp)
+        if (bCommandLineSetNthreadsOMP && nt_omp != *nthreads_omp)
         {
-            gmx_fatal(FARGS,"OMP_NUM_THREADS (%d) and the number of threads requested on the command line (%d) have different values",nt_omp,*nthreads_omp);
+            gmx_fatal(FARGS,"Environment variable OMP_NUM_THREADS (%d) and the number of threads requested on the command line (%d) have different values. Either omit one, or set them both to the same value.",nt_omp,*nthreads_omp);
         }
 
-        /* Setting the number of OpenMP threads.
-         * NOTE: with tMPI this function is only called on the master node,
-         * but with MPI on all nodes which means lots of messages on stderr.
-         */
-        fprintf(stderr,"Getting the number of OpenMP threads from OMP_NUM_THREADS: %d\n",nt_omp);
+        /* Setting the number of OpenMP threads. */
         *nthreads_omp = nt_omp;
+
+        /* Output the results */
+        sprintf(buffer,
+                "The number of OpenMP threads was set by environment variable OMP_NUM_THREADS to %d%s\n",
+                nt_omp,
+                bCommandLineSetNthreadsOMP ? " (and the command-line setting agreed with that)" : "");
+        if (bIsSimMaster)
+        {
+            /* This prints once per simulation for multi-simulations,
+             * which might help diagnose issues with inhomogenous
+             * cluster setups. */
+            fputs(buffer, stderr);
+        }
+        if (debug)
+        {
+            /* This prints once per process for real MPI (i.e. once
+             * per debug file), and once per simulation for thread MPI
+             * (because of logic in the calling function). */
+            fputs(buffer, debug);
+        }
     }
 }
 
index 6b8c27b0245ca2af909bbcd92bd558790ca3eeaa..2d583ce8546b5fb6b035358c88ef558d54c69d61 100644 (file)
@@ -127,12 +127,13 @@ static void par_fn(char *base,int ftp,const t_commrec *cr,
 }
 
 void check_multi_int(FILE *log,const gmx_multisim_t *ms,int val,
-                     const char *name)
+                     const char *name,
+                     gmx_bool bQuiet)
 {
   int  *ibuf,p;
   gmx_bool bCompatible;
 
-  if (NULL != log)
+  if (NULL != log && !bQuiet)
       fprintf(log,"Multi-checking %s ... ",name);
   
   if (ms == NULL)
@@ -149,7 +150,7 @@ void check_multi_int(FILE *log,const gmx_multisim_t *ms,int val,
   
   if (bCompatible) 
   {
-      if (NULL != log)
+      if (NULL != log && !bQuiet)
           fprintf(log,"OK\n");
   }
   else 
@@ -167,13 +168,14 @@ void check_multi_int(FILE *log,const gmx_multisim_t *ms,int val,
 }
 
 void check_multi_large_int(FILE *log,const gmx_multisim_t *ms,
-                           gmx_large_int_t val, const char *name)
+                           gmx_large_int_t val, const char *name,
+                           gmx_bool bQuiet)
 {
   gmx_large_int_t  *ibuf;
   int p;
   gmx_bool bCompatible;
 
-  if (NULL != log)
+  if (NULL != log && !bQuiet)
       fprintf(log,"Multi-checking %s ... ",name);
   
   if (ms == NULL)
@@ -190,7 +192,7 @@ void check_multi_large_int(FILE *log,const gmx_multisim_t *ms,
   
   if (bCompatible) 
   {
-      if (NULL != log)
+      if (NULL != log && !bQuiet)
           fprintf(log,"OK\n");
   }
   else 
@@ -350,7 +352,10 @@ static void comm_args(const t_commrec *cr,int *argc,char ***argv)
   
   if (!MASTER(cr))
     snew(*argv,*argc+1);
-  fprintf(stderr,"NODEID=%d argc=%d\n",cr->nodeid,*argc);
+  if (debug)
+  {
+      fprintf(debug,"NODEID=%d argc=%d\n",cr->nodeid,*argc);
+  }
   for(i=0; (i<*argc); i++) {
     if (MASTER(cr))
       len = strlen((*argv)[i])+1;
index 71ceb9862264dd888555a2e62b6cd6c7b885f0b1..c792fdda391fcc8e8549d55beb38695715de8d7c 100644 (file)
@@ -225,8 +225,11 @@ int gmx_setup(int *argc,char **argv,int *nnodes)
 #endif
  
 #ifdef GMX_LIB_MPI 
-  fprintf(stderr,"NNODES=%d, MYRANK=%d, HOSTNAME=%s\n",
-         mpi_num_nodes,mpi_my_rank,mpi_hostname);
+  if (debug)
+  {
+      fprintf(debug,"NNODES=%d, MYRANK=%d, HOSTNAME=%s\n",
+              mpi_num_nodes,mpi_my_rank,mpi_hostname);
+  }
 #endif
   
   *nnodes=mpi_num_nodes;
index 6b011a1dd7c8a836864d268c6f69f794d21a1300..f1c22f76237152c124fade6a1c2d7982cce0e621 100644 (file)
@@ -226,10 +226,12 @@ void init_orires(FILE *fplog,const gmx_mtop_t *mtop,
         fprintf(fplog,"  the orientation restraints are ensemble averaged over %d systems\n",ms->nsim);
         
         check_multi_int(fplog,ms,od->nr,
-                        "the number of orientation restraints");
+                        "the number of orientation restraints",
+                        FALSE);
         check_multi_int(fplog,ms,od->nref,
-                        "the number of fit atoms for orientation restraining");
-        check_multi_int(fplog,ms,ir->nsteps,"nsteps");
+                        "the number of fit atoms for orientation restraining",
+                        FALSE);
+        check_multi_int(fplog,ms,ir->nsteps,"nsteps",FALSE);
         /* Copy the reference coordinates from the master to the other nodes */
         gmx_sum_sim(DIM*od->nref,od->xref[0],ms);
     }
index 9d9ab37d492556d39631810ed4aa2f5eb429074a..5d80e8422ba3c27f1c27e40fb201f9fe47b84f84 100644 (file)
@@ -619,7 +619,7 @@ int cmain(int argc,char *argv[])
   ivec     ddxyz;
   int      dd_node_order;
   gmx_bool     bAddPart;
-  FILE     *fplog,*fptest;
+  FILE     *fplog,*fpmulti;
   int      sim_part,sim_part_fn;
   const char *part_suffix=".part";
   char     suffix[STRLEN];
@@ -703,7 +703,7 @@ int cmain(int argc,char *argv[])
                                                 &sim_part_fn,NULL,cr,
                                                 bAppendFiles,NFILE,fnm,
                                                 part_suffix,&bAddPart);
-      if (sim_part_fn==0 && MASTER(cr))
+      if (sim_part_fn==0 && MULTIMASTER(cr))
       {
           fprintf(stdout,"No previous checkpoint file present, assuming this is a new run.\n");
       }
@@ -714,7 +714,17 @@ int cmain(int argc,char *argv[])
 
       if (MULTISIM(cr) && MASTER(cr))
       {
-          check_multi_int(stdout,cr->ms,sim_part,"simulation part");
+          if (MULTIMASTER(cr))
+          {
+              /* Log file is not yet available, so if there's a
+               * problem we can only write to stderr. */
+              fpmulti = stderr;
+          }
+          else
+          {
+              fpmulti = NULL;
+          }
+          check_multi_int(fpmulti,cr->ms,sim_part,"simulation part",TRUE);
       }
   } 
   else
@@ -734,7 +744,7 @@ int cmain(int argc,char *argv[])
       sprintf(suffix,"%s%04d",part_suffix,sim_part_fn);
 
       add_suffix_to_output_names(fnm,NFILE,suffix);
-      if (MASTER(cr))
+      if (MULTIMASTER(cr))
       {
           fprintf(stdout,"Checkpoint file is from part %d, new output files will be suffixed '%s'.\n",sim_part-1,suffix);
       }
index 2b3a7c2b93ff81c31d2b1b8c302bef395a5cf062..4b1964aa32664708691f34469dd9335e020dc958 100644 (file)
@@ -162,17 +162,17 @@ gmx_repl_ex_t init_replica_exchange(FILE *fplog,
 
     fprintf(fplog,"Repl  There are %d replicas:\n",re->nrepl);
 
-    check_multi_int(fplog,ms,state->natoms,"the number of atoms");
-    check_multi_int(fplog,ms,ir->eI,"the integrator");
-    check_multi_large_int(fplog,ms,ir->init_step+ir->nsteps,"init_step+nsteps");
+    check_multi_int(fplog,ms,state->natoms,"the number of atoms",FALSE);
+    check_multi_int(fplog,ms,ir->eI,"the integrator",FALSE);
+    check_multi_large_int(fplog,ms,ir->init_step+ir->nsteps,"init_step+nsteps",FALSE);
     check_multi_large_int(fplog,ms,(ir->init_step+nst-1)/nst,
-                          "first exchange step: init_step/-replex");
-    check_multi_int(fplog,ms,ir->etc,"the temperature coupling");
+                          "first exchange step: init_step/-replex",FALSE);
+    check_multi_int(fplog,ms,ir->etc,"the temperature coupling",FALSE);
     check_multi_int(fplog,ms,ir->opts.ngtc,
-                    "the number of temperature coupling groups");
-    check_multi_int(fplog,ms,ir->epc,"the pressure coupling");
-    check_multi_int(fplog,ms,ir->efep,"free energy");
-    check_multi_int(fplog,ms,ir->fepvals->n_lambda,"number of lambda states");
+                    "the number of temperature coupling groups",FALSE);
+    check_multi_int(fplog,ms,ir->epc,"the pressure coupling",FALSE);
+    check_multi_int(fplog,ms,ir->efep,"free energy",FALSE);
+    check_multi_int(fplog,ms,ir->fepvals->n_lambda,"number of lambda states",FALSE);
 
     re->temp = ir->opts.ref_t[0];
     for(i=1; (i<ir->opts.ngtc); i++)
index 773b9160e701ae559a0ce93a725498854ba288be..d6cfe63e22f94524c067289d377ac39bd97ac8a9 100644 (file)
@@ -1079,9 +1079,10 @@ static void set_cpu_affinity(FILE *fplog,
 
 
 static void check_and_update_hw_opt(gmx_hw_opt_t *hw_opt,
-                                    int cutoff_scheme)
+                                    int cutoff_scheme,
+                                    gmx_bool bIsSimMaster)
 {
-    gmx_omp_nthreads_read_env(&hw_opt->nthreads_omp);
+    gmx_omp_nthreads_read_env(&hw_opt->nthreads_omp, bIsSimMaster);
 
 #ifndef GMX_THREAD_MPI
     if (hw_opt->nthreads_tot > 0)
@@ -1342,7 +1343,7 @@ int mdrunner(gmx_hw_opt_t *hw_opt,
     if (SIMMASTER(cr))
 #endif
     {
-        check_and_update_hw_opt(hw_opt,minf.cutoff_scheme);
+        check_and_update_hw_opt(hw_opt,minf.cutoff_scheme,SIMMASTER(cr));
 
 #ifdef GMX_THREAD_MPI
         /* Early check for externally set process affinity. Can't do over all
index a441edf36f0fd28da28662ed926a990e9c67a6e2..09b23e556d336bf5d73902bfd560421e526c1de1 100644 (file)
@@ -5574,11 +5574,6 @@ void setup_dd_grid(FILE *fplog,gmx_domdec_t *dd)
         }
     }
     
-    if (DDMASTER(dd))
-    {
-        fprintf(stderr,"Making %dD domain decomposition %d x %d x %d\n",
-           dd->ndim,dd->nc[XX],dd->nc[YY],dd->nc[ZZ]);
-    }
     if (fplog)
     {
         fprintf(fplog,"\nMaking %dD domain decomposition grid %d x %d x %d, home cell index %d %d %d\n\n",