Use MDLogger in grompp
authorPaul Bauer <paul.bauer.q@gmail.com>
Wed, 26 Jun 2019 15:17:01 +0000 (17:17 +0200)
committerPaul Bauer <paul.bauer.q@gmail.com>
Fri, 31 Jan 2020 20:11:53 +0000 (21:11 +0100)
Start using MDLogger in main grompp routine.

Refs #2999, #3005

Change-Id: Ic88e9a2ec271a4bf9502da5a268f532d422a1f1b

src/gromacs/gmxpreprocess/grompp.cpp

index 28562ffeb4904b5d305e5decd32a893e363fb8c0..a5a830c4d869879ca9df848910a86bd2c5df58f8 100644 (file)
 #include "gromacs/utility/cstringutil.h"
 #include "gromacs/utility/exceptions.h"
 #include "gromacs/utility/fatalerror.h"
+#include "gromacs/utility/filestream.h"
 #include "gromacs/utility/futil.h"
 #include "gromacs/utility/gmxassert.h"
 #include "gromacs/utility/keyvaluetreebuilder.h"
 #include "gromacs/utility/listoflists.h"
+#include "gromacs/utility/logger.h"
+#include "gromacs/utility/loggerbuilder.h"
 #include "gromacs/utility/mdmodulenotification.h"
 #include "gromacs/utility/smalloc.h"
 #include "gromacs/utility/snprintf.h"
@@ -263,11 +266,16 @@ static int rm_interactions(int ifunc, gmx::ArrayRef<MoleculeInformation> mols)
     return n;
 }
 
