Generic gmx dump implementation for KeyValueTree
authorTeemu Murtola <teemu.murtola@gmail.com>
Sat, 22 Oct 2016 03:45:48 +0000 (06:45 +0300)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 8 Feb 2017 21:32:29 +0000 (22:32 +0100)
Instead of calling IInputRecExtension, make gmx dump directly write out
the KeyValueTree structure.  This removes the need for any special code
in the modules to support gmx dump.

This does not work with the mdp output from gmx dump, but neither did
the old code.  This should be relatively easy to put back once we have a
structured mdp format.

Change-Id: Ie1420ce6e1ae14e73da1be09ce48faf1b5ed9daf

src/gromacs/applied-forces/electricfield.cpp
src/gromacs/mdlib/broadcaststructs.cpp
src/gromacs/mdrunutility/mdmodules.cpp
src/gromacs/mdrunutility/mdmodules.h
src/gromacs/mdtypes/inputrec.cpp
src/gromacs/mdtypes/inputrec.h
src/gromacs/tools/dump.cpp
src/gromacs/utility/keyvaluetree.cpp
src/gromacs/utility/keyvaluetree.h
src/gromacs/utility/textwriter.cpp

index 964ad3411a78184eb6ec23aaf0e3714fab553d9d..bcfe74357d3c2aaab27099214754ae8823f8dfe7 100644 (file)
@@ -66,7 +66,6 @@
 #include "gromacs/utility/pleasecite.h"
 #include "gromacs/utility/strconvert.h"
 #include "gromacs/utility/stringutil.h"
-#include "gromacs/utility/txtdump.h"
 
 namespace gmx
 {
@@ -169,7 +168,6 @@ class ElectricField : public IInputRecExtension, public IForceProvider
                              const IInputRecExtension *field2,
                              real                      reltol,
                              real                      abstol);
-        virtual void printParameters(FILE *fp, int indent);
         virtual void initOutput(FILE *fplog, int nfile, const t_filenm fnm[],
                                 bool bAppendFiles, const gmx_output_env_t *oenv);
         virtual void finishOutput();
@@ -406,18 +404,6 @@ void ElectricField::initForcerec(t_forcerec *fr)
     }
 }
 
-void ElectricField::printParameters(FILE *fp, int indent)
-{
-    const char *const dimension[DIM] = { "X", "Y", "Z" };
-    indent = pr_title(fp, indent, "ElectricField");
-    for (int m = 0; m < DIM; m++)
-    {
-        pr_indent(fp, indent);
-        fprintf(fp, "-%s E0 = %g omega = %g t0 = %g sigma = %g\n",
-                dimension[m], a(m), omega(m), t0(m), sigma(m));
-    }
-}
-
 void ElectricField::setFieldTerm(int dim, real a, real omega, real t0, real sigma)
 {
     range_check(dim, 0, DIM);
index 35f0e9e6da4edab3cbe17b27805cd0901119d41e..5a27ae4dc425ece145527b262c9670008423145d 100644 (file)
@@ -668,7 +668,10 @@ static void bc_inputrec(const t_commrec *cr, t_inputrec *inputrec)
     gmx::IInputRecExtension *eptr = inputrec->efield;
     block_bc(cr, *inputrec);
     inputrec->efield = eptr;
-    inputrec->params = nullptr;
+    if (!SIMMASTER(cr))
+    {
+        inputrec->params = nullptr;
+    }
 
     bc_grpopts(cr, &(inputrec->opts));
 
index c7bcd9a5be4d41c784affcf2f92c93418eca6b0e..a97be4dc9b13a4497919601f51ecff8926712c94 100644 (file)
@@ -144,4 +144,15 @@ void MDModules::assignOptionsToModulesFromTpr()
     assignOptionsFromKeyValueTree(&moduleOptions, *impl_->ir_->params, nullptr);
 }
 
+void MDModules::adjustInputrecBasedOnModules()
+{
+    gmx::Options                        options;
+    impl_->field_->initMdpOptions(&options);
+    std::unique_ptr<KeyValueTreeObject> params(impl_->ir_->params);
+    // Avoid double freeing if the next operation throws.
+    impl_->ir_->params = nullptr;
+    impl_->ir_->params = new KeyValueTreeObject(
+                gmx::adjustKeyValueTreeFromOptions(*params, options));
+}
+
 } // namespace gmx
index 916caa59d3b95c653d6ed8e1872fe685f1a3e4e3..d221556781c36710a171833d8907e1f7b4d68b63 100644 (file)
@@ -137,6 +137,14 @@ class MDModules
          */
         void assignOptionsToModulesFromTpr();
 
+        /*! \brief
+         * Normalizes inputrec parameters to match current code version.
+         *
+         * This orders the parameters in inputrec to match the current code and
+         * adds any missing defaults.
+         */
+        void adjustInputrecBasedOnModules();
+
     private:
         class Impl;
 
