Generic KeyValueTree broadcast
authorTeemu Murtola <teemu.murtola@gmail.com>
Sun, 23 Oct 2016 07:14:42 +0000 (10:14 +0300)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 8 Feb 2017 21:33:48 +0000 (22:33 +0100)
Broadcast the KeyValueTree as part of inputrec, instead of having
special code in each module and a method in IInputRecExtension for this
purpose.  Now the only thing the electric field code needs to do for mdp
input/inputrec is to declare the options it accepts (and how they map
from mdp options).

Change-Id: Ia40120d4dd887c533cb6fcb60535852436dbd781

src/gromacs/applied-forces/electricfield.cpp
src/gromacs/mdlib/broadcaststructs.cpp
src/gromacs/mdtypes/inputrec.h
src/gromacs/utility/inmemoryserializer.cpp [new file with mode: 0644]
src/gromacs/utility/inmemoryserializer.h [new file with mode: 0644]
src/programs/mdrun/runner.cpp

index 5ee78d809e9f2d273325f1257070e529ae883383..7dbe79da07d303623b5ebb6cbdb02f6727511980 100644 (file)
@@ -48,7 +48,6 @@
 #include "gromacs/commandline/filenm.h"
 #include "gromacs/fileio/gmxfio.h"
 #include "gromacs/fileio/xvgr.h"
-#include "gromacs/gmxlib/network.h"
 #include "gromacs/math/units.h"
 #include "gromacs/math/vec.h"
 #include "gromacs/mdtypes/commrec.h"
@@ -59,7 +58,6 @@
 #include "gromacs/options/ioptionscontainerwithsections.h"
 #include "gromacs/options/optionsection.h"
 #include "gromacs/utility/exceptions.h"
-#include "gromacs/utility/fatalerror.h"
 #include "gromacs/utility/keyvaluetreebuilder.h"
 #include "gromacs/utility/keyvaluetreetransform.h"
 #include "gromacs/utility/pleasecite.h"
@@ -162,7 +160,7 @@ class ElectricField : public IInputRecExtension, public IForceProvider
         // From IInputRecExtension
         virtual void initMdpTransform(IKeyValueTreeTransformRules *transform);
         virtual void initMdpOptions(IOptionsContainerWithSections *options);
-        virtual void broadCast(const t_commrec *cr);
+
         virtual void initOutput(FILE *fplog, int nfile, const t_filenm fnm[],
                                 bool bAppendFiles, const gmx_output_env_t *oenv);
         virtual void finishOutput();
@@ -178,18 +176,6 @@ class ElectricField : public IInputRecExtension, public IForceProvider
         //! Return whether or not to apply a field
         bool isActive() const;
 