-static int check_atom_names(const char* fn1, const char* fn2, gmx_mtop_t* mtop, const t_atoms* at)
+static int check_atom_names(const char*          fn1,
+                            const char*          fn2,
+                            gmx_mtop_t*          mtop,
+                            const t_atoms*       at,
+                            const gmx::MDLogger& logger)
 {
     int      m, i, j, nmismatch;
     t_atoms* tat;
-#define MAXMISMATCH 20
+
+    constexpr int c_maxNumberOfMismatches = 20;
 
     if (mtop->natoms != at->nr)
     {
@@ -285,15 +293,20 @@ static int check_atom_names(const char* fn1, const char* fn2, gmx_mtop_t* mtop,
             {
                 if (strcmp(*(tat->atomname[j]), *(at->atomname[i])) != 0)
                 {
-                    if (nmismatch < MAXMISMATCH)
+                    if (nmismatch < c_maxNumberOfMismatches)
                     {
-                        fprintf(stderr,
-                                "Warning: atom name %d in %s and %s does not match (%s - %s)\n",
-                                i + 1, fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
+                        GMX_LOG(logger.warning)
+                                .asParagraph()
+                                .appendTextFormatted(
+                                        "atom name %d in %s and %s does not match (%s - %s)", i + 1,
+                                        fn1, fn2, *(tat->atomname[j]), *(at->atomname[i]));
                     }
-                    else if (nmismatch == MAXMISMATCH)
+                    else if (nmismatch == c_maxNumberOfMismatches)
                     {
-                        fprintf(stderr, "(more than %d non-matching atom names)\n", MAXMISMATCH);
+                        GMX_LOG(logger.warning)
+                                .asParagraph()
+                                .appendTextFormatted("(more than %d non-matching atom names)",
+                                                     c_maxNumberOfMismatches);
                     }
                     nmismatch++;
                 }
@@ -331,7 +344,6 @@ static void check_bonds_timestep(const gmx_mtop_t* mtop, double dt, warninp* wi)
     int  i, a1, a2, w_a1, w_a2, j;
     real twopi2, limit2, fc, re, m1, m2, period2, w_period2;
     bool bFound, bWater, bWarn;
-    char warn_buf[STRLEN];
 
     /* Get the interaction parameters */
     gmx::ArrayRef<const t_iparams> ip = mtop->ffparams.iparams;
@@ -420,7 +432,7 @@ static void check_bonds_timestep(const gmx_mtop_t* mtop, double dt, warninp* wi)
         bWarn = (w_period2 < gmx::square(min_steps_warn * dt));
         /* A check that would recognize most water models */
         bWater = ((*w_moltype->atoms.atomname[0])[0] == 'O' && w_moltype->atoms.nr <= 5);
-        sprintf(warn_buf,
+        std::string warningMessage = gmx::formatString(
                 "The bond in molecule-type %s between atoms %d %s and %d %s has an estimated "
                 "oscillational period of %.1e ps, which is less than %d times the time step of "
                 "%.1e ps.\n"
@@ -432,11 +444,11 @@ static void check_bonds_timestep(const gmx_mtop_t* mtop, double dt, warninp* wi)
                        : "Maybe you forgot to change the constraints mdp option.");
         if (bWarn)
         {
-            warning(wi, warn_buf);
+            warning(wi, warningMessage.c_str());
         }
         else
         {
-            warning_note(wi, warn_buf);
+            warning_note(wi, warningMessage.c_str());
         }
     }
 }
@@ -456,8 +468,7 @@ static void check_vel(gmx_mtop_t* mtop, rvec v[])
 
 static void check_shells_inputrec(gmx_mtop_t* mtop, t_inputrec* ir, warninp* wi)
 {
-    int  nshells = 0;
-    char warn_buf[STRLEN];
+    int nshells = 0;
 
     for (const AtomProxy atomP : AtomRange(*mtop))
     {
@@ -470,10 +481,10 @@ static void check_shells_inputrec(gmx_mtop_t* mtop, t_inputrec* ir, warninp* wi)
     if ((nshells > 0) && (ir->nstcalcenergy != 1))
     {
         set_warning_line(wi, "unknown", -1);
-        snprintf(warn_buf, STRLEN, "There are %d shells, changing nstcalcenergy from %d to 1",
-                 nshells, ir->nstcalcenergy);
+        std::string warningMessage = gmx::formatString(
+                "There are %d shells, changing nstcalcenergy from %d to 1", nshells, ir->nstcalcenergy);
         ir->nstcalcenergy = 1;
-        warning(wi, warn_buf);
+        warning(wi, warningMessage.c_str());
     }
 }
 
@@ -570,13 +581,12 @@ static void new_status(const char*                           topfile,
                        double*                               reppow,
                        real*                                 fudgeQQ,
                        gmx_bool                              bMorse,
-                       warninp*                              wi)
+                       warninp*                              wi,
+                       const gmx::MDLogger&                  logger)
 {
     std::vector<gmx_molblock_t> molblock;
     int                         i, nmismatch;
     bool                        ffParametrizedWithHBondConstraints;
-    char                        buf[STRLEN];
-    char                        warn_buf[STRLEN];
 
     /* TOPOLOGY processing */
     sys->name = do_top(bVerbose, topfile, topppfile, opts, bZero, &(sys->symtab), interactions,
@@ -618,8 +628,9 @@ static void new_status(const char*                           topfile,
         if (i > 0)
         {
             set_warning_line(wi, "unknown", -1);
-            sprintf(warn_buf, "disre = no, removed %d distance restraints", i);
-            warning_note(wi, warn_buf);
+            std::string warningMessage =
+                    gmx::formatString("disre = no, removed %d distance restraints", i);
+            warning_note(wi, warningMessage.c_str());
         }
     }
     if (!opts->bOrire)
@@ -628,8 +639,9 @@ static void new_status(const char*                           topfile,
         if (i > 0)
         {
             set_warning_line(wi, "unknown", -1);
-            sprintf(warn_buf, "orire = no, removed %d orientation restraints", i);
-            warning_note(wi, warn_buf);
+            std::string warningMessage =
+                    gmx::formatString("orire = no, removed %d orientation restraints", i);
+            warning_note(wi, warningMessage.c_str());
         }
     }
 
@@ -656,7 +668,7 @@ static void new_status(const char*                           topfile,
     /* COORDINATE file processing */
     if (bVerbose)
     {
-        fprintf(stderr, "processing coordinates...\n");
+        GMX_LOG(logger.info).asParagraph().appendTextFormatted("processing coordinates...");
     }
 
     t_topology* conftop;
@@ -691,24 +703,26 @@ static void new_status(const char*                           topfile,
     /* This call fixes the box shape for runs with pressure scaling */
     set_box_rel(ir, state);
 
-    nmismatch = check_atom_names(topfile, confin, sys, &conftop->atoms);
+    nmismatch = check_atom_names(topfile, confin, sys, &conftop->atoms, logger);
     done_top(conftop);
     sfree(conftop);
 
     if (nmismatch)
     {
-        sprintf(buf,
+        std::string warningMessage = gmx::formatString(
                 "%d non-matching atom name%s\n"
                 "atom names from %s will be used\n"
                 "atom names from %s will be ignored\n",
                 nmismatch, (nmismatch == 1) ? "" : "s", topfile, confin);
-        warning(wi, buf);
+        warning(wi, warningMessage.c_str());
     }
 
     /* Do more checks, mostly related to constraints */
     if (bVerbose)
     {
-        fprintf(stderr, "double-checking input for internal consistency...\n");
+        GMX_LOG(logger.info)
+                .asParagraph()
+                .appendTextFormatted("double-checking input for internal consistency...");
     }
     {
         bool bHasNormalConstraints =
@@ -732,7 +746,7 @@ static void new_status(const char*                           topfile,
         if (opts->seed == -1)
         {
             opts->seed = static_cast<int>(gmx::makeRandomSeed());
-            fprintf(stderr, "Setting gen_seed to %d\n", opts->seed);
+            GMX_LOG(logger.info).asParagraph().appendTextFormatted("Setting gen_seed to %d", opts->seed);
         }
         state->flags |= (1 << estV);
         maxwell_speed(opts->tempi, opts->seed, sys, state->v.rvec_array());
@@ -778,7 +792,8 @@ static void cont_status(const char*             slog,
                         t_inputrec*             ir,
                         t_state*                state,
                         gmx_mtop_t*             sys,
-                        const gmx_output_env_t* oenv)
+                        const gmx_output_env_t* oenv,
+                        const gmx::MDLogger&    logger)
 /* If fr_time == -1 read the last frame available which is complete */
 {
     bool         bReadVel;
@@ -788,23 +803,27 @@ static void cont_status(const char*             slog,
 
     bReadVel = (bNeedVel && !bGenVel);
 
-    fprintf(stderr, "Reading Coordinates%s and Box size from old trajectory\n",
-            bReadVel ? ", Velocities" : "");
+    GMX_LOG(logger.info)
+            .asParagraph()
+            .appendTextFormatted("Reading Coordinates%s and Box size from old trajectory",
+                                 bReadVel ? ", Velocities" : "");
     if (fr_time == -1)
     {
-        fprintf(stderr, "Will read whole trajectory\n");
+        GMX_LOG(logger.info).asParagraph().appendTextFormatted("Will read whole trajectory");
     }
     else
     {
-        fprintf(stderr, "Will read till time %g\n", fr_time);
+        GMX_LOG(logger.info).asParagraph().appendTextFormatted("Will read till time %g", fr_time);
     }
     if (!bReadVel)
     {
         if (bGenVel)
         {
-            fprintf(stderr,
-                    "Velocities generated: "
-                    "ignoring velocities in input trajectory\n");
+            GMX_LOG(logger.info)
+                    .asParagraph()
+                    .appendTextFormatted(
+                            "Velocities generated: "
+                            "ignoring velocities in input trajectory");
         }
         read_first_frame(oenv, &fp, slog, &fr, TRX_NEED_X);
     }
@@ -814,11 +833,12 @@ static void cont_status(const char*             slog,
 
         if (!fr.bV)
         {
-            fprintf(stderr,
-                    "\n"
-                    "WARNING: Did not find a frame with velocities in file %s,\n"
-                    "         all velocities will be set to zero!\n\n",
-                    slog);
+            GMX_LOG(logger.warning)
+                    .asParagraph()
+                    .appendTextFormatted(
+                            "WARNING: Did not find a frame with velocities in file %s,\n"
+                            "         all velocities will be set to zero!",
+                            slog);
             for (auto& vi : makeArrayRef(state->v))
             {
                 vi = { 0, 0, 0 };
@@ -854,8 +874,8 @@ static void cont_status(const char*             slog,
      */
     set_box_rel(ir, state);
 
-    fprintf(stderr, "Using frame at t = %g ps\n", use_time);
-    fprintf(stderr, "Starting time for run is %g ps\n", ir->init_t);
+    GMX_LOG(logger.info).asParagraph().appendTextFormatted("Using frame at t = %g ps", use_time);
+    GMX_LOG(logger.info).asParagraph().appendTextFormatted("Starting time for run is %g ps", ir->init_t);
 
     if ((ir->epc != epcNO || ir->etc == etcNOSEHOOVER) && ener)
     {
@@ -871,7 +891,8 @@ static void read_posres(gmx_mtop_t*                              mtop,
                         int                                      rc_scaling,
                         PbcType                                  pbcType,
                         rvec                                     com,
-                        warninp*                                 wi)
+                        warninp*                                 wi,
+                        const gmx::MDLogger&                     logger)
 {
     gmx_bool*   hadAtom;
     rvec *      x, *v;
@@ -880,7 +901,6 @@ static void read_posres(gmx_mtop_t*                              mtop,
     t_topology* top;
     matrix      box, invbox;
     int         natoms, npbcdim = 0;
-    char        warn_buf[STRLEN];
     int         a, nat_molb;
     t_atom*     atom;
 
@@ -891,11 +911,11 @@ static void read_posres(gmx_mtop_t*                              mtop,
     sfree(top);
     if (natoms != mtop->natoms)
     {
-        sprintf(warn_buf,
+        std::string warningMessage = gmx::formatString(
                 "The number of atoms in %s (%d) does not match the number of atoms in the topology "
                 "(%d). Will assume that the first %d atoms in the topology and %s match.",
                 fn, natoms, mtop->natoms, std::min(mtop->natoms, natoms), fn);
-        warning(wi, warn_buf);
+        warning(wi, warningMessage.c_str());
     }
 
     npbcdim = numPbcDimensions(pbcType);
@@ -996,9 +1016,11 @@ static void read_posres(gmx_mtop_t*                              mtop,
         {
             com[j] = sum[j] / totmass;
         }
-        fprintf(stderr,
-                "The center of mass of the position restraint coord's is %6.3f %6.3f %6.3f\n",
-                com[XX], com[YY], com[ZZ]);
+        GMX_LOG(logger.info)
+                .asParagraph()
+                .appendTextFormatted(
+                        "The center of mass of the position restraint coord's is %6.3f %6.3f %6.3f",
+                        com[XX], com[YY], com[ZZ]);
     }
 
     if (rc_scaling != erscNO)
@@ -1061,31 +1083,36 @@ static void gen_posres(gmx_mtop_t*                              mtop,
                        PbcType                                  pbcType,
                        rvec                                     com,
                        rvec                                     comB,
-                       warninp*                                 wi)
+                       warninp*                                 wi,
+                       const gmx::MDLogger&                     logger)
 {
-    read_posres(mtop, mi, FALSE, fnA, rc_scaling, pbcType, com, wi);
+    read_posres(mtop, mi, FALSE, fnA, rc_scaling, pbcType, com, wi, logger);
     /* It is safer to simply read the b-state posres rather than trying
      * to be smart and copy the positions.
      */
-    read_posres(mtop, mi, TRUE, fnB, rc_scaling, pbcType, comB, wi);
+    read_posres(mtop, mi, TRUE, fnB, rc_scaling, pbcType, comB, wi, logger);
 }
 
-static void set_wall_atomtype(PreprocessingAtomTypes* at, t_gromppopts* opts, t_inputrec* ir, warninp* wi)
+static void set_wall_atomtype(PreprocessingAtomTypes* at,
+                              t_gromppopts*           opts,
+                              t_inputrec*             ir,
+                              warninp*                wi,
+                              const gmx::MDLogger&    logger)
 {
-    int  i;
-    char warn_buf[STRLEN];
+    int i;
 
     if (ir->nwall > 0)
     {
-        fprintf(stderr, "Searching the wall atom type(s)\n");
+        GMX_LOG(logger.info).asParagraph().appendTextFormatted("Searching the wall atom type(s)");
     }
     for (i = 0; i < ir->nwall; i++)
     {
         ir->wall_atomtype[i] = at->atomTypeFromName(opts->wall_atomtype[i]);
         if (ir->wall_atomtype[i] == NOTSET)
         {
-            sprintf(warn_buf, "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
-            warning_error(wi, warn_buf);
+            std::string warningMessage = gmx::formatString(
+                    "Specified wall atom type %s is not defined", opts->wall_atomtype[i]);
+            warning_error(wi, warningMessage.c_str());
         }
     }
 }
@@ -1248,8 +1275,7 @@ static void init_cmap_grid(gmx_cmap_t* cmap_grid, int ngrid, int grid_spacing)
 
 static int count_constraints(const gmx_mtop_t* mtop, gmx::ArrayRef<const MoleculeInformation> mi, warninp* wi)
 {
-    int  count, count_mol;
-    char buf[STRLEN];
+    int count, count_mol;
 
     count = 0;
     for (const gmx_molblock_t& molb : mtop->molblock)
@@ -1271,12 +1297,12 @@ static int count_constraints(const gmx_mtop_t* mtop, gmx::ArrayRef<const Molecul
 
         if (count_mol > nrdf_internal(&mi[molb.type].atoms))
         {
-            sprintf(buf,
+            std::string warningMessage = gmx::formatString(
                     "Molecule type '%s' has %d constraints.\n"
                     "For stability and efficiency there should not be more constraints than "
                     "internal number of degrees of freedom: %d.\n",
                     *mi[molb.type].name, count_mol, nrdf_internal(&mi[molb.type].atoms));
-            warning(wi, buf);
+            warning(wi, warningMessage.c_str());
         }
         count += molb.nmol * count_mol;
     }
@@ -1325,14 +1351,12 @@ static real get_max_reference_temp(const t_inputrec* ir, warninp* wi)
 
     if (bNoCoupl)
     {
-        char buf[STRLEN];
-
-        sprintf(buf,
+        std::string warningMessage = gmx::formatString(
                 "Some temperature coupling groups do not use temperature coupling. We will assume "
                 "their temperature is not more than %.3f K. If their temperature is higher, the "
                 "energy error and the Verlet buffer might be underestimated.",
                 ref_t);
-        warning(wi, buf);
+        warning(wi, warningMessage.c_str());
     }
 
     return ref_t;
@@ -1341,7 +1365,7 @@ static real get_max_reference_temp(const t_inputrec* ir, warninp* wi)
 /* Checks if there are unbound atoms in moleculetype molt.
  * Prints a note for each unbound atoms and a warning if any is present.
  */
-static void checkForUnboundAtoms(const gmx_moltype_t* molt, gmx_bool bVerbose, warninp* wi)
+static void checkForUnboundAtoms(const gmx_moltype_t* molt, gmx_bool bVerbose, warninp* wi, const gmx::MDLogger& logger)
 {
     const t_atoms* atoms = &molt->atoms;
 
@@ -1378,10 +1402,12 @@ static void checkForUnboundAtoms(const gmx_moltype_t* molt, gmx_bool bVerbose, w
         {
             if (bVerbose)
             {
-                fprintf(stderr,
-                        "\nAtom %d '%s' in moleculetype '%s' is not bound by a potential or "
-                        "constraint to any other atom in the same moleculetype.\n",
-                        a + 1, *atoms->atomname[a], *molt->name);
+                GMX_LOG(logger.warning)
+                        .asParagraph()
+                        .appendTextFormatted(
+                                "Atom %d '%s' in moleculetype '%s' is not bound by a potential or "
+                                "constraint to any other atom in the same moleculetype.",
+                                a + 1, *atoms->atomname[a], *molt->name);
             }
             numDanglingAtoms++;
         }
@@ -1389,24 +1415,23 @@ static void checkForUnboundAtoms(const gmx_moltype_t* molt, gmx_bool bVerbose, w
 
     if (numDanglingAtoms > 0)
     {
-        char buf[STRLEN];
-        sprintf(buf,
+        std::string warningMessage = gmx::formatString(
                 "In moleculetype '%s' %d atoms are not bound by a potential or constraint to any "
                 "other atom in the same moleculetype. Although technically this might not cause "
                 "issues in a simulation, this often means that the user forgot to add a "
                 "bond/potential/constraint or put multiple molecules in the same moleculetype "
                 "definition by mistake. Run with -v to get information for each atom.",
                 *molt->name, numDanglingAtoms);
-        warning_note(wi, buf);
+        warning_note(wi, warningMessage.c_str());
     }
 }
 
 /* Checks all moleculetypes for unbound atoms */
-static void checkForUnboundAtoms(const gmx_mtop_t* mtop, gmx_bool bVerbose, warninp* wi)
+static void checkForUnboundAtoms(const gmx_mtop_t* mtop, gmx_bool bVerbose, warninp* wi, const gmx::MDLogger& logger)
 {
     for (const gmx_moltype_t& molt : mtop->moltype)
     {
-        checkForUnboundAtoms(&molt, bVerbose, wi);
+        checkForUnboundAtoms(&molt, bVerbose, wi, logger);
     }
 }
 
@@ -1574,12 +1599,18 @@ static void checkDecoupledModeAccuracy(const gmx_mtop_t* mtop, const t_inputrec*
     }
 }
 
-static void set_verlet_buffer(const gmx_mtop_t* mtop, t_inputrec* ir, real buffer_temp, matrix box, warninp* wi)
+static void set_verlet_buffer(const gmx_mtop_t*    mtop,
+                              t_inputrec*          ir,
+                              real                 buffer_temp,
+                              matrix               box,
+                              warninp*             wi,
+                              const gmx::MDLogger& logger)
 {
-    char warn_buf[STRLEN];
-
-    printf("Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K\n", ir->verletbuf_tol,
-           buffer_temp);
+    GMX_LOG(logger.info)
+            .asParagraph()
+            .appendTextFormatted(
+                    "Determining Verlet buffer for a tolerance of %g kJ/mol/ps at %g K",
+                    ir->verletbuf_tol, buffer_temp);
 
     /* Calculate the buffer size for simple atom vs atoms list */
     VerletbufListSetup listSetup1x1;
@@ -1596,22 +1627,31 @@ static void set_verlet_buffer(const gmx_mtop_t* mtop, t_inputrec* ir, real buffe
     const int n_nonlin_vsite = countNonlinearVsites(*mtop);
     if (n_nonlin_vsite > 0)
     {
-        sprintf(warn_buf,
+        std::string warningMessage = gmx::formatString(
                 "There are %d non-linear virtual site constructions. Their contribution to the "
                 "energy error is approximated. In most cases this does not affect the error "
                 "significantly.",
                 n_nonlin_vsite);
-        warning_note(wi, warn_buf);
+        warning_note(wi, warningMessage);
     }
 
-    printf("Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm\n", 1, 1,
-           rlist_1x1, rlist_1x1 - std::max(ir->rvdw, ir->rcoulomb));
+    GMX_LOG(logger.info)
+            .asParagraph()
+            .appendTextFormatted(
+                    "Calculated rlist for %dx%d atom pair-list as %.3f nm, buffer size %.3f nm", 1,
+                    1, rlist_1x1, rlist_1x1 - std::max(ir->rvdw, ir->rcoulomb));
 
-    printf("Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm\n",
-           listSetup4x4.cluster_size_i, listSetup4x4.cluster_size_j, ir->rlist,
-           ir->rlist - std::max(ir->rvdw, ir->rcoulomb));
+    GMX_LOG(logger.info)
+            .asParagraph()
+            .appendTextFormatted(
+                    "Set rlist, assuming %dx%d atom pair-list, to %.3f nm, buffer size %.3f nm",
+                    listSetup4x4.cluster_size_i, listSetup4x4.cluster_size_j, ir->rlist,
+                    ir->rlist - std::max(ir->rvdw, ir->rcoulomb));
 
-    printf("Note that mdrun will redetermine rlist based on the actual pair-list setup\n");
+    GMX_LOG(logger.info)
+            .asParagraph()
+            .appendTextFormatted(
+                    "Note that mdrun will redetermine rlist based on the actual pair-list setup");
 
     if (gmx::square(ir->rlist) >= max_cutoff2(ir->pbcType, box))
     {
@@ -1735,7 +1775,6 @@ int gmx_grompp(int argc, char* argv[])
     gmx_output_env_t*                    oenv;
     gmx_bool                             bVerbose = FALSE;
     warninp*                             wi;
-    char                                 warn_buf[STRLEN];
 
     t_filenm fnm[] = { { efMDP, nullptr, nullptr, ffREAD },
                        { efMDP, "-po", "mdout", ffWRITE },
@@ -1799,6 +1838,13 @@ int gmx_grompp(int argc, char* argv[])
     snew(opts->include, STRLEN);
     snew(opts->define, STRLEN);
 
+    gmx::LoggerBuilder builder;
+    builder.addTargetStream(gmx::MDLogger::LogLevel::Info, &gmx::TextOutputFile::standardOutput());
+    builder.addTargetStream(gmx::MDLogger::LogLevel::Warning, &gmx::TextOutputFile::standardError());
+    gmx::LoggerOwner    logOwner(builder.build());
+    const gmx::MDLogger logger(logOwner.logger());
+
+
     wi = init_warning(TRUE, maxwarn);
 
     /* PARAMETER file processing */
@@ -1812,31 +1858,37 @@ int gmx_grompp(int argc, char* argv[])
 
     if (bVerbose)
     {
-        fprintf(stderr, "checking input for internal consistency...\n");
+        GMX_LOG(logger.info)
+                .asParagraph()
+                .appendTextFormatted("checking input for internal consistency...");
     }
     check_ir(mdparin, mdModules.notifier(), ir, opts, wi);
 
     if (ir->ld_seed == -1)
     {
         ir->ld_seed = static_cast<int>(gmx::makeRandomSeed());
-        fprintf(stderr, "Setting the LD random seed to %" PRId64 "\n", ir->ld_seed);
+        GMX_LOG(logger.info)
+                .asParagraph()
+                .appendTextFormatted("Setting the LD random seed to %" PRId64 "", ir->ld_seed);
     }
 
     if (ir->expandedvals->lmc_seed == -1)
     {
         ir->expandedvals->lmc_seed = static_cast<int>(gmx::makeRandomSeed());
-        fprintf(stderr, "Setting the lambda MC random seed to %d\n", ir->expandedvals->lmc_seed);
+        GMX_LOG(logger.info)
+                .asParagraph()
+                .appendTextFormatted("Setting the lambda MC random seed to %d", ir->expandedvals->lmc_seed);
     }
 
     bNeedVel = EI_STATE_VELOCITY(ir->eI);
     bGenVel  = (bNeedVel && opts->bGenVel);
     if (bGenVel && ir->bContinuation)
     {
-        sprintf(warn_buf,
+        std::string warningMessage = gmx::formatString(
                 "Generating velocities is inconsistent with attempting "
                 "to continue a previous run. Choose only one of "
                 "gen-vel = yes and continuation = yes.");
-        warning_error(wi, warn_buf);
+        warning_error(wi, warningMessage);
     }
 
     std::array<InteractionsOfType, F_NRE> interactions;
@@ -1856,7 +1908,7 @@ int gmx_grompp(int argc, char* argv[])
     t_state state;
     new_status(fn, opt2fn_null("-pp", NFILE, fnm), opt2fn("-c", NFILE, fnm), opts, ir, bZero,
                bGenVel, bVerbose, &state, &atypes, &sys, &mi, &intermolecular_interactions,
-               interactions, &comb, &reppow, &fudgeQQ, opts->bMorse, wi);
+               interactions, &comb, &reppow, &fudgeQQ, opts->bMorse, wi, logger);
 
     if (debug)
     {
@@ -1883,15 +1935,17 @@ int gmx_grompp(int argc, char* argv[])
     {
         if (ir->eI == eiCG || ir->eI == eiLBFGS)
         {
-            sprintf(warn_buf, "Can not do %s with %s, use %s", EI(ir->eI),
-                    econstr_names[econtSHAKE], econstr_names[econtLINCS]);
-            warning_error(wi, warn_buf);
+            std::string warningMessage =
+                    gmx::formatString("Can not do %s with %s, use %s", EI(ir->eI),
+                                      econstr_names[econtSHAKE], econstr_names[econtLINCS]);
+            warning_error(wi, warningMessage);
         }
         if (ir->bPeriodicMols)
         {
-            sprintf(warn_buf, "Can not do periodic molecules with %s, use %s",
-                    econstr_names[econtSHAKE], econstr_names[econtLINCS]);
-            warning_error(wi, warn_buf);
+            std::string warningMessage =
+                    gmx::formatString("Can not do periodic molecules with %s, use %s",
+                                      econstr_names[econtSHAKE], econstr_names[econtLINCS]);
+            warning_error(wi, warningMessage);
         }
     }
 
@@ -1928,12 +1982,12 @@ int gmx_grompp(int argc, char* argv[])
     {
         if (ir->epc == epcPARRINELLORAHMAN || ir->epc == epcMTTK)
         {
-            sprintf(warn_buf,
+            std::string warningMessage = gmx::formatString(
                     "You are combining position restraints with %s pressure coupling, which can "
                     "lead to instabilities. If you really want to combine position restraints with "
                     "pressure coupling, we suggest to use %s pressure coupling instead.",
                     EPCOUPLTYPE(ir->epc), EPCOUPLTYPE(epcBERENDSEN));
-            warning_note(wi, warn_buf);
+            warning_note(wi, warningMessage);
         }
 
         const char* fn = opt2fn("-r", NFILE, fnm);
@@ -1969,18 +2023,15 @@ int gmx_grompp(int argc, char* argv[])
 
         if (bVerbose)
         {
-            fprintf(stderr, "Reading position restraint coords from %s", fn);
-            if (strcmp(fn, fnB) == 0)
-            {
-                fprintf(stderr, "\n");
-            }
-            else
+            std::string message = gmx::formatString("Reading position restraint coords from %s", fn);
+            if (strcmp(fn, fnB) != 0)
             {
-                fprintf(stderr, " and %s\n", fnB);
+                message += gmx::formatString(" and %s", fnB);
             }
+            GMX_LOG(logger.info).asParagraph().appendText(message);
         }
         gen_posres(&sys, mi, fn, fnB, ir->refcoord_scaling, ir->pbcType, ir->posres_com,
-                   ir->posres_comB, wi);
+                   ir->posres_comB, wi, logger);
     }
 
     /* If we are using CMAP, setup the pre-interpolation grid */
@@ -1992,7 +2043,7 @@ int gmx_grompp(int argc, char* argv[])
                    interactions[F_CMAP].cmap, &sys.ffparams.cmap_grid);
     }
 
-    set_wall_atomtype(&atypes, opts, ir, wi);
+    set_wall_atomtype(&atypes, opts, ir, wi, logger);
     if (bRenum)
     {
         atypes.renumberTypes(interactions, &sys, ir->wall_atomtype, bVerbose);
@@ -2013,7 +2064,7 @@ int gmx_grompp(int argc, char* argv[])
 
     if (bVerbose)
     {
-        fprintf(stderr, "converting bonded parameters...\n");
+        GMX_LOG(logger.info).asParagraph().appendTextFormatted("converting bonded parameters...");
     }
 
     const int ntype = atypes.size();
@@ -2046,7 +2097,7 @@ int gmx_grompp(int argc, char* argv[])
     /* check masses */
     check_mol(&sys, wi);
 
-    checkForUnboundAtoms(&sys, bVerbose, wi);
+    checkForUnboundAtoms(&sys, bVerbose, wi, logger);
 
     if (EI_DYNAMICS(ir->eI) && ir->eI != eiBD)
     {
@@ -2069,7 +2120,7 @@ int gmx_grompp(int argc, char* argv[])
 
     if (bVerbose)
     {
-        fprintf(stderr, "initialising group options...\n");
+        GMX_LOG(logger.info).asParagraph().appendTextFormatted("initialising group options...");
     }
     do_index(mdparin, ftp2fn_null(efNDX, NFILE, fnm), &sys, bVerbose, mdModules.notifier(), ir, wi);
 
@@ -2091,19 +2142,19 @@ int gmx_grompp(int argc, char* argv[])
                 }
                 if (buffer_temp > 0)
                 {
-                    sprintf(warn_buf,
+                    std::string warningMessage = gmx::formatString(
                             "NVE simulation: will use the initial temperature of %.3f K for "
                             "determining the Verlet buffer size",
                             buffer_temp);
-                    warning_note(wi, warn_buf);
+                    warning_note(wi, warningMessage);
                 }
                 else
                 {
-                    sprintf(warn_buf,
+                    std::string warningMessage = gmx::formatString(
                             "NVE simulation with an initial temperature of zero: will use a Verlet "
                             "buffer of %d%%. Check your energy drift!",
                             gmx::roundToInt(verlet_buffer_ratio_NVE_T0 * 100));
-                    warning_note(wi, warn_buf);
+                    warning_note(wi, warningMessage);
                 }
             }
             else
@@ -2138,7 +2189,7 @@ int gmx_grompp(int argc, char* argv[])
 
                     if (ir->verletbuf_tol > 1.1 * driftTolerance * totalEnergyDriftPerAtomPerPicosecond)
                     {
-                        sprintf(warn_buf,
+                        std::string warningMessage = gmx::formatString(
                                 "You are using a Verlet buffer tolerance of %g kJ/mol/ps for an "
                                 "NVE simulation of length %g ps, which can give a final drift of "
                                 "%d%%. For conserving energy to %d%% when using constraints, you "
@@ -2147,11 +2198,11 @@ int gmx_grompp(int argc, char* argv[])
                                 gmx::roundToInt(ir->verletbuf_tol / totalEnergyDriftPerAtomPerPicosecond * 100),
                                 gmx::roundToInt(100 * driftTolerance),
                                 driftTolerance * totalEnergyDriftPerAtomPerPicosecond);
-                        warning_note(wi, warn_buf);
+                        warning_note(wi, warningMessage);
                     }
                 }
 
-                set_verlet_buffer(&sys, ir, buffer_temp, state.box, wi);
+                set_verlet_buffer(&sys, ir, buffer_temp, state.box, wi, logger);
             }
         }
     }
@@ -2191,10 +2242,12 @@ int gmx_grompp(int argc, char* argv[])
     {
         if (bVerbose)
         {
-            fprintf(stderr, "getting data from old trajectory ...\n");
+            GMX_LOG(logger.info)
+                    .asParagraph()
+                    .appendTextFormatted("getting data from old trajectory ...");
         }
         cont_status(ftp2fn(efTRN, NFILE, fnm), ftp2fn_null(efEDR, NFILE, fnm), bNeedVel, bGenVel,
-                    fr_time, ir, &state, &sys, oenv);
+                    fr_time, ir, &state, &sys, oenv, logger);
     }
 
     if (ir->pbcType == PbcType::XY && ir->nwall != 2)
@@ -2297,7 +2350,10 @@ int gmx_grompp(int argc, char* argv[])
     if (EEL_PME(ir->coulombtype))
     {
         float ratio = pme_load_estimate(sys, *ir, state.box);
-        fprintf(stderr, "Estimate for the relative computational load of the PME mesh part: %.2f\n", ratio);
+        GMX_LOG(logger.info)
+                .asParagraph()
+                .appendTextFormatted(
+                        "Estimate for the relative computational load of the PME mesh part: %.2f", ratio);
         /* With free energy we might need to do PME both for the A and B state
          * charges. This will double the cost, but the optimal performance will
          * then probably be at a slightly larger cut-off and grid spacing.
@@ -2319,17 +2375,17 @@ int gmx_grompp(int argc, char* argv[])
     }
 
     {
-        char   warn_buf[STRLEN];
-        double cio = compute_io(ir, sys.natoms, sys.groups, F_NRE, 1);
-        sprintf(warn_buf, "This run will generate roughly %.0f Mb of data", cio);
+        double      cio = compute_io(ir, sys.natoms, sys.groups, F_NRE, 1);
+        std::string warningMessage =
+                gmx::formatString("This run will generate roughly %.0f Mb of data", cio);
         if (cio > 2000)
         {
             set_warning_line(wi, mdparin, -1);
-            warning_note(wi, warn_buf);
+            warning_note(wi, warningMessage);
         }
         else
         {
-            printf("%s\n", warn_buf);
+            GMX_LOG(logger.info).asParagraph().appendText(warningMessage);
         }
     }
 
@@ -2345,7 +2401,7 @@ int gmx_grompp(int argc, char* argv[])
 
     if (bVerbose)
     {
-        fprintf(stderr, "writing run input file...\n");
+        GMX_LOG(logger.info).asParagraph().appendTextFormatted("writing run input file...");
     }
 
     done_warning(wi, FARGS);