index c73c2fa255bb2494759988b83033e126b2ce01d5..96e91bd923ef63eb8886eaf9eaafe7606fc6a726 100644 (file)
@@ -55,6 +55,7 @@
 #include "gromacs/utility/smalloc.h"
 #include "gromacs/utility/snprintf.h"
 #include "gromacs/utility/stringutil.h"
+#include "gromacs/utility/textwriter.h"
 #include "gromacs/utility/txtdump.h"
 
 //! Macro to select a bool name
@@ -951,9 +952,6 @@ void pr_inputrec(FILE *fp, int indent, const char *title, const t_inputrec *ir,
             pr_simtempvals(fp, indent, ir->simtempvals, ir->fepvals->n_lambda);
         }
 
-        /* ELECTRIC FIELDS */
-        ir->efield->printParameters(fp, indent);
-
         /* ION/WATER SWAPPING FOR COMPUTATIONAL ELECTROPHYSIOLOGY */
         PS("swapcoords", ESWAPTYPE(ir->eSwapCoords));
         if (ir->eSwapCoords != eswapNO)
@@ -971,6 +969,13 @@ void pr_inputrec(FILE *fp, int indent, const char *title, const t_inputrec *ir,
         PR("userreal3", ir->userreal3);
         PR("userreal4", ir->userreal4);
 
+        if (!bMDPformat)
+        {
+            gmx::TextWriter writer(fp);
+            writer.wrapperSettings().setIndent(indent);
+            ir->params->writeUsing(&writer);
+        }
+
         pr_grp_opts(fp, indent, "grpopts", &(ir->opts), bMDPformat);
     }
 }
index 3a43a39eecc244addfa2bb7360ab77d3898aabe5..e74fc4aee9e12c7a3d653d15133d2d5642a137fa 100644 (file)
@@ -124,13 +124,6 @@ class IInputRecExtension
                              real                           reltol,
                              real                           abstol)  = 0;
 
-        /*! \brief Print parameters belonging to this class to a file
-         *
-         * \param[in] fp     File pointer
-         * \param[in] indent Initial indentation level for printing
-         */
-        virtual void printParameters(FILE *fp, int indent)           = 0;
-
         /*! \brief Initiate output parameters
          *
          * \param[in] fplog File pointer for log messages
index 4cc89d9bf9b270964968c49debb3c482c88e0c51..e2af9cb1f50a032e195f53bec67fdd51effdc01e 100644 (file)
@@ -74,7 +74,8 @@ static void list_tpx(const char *fn,
                      gmx_bool    bShowNumbers,
                      gmx_bool    bShowParameters,
                      const char *mdpfn,
-                     gmx_bool    bSysTop)
+                     gmx_bool    bSysTop,
+                     gmx_bool    bOriginalInputrec)
 {
     FILE         *gp;
     int           indent, i, j, **gcount, atot;
@@ -95,9 +96,9 @@ static void list_tpx(const char *fn,
                    ir,
                    &state,
                    tpx.bTop ? &mtop : nullptr);
-    if (tpx.bIr)
+    if (tpx.bIr && !bOriginalInputrec)
     {
-        mdModules.assignOptionsToModulesFromTpr();
+        mdModules.adjustInputrecBasedOnModules();
     }
 
     if (mdpfn && tpx.bIr)
@@ -629,13 +630,15 @@ int gmx_dump(int argc, char *argv[])
 
     gmx_output_env_t *oenv;
     /* Command line options */