-        /*! \brief Add a component to the electric field
-         *
-         * The electric field has three spatial dimensions that are
-         * added to the data structure one at a time.
-         * \param[in] dim   Dimension, XX, YY, ZZ (0, 1, 2)
-         * \param[in] a     Amplitude of the field in V/nm
-         * \param[in] omega Frequency (1/ps)
-         * \param[in] t0    Time of pulse peak (ps)
-         * \param[in] sigma Width of peak (ps)
-         */
-        void setFieldTerm(int dim, real a, real omega, real t0, real sigma);
-
         /*! \brief Return the field strength
          *
          * \param[in] dim The spatial direction
@@ -324,37 +310,6 @@ void ElectricField::initMdpOptions(IOptionsContainerWithSections *options)
     efield_[ZZ].initMdpOptions(&section, "z");
 }
 
-void ElectricField::broadCast(const t_commrec *cr)
-{
-    rvec a1, omega1, sigma1, t01;
-
-    if (MASTER(cr))
-    {
-        // Load the parameters read from tpr into temp vectors
-        for (int m = 0; m < DIM; m++)
-        {
-            a1[m]     = a(m);
-            omega1[m] = omega(m);
-            sigma1[m] = sigma(m);
-            t01[m]    = t0(m);
-        }
-    }
-    // Broadcasting the parameters
-    gmx_bcast(DIM*sizeof(a1[0]), a1, cr);
-    gmx_bcast(DIM*sizeof(omega1[0]), omega1, cr);
-    gmx_bcast(DIM*sizeof(t01[0]), t01, cr);
-    gmx_bcast(DIM*sizeof(sigma1[0]), sigma1, cr);
-
-    // And storing them locally
-    if (!MASTER(cr))
-    {
-        for (int m = 0; m < DIM; m++)
-        {
-            setFieldTerm(m, a1[m], omega1[m], t01[m], sigma1[m]);
-        }
-    }
-}
-
 void ElectricField::initOutput(FILE *fplog, int nfile, const t_filenm fnm[],
                                bool bAppendFiles, const gmx_output_env_t *oenv)
 {
@@ -399,12 +354,6 @@ void ElectricField::initForcerec(t_forcerec *fr)
     }
 }
 
-void ElectricField::setFieldTerm(int dim, real a, real omega, real t0, real sigma)
-{
-    range_check(dim, 0, DIM);
-    efield_[dim].setField(a, omega, t0, sigma);
-}
-
 real ElectricField::field(int dim, real t) const
 {
     return efield_[dim].evaluate(t);
index 5a27ae4dc425ece145527b262c9670008423145d..adb857578bf492fd03bb11466571c5de4014bd16 100644 (file)
@@ -53,6 +53,9 @@
 #include "gromacs/topology/symtab.h"
 #include "gromacs/topology/topology.h"
 #include "gromacs/utility/fatalerror.h"
+#include "gromacs/utility/inmemoryserializer.h"
+#include "gromacs/utility/keyvaluetree.h"
+#include "gromacs/utility/keyvaluetreeserializer.h"
 #include "gromacs/utility/smalloc.h"
 
 static void bc_cstring(const t_commrec *cr, char **s)
@@ -668,9 +671,29 @@ static void bc_inputrec(const t_commrec *cr, t_inputrec *inputrec)
     gmx::IInputRecExtension *eptr = inputrec->efield;
     block_bc(cr, *inputrec);
     inputrec->efield = eptr;
-    if (!SIMMASTER(cr))
+    if (SIMMASTER(cr))
     {
+        gmx::InMemorySerializer serializer;
+        gmx::serializeKeyValueTree(*inputrec->params, &serializer);
+        std::vector<char>       buffer = serializer.finishAndGetBuffer();
+        size_t                  size   = buffer.size();
+        block_bc(cr, size);
+        nblock_bc(cr, size, buffer.data());
+    }
+    else
+    {
+        // block_bc() above overwrites the old pointer, so set it to a
+        // reasonable value in case code below throws.
+        // cppcheck-suppress redundantAssignment
         inputrec->params = nullptr;
+        std::vector<char> buffer;
+        size_t            size;
+        block_bc(cr, size);
+        nblock_abc(cr, size, &buffer);
+        gmx::InMemoryDeserializer serializer(buffer);
+        // cppcheck-suppress redundantAssignment
+        inputrec->params = new gmx::KeyValueTreeObject(
+                    gmx::deserializeKeyValueTree(&serializer));
     }
 
     bc_grpopts(cr, &(inputrec->opts));
@@ -709,7 +732,6 @@ static void bc_inputrec(const t_commrec *cr, t_inputrec *inputrec)
         snew_bc(cr, inputrec->imd, 1);
         bc_imd(cr, inputrec->imd);
     }
-    inputrec->efield->broadCast(cr);
     if (inputrec->eSwapCoords != eswapNO)
     {
         snew_bc(cr, inputrec->swap, 1);
index a7ef04436ac58f253178c5d2ae0756f7fa2614d8..b19775d8936e5717b49aaef7f7dd7ad80f809fad 100644 (file)
 
 struct gmx_output_env_t;
 struct pull_params_t;
-struct t_commrec;
-struct t_fileio;
 struct t_filenm;
 struct t_forcerec;
-struct t_inpfile;
-struct warninp;
 
 namespace gmx
 {
@@ -105,12 +101,6 @@ class IInputRecExtension
          */
         virtual void initMdpOptions(IOptionsContainerWithSections *options) = 0;
 
