Warn for type mismatch for gmx printf like functions 2/3
authorRoland Schulz <roland.schulz@intel.com>
Fri, 13 Jul 2018 04:18:58 +0000 (21:18 -0700)
committerRoland Schulz <roland.schulz@intel.com>
Wed, 12 Sep 2018 22:01:39 +0000 (15:01 -0700)
Most manual fixes. Inclduing those with missing/extra
arguments.

Related #2570

Change-Id: I2046d323dad6420e61f47918b12aa538f4a6d4a0

19 files changed:
src/gromacs/awh/biasstate.cpp
src/gromacs/awh/read-params.cpp
src/gromacs/fileio/enxio.cpp
src/gromacs/gmxana/dlist.cpp
src/gromacs/gmxana/gmx_bar.cpp
src/gromacs/gmxana/gmx_densorder.cpp
src/gromacs/gmxana/gmx_do_dssp.cpp
src/gromacs/gmxana/gmx_nmr.cpp
src/gromacs/gmxana/gmx_trjcat.cpp
src/gromacs/gmxana/gmx_trjorder.cpp
src/gromacs/gmxana/gmx_wham.cpp
src/gromacs/gmxpreprocess/pdb2gmx.cpp
src/gromacs/hardware/detecthardware.cpp
src/gromacs/imd/imd.cpp
src/gromacs/mdlib/mdebin.cpp
src/gromacs/mdlib/sim_util.cpp
src/gromacs/mdrun/runner.cpp
src/gromacs/pulling/pull_rotation.cpp
src/gromacs/utility/smalloc.cpp