-    static gmx_bool   bShowNumbers = TRUE;
-    static gmx_bool   bShowParams  = FALSE;
-    static gmx_bool   bSysTop      = FALSE;
-    t_pargs           pa[]         = {
+    gmx_bool          bShowNumbers      = TRUE;
+    gmx_bool          bShowParams       = FALSE;
+    gmx_bool          bSysTop           = FALSE;
+    gmx_bool          bOriginalInputrec = FALSE;
+    t_pargs           pa[]              = {
         { "-nr", FALSE, etBOOL, {&bShowNumbers}, "Show index numbers in output (leaving them out makes comparison easier, but creates a useless topology)" },
         { "-param", FALSE, etBOOL, {&bShowParams}, "Show parameters for each bonded interaction (for comparing dumps, it is useful to combine this with -nonr)" },
-        { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" }
+        { "-sys", FALSE, etBOOL, {&bSysTop}, "List the atoms and bonded interactions for the whole system instead of for each molecule type" },
+        { "-orgir", FALSE, etBOOL, {&bOriginalInputrec}, "Show input parameters from tpr as they were written by the version that produced the file, instead of how the current version reads them" }
     };
 
     if (!parse_common_args(&argc, argv, 0, NFILE, fnm, asize(pa), pa,
@@ -648,7 +651,7 @@ int gmx_dump(int argc, char *argv[])
     if (ftp2bSet(efTPR, NFILE, fnm))
     {
         list_tpx(ftp2fn(efTPR, NFILE, fnm), bShowNumbers, bShowParams,
-                 ftp2fn_null(efMDP, NFILE, fnm), bSysTop);
+                 ftp2fn_null(efMDP, NFILE, fnm), bSysTop, bOriginalInputrec);
     }
     else if (ftp2bSet(efTRX, NFILE, fnm))
     {
index a71c1dab6944718f7bf8a1e27535822d9a1f9467..5b752fd66aa5857c651efe1657de874fd6986835 100644 (file)
@@ -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,2017, 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.
 
 #include "gromacs/utility/gmxassert.h"
 #include "gromacs/utility/stringutil.h"
+#include "gromacs/utility/textwriter.h"
 
 namespace gmx
 {
 
-/********************************************************************
- * KeyValueTreePath
- */
-
 namespace
 {
 
@@ -60,8 +57,27 @@ std::vector<std::string> splitPathElements(const std::string &path)
     return splitDelimitedString(path.substr(1), '/');
 }
 
+//! Helper function to format a simple KeyValueTreeValue.
+std::string formatSingleValue(const KeyValueTreeValue &value)
+{
+    if (value.isType<float>())
+    {
+        return formatString("%g", value.cast<float>());
+    }
+    else if (value.isType<double>())
+    {
+        return formatString("%g", value.cast<double>());
+    }
+    GMX_RELEASE_ASSERT(false, "Unknown value type");
+    return std::string();
+}
+
 }   // namespace
 
+/********************************************************************
+ * KeyValueTreePath
+ */
+
 KeyValueTreePath::KeyValueTreePath(const std::string &path)
     : path_(splitPathElements(path))
 {
@@ -72,4 +88,48 @@ std::string KeyValueTreePath::toString() const
     return "/" + joinStrings(path_, "/");
 }
 
+/********************************************************************
+ * KeyValueTreeObject
+ */
+
+void KeyValueTreeObject::writeUsing(TextWriter *writer) const
+{
+    for (const auto &prop : properties())
+    {
+        const auto &value = prop.value();
+        if (value.isObject())
+        {
+            writer->writeString(prop.key());
+            writer->writeLine(":");
+            int oldIndent = writer->wrapperSettings().indent();
+            writer->wrapperSettings().setIndent(oldIndent + 2);
+            value.asObject().writeUsing(writer);
+            writer->wrapperSettings().setIndent(oldIndent);
+        }
+        else
+        {
+            int indent = writer->wrapperSettings().indent();
+            writer->writeString(formatString("%*s", -(33-indent), prop.key().c_str()));
+            writer->writeString(" = ");
+            if (value.isArray())
+            {
+                writer->writeString("[");
+                for (const auto &elem : value.asArray().values())
+                {
+                    GMX_RELEASE_ASSERT(!elem.isObject() && !elem.isArray(),
+                                       "Arrays of objects not currently implemented");
+                    writer->writeString(" ");
+                    writer->writeString(formatSingleValue(elem));
+                }
+                writer->writeString(" ]");
+            }
+            else
+            {
+                writer->writeString(formatSingleValue(value));
+            }
+            writer->writeLine();
+        }
+    }
+}
+
 } // namespace gmx
index f0ad62674e2440a53c6df1aa4fe25ca57b5823fa..88b5b29e2459dd1a3cae1f01eae2a2ccdd1d21c8 100644 (file)
@@ -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,2017, 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.
@@ -57,6 +57,7 @@ namespace gmx
 
 class KeyValueTreeArray;
 class KeyValueTreeObject;
+class TextWriter;
 
 class KeyValueTreePath
 {
@@ -180,6 +181,8 @@ class KeyValueTreeObject
             return valueMap_.at(key);
         }
 
+        void writeUsing(TextWriter *writer) const;
+
     private:
         KeyValueTreeValue &operator[](const std::string &key)
         {
index eaab84f0c52e7d779104a46d00c8e4dce47cc5cd..7701950eaf6fdc7b50a3185d5bd167298b5c7472 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2015, by the GROMACS development team, led by
+ * Copyright (c) 2015,2016,2017, 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.
@@ -105,7 +105,14 @@ class TextWriter::Impl
 
         void writeWrappedString(const std::string &str)
         {
-            writeRawString(wrapper_.wrapToString(str));
+            if (newLineCount_ > 0)
+            {
+                writeRawString(wrapper_.wrapToString(str));
+            }
+            else
+            {
+                writeRawString(str);
+            }
         }
 
         TextOutputStreamPointer stream_;