Basic support for transforming KeyValueTrees
authorTeemu Murtola <teemu.murtola@gmail.com>
Tue, 26 Jul 2016 13:05:13 +0000 (16:05 +0300)
committerMark Abraham <mark.j.abraham@gmail.com>
Tue, 30 Aug 2016 14:41:16 +0000 (16:41 +0200)
Add support for writing transforms from one KeyValueTree to another.
For now, only a small subset of simple transforms is supported, but this
should be mostly sufficient for a proof-of-concept for MDP processing on
the electric field code, and new cases can be relatively easily added.

These transforms can support two use cases: specifying how a legacy MDP
file maps to a more structured format, and for supporting changes from
version to version in the structured format.  However, versioning for
the structured format (in particular, if/when it will be necessary to
recognize the version from input) likely still needs some work.

There is not much in terms of consistency checks, error handling, or
documentation, but that is better to add once using code exists.

Change-Id: Ie2570eee8a62c2d64b79e00925bbf1e40b62b261

src/gromacs/utility/keyvaluetree.h
src/gromacs/utility/keyvaluetreebuilder.h
src/gromacs/utility/keyvaluetreetransform.cpp [new file with mode: 0644]
src/gromacs/utility/keyvaluetreetransform.h [new file with mode: 0644]
src/gromacs/utility/tests/CMakeLists.txt
src/gromacs/utility/tests/keyvaluetreetransform.cpp [new file with mode: 0644]
src/gromacs/utility/tests/refdata/TreeValueTransformTest_ObjectFromMultipleStrings.xml [new file with mode: 0644]
src/gromacs/utility/tests/refdata/TreeValueTransformTest_ObjectFromString.xml [new file with mode: 0644]
src/gromacs/utility/tests/refdata/TreeValueTransformTest_SimpleTransforms.xml [new file with mode: 0644]
src/gromacs/utility/tests/refdata/TreeValueTransformTest_SimpleTransformsToObject.xml [new file with mode: 0644]

index 6c4cbbd0eba70c84d688ef6f9824765880ed0c3e..6dd0c0af4b73b1dbcd23e59b57753fbeab6c33d4 100644 (file)
@@ -79,6 +79,7 @@ class KeyValueTreeValue
 
         friend class KeyValueTreeBuilder;
         friend class KeyValueTreeObjectBuilder;
+        friend class KeyValueTreeValueBuilder;
 };
 
 class KeyValueTreeArray
@@ -139,7 +140,20 @@ class KeyValueTreeObject
 
         const std::vector<KeyValueTreeProperty> &properties() const { return values_; }
 
+        bool keyExists(const std::string &key) const
+        {
+            return valueMap_.find(key) != valueMap_.end();
+        }
+        const KeyValueTreeValue &operator[](const std::string &key) const
+        {
+            return valueMap_.at(key);
+        }
+
     private:
+        KeyValueTreeValue &operator[](const std::string &key)
+        {
+            return valueMap_.at(key);
+        }
         std::map<std::string, KeyValueTreeValue>::iterator
         addProperty(const std::string &key, KeyValueTreeValue &&value)
         {
index 1a01dc04e3d0b8d1c3a3a7a104fc9aabb63be5a9..94cd6606b731145f01c76eeb50e98785e828d4c8 100644 (file)
@@ -82,6 +82,26 @@ class KeyValueTreeBuilder
         friend class KeyValueTreeUniformArrayBuilder;
 };
 
