From: Roland Schulz Date: Sat, 14 Jul 2018 05:57:53 +0000 (-0700) Subject: Warn for type mismatch for gmx printf like functions 3/3 X-Git-Url: http://biod.pnpi.spb.ru/gitweb/?a=commitdiff_plain;h=e5a941a4649d73e4a4d19104d5f4ffa19f5d8fd6;p=alexxy%2Fgromacs.git Warn for type mismatch for gmx printf like functions 3/3 Add attributes for gcc (and compatible) and msvc. And fix remaining issues. Related #2570 Change-Id: I9b46559c2554404b309aa8710eb34ffa128d8cae --- diff --git a/src/gromacs/awh/bias.cpp b/src/gromacs/awh/bias.cpp index bde1a68fbe..437784e24d 100644 --- a/src/gromacs/awh/bias.cpp +++ b/src/gromacs/awh/bias.cpp @@ -47,6 +47,7 @@ #include "bias.h" #include +#include #include #include #include @@ -242,7 +243,7 @@ static void ensureStateAndRunConsistency(const BiasParams ¶ms, int64_t numUpdatesExpected = state.histogramSize().numUpdates(); if (numUpdatesFromSamples != numUpdatesExpected) { - std::string mesg = gmx::formatString("The number of AWH updates in the checkpoint file (%ld) does not match the total number of AWH samples divided by the number of samples per update for %d sharing AWH bias(es) (%ld/%d=%ld)", + std::string mesg = gmx::formatString("The number of AWH updates in the checkpoint file (%" PRId64 ") does not match the total number of AWH samples divided by the number of samples per update for %d sharing AWH bias(es) (%" PRId64 "/%d=%" PRId64 ")", numUpdatesExpected, params.numSharedUpdate, numSamples, @@ -255,7 +256,7 @@ static void ensureStateAndRunConsistency(const BiasParams ¶ms, */ if (numUpdatesFromSamples % state.histogramSize().numUpdates() == 0) { - mesg += gmx::formatString(" Or the run you continued from used %ld sharing simulations, whereas you now specified %d sharing simulations.", + mesg += gmx::formatString(" Or the run you continued from used %" PRId64 " sharing simulations, whereas you now specified %d sharing simulations.", numUpdatesFromSamples/state.histogramSize().numUpdates(), params.numSharedUpdate); } diff --git a/src/gromacs/domdec/domdec.cpp b/src/gromacs/domdec/domdec.cpp index 3e1bd81c8c..d9f565e980 100644 --- a/src/gromacs/domdec/domdec.cpp +++ b/src/gromacs/domdec/domdec.cpp @@ -40,6 +40,7 @@ #include "config.h" #include +#include #include #include #include @@ -3399,7 +3400,7 @@ static DlbState determineInitialDlbState(const gmx::MDLogger &mdlog, case DlbState::onUser: return forceDlbOffOrBail(dlbState, reasonStr + " In load balanced runs binary reproducibility cannot be ensured.", mdlog); default: - gmx_fatal(FARGS, "Death horror: undefined case (%d) for load balancing choice", dlbState); + gmx_fatal(FARGS, "Death horror: undefined case (%d) for load balancing choice", static_cast(dlbState)); } } @@ -6451,7 +6452,7 @@ void dd_partition_system(FILE *fplog, { if (state_local->ddp_count > dd->ddp_count) { - gmx_fatal(FARGS, "Internal inconsistency state_local->ddp_count (%d) > dd->ddp_count (%ld)", state_local->ddp_count, dd->ddp_count); + gmx_fatal(FARGS, "Internal inconsistency state_local->ddp_count (%d) > dd->ddp_count (%" PRId64 ")", state_local->ddp_count, dd->ddp_count); } if (state_local->ddp_count_cg_gl != state_local->ddp_count) diff --git a/src/gromacs/domdec/domdec_setup.cpp b/src/gromacs/domdec/domdec_setup.cpp index 568c948d65..8df4070901 100644 --- a/src/gromacs/domdec/domdec_setup.cpp +++ b/src/gromacs/domdec/domdec_setup.cpp @@ -45,6 +45,7 @@ #include "gmxpre.h" #include +#include #include #include @@ -760,7 +761,7 @@ real dd_choose_grid(const gmx::MDLogger &mdlog, /* Check if the largest divisor is more than nnodes^2/3 */ if (ldiv*ldiv*ldiv > nnodes_div*nnodes_div) { - gmx_fatal(FARGS, "The number of ranks you selected (%ld) contains a large prime factor %ld. In most cases this will lead to bad performance. Choose a number with smaller prime factors or set the decomposition (option -dd) manually.", + gmx_fatal(FARGS, "The number of ranks you selected (%" PRId64 ") contains a large prime factor %" PRId64 ". In most cases this will lead to bad performance. Choose a number with smaller prime factors or set the decomposition (option -dd) manually.", nnodes_div, ldiv); } } diff --git a/src/gromacs/essentialdynamics/edsam.cpp b/src/gromacs/essentialdynamics/edsam.cpp index 4e881bfcdc..db1b4d12de 100644 --- a/src/gromacs/essentialdynamics/edsam.cpp +++ b/src/gromacs/essentialdynamics/edsam.cpp @@ -2363,7 +2363,7 @@ static void crosscheck_edi_file_vs_checkpoint(const gmx_edsam &ed, edsamhistory_ if (static_cast(ed.edpar.size()) != EDstate->nED) { gmx_fatal(FARGS, "The number of essential dynamics / flooding groups is not consistent.\n" - "There are %d ED groups in the .cpt file, but %d in the .edi file!\n" + "There are %d ED groups in the .cpt file, but %zu in the .edi file!\n" "Are you sure this is the correct .edi file?\n", EDstate->nED, ed.edpar.size()); } } diff --git a/src/gromacs/ewald/tests/pmesplinespreadtest.cpp b/src/gromacs/ewald/tests/pmesplinespreadtest.cpp index 43c4b02852..3e593791bf 100644 --- a/src/gromacs/ewald/tests/pmesplinespreadtest.cpp +++ b/src/gromacs/ewald/tests/pmesplinespreadtest.cpp @@ -179,7 +179,7 @@ class PmeSplineAndSpreadTest : public ::testing::TestWithParam(ceil(sqrt(atomCount))); /* 2 is empiric; sqrt(atomCount) assumes all the input charges may spread onto the same cell */ - SCOPED_TRACE(formatString("Testing grid values with tolerance of %ld", ulpToleranceGrid)); + SCOPED_TRACE(formatString("Testing grid values with tolerance of %d", ulpToleranceGrid)); if (!gridValuesSizeAssigned) { previousGridValuesSize = nonZeroGridValues.size(); diff --git a/src/gromacs/fileio/tpxio.cpp b/src/gromacs/fileio/tpxio.cpp index fa44d4ed05..364bddea21 100644 --- a/src/gromacs/fileio/tpxio.cpp +++ b/src/gromacs/fileio/tpxio.cpp @@ -2599,7 +2599,7 @@ static void do_tpxheader(t_fileio *fio, gmx_bool bRead, t_tpxheader *tpx, if ((precision != sizeof(float)) && !bDouble) { gmx_fatal(FARGS, "Unknown precision in file %s: real is %d bytes " - "instead of %lu or %lu", + "instead of %zu or %zu", gmx_fio_getname(fio), precision, sizeof(float), sizeof(double)); } gmx_fio_setprecision(fio, bDouble); diff --git a/src/gromacs/gmxana/gmx_hydorder.cpp b/src/gromacs/gmxana/gmx_hydorder.cpp index 4fe52c384e..403fd5ae95 100644 --- a/src/gromacs/gmxana/gmx_hydorder.cpp +++ b/src/gromacs/gmxana/gmx_hydorder.cpp @@ -665,7 +665,7 @@ int gmx_hydorder(int argc, char *argv[]) gmx::ArrayRef intfn = opt2fns("-o", NFILE, fnm); if (intfn.size() != 2) { - gmx_fatal(FARGS, "No or not correct number (2) of output-files: %ld", intfn.size()); + gmx_fatal(FARGS, "No or not correct number (2) of output-files: %td", intfn.size()); } calc_tetra_order_interface(ndxfnm, tpsfnm, trxfnm, binwidth, nsttblock, &frames, &xslices, &yslices, sg1, sg2, &intfpos, oenv); writesurftoxpms(intfpos, frames, xslices, yslices, binwidth, intfn, nlevels); @@ -675,7 +675,7 @@ int gmx_hydorder(int argc, char *argv[]) gmx::ArrayRef spectra = opt2fns("-Spect", NFILE, fnm); if (spectra.size() != 2) { - gmx_fatal(FARGS, "No or not correct number (2) of output-files: %ld", spectra.size()); + gmx_fatal(FARGS, "No or not correct number (2) of output-files: %td", spectra.size()); } powerspectavg(intfpos, frames, xslices, yslices, spectra); } @@ -685,7 +685,7 @@ int gmx_hydorder(int argc, char *argv[]) gmx::ArrayRef raw = opt2fns("-or", NFILE, fnm); if (raw.size() != 2) { - gmx_fatal(FARGS, "No or not correct number (2) of output-files: %ld", raw.size()); + gmx_fatal(FARGS, "No or not correct number (2) of output-files: %td", raw.size()); } writeraw(intfpos, frames, xslices, yslices, raw); } diff --git a/src/gromacs/gmxana/gmx_make_edi.cpp b/src/gromacs/gmxana/gmx_make_edi.cpp index 77d361396d..1e9608e3ec 100644 --- a/src/gromacs/gmxana/gmx_make_edi.cpp +++ b/src/gromacs/gmxana/gmx_make_edi.cpp @@ -291,12 +291,12 @@ static int sscan_list(int *list[], const char *str, const char *listname) /* format error occured */ case sError: - gmx_fatal(FARGS, "Error in the list of eigenvectors for %s at pos %ld with char %c", listname, pos-startpos, *(pos-1)); + gmx_fatal(FARGS, "Error in the list of eigenvectors for %s at pos %td with char %c", listname, pos-startpos, *(pos-1)); /* logical error occured */ case sZero: - gmx_fatal(FARGS, "Error in the list of eigenvectors for %s at pos %ld: eigenvector 0 is not valid", listname, pos-startpos); + gmx_fatal(FARGS, "Error in the list of eigenvectors for %s at pos %td: eigenvector 0 is not valid", listname, pos-startpos); case sSmaller: - gmx_fatal(FARGS, "Error in the list of eigenvectors for %s at pos %ld: second index %d is not bigger than %d", listname, pos-startpos, end_number, number); + gmx_fatal(FARGS, "Error in the list of eigenvectors for %s at pos %td: second index %d is not bigger than %d", listname, pos-startpos, end_number, number); } ++pos; /* read next character */ } /*scanner has finished */ diff --git a/src/gromacs/gmxlib/network.cpp b/src/gromacs/gmxlib/network.cpp index 41f6826ba2..9276f5a444 100644 --- a/src/gromacs/gmxlib/network.cpp +++ b/src/gromacs/gmxlib/network.cpp @@ -756,7 +756,7 @@ const char *opt2fn_master(const char *opt, int nfile, const t_filenm fnm[], void gmx_fatal_collective(int f_errno, const char *file, int line, MPI_Comm comm, gmx_bool bMaster, - const char *fmt, ...) + gmx_fmtstr const char *fmt, ...) { va_list ap; gmx_bool bFinalize; diff --git a/src/gromacs/gmxlib/network.h b/src/gromacs/gmxlib/network.h index 233d64c469..5919aff437 100644 --- a/src/gromacs/gmxlib/network.h +++ b/src/gromacs/gmxlib/network.h @@ -45,6 +45,7 @@ #include "gromacs/utility/basedefinitions.h" #include "gromacs/utility/gmxmpi.h" +#include "gromacs/utility/stringutil.h" struct gmx_multisim_t; struct t_commrec; @@ -131,7 +132,7 @@ const char *opt2fn_master(const char *opt, int nfile, [[ noreturn ]] void gmx_fatal_collective(int f_errno, const char *file, int line, MPI_Comm comm, gmx_bool bMaster, - const char *fmt, ...); + gmx_fmtstr const char *fmt, ...) gmx_format(printf, 6, 7); /* As gmx_fatal declared in utility/fatalerror.h, * but only the master process prints the error message. * This should only be called one of the following two situations: diff --git a/src/gromacs/gmxpreprocess/readir.cpp b/src/gromacs/gmxpreprocess/readir.cpp index b27b707f1b..355b112447 100644 --- a/src/gromacs/gmxpreprocess/readir.cpp +++ b/src/gromacs/gmxpreprocess/readir.cpp @@ -3493,7 +3493,7 @@ void do_index(const char* mdparin, const char *ndx, auto accelerationGroupNames = gmx::splitString(is->accgrps); if (accelerationGroupNames.size() * DIM != accelerations.size()) { - gmx_fatal(FARGS, "Invalid Acceleration input: %d groups and %d acc. values", + gmx_fatal(FARGS, "Invalid Acceleration input: %zu groups and %zu acc. values", accelerationGroupNames.size(), accelerations.size()); } do_numbering(natoms, groups, accelerationGroupNames, grps, gnames, egcACC, @@ -3508,7 +3508,7 @@ void do_index(const char* mdparin, const char *ndx, auto freezeGroupNames = gmx::splitString(is->freeze); if (freezeDims.size() != DIM * freezeGroupNames.size()) { - gmx_fatal(FARGS, "Invalid Freezing input: %d groups and %d freeze values", + gmx_fatal(FARGS, "Invalid Freezing input: %zu groups and %zu freeze values", freezeGroupNames.size(), freezeDims.size()); } do_numbering(natoms, groups, freezeGroupNames, grps, gnames, egcFREEZE, diff --git a/src/gromacs/gmxpreprocess/toputil.cpp b/src/gromacs/gmxpreprocess/toputil.cpp index 8a66e7bf73..96a183e863 100644 --- a/src/gromacs/gmxpreprocess/toputil.cpp +++ b/src/gromacs/gmxpreprocess/toputil.cpp @@ -66,7 +66,7 @@ void set_p_string(t_param *p, const char *s) } else { - gmx_fatal(FARGS, "Increase MAXSLEN in the grompp code to at least %lu," + gmx_fatal(FARGS, "Increase MAXSLEN in the grompp code to at least %zu," " or shorten your definition of bonds like %s to at most %d", strlen(s)+1, s, MAXSLEN-1); } diff --git a/src/gromacs/gpu_utils/cudautils.cuh b/src/gromacs/gpu_utils/cudautils.cuh index 6eabc3ce14..dac92fa6fe 100644 --- a/src/gromacs/gpu_utils/cudautils.cuh +++ b/src/gromacs/gpu_utils/cudautils.cuh @@ -85,7 +85,7 @@ static inline void ensureNoPendingCudaError(const char *errorMessage) GMX_ASSERT(stat == cudaSuccess, fullMessage.c_str()); // TODO When we evolve a better logging framework, use that // for release-build error reporting. - gmx_warning(fullMessage.c_str()); + gmx_warning("%s", fullMessage.c_str()); } } // namespace diff --git a/src/gromacs/gpu_utils/gpu_utils.cu b/src/gromacs/gpu_utils/gpu_utils.cu index 107d90fb31..15e467308f 100644 --- a/src/gromacs/gpu_utils/gpu_utils.cu +++ b/src/gromacs/gpu_utils/gpu_utils.cu @@ -117,8 +117,8 @@ static void checkCompiledTargetCompatibility(const gmx_device_info_t *devInfo) "The %s binary does not include support for the CUDA architecture " "of the selected GPU (device ID #%d, compute capability %d.%d). " "By default, GROMACS supports all common architectures, so your GPU " - "might be rare, or some architectures were disabled in the build. ", - "Consult the install guide for how to use the GMX_CUDA_TARGET_SM and ", + "might be rare, or some architectures were disabled in the build. " + "Consult the install guide for how to use the GMX_CUDA_TARGET_SM and " "GMX_CUDA_TARGET_COMPUTE CMake variables to add this architecture.", gmx::getProgramContext().displayName(), devInfo->id, devInfo->prop.major, devInfo->prop.minor); @@ -757,8 +757,8 @@ void findGpus(gmx_gpu_info_t *gpu_info) // errors during sanity checks don't propagate. if ((stat = cudaGetLastError()) != cudaSuccess) { - gmx_warning(gmx::formatString("An error occurred while sanity checking device #%d; %s: %s", - devs[i].id, cudaGetErrorName(stat), cudaGetErrorString(stat)).c_str()); + gmx_warning("An error occurred while sanity checking device #%d; %s: %s", + devs[i].id, cudaGetErrorName(stat), cudaGetErrorString(stat)); } } } diff --git a/src/gromacs/gpu_utils/gpu_utils_ocl.cpp b/src/gromacs/gpu_utils/gpu_utils_ocl.cpp index ea8451fb2e..1936cc5937 100644 --- a/src/gromacs/gpu_utils/gpu_utils_ocl.cpp +++ b/src/gromacs/gpu_utils/gpu_utils_ocl.cpp @@ -167,7 +167,7 @@ bool canDetectGpus(std::string *errorMessage) } #endif GMX_RELEASE_ASSERT(status == CL_SUCCESS, - gmx::formatString("An unexpected value was returned from clGetPlatformIDs %u: %s", + gmx::formatString("An unexpected value was returned from clGetPlatformIDs %d: %s", status, ocl_get_error_string(status).c_str()).c_str()); bool foundPlatform = (numPlatforms > 0); if (!foundPlatform && errorMessage != nullptr) @@ -196,7 +196,7 @@ void findGpus(gmx_gpu_info_t *gpu_info) cl_int status = clGetPlatformIDs(0, nullptr, &ocl_platform_count); if (CL_SUCCESS != status) { - GMX_THROW(gmx::InternalError(gmx::formatString("An unexpected value %u was returned from clGetPlatformIDs: ", + GMX_THROW(gmx::InternalError(gmx::formatString("An unexpected value %d was returned from clGetPlatformIDs: ", status) + ocl_get_error_string(status))); } @@ -211,7 +211,7 @@ void findGpus(gmx_gpu_info_t *gpu_info) status = clGetPlatformIDs(ocl_platform_count, ocl_platform_ids, nullptr); if (CL_SUCCESS != status) { - GMX_THROW(gmx::InternalError(gmx::formatString("An unexpected value %u was returned from clGetPlatformIDs: ", + GMX_THROW(gmx::InternalError(gmx::formatString("An unexpected value %d was returned from clGetPlatformIDs: ", status) + ocl_get_error_string(status))); } diff --git a/src/gromacs/hardware/printhardware.cpp b/src/gromacs/hardware/printhardware.cpp index 9748c8e4a0..dc58b6025b 100644 --- a/src/gromacs/hardware/printhardware.cpp +++ b/src/gromacs/hardware/printhardware.cpp @@ -302,7 +302,7 @@ static std::string detected_hardware_string(const gmx_hw_info_t *hwinfo, s += gmx::formatString(" Numa nodes:\n"); for (auto &n : hwTop.machine().numa.nodes) { - s += gmx::formatString(" Node %2d (%" PRIu64 " bytes mem):", n.id, n.memory); + s += gmx::formatString(" Node %2d (%zu bytes mem):", n.id, n.memory); for (auto &l : n.logicalProcessorId) { s += gmx::formatString(" %3d", l); @@ -312,12 +312,12 @@ static std::string detected_hardware_string(const gmx_hw_info_t *hwinfo, s += gmx::formatString(" Latency:\n "); for (std::size_t j = 0; j < hwTop.machine().numa.nodes.size(); j++) { - s += gmx::formatString(" %5lu", j); + s += gmx::formatString(" %5zu", j); } s += gmx::formatString("\n"); for (std::size_t i = 0; i < hwTop.machine().numa.nodes.size(); i++) { - s += gmx::formatString(" %5lu", i); + s += gmx::formatString(" %5zu", i); for (std::size_t j = 0; j < hwTop.machine().numa.nodes.size(); j++) { s += gmx::formatString(" %5.2f", hwTop.machine().numa.relativeLatency[i][j]); @@ -329,7 +329,7 @@ static std::string detected_hardware_string(const gmx_hw_info_t *hwinfo, s += gmx::formatString(" Caches:\n"); for (auto &c : hwTop.machine().caches) { - s += gmx::formatString(" L%d: %" PRIu64 " bytes, linesize %d bytes, assoc. %d, shared %d ways\n", + s += gmx::formatString(" L%d: %zu bytes, linesize %d bytes, assoc. %d, shared %d ways\n", c.level, c.size, c.linesize, c.associativity, c.shared); } } diff --git a/src/gromacs/listed-forces/manage-threading.cpp b/src/gromacs/listed-forces/manage-threading.cpp index e9d7706127..af96415370 100644 --- a/src/gromacs/listed-forces/manage-threading.cpp +++ b/src/gromacs/listed-forces/manage-threading.cpp @@ -48,6 +48,7 @@ #include "config.h" #include +#include #include #include @@ -450,7 +451,7 @@ void setup_bonded_threading(bonded_threading_t *bt, if (gmx_debug_at) { #if BITMASK_SIZE <= 64 //move into bitmask when it is C++ - std::string flags = gmx::formatString("%lx", *mask); + std::string flags = gmx::formatString("%" PRIx64, *mask); #else std::string flags = gmx::formatAndJoin(*mask, "", gmx::StringFormatter("%x")); diff --git a/src/gromacs/mdlib/lincs.cpp b/src/gromacs/mdlib/lincs.cpp index e37b56ae78..85427ce162 100644 --- a/src/gromacs/mdlib/lincs.cpp +++ b/src/gromacs/mdlib/lincs.cpp @@ -1318,7 +1318,7 @@ static void set_lincs_matrix_task(Lincs *li, li_task->ntriangle++; if (li->blnr[i+1] - li->blnr[i] > static_cast(sizeof(li_task->tri_bits[0])*8 - 1)) { - gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %lu allowed for constraints participating in triangles", + gmx_fatal(FARGS, "A constraint is connected to %d constraints, this is more than the %zu allowed for constraints participating in triangles", li->blnr[i+1] - li->blnr[i], sizeof(li_task->tri_bits[0])*8-1); } diff --git a/src/gromacs/mdlib/nbnxn_ocl/nbnxn_ocl.cpp b/src/gromacs/mdlib/nbnxn_ocl/nbnxn_ocl.cpp index 3ef88eaea1..443d656366 100644 --- a/src/gromacs/mdlib/nbnxn_ocl/nbnxn_ocl.cpp +++ b/src/gromacs/mdlib/nbnxn_ocl/nbnxn_ocl.cpp @@ -139,7 +139,7 @@ static inline void validate_global_work_size(const KernelLaunchConfig &config, i { gmx_fatal(FARGS, "Watch out, the input system is too large to simulate!\n" "The number of nonbonded work units (=number of super-clusters) exceeds the" - "device capabilities. Global work size limit exceeded (%d > %d)!", + "device capabilities. Global work size limit exceeded (%zu > %zu)!", global_work_size[i], device_limit); } } diff --git a/src/gromacs/mdlib/nbnxn_search.cpp b/src/gromacs/mdlib/nbnxn_search.cpp index eae0873ae6..fbb2ab1bb9 100644 --- a/src/gromacs/mdlib/nbnxn_search.cpp +++ b/src/gromacs/mdlib/nbnxn_search.cpp @@ -1726,7 +1726,7 @@ static void make_fep_list(const nbnxn_search *nbs, if (ngid*gridj->na_cj > gmx::index(sizeof(gid_cj)*8)) { - gmx_fatal(FARGS, "The Verlet scheme with %dx%d kernels and free-energy only supports up to %lu energy groups", + gmx_fatal(FARGS, "The Verlet scheme with %dx%d kernels and free-energy only supports up to %zu energy groups", gridi->na_c, gridj->na_cj, (sizeof(gid_cj)*8)/gridj->na_cj); } diff --git a/src/gromacs/mdlib/nbnxn_tuning.cpp b/src/gromacs/mdlib/nbnxn_tuning.cpp index a7926bce9f..b34c2544f2 100644 --- a/src/gromacs/mdlib/nbnxn_tuning.cpp +++ b/src/gromacs/mdlib/nbnxn_tuning.cpp @@ -471,7 +471,7 @@ static std::string formatListSetup(const std::string &listName, } listSetup += "updated every "; // Make the shortest int format string that fits nstListForSpacing - std::string nstListFormat = "%" + gmx::formatString("%lu", gmx::formatString("%d", nstListForSpacing).size()) + "d"; + std::string nstListFormat = "%" + gmx::formatString("%zu", gmx::formatString("%d", nstListForSpacing).size()) + "d"; listSetup += gmx::formatString(nstListFormat.c_str(), nstList); listSetup += gmx::formatString(" steps, buffer %.3f nm, rlist %.3f nm\n", rList - interactionCutoff, rList); diff --git a/src/gromacs/mdlib/shellfc.cpp b/src/gromacs/mdlib/shellfc.cpp index e0b3fefb67..51d40ab012 100644 --- a/src/gromacs/mdlib/shellfc.cpp +++ b/src/gromacs/mdlib/shellfc.cpp @@ -501,7 +501,7 @@ gmx_shellfc_t *init_shell_flexcon(FILE *fplog, case F_ANHARM_POL: if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10)) { - gmx_fatal(FARGS, "polarize can not be used with qA(%e) != qB(%e) for atom %d of molecule block %lu", qS, atom[aS].qB, aS+1, mb+1); + gmx_fatal(FARGS, "polarize can not be used with qA(%e) != qB(%e) for atom %d of molecule block %zu", qS, atom[aS].qB, aS+1, mb+1); } shell[nsi].k += gmx::square(qS)*ONE_4PI_EPS0/ ffparams->iparams[type].polarize.alpha; @@ -509,7 +509,7 @@ gmx_shellfc_t *init_shell_flexcon(FILE *fplog, case F_WATER_POL: if (!gmx_within_tol(qS, atom[aS].qB, GMX_REAL_EPS*10)) { - gmx_fatal(FARGS, "water_pol can not be used with qA(%e) != qB(%e) for atom %d of molecule block %lu", qS, atom[aS].qB, aS+1, mb+1); + gmx_fatal(FARGS, "water_pol can not be used with qA(%e) != qB(%e) for atom %d of molecule block %zu", qS, atom[aS].qB, aS+1, mb+1); } alpha = (ffparams->iparams[type].wpol.al_x+ ffparams->iparams[type].wpol.al_y+ diff --git a/src/gromacs/mdrun/md.cpp b/src/gromacs/mdrun/md.cpp index 05daa4b60a..68d8b0bfc4 100644 --- a/src/gromacs/mdrun/md.cpp +++ b/src/gromacs/mdrun/md.cpp @@ -45,6 +45,7 @@ #include "config.h" +#include #include #include #include @@ -772,11 +773,11 @@ void gmx::Integrator::do_md() { if (!rerun_fr.bBox) { - gmx_fatal(FARGS, "Rerun trajectory frame step %ld time %f does not contain a box, while pbc is used", rerun_fr.step, rerun_fr.time); + gmx_fatal(FARGS, "Rerun trajectory frame step %" PRId64 " time %f does not contain a box, while pbc is used", rerun_fr.step, rerun_fr.time); } if (max_cutoff2(ir->ePBC, rerun_fr.box) < gmx::square(fr->rlist)) { - gmx_fatal(FARGS, "Rerun trajectory frame step %ld time %f has too small box dimensions", rerun_fr.step, rerun_fr.time); + gmx_fatal(FARGS, "Rerun trajectory frame step %" PRId64 " time %f has too small box dimensions", rerun_fr.step, rerun_fr.time); } } } diff --git a/src/gromacs/mdrun/runner.cpp b/src/gromacs/mdrun/runner.cpp index 725857f909..79af995e5b 100644 --- a/src/gromacs/mdrun/runner.cpp +++ b/src/gromacs/mdrun/runner.cpp @@ -48,6 +48,7 @@ #include "config.h" #include +#include #include #include #include @@ -303,7 +304,7 @@ static void override_nsteps_cmdline(const gmx::MDLogger &mdlog, } else if (nsteps_cmdline < -2) { - gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %ld", + gmx_fatal(FARGS, "Invalid nsteps value passed on the command line: %" PRId64, nsteps_cmdline); } /* Do nothing if nsteps_cmdline == -2 */ diff --git a/src/gromacs/pulling/pull.cpp b/src/gromacs/pulling/pull.cpp index 792026f3ed..8c8fa5527e 100644 --- a/src/gromacs/pulling/pull.cpp +++ b/src/gromacs/pulling/pull.cpp @@ -41,6 +41,7 @@ #include "config.h" #include +#include #include #include #include @@ -1112,7 +1113,7 @@ static void do_constraint(struct pull_t *pull, t_pbc *pbc, if (dnorm2(r_ij[c]) == 0) { - gmx_fatal(FARGS, "Distance for pull coordinate %lu is zero with constraint pulling, which is not allowed.", c + 1); + gmx_fatal(FARGS, "Distance for pull coordinate %zu is zero with constraint pulling, which is not allowed.", c + 1); } } diff --git a/src/gromacs/selection/selmethod.cpp b/src/gromacs/selection/selmethod.cpp index 6ac6d36ad1..bf9b2ee572 100644 --- a/src/gromacs/selection/selmethod.cpp +++ b/src/gromacs/selection/selmethod.cpp @@ -49,6 +49,7 @@ #include "gromacs/utility/arraysize.h" #include "gromacs/utility/cstringutil.h" #include "gromacs/utility/exceptions.h" +#include "gromacs/utility/stringutil.h" #include "selmethod-impl.h" #include "symrec.h" @@ -72,7 +73,10 @@ typedef struct { * Convenience function for reporting errors found in selection methods. */ static void -report_error(FILE *fp, const char *name, const char *fmt, ...) +report_error(FILE *fp, const char *name, gmx_fmtstr const char *fmt, ...) gmx_format(printf, 3, 4); + +static void +report_error(FILE *fp, const char *name, gmx_fmtstr const char *fmt, ...) { va_list ap; va_start(ap, fmt); @@ -90,7 +94,10 @@ report_error(FILE *fp, const char *name, const char *fmt, ...) */ static void report_param_error(FILE *fp, const char *mname, const char *pname, - const char *fmt, ...) + gmx_fmtstr const char *fmt, ...) gmx_format(printf, 4, 5); +static void +report_param_error(FILE *fp, const char *mname, const char *pname, + gmx_fmtstr const char *fmt, ...) { va_list ap; va_start(ap, fmt); diff --git a/src/gromacs/tables/splineutil.cpp b/src/gromacs/tables/splineutil.cpp index fc1f4941b7..9a21fe2e03 100644 --- a/src/gromacs/tables/splineutil.cpp +++ b/src/gromacs/tables/splineutil.cpp @@ -143,7 +143,7 @@ throwUnlessDerivativeIsConsistentWithFunction(ArrayRef func } if (!isConsistent) { - GMX_THROW(InconsistentInputError(formatString("Derivative inconsistent with numerical vector for elements %lu-%lu", minFail+1, maxFail+1))); + GMX_THROW(InconsistentInputError(formatString("Derivative inconsistent with numerical vector for elements %zu-%zu", minFail+1, maxFail+1))); } } diff --git a/src/gromacs/utility/fatalerror.cpp b/src/gromacs/utility/fatalerror.cpp index c3dfac65f8..bfb27986b8 100644 --- a/src/gromacs/utility/fatalerror.cpp +++ b/src/gromacs/utility/fatalerror.cpp @@ -226,7 +226,7 @@ void gmx_fatal_mpi_va(int /*f_errno*/, const char *file, int line, gmx_exit_on_fatal_error(exitType, 1); } -void gmx_fatal(int f_errno, const char *file, int line, const char *fmt, ...) +void gmx_fatal(int f_errno, const char *file, int line, gmx_fmtstr const char *fmt, ...) { va_list ap; va_start(ap, fmt); @@ -264,7 +264,7 @@ void _range_check(int n, int n_min, int n_max, const char *warn_str, } } -void gmx_warning(const char *fmt, ...) +void gmx_warning(gmx_fmtstr const char *fmt, ...) { va_list ap; char msg[STRLEN]; diff --git a/src/gromacs/utility/fatalerror.h b/src/gromacs/utility/fatalerror.h index 329de245dd..0f5ce2b507 100644 --- a/src/gromacs/utility/fatalerror.h +++ b/src/gromacs/utility/fatalerror.h @@ -48,6 +48,7 @@ #include #include "gromacs/utility/basedefinitions.h" +#include "gromacs/utility/stringutil.h" /*! \brief * Debug log file. @@ -173,7 +174,7 @@ gmx_fatal_mpi_va(int fatal_errno, const char *file, int line, \endcode */ [[noreturn]] void -gmx_fatal(int fatal_errno, const char *file, int line, const char *fmt, ...); +gmx_fatal(int fatal_errno, const char *file, int line, gmx_fmtstr const char *fmt, ...) gmx_format(printf, 4, 5); /** Helper macro to pass first three parameters to gmx_fatal(). */ #define FARGS 0, __FILE__, __LINE__ @@ -235,6 +236,6 @@ void _range_check(int n, int n_min, int n_max, const char *warn_str, * The message string should NOT start with "WARNING" * and should NOT end with a newline. */ -void gmx_warning(const char *fmt, ...); +void gmx_warning(gmx_fmtstr const char *fmt, ...) gmx_format(printf, 1, 2); #endif diff --git a/src/gromacs/utility/logger.cpp b/src/gromacs/utility/logger.cpp index 7be58088b4..e334634a7e 100644 --- a/src/gromacs/utility/logger.cpp +++ b/src/gromacs/utility/logger.cpp @@ -1,7 +1,7 @@ /* * This file is part of the GROMACS molecular simulation package. * - * Copyright (c) 2016, by the GROMACS development team, led by + * Copyright (c) 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. @@ -60,7 +60,7 @@ ILogTarget::~ILogTarget() } -LogEntryWriter &LogEntryWriter::appendTextFormatted(const char *fmt, ...) +LogEntryWriter &LogEntryWriter::appendTextFormatted(gmx_fmtstr const char *fmt, ...) { va_list ap; diff --git a/src/gromacs/utility/logger.h b/src/gromacs/utility/logger.h index f04af827fa..dc09a8d3fa 100644 --- a/src/gromacs/utility/logger.h +++ b/src/gromacs/utility/logger.h @@ -46,6 +46,8 @@ #include +#include "gromacs/utility/stringutil.h" + namespace gmx { @@ -92,7 +94,7 @@ class LogEntryWriter return *this; } //! Appends given text as a line in the log entry, with printf-style formatting. - LogEntryWriter &appendTextFormatted(const char *fmt, ...); + LogEntryWriter &appendTextFormatted(gmx_fmtstr const char *fmt, ...) gmx_format(printf, 2, 3); //! Writes the log entry with empty lines before and after. LogEntryWriter &asParagraph() { diff --git a/src/gromacs/utility/stringutil.cpp b/src/gromacs/utility/stringutil.cpp index ad0f963bf3..f42a3f16d9 100644 --- a/src/gromacs/utility/stringutil.cpp +++ b/src/gromacs/utility/stringutil.cpp @@ -133,7 +133,7 @@ std::string stripString(const std::string &str) return std::string(start, end); } -std::string formatString(const char *fmt, ...) +std::string formatString(gmx_fmtstr const char *fmt, ...) { va_list ap; va_start(ap, fmt); diff --git a/src/gromacs/utility/stringutil.h b/src/gromacs/utility/stringutil.h index 2c0be63baa..b3cca7e3e3 100644 --- a/src/gromacs/utility/stringutil.h +++ b/src/gromacs/utility/stringutil.h @@ -1,7 +1,7 @@ /* * This file is part of the GROMACS molecular simulation package. * - * Copyright (c) 2011,2012,2013,2014,2015,2016,2017, by the GROMACS development team, led by + * Copyright (c) 2011,2012,2013,2014,2015,2016,2017,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. @@ -168,6 +168,31 @@ std::string stripSuffixIfPresent(const std::string &str, const char *suffix); * \throws std::bad_alloc if out of memory. */ std::string stripString(const std::string &str); +#ifdef __GNUC__ +#define gmx_format(archetype, string_index, first_to_check) \ + __attribute__ ((format (archetype, string_index, first_to_check))) +#else +/*! \brief GCC like function format attribute + * + * The format attribute specifies that a function takes printf, scanf, ... + * style arguments that should be type-checked against a format string. + * The attribute has to be placed after the function. + * This attribute is only valid for function declarations and not function + * definitions (GCC limitation). For member functions the implicit `this` + * pointer is included in the argument count. + */ +#define gmx_format(archetype, string_index, first_to_check) +#endif +#ifdef _MSC_VER +#define gmx_fmtstr _In_ _Printf_format_string_ +#else +/*! \brief MSVC like function format attribute + * + * Does type checking for printf like format strings in MSVC style. + * Attribute has to be placed before format string. + */ +#define gmx_fmtstr +#endif /*! \brief * Formats a string (snprintf() wrapper). * @@ -177,7 +202,8 @@ std::string stripString(const std::string &str); * instead of requiring a preallocated buffer. Arbitrary length output is * supported. */ -std::string formatString(const char *fmt, ...); +std::string formatString(gmx_fmtstr const char *fmt, ...) gmx_format(printf, 1, 2); + /*! \brief * Formats a string (vsnprintf() wrapper). *