-        /*! \brief Broadcast input parameters to all ranks
-         *
-         * \param[in] cr  Communication record, gromacs structure
-         */
-        virtual void broadCast(const t_commrec *cr)     = 0;
-
         /*! \brief Initiate output parameters
          *
          * \param[in] fplog File pointer for log messages
diff --git a/src/gromacs/utility/inmemoryserializer.cpp b/src/gromacs/utility/inmemoryserializer.cpp
new file mode 100644 (file)
index 0000000..fb7d4b7
--- /dev/null
@@ -0,0 +1,190 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * 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.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+#include "gmxpre.h"
+
+#include "inmemoryserializer.h"
+
+#include <algorithm>
+#include <vector>
+
+#include "gromacs/utility/gmxassert.h"
+
+namespace gmx
+{
+
+namespace
+{
+
+//! Returns offset to add to `pos` to get it aligned at `alignment` bytes.
+size_t alignedOffset(size_t pos, size_t alignment)
+{
+    return (alignment - pos % alignment) % alignment;
+}
+
+}   // namespace
+
+/********************************************************************
+ * InMemorySerializer
+ */
+
+class InMemorySerializer::Impl
+{
+    public:
+        template <typename T>
+        void doValue(T value)
+        {
+            // Here, we assume that the vector memory buffer start is aligned,
+            // similar to what malloc() guarantees.
+            const size_t size = buffer_.size();
+            const size_t pos  = size + alignedOffset(size, alignof(T));
+            buffer_.resize(pos + sizeof(T));
+            *reinterpret_cast<T *>(&buffer_[pos]) = value;
+        }
+        void doString(const std::string &value)
+        {
+            doValue<size_t>(value.size());
+            const size_t pos = buffer_.size();
+            buffer_.resize(pos + value.size());
+            std::copy(value.begin(), value.end(), &buffer_[pos]);
+        }
+
+        std::vector<char> buffer_;
+};
+
+InMemorySerializer::InMemorySerializer()
+    : impl_(new Impl)
+{
+}
+
+InMemorySerializer::~InMemorySerializer()
+{
+}
+
+std::vector<char> InMemorySerializer::finishAndGetBuffer()
+{
+    return std::move(impl_->buffer_);
+}
+
+void InMemorySerializer::doUChar(unsigned char *value)
+{
+    impl_->doValue(*value);
+}
+
+void InMemorySerializer::doInt(int *value)
+{
+    impl_->doValue(*value);
+}
+
+void InMemorySerializer::doFloat(float *value)
+{
+    impl_->doValue(*value);
+}
+
+void InMemorySerializer::doDouble(double *value)
+{
+    impl_->doValue(*value);
+}
+
+void InMemorySerializer::doString(std::string *value)
+{
+    impl_->doString(*value);
+}
+
+/********************************************************************
+ * InMemoryDeserializer
+ */
+
+class InMemoryDeserializer::Impl
+{
+    public:
+        explicit Impl(const std::vector<char> &buffer)
+            : buffer_(buffer), pos_(0)
+        {
+        }
+
+        template <typename T>
+        void doValue(T *value)
+        {
+            pos_  += alignedOffset(pos_, alignof(T));
+            *value = *reinterpret_cast<const T *>(&buffer_[pos_]);
+            pos_  += sizeof(T);
+        }
+        void doString(std::string *value)
+        {
+            size_t size;
+            doValue<size_t>(&size);
+            *value = std::string(&buffer_[pos_], &buffer_[pos_ + size]);
+            pos_  += size;
+        }
+
+        const std::vector<char> &buffer_;
+        size_t                   pos_;
+};
+
+InMemoryDeserializer::InMemoryDeserializer(const std::vector<char> &buffer)
+    : impl_(new Impl(buffer))
+{
+}
+
+InMemoryDeserializer::~InMemoryDeserializer()
+{
+}
+
+void InMemoryDeserializer::doUChar(unsigned char *value)
+{
+    impl_->doValue(value);
+}
+
+void InMemoryDeserializer::doInt(int *value)
+{
+    impl_->doValue(value);
+}
+
+void InMemoryDeserializer::doFloat(float *value)
+{
+    impl_->doValue(value);
+}
+
+void InMemoryDeserializer::doDouble(double *value)
+{
+    impl_->doValue(value);
+}
+
+void InMemoryDeserializer::doString(std::string *value)
+{
+    impl_->doString(value);
+}
+
+} // namespace gmx
diff --git a/src/gromacs/utility/inmemoryserializer.h b/src/gromacs/utility/inmemoryserializer.h
new file mode 100644 (file)
index 0000000..55d315c
--- /dev/null
@@ -0,0 +1,98 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * 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.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \libinternal \file
+ * \brief
+ * Declares gmx::ISerializer implementation for in-memory serialization.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \inlibraryapi
+ * \ingroup module_utility
+ */
+#ifndef GMX_UTILITY_INMEMORYSERIALIZER_H
+#define GMX_UTILITY_INMEMORYSERIALIZER_H
+
+#include <vector>
+
+#include "gromacs/utility/classhelpers.h"
+#include "gromacs/utility/iserializer.h"
+
+namespace gmx
+{
+
+class InMemorySerializer : public ISerializer
+{
+    public:
+        InMemorySerializer();
+        virtual ~InMemorySerializer();
+
+        std::vector<char> finishAndGetBuffer();
+
+        // From ISerializer
+        virtual bool reading() const { return false; }
+        virtual void doUChar(unsigned char *value);
+        virtual void doInt(int *value);
+        virtual void doFloat(float *value);
+        virtual void doDouble(double *value);
+        virtual void doString(std::string *value);
+
+    private:
+        class Impl;
+
+        PrivateImplPointer<Impl> impl_;
+};
+
+class InMemoryDeserializer : public ISerializer
+{
+    public:
+        explicit InMemoryDeserializer(const std::vector<char> &buffer);
+        virtual ~InMemoryDeserializer();
+
+        // From ISerializer
+        virtual bool reading() const { return true; }
+        virtual void doUChar(unsigned char *value);
+        virtual void doInt(int *value);
+        virtual void doFloat(float *value);
+        virtual void doDouble(double *value);
+        virtual void doString(std::string *value);
+
+    private:
+        class Impl;
+
+        PrivateImplPointer<Impl> impl_;
+};
+
+} // namespace gmx
+
+#endif
index d17c5f9f4a16e4a661bff574e1ba9aba605bad16..78ab6dc34abb9728024d44d24e668aaca0b18886 100644 (file)
@@ -779,7 +779,6 @@ int mdrunner(gmx_hw_opt_t *hw_opt,
     {
         /* Read (nearly) all data required for the simulation */
         read_tpx_state(ftp2fn(efTPR, nfile, fnm), inputrec, state, mtop);
-        mdModules.assignOptionsToModulesFromTpr();
 
         if (inputrec->cutoff_scheme == ecutsVERLET)
         {
@@ -897,6 +896,7 @@ int mdrunner(gmx_hw_opt_t *hw_opt,
          */
         gmx_bcast_sim(sizeof(bUseGPU), &bUseGPU, cr);
     }
+    mdModules.assignOptionsToModulesFromTpr();
 
     if (fplog != nullptr)
     {