+class KeyValueTreeValueBuilder
+{
+    public:
+        template <typename T>
+        void setValue(const T &value)
+        {
+            value_ = Variant::create<T>(value);
+        }
+        void setVariantValue(Variant &&value)
+        {
+            value_ = std::move(value);
+        }
+        KeyValueTreeObjectBuilder createObject();
+
+        KeyValueTreeValue build() { return KeyValueTreeValue(std::move(value_)); }
+
+    private:
+        Variant value_;
+};
+
 class KeyValueTreeArrayBuilderBase
 {
     protected:
@@ -135,10 +155,14 @@ class KeyValueTreeObjectArrayBuilder : public KeyValueTreeArrayBuilderBase
 class KeyValueTreeObjectBuilder
 {
     public:
+        void addRawValue(const std::string &key, KeyValueTreeValue &&value)
+        {
+            object_->addProperty(key, std::move(value));
+        }
         template <typename T>
         void addValue(const std::string &key, const T &value)
         {
-            object_->addProperty(key, KeyValueTreeBuilder::createValue<T>(value));
+            addRawValue(key, KeyValueTreeBuilder::createValue<T>(value));
         }
         KeyValueTreeObjectBuilder addObject(const std::string &key)
         {
@@ -156,6 +180,20 @@ class KeyValueTreeObjectBuilder
             auto iter = object_->addProperty(key, KeyValueTreeBuilder::createValue<KeyValueTreeArray>());
             return KeyValueTreeObjectArrayBuilder(&iter->second.asArray());
         }
+        void mergeObject(KeyValueTreeValue &&value)
+        {
+            KeyValueTreeObject &obj = value.asObject();
+            for (auto &prop : obj.valueMap_)
+            {
+                addRawValue(prop.first, std::move(prop.second));
+            }
+        }
+
+        bool keyExists(const std::string &key) const { return object_->keyExists(key); }
+        KeyValueTreeObjectBuilder getObject(const std::string &key) const
+        {
+            return KeyValueTreeObjectBuilder(&(*object_)[key].asObject());
+        }
 
     private:
         explicit KeyValueTreeObjectBuilder(KeyValueTreeObject *object)
@@ -170,6 +208,7 @@ class KeyValueTreeObjectBuilder
         KeyValueTreeObject *object_;
 
         friend class KeyValueTreeBuilder;
+        friend class KeyValueTreeValueBuilder;
         friend class KeyValueTreeObjectArrayBuilder;
 };
 
@@ -182,6 +221,12 @@ inline KeyValueTreeObjectBuilder KeyValueTreeBuilder::rootObject()
     return KeyValueTreeObjectBuilder(&root_);
 }
 