index 4bf7bf3263122fa654f4613b9013b6750a17a08b..c45f8f19b5e01827195cd768beb1602bbb2c65cd 100644 (file)
@@ -1523,9 +1523,9 @@ static void readUserPmfAndTargetDistribution(const std::vector<DimParams> &dimPa
     if (numColumns < numColumnsMin)
     {
         std::string mesg =
-            gmx::formatString("The number of columns in %s (%s) should be at least %d."
+            gmx::formatString("The number of columns in %s should be at least %d."
                               "\n\n%s",
-                              filename.c_str(), correctFormatMessage.c_str());
+                              filename.c_str(), numColumnsMin, correctFormatMessage.c_str());
         GMX_THROW(InvalidInputError(mesg));
     }
 
index 63e4fc0696d0e10227df365d6051ea3f3b4d91ab..353ec243505ccf3763bb49db1fd693ad172921ae 100644 (file)
@@ -271,7 +271,7 @@ static void read_bias_params(std::vector<t_inpfile> *inp, AwhBiasParams *awhBias
     awhBiasParams->errorInitial = get_ereal(inp, opt, 10, wi);
     if (awhBiasParams->errorInitial <= 0)
     {
-        gmx_fatal(FARGS, "%s (%d) needs to be > 0.", opt);
+        gmx_fatal(FARGS, "%s needs to be > 0.", opt);
     }
 
     if (bComment)
@@ -703,7 +703,7 @@ static void checkInputConsistencyInterval(const AwhParams *awhParams, warninp_t
 
             if ((period == 0) && (origin > end))
             {
-                gmx_fatal(FARGS, "For the non-periodic pull coordinates awh%d-dim%d-start cannot be larger than awh%f-dim%d-end",
+                gmx_fatal(FARGS, "For the non-periodic pull coordinates awh%d-dim%d-start (%f) cannot be larger than awh%d-dim%d-end (%f)",
                           k + 1, d + 1, origin, k + 1, d + 1, end);
             }
 
index 78aaf2ae18cc5215e14de805e68c07a6d74e0094..dc0a177c9ff4609568915a5d4946297e06b0a8c6 100644 (file)
@@ -505,7 +505,7 @@ static gmx_bool do_eheader(ener_file_t ef, int *file_version, t_enxframe *fr,
         }
         if (*bOK && *file_version > enx_version)
         {
-            gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program", gmx_fio_getname(ef->fio), file_version, enx_version);
+            gmx_fatal(FARGS, "reading tpx file (%s) version %d with version %d program", gmx_fio_getname(ef->fio), *file_version, enx_version);
         }
         if (!gmx_fio_do_double(ef->fio, fr->t))
         {
index 0b09a2842bf43bbcf82266fdf6004aab70c85ede..6661f757f26ba31a482cf54cc0572de4ee27b06e 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013,2014,2015,2016, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014,2015,2016,2018, 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.
@@ -229,7 +229,7 @@ t_dlist *mk_dlist(FILE *log,
             {
                 gmx_fatal(FARGS, "Unknown residue %s when searching for residue type.\n"
                           "Maybe you need to add a custom residue in residuetypes.dat.",
-                          thisres, __FILE__, __LINE__);
+                          thisres);
             }
 
             sprintf(dl[nl].name, "%s%d", thisres, ires+r0);
index 76bd921c055289747464bfc57b94518f2a361a67..446d2655225aa68d3cc40052f1036c4b8d07e004 100644 (file)
@@ -3316,7 +3316,7 @@ static void read_barsim_edr(const char *fn, real *temp, sim_data_t *sd)
             if (!lambda_vec_same(&start_lambda, native_lambda) )
             {
                 gmx_fatal(FARGS, "Native lambda not constant in file %s: started at %f, and becomes %f at time %f",
-                          fn, native_lambda->val, start_lambda.val, start_time);
+                          fn, native_lambda->val[0], start_lambda.val[0], start_time);
             }
             /* check the number of samples against the previous number */
             if ( ((nblocks_raw+nblocks_hist) != nsamples) || (nlam != 1 ) )
index 94997770cc7ee988d06940d30c05f295bc324548..cdc9d9b1dc3b9b095662d16fdd74651adc8b6794 100644 (file)
@@ -775,7 +775,7 @@ int gmx_densorder(int argc, char *argv[])
         gmx::ArrayRef<const std::string> graphFiles = opt2fns("-og", NFILE, fnm);
         if (graphFiles.size() != 2)
         {
-            gmx_fatal(FARGS, "No or not correct number (2) of output-files: %zu", graphFiles.size());
+            gmx_fatal(FARGS, "No or not correct number (2) of output-files: %td", graphFiles.size());
         }
         writesurftoxpms(surf1, surf2, tblock, xslices, yslices, zslices, binw, binwz, graphFiles, zslices);
     }
@@ -790,7 +790,7 @@ int gmx_densorder(int argc, char *argv[])
         gmx::ArrayRef<const std::string> rawFiles = opt2fns("-or", NFILE, fnm);
         if (rawFiles.size() != 2)
         {
-            gmx_fatal(FARGS, "No or not correct number (2) of output-files: %zu", rawFiles.size());
+            gmx_fatal(FARGS, "No or not correct number (2) of output-files: %td", rawFiles.size());
         }
         writeraw(surf1, surf2, tblock, xslices, yslices, rawFiles, oenv);
     }
@@ -802,7 +802,7 @@ int gmx_densorder(int argc, char *argv[])
         gmx::ArrayRef<const std::string> spectra = opt2fns("-Spect", NFILE, fnm);
         if (spectra.size() != 2)
         {
-            gmx_fatal(FARGS, "No or not correct number (2) of output-file-series: %zu",
+            gmx_fatal(FARGS, "No or not correct number (2) of output-file-series: %td",
                       spectra.size());
         }
         powerspectavg_intf(surf1, surf2, tblock, xslices, yslices, spectra);
index cdf208521243115a37a2385aed0d82e0453d6c6b..cc8c87c6925a13e911f7761c991c4620aa3b89c3 100644 (file)
@@ -665,7 +665,7 @@ int gmx_do_dssp(int argc, char *argv[])
 
         if (0 != system(dssp))
         {
-            gmx_fatal(FARGS, "Failed to execute command: %s\n",
+            gmx_fatal(FARGS, "Failed to execute command: %s\n"
                       "Try specifying your dssp version with the -ver option.", dssp);
         }
 
index 4bf48a4eef90132b02470b2fe0820c1eac454d82..b8f15180e4bf6415862817ff1b8209b509e2bc4a 100644 (file)
@@ -786,7 +786,7 @@ int gmx_nmr(int argc, char *argv[])
 
                         if (blk->sub[0].nr != nor)
                         {
-                            gmx_fatal(FARGS, "Number of orientation restraints in energy file (%d) does not match with the topology (%d)", blk->sub[0].nr);
+                            gmx_fatal(FARGS, "Number of orientation restraints in energy file (%d) does not match with the topology (%d)", blk->sub[0].nr, nor);
                         }
                         if (bORA || bODA)
                         {
index 8521164281f40ca5e2806ac27bba790b289e484f..9d5ed554e24523ca991da40de3985aba85384aad 100644 (file)
@@ -552,7 +552,7 @@ int gmx_trjcat(int argc, char *argv[])
 
     if (bDeMux && static_cast<int>(inFiles.size()) != nset)
     {
-        gmx_fatal(FARGS, "You have specified %zu files and %d entries in the demux table", inFiles.size(), nset);
+        gmx_fatal(FARGS, "You have specified %td files and %d entries in the demux table", inFiles.size(), nset);
     }
 
     ftpin = fn2ftp(inFiles[0].c_str());
@@ -581,7 +581,7 @@ int gmx_trjcat(int argc, char *argv[])
     }
     else if (bDeMux && static_cast<int>(outFiles.size()) != nset && outFiles.size() != 1)
     {
-        gmx_fatal(FARGS, "Number of output files should be 1 or %d (#input files), not %zu", nset, outFiles.size());
+        gmx_fatal(FARGS, "Number of output files should be 1 or %d (#input files), not %td", nset, outFiles.size());
     }
     if (bDeMux)
     {
index f34e6077c603c141f49ba16aaf566495f936ec2b..698107ca1218e29674996d6c8d7900d4f5128287 100644 (file)
@@ -189,7 +189,7 @@ int gmx_trjorder(int argc, char *argv[])
         {
             if (index[i][j] > natoms)
             {
-                gmx_fatal(FARGS, "An atom number in group %s is larger than the number of atoms in the trajectory");
+                gmx_fatal(FARGS, "An atom number in group %s is larger than the number of atoms in the trajectory", grpname[i]);
             }
         }
     }
index 395d9e96deb51234b9e9de12af56ffc91858d46a..a47e10c4e3b5d6ac8b86a17c8e78db574eb2fc3a 100644 (file)
@@ -417,7 +417,7 @@ static void setup_tab(const char *fn, t_UmbrellaOptions *opt)
     {
         if  (std::abs(y[0][i+1]-y[0][i]-opt->tabDz) > opt->tabDz/1e6)
         {
-            gmx_fatal(FARGS, "z-values in %d are not equally spaced.\n", ny, fn);
+            gmx_fatal(FARGS, "z-values in %s are not equally spaced.\n", fn);
         }
     }
     snew(opt->tabY, nl);
@@ -1871,7 +1871,7 @@ static FILE *open_pdo_pipe(const char *fn, t_UmbrellaOptions *opt, gmx_bool *bPi
             {
                 gmx_fatal(FARGS, "Cannot find executable gunzip in /bin or /usr/bin.\n"
                           "You may want to define the path to gunzip "
-                          "with the environment variable GMX_PATH_GZIP.", gunzip);
+                          "with the environment variable GMX_PATH_GZIP.");
             }
         }
         else
@@ -2126,7 +2126,7 @@ static void read_tpr_header(const char *fn, t_UmbrellaHeader* header, t_Umbrella
             }
             if (header->pcrd[i].k <= 0.0)
             {
-                gmx_fatal(FARGS, "%s: Pull coordinate %d has force constant of of %g in %s.\n"
+                gmx_fatal(FARGS, "%s: Pull coordinate %d has force constant of of %g.\n"
                           "That doesn't seem to be an Umbrella tpr.\n",
                           fn, i+1, header->pcrd[i].k);
             }
index d3310366284802004a65cd2dd40a7e9891f63d51..87ec06ae800bb2f2806536b2ed1bd20a2491c377 100644 (file)
@@ -249,7 +249,7 @@ static void read_rtprename(const char *fname, FILE *fp,
         {
             if (nc != 2 && nc != 5)
             {
-                gmx_fatal(FARGS, "Residue renaming database '%s' has %d columns instead of %d, %d or %d", fname, ncol, 2, 5);
+                gmx_fatal(FARGS, "Residue renaming database '%s' has %d columns instead of %d or %d", fname, ncol, 2, 5);
             }
             ncol = nc;
         }
@@ -1586,8 +1586,7 @@ int gmx_pdb2gmx(int argc, char *argv[])
             bVsiteAromatics = TRUE;
             break;
         default:
-            gmx_fatal(FARGS, "DEATH HORROR in $s (%s): vsitestr[0]='%d'",
-                      __FILE__, __LINE__, vsitestr[0]);
+            gmx_fatal(FARGS, "Internal inconsistency: vsitestr[0]='%s'", vsitestr[0]);
     } /* end switch */
 
     /* Open the symbol table */
index 81a8f98d73c82c3981545436fed474fe7bacd932..aea9a1a3e6fef03dc48374eac44285121de5b85a 100644 (file)
@@ -144,7 +144,7 @@ static void gmx_detect_gpus(const gmx::MDLogger            &mdlog,
         {
             GMX_LOG(mdlog.info).asParagraph().appendTextFormatted(
                     "NOTE: Detection of GPUs failed. The API reported:\n"
-                    "      %s\n",
+                    "      %s\n"
                     "      GROMACS cannot run tasks on a GPU.",
                     errorMessage.c_str());
         }
index 8c2af328dd605f4ce057bb07a70c7ac3b42fa861..13d54d66bfce908e208d011837a5283673c775af 100644 (file)
@@ -534,7 +534,7 @@ static void imd_prepare_master_socket(t_gmx_IMD_setup *IMDsetup)
 
     if (imdsock_getport(IMDsetup->socket, &IMDsetup->port))
     {
-        gmx_fatal(FARGS, "%s Could not determine port number.\n", IMDstr, ret);
+        gmx_fatal(FARGS, "%s Could not determine port number.\n", IMDstr);
     }
 
     fprintf(stderr, "%s Listening for IMD connection on port %d.\n", IMDstr, IMDsetup->port);
index ddd39c3925e76981e03f9aef2079a7d7e0746d17..b8efb73d5862cb93c64511969a394a5bbd849881 100644 (file)
@@ -1567,7 +1567,7 @@ void restore_energyhistory_from_state(t_mdebin              * mdebin,
     if ((enerhist->nsum     > 0 && nener != enerhist->ener_sum.size()) ||
         (enerhist->nsum_sim > 0 && nener != enerhist->ener_sum_sim.size()))
     {
-        gmx_fatal(FARGS, "Mismatch between number of energies in run input (%d) and checkpoint file (%zu or %zu).",
+        gmx_fatal(FARGS, "Mismatch between number of energies in run input (%u) and checkpoint file (%zu or %zu).",
                   nener, enerhist->ener_sum.size(), enerhist->ener_sum_sim.size());
     }
 
index 03913bf416d274118882f00927053efb99fa811e..042bb6ceabf56ab6b12bacc5e634104976daa3b0 100644 (file)
@@ -333,7 +333,7 @@ static void print_large_forces(FILE            *fp,
                                const rvec      *f)
 {
     real           force2Tolerance = gmx::square(forceTolerance);
-    std::uintmax_t numNonFinite    = 0;
+    gmx::index     numNonFinite    = 0;
     for (int i = 0; i < md->homenr; i++)
     {
         real force2    = norm2(f[i]);
@@ -355,7 +355,7 @@ static void print_large_forces(FILE            *fp,
          * the printing on other ranks. But we can only avoid that with
          * an expensive MPI barrier that we would need at each step.
          */
-        gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %ju atoms", step, numNonFinite);
+        gmx_fatal(FARGS, "At step %" PRId64 " detected non-finite forces on %td atoms", step, numNonFinite);
     }
 }
 
index 82c96097f4250bb3071738d61a90d3822d985be7..7337df0886251b4c24d86a565cd26ba2edae96a4 100644 (file)
@@ -679,16 +679,14 @@ int Mdrunner::mdrunner()
         gmx_fatal(FARGS,
                   "The -dd or -npme option request a parallel simulation, "
 #if !GMX_MPI
-                  "but %s was compiled without threads or MPI enabled"
+                  "but %s was compiled without threads or MPI enabled", output_env_get_program_display_name(oenv));
 #else
 #if GMX_THREAD_MPI
-                  "but the number of MPI-threads (option -ntmpi) is not set or is 1"
+                  "but the number of MPI-threads (option -ntmpi) is not set or is 1");
 #else
-                  "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec"
+                  "but %s was not started through mpirun/mpiexec or only one rank was requested through mpirun/mpiexec", output_env_get_program_display_name(oenv));
 #endif
 #endif
-                  , output_env_get_program_display_name(oenv)
-                  );
     }
 
     if (doRerun &&
index 76badaeda53aed597e6b1f703bfdee899870d759..f9fd3abeeffe00d555a6d2a17d45023930784fea 100644 (file)
@@ -644,7 +644,7 @@ static double calc_beta_max(real min_gaussian, real slab_dist)
     }
     if (min_gaussian <= 0)
     {
-        gmx_fatal(FARGS, "Cutoff value for Gaussian must be > 0. (You requested %f)");
+        gmx_fatal(FARGS, "Cutoff value for Gaussian must be > 0. (You requested %f)", min_gaussian);
     }
 
     /* Define the sigma value */
index f329a2c8ec60e0bf459d16ae43dc5ba43face234..bf5015c006f71ad685e9fb65933153de8b71e769 100644 (file)
@@ -165,9 +165,9 @@ void *save_realloc(const char *name, const char *file, int line, void *ptr,
         if (p == nullptr)
         {
             gmx_fatal(errno, __FILE__, __LINE__,
-                      "Not enough memory. Failed to realloc %" PRId64 " bytes for %s, %s=%x\n"
+                      "Not enough memory. Failed to realloc %zu bytes for %s, %s=%p\n"
                       "(called from file %s, line %d)",
-                      static_cast<int64_t>(size), name, name, ptr, file, line);
+                      size, name, name, ptr, file, line);
         }
     }
     return p;