+inline KeyValueTreeObjectBuilder KeyValueTreeValueBuilder::createObject()
+{
+    value_ = Variant::create<KeyValueTreeObject>(KeyValueTreeObject());
+    return KeyValueTreeObjectBuilder(&value_.castRef<KeyValueTreeObject>());
+}
+
 inline KeyValueTreeObjectBuilder KeyValueTreeObjectArrayBuilder::addObject()
 {
     auto &value = addRawValue(KeyValueTreeBuilder::createValue<KeyValueTreeObject>());
diff --git a/src/gromacs/utility/keyvaluetreetransform.cpp b/src/gromacs/utility/keyvaluetreetransform.cpp
new file mode 100644 (file)
index 0000000..00eb73c
--- /dev/null
@@ -0,0 +1,283 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2016, 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 "keyvaluetreetransform.h"
+
+#include <exception>
+#include <vector>
+
+#include "gromacs/utility/keyvaluetreebuilder.h"
+#include "gromacs/utility/stringutil.h"
+
+namespace gmx
+{
+
+/********************************************************************
+ * IKeyValueTreeTransformRules
+ */
+
+IKeyValueTreeTransformRules::~IKeyValueTreeTransformRules()
+{
+}
+
+/********************************************************************
+ * KeyValueTreeTransformRule
+ */
+
+namespace internal
+{
+
+class KeyValueTreeTransformRule
+{
+    public:
+
+    private:
+        std::function<KeyValueTreeValue(const KeyValueTreeValue &)> transform_;
+};
+
+/********************************************************************
+ * KeyValueTreeTransformer::Impl
+ */
+
+class KeyValueTreeTransformerImpl : public IKeyValueTreeTransformRules
+{
+    public:
+        class Rule
+        {
+            public:
+                typedef std::function<void(KeyValueTreeValueBuilder *, const KeyValueTreeValue &)>
+                    TransformFunction;
+                void doTransform(KeyValueTreeBuilder     *builder,
+                                 const KeyValueTreeValue &value) const;
+                void doChildTransforms(KeyValueTreeBuilder      *builder,
+                                       const KeyValueTreeObject &object) const;
+                void applyTransformedValue(KeyValueTreeBuilder  *builder,
+                                           KeyValueTreeValue   &&value) const;
+
+                const Rule *findMatchingChildRule(const std::string &key) const
+                {
+                    auto iter = childRules_.find(key);
+                    if (iter == childRules_.end())
+                    {
+                        return nullptr;
+                    }
+                    return &iter->second;
+                }
+                Rule *getOrCreateChildRule(const std::string &key)
+                {
+                    auto iter = childRules_.find(key);
+                    if (iter == childRules_.end())
+                    {
+                        iter = childRules_.insert(std::make_pair(key, Rule())).first;
+                    }
+                    return &iter->second;
+                }
+
+                std::vector<std::string>    targetPath_;
+                std::string                 targetKey_;
+                TransformFunction           transform_;
+                std::map<std::string, Rule> childRules_;
+        };
+
+        virtual KeyValueTreeTransformRuleBuilder addRule()
+        {
+            return KeyValueTreeTransformRuleBuilder(this);
+        }
+
+        Rule  rootRule_;
+};
+
+void KeyValueTreeTransformerImpl::Rule::doTransform(
+        KeyValueTreeBuilder *builder, const KeyValueTreeValue &value) const
+{
+    if (transform_)
+    {
+        KeyValueTreeValueBuilder valueBuilder;
+        transform_(&valueBuilder, value);
+        applyTransformedValue(builder, valueBuilder.build());
+        return;
+    }
+    if (!childRules_.empty())
+    {
+        doChildTransforms(builder, value.asObject());
+    }
+}
+
+void KeyValueTreeTransformerImpl::Rule::doChildTransforms(
+        KeyValueTreeBuilder *builder, const KeyValueTreeObject &object) const
+{
+    for (const auto &prop : object.properties())
+    {
+        const Rule *childRule = findMatchingChildRule(prop.key());
+        if (childRule != nullptr)
+        {
+            childRule->doTransform(builder, prop.value());
+        }
+    }
+}
+
+void KeyValueTreeTransformerImpl::Rule::applyTransformedValue(
+        KeyValueTreeBuilder *builder, KeyValueTreeValue &&value) const
+{
+    KeyValueTreeObjectBuilder objBuilder = builder->rootObject();
+    for (const std::string &key : targetPath_)
+    {
+        if (objBuilder.keyExists(key))
+        {
+            objBuilder = objBuilder.getObject(key);
+        }
+        else
+        {
+            objBuilder = objBuilder.addObject(key);
+        }
+    }
+    if (objBuilder.keyExists(targetKey_))
+    {
+        objBuilder.getObject(targetKey_).mergeObject(std::move(value));
+    }
+    else
+    {
+        objBuilder.addRawValue(targetKey_, std::move(value));
+    }
+}
+
+}   // namespace internal
+
+/********************************************************************
+ * KeyValueTreeTransformer
+ */
+
+KeyValueTreeTransformer::KeyValueTreeTransformer()
+    : impl_(new internal::KeyValueTreeTransformerImpl)
+{
+}
+
+KeyValueTreeTransformer::~KeyValueTreeTransformer()
+{
+}
+
+IKeyValueTreeTransformRules *KeyValueTreeTransformer::rules()
+{
+    return impl_.get();
+}
+
+KeyValueTreeObject KeyValueTreeTransformer::transform(const KeyValueTreeObject &tree) const
+{
+    gmx::KeyValueTreeBuilder builder;
+    impl_->rootRule_.doChildTransforms(&builder, tree);
+    return builder.build();
+}
+
+/********************************************************************
+ * KeyValueTreeTransformRuleBuilder::Data
+ */
+
+class KeyValueTreeTransformRuleBuilder::Data
+{
+    public:
+        typedef internal::KeyValueTreeTransformerImpl::Rule::TransformFunction
+            TransformFunction;
+
+        void createRule(internal::KeyValueTreeTransformerImpl *impl)
+        {
+            std::vector<std::string>                     from = splitDelimitedString(fromPath_.substr(1), '/');
+            std::vector<std::string>                     to   = splitDelimitedString(toPath_.substr(1), '/');
+            internal::KeyValueTreeTransformerImpl::Rule *rule = &impl->rootRule_;
+            for (const std::string &key : from)
+            {
+                rule = rule->getOrCreateChildRule(key);
+            }
+            rule->targetKey_  = to.back();
+            to.pop_back();
+            rule->targetPath_ = std::move(to);
+            rule->transform_  = transform_;
+        }
+
+        std::string       fromPath_;
+        std::string       toPath_;
+        TransformFunction transform_;
+};
+
+/********************************************************************
+ * KeyValueTreeTransformRuleBuilder
+ */
+
+KeyValueTreeTransformRuleBuilder::KeyValueTreeTransformRuleBuilder(
+        internal::KeyValueTreeTransformerImpl *impl)
+    : impl_(impl), data_(new Data)
+{
+}
+
+KeyValueTreeTransformRuleBuilder::~KeyValueTreeTransformRuleBuilder()
+{
+    if (!std::uncaught_exception())
+    {
+        data_->createRule(impl_);
+    }
+}
+
+void KeyValueTreeTransformRuleBuilder::setFromPath(const std::string &path)
+{
+    data_->fromPath_ = path;
+}
+
+void KeyValueTreeTransformRuleBuilder::setToPath(const std::string &path)
+{
+    data_->toPath_ = path;
+}
+
+void KeyValueTreeTransformRuleBuilder::addTransformToVariant(
+        std::function<Variant(const Variant &)> transform)
+{
+    data_->transform_ =
+        [transform] (KeyValueTreeValueBuilder *builder, const KeyValueTreeValue &value)
+        {
+            builder->setVariantValue(transform(value.asVariant()));
+        };
+}
+
+void KeyValueTreeTransformRuleBuilder::addTransformToObject(
+        std::function<void(KeyValueTreeObjectBuilder *, const Variant &)> transform)
+{
+    data_->transform_ =
+        [transform] (KeyValueTreeValueBuilder *builder, const KeyValueTreeValue &value)
+        {
+            KeyValueTreeObjectBuilder obj = builder->createObject();
+            transform(&obj, value.asVariant());
+        };
+}
+
+} // namespace gmx
diff --git a/src/gromacs/utility/keyvaluetreetransform.h b/src/gromacs/utility/keyvaluetreetransform.h
new file mode 100644 (file)
index 0000000..102df2b
--- /dev/null
@@ -0,0 +1,189 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2016, 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 utilities for transforming key-value trees.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \inlibraryapi
+ * \ingroup module_utility
+ */
+#ifndef GMX_UTILITY_KEYVALUETREETRANSFORM_H
+#define GMX_UTILITY_KEYVALUETREETRANSFORM_H
+
+#include <functional>
+#include <string>
+
+#include "gromacs/utility/classhelpers.h"
+#include "gromacs/utility/variant.h"
+
+namespace gmx
+{
+
+class KeyValueTreeObject;
+class KeyValueTreeObjectBuilder;
+
+class KeyValueTreeTransformRuleBuilder;
+
+namespace internal
+{
+class KeyValueTreeTransformerImpl;
+}
+
+class IKeyValueTreeTransformRules
+{
+    public:
+        virtual KeyValueTreeTransformRuleBuilder addRule() = 0;
+
+    protected:
+        ~IKeyValueTreeTransformRules();
+};
+
+class KeyValueTreeTransformer
+{
+    public:
+        KeyValueTreeTransformer();
+        ~KeyValueTreeTransformer();
+
+        IKeyValueTreeTransformRules *rules();
+
+        KeyValueTreeObject transform(const KeyValueTreeObject &tree) const;
+
+    private:
+        PrivateImplPointer<internal::KeyValueTreeTransformerImpl> impl_;
+};
+
+class KeyValueTreeTransformRuleBuilder
+{
+    public:
+        class Base
+        {
+            protected:
+                explicit Base(KeyValueTreeTransformRuleBuilder *builder)
+                    : builder_(builder)
+                {
+                }
+
+                KeyValueTreeTransformRuleBuilder *builder_;
+        };
+
+        template <typename FromType, typename ToType>
+        class ToValue : public Base
+        {
+            public:
+                explicit ToValue(KeyValueTreeTransformRuleBuilder *builder)
+                    : Base(builder)
+                {
+                }
+
+                void transformWith(std::function<ToType(const FromType &)> transform)
+                {
+                    builder_->addTransformToVariant(
+                            [transform] (const Variant &value)
+                            {
+                                return Variant::create<ToType>(transform(value.cast<FromType>()));
+                            });
+                }
+        };
+
+        template <typename FromType>
+        class ToObject : public Base
+        {
+            public:
+                explicit ToObject(KeyValueTreeTransformRuleBuilder *builder)
+                    : Base(builder)
+                {
+                }
+
+                void transformWith(std::function<void(KeyValueTreeObjectBuilder *, const FromType &)> transform)
+                {
+                    builder_->addTransformToObject(
+                            [transform] (KeyValueTreeObjectBuilder *builder, const Variant &value)
+                            {
+                                transform(builder, value.cast<FromType>());
+                            });
+                }
+        };
+
+        template <typename FromType>
+        class AfterFrom : public Base
+        {
+            public:
+                explicit AfterFrom(KeyValueTreeTransformRuleBuilder *builder)
+                    : Base(builder)
+                {
+                }
+
+                template <typename ToType>
+                ToValue<FromType, ToType> to(const std::string &path)
+                {
+                    builder_->setToPath(path);
+                    return ToValue<FromType, ToType>(builder_);
+                }
+
+                ToObject<FromType> toObject(const std::string &path)
+                {
+                    builder_->setToPath(path);
+                    return ToObject<FromType>(builder_);
+                }
+        };
+
+        explicit KeyValueTreeTransformRuleBuilder(internal::KeyValueTreeTransformerImpl *impl);
+        KeyValueTreeTransformRuleBuilder(KeyValueTreeTransformRuleBuilder &&)            = default;
+        KeyValueTreeTransformRuleBuilder &operator=(KeyValueTreeTransformRuleBuilder &&) = default;
+        ~KeyValueTreeTransformRuleBuilder();
+
+        template <typename FromType>
+        AfterFrom<FromType> from(const std::string &path)
+        {
+            setFromPath(path);
+            return AfterFrom<FromType>(this);
+        }
+
+    private:
+        void setFromPath(const std::string &path);
+        void setToPath(const std::string &path);
+        void addTransformToVariant(std::function<Variant(const Variant &)> transform);
+        void addTransformToObject(std::function<void(KeyValueTreeObjectBuilder *, const Variant &)> transform);
+
+        class Data;
+
+        internal::KeyValueTreeTransformerImpl *impl_;
+        std::unique_ptr<Data>                  data_;
+};
+
+} // namespace gmx
+
+#endif
index 6cb80ec965842cec82943e94f2e71ed847994958..dd562a6bb295d6e746ed589f3d389ea861d46e14 100644 (file)
@@ -37,6 +37,7 @@ gmx_add_unit_test(UtilityUnitTests utility-test
                   arrayref.cpp
                   basedefinitions.cpp
                   bitmask32.cpp bitmask64.cpp bitmask128.cpp
+                  keyvaluetreetransform.cpp
                   logger.cpp
                   path.cpp
                   stringutil.cpp
diff --git a/src/gromacs/utility/tests/keyvaluetreetransform.cpp b/src/gromacs/utility/tests/keyvaluetreetransform.cpp
new file mode 100644 (file)
index 0000000..11834d8
--- /dev/null
@@ -0,0 +1,143 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2016, 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 "gromacs/utility/keyvaluetreetransform.h"
+
+#include <string>
+#include <vector>
+
+#include <gtest/gtest.h>
+
+#include "gromacs/utility/keyvaluetree.h"
+#include "gromacs/utility/keyvaluetreebuilder.h"
+#include "gromacs/utility/strconvert.h"
+#include "gromacs/utility/stringutil.h"
+
+#include "testutils/refdata.h"
+
+namespace
+{
+
+class TreeValueTransformTest : public ::testing::Test
+{
+    public:
+        void testTransform(const gmx::KeyValueTreeObject      &input,
+                           const gmx::KeyValueTreeTransformer &transform)
+        {
+            gmx::KeyValueTreeObject         result = transform.transform(input);
+
+            gmx::test::TestReferenceData    data;
+            gmx::test::TestReferenceChecker checker(data.rootChecker());
+            checker.checkKeyValueTreeObject(input, "Input");
+            checker.checkKeyValueTreeObject(result, "Tree");
+        }
+};
+
+TEST_F(TreeValueTransformTest, SimpleTransforms)
+{
+    gmx::KeyValueTreeBuilder     builder;
+    builder.rootObject().addValue<std::string>("a", "1");
+    builder.rootObject().addValue<std::string>("b", "2");
+    gmx::KeyValueTreeObject      input = builder.build();
+
+    gmx::KeyValueTreeTransformer transform;
+    transform.rules()->addRule()
+        .from<std::string>("/a").to<int>("/i").transformWith(&gmx::fromStdString<int>);
+    transform.rules()->addRule()
+        .from<std::string>("/b").to<int>("/j").transformWith(&gmx::fromStdString<int>);
+
+    testTransform(input, transform);
+}
+
+TEST_F(TreeValueTransformTest, SimpleTransformsToObject)
+{
+    gmx::KeyValueTreeBuilder     builder;
+    builder.rootObject().addValue<std::string>("a", "1");
+    builder.rootObject().addValue<std::string>("b", "2");
+    gmx::KeyValueTreeObject      input = builder.build();
+
+    gmx::KeyValueTreeTransformer transform;
+    transform.rules()->addRule()
+        .from<std::string>("/a").to<int>("/foo/i").transformWith(&gmx::fromStdString<int>);
+    transform.rules()->addRule()
+        .from<std::string>("/b").to<int>("/foo/j").transformWith(&gmx::fromStdString<int>);
+
+    testTransform(input, transform);
+}
+
+
+TEST_F(TreeValueTransformTest, ObjectFromString)
+{
+    gmx::KeyValueTreeBuilder     builder;
+    builder.rootObject().addValue<std::string>("a", "1 2");
+    gmx::KeyValueTreeObject      input = builder.build();
+
+    gmx::KeyValueTreeTransformer transform;
+    transform.rules()->addRule()
+        .from<std::string>("/a").toObject("/foo").transformWith(
+            [] (gmx::KeyValueTreeObjectBuilder *builder, const std::string &value)
+            {
+                std::vector<std::string> values = gmx::splitString(value);
+                builder->addValue<int>("a", gmx::fromString<int>(values[0]));
+                builder->addValue<int>("b", gmx::fromString<int>(values[1]));
+            });
+
+    testTransform(input, transform);
+}
+
+TEST_F(TreeValueTransformTest, ObjectFromMultipleStrings)
+{
+    gmx::KeyValueTreeBuilder     builder;
+    builder.rootObject().addValue<std::string>("a", "1");
+    builder.rootObject().addValue<std::string>("b", "2 3");
+    gmx::KeyValueTreeObject      input = builder.build();
+
+    gmx::KeyValueTreeTransformer transform;
+    transform.rules()->addRule()
+        .from<std::string>("/a").to<int>("/foo/a").transformWith(&gmx::fromStdString<int>);
+    transform.rules()->addRule()
+        .from<std::string>("/b").toObject("/foo").transformWith(
+            [] (gmx::KeyValueTreeObjectBuilder *builder, const std::string &value)
+            {
+                std::vector<std::string> values = gmx::splitString(value);
+                builder->addValue<int>("b", gmx::fromString<int>(values[0]));
+                builder->addValue<int>("c", gmx::fromString<int>(values[1]));
+            });
+
+    testTransform(input, transform);
+}
+
+} // namespace
diff --git a/src/gromacs/utility/tests/refdata/TreeValueTransformTest_ObjectFromMultipleStrings.xml b/src/gromacs/utility/tests/refdata/TreeValueTransformTest_ObjectFromMultipleStrings.xml
new file mode 100644 (file)
index 0000000..42e41f4
--- /dev/null
@@ -0,0 +1,15 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <Object Name="Input">
+    <String Name="a">1</String>
+    <String Name="b">2 3</String>
+  </Object>
+  <Object Name="Tree">
+    <Object Name="foo">
+      <Int Name="a">1</Int>
+      <Int Name="b">2</Int>
+      <Int Name="c">3</Int>
+    </Object>
+  </Object>
+</ReferenceData>
diff --git a/src/gromacs/utility/tests/refdata/TreeValueTransformTest_ObjectFromString.xml b/src/gromacs/utility/tests/refdata/TreeValueTransformTest_ObjectFromString.xml
new file mode 100644 (file)
index 0000000..a47a380
--- /dev/null
@@ -0,0 +1,13 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <Object Name="Input">
+    <String Name="a">1 2</String>
+  </Object>
+  <Object Name="Tree">
+    <Object Name="foo">
+      <Int Name="a">1</Int>
+      <Int Name="b">2</Int>
+    </Object>
+  </Object>
+</ReferenceData>
diff --git a/src/gromacs/utility/tests/refdata/TreeValueTransformTest_SimpleTransforms.xml b/src/gromacs/utility/tests/refdata/TreeValueTransformTest_SimpleTransforms.xml
new file mode 100644 (file)
index 0000000..eb275c7
--- /dev/null
@@ -0,0 +1,12 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <Object Name="Input">
+    <String Name="a">1</String>
+    <String Name="b">2</String>
+  </Object>
+  <Object Name="Tree">
+    <Int Name="i">1</Int>
+    <Int Name="j">2</Int>
+  </Object>
+</ReferenceData>
diff --git a/src/gromacs/utility/tests/refdata/TreeValueTransformTest_SimpleTransformsToObject.xml b/src/gromacs/utility/tests/refdata/TreeValueTransformTest_SimpleTransformsToObject.xml
new file mode 100644 (file)
index 0000000..72e2c38
--- /dev/null
@@ -0,0 +1,14 @@
+<?xml version="1.0"?>
+<?xml-stylesheet type="text/xsl" href="referencedata.xsl"?>
+<ReferenceData>
+  <Object Name="Input">
+    <String Name="a">1</String>
+    <String Name="b">2</String>
+  </Object>
+  <Object Name="Tree">
+    <Object Name="foo">
+      <Int Name="i">1</Int>
+      <Int Name="j">2</Int>
+    </Object>
+  </Object>
+</ReferenceData>