Generalize option section internals
authorTeemu Murtola <teemu.murtola@gmail.com>
Sat, 2 Jul 2016 18:14:13 +0000 (21:14 +0300)
committerMark Abraham <mark.j.abraham@gmail.com>
Wed, 17 Aug 2016 16:30:41 +0000 (18:30 +0200)
Use a similar mechanism as for options to allow multiple different
section types to be handled through addSection().  Prerequisite for
adding support for repeating sections.  Support for iterating over
option sections is not of particular focus at this point; the initial
focus is to get the API for declaring sections and options working.

Also change the approach in IOptionsContainer to remove the need for
using directives in classes implementing it, and make the build system
list the source files explicitly.

Change-Id: I3f70f6937cb6b594c5bff120d5f2f87e42841ac4

src/gromacs/options/CMakeLists.txt
src/gromacs/options/abstractsection.cpp [new file with mode: 0644]
src/gromacs/options/abstractsection.h [new file with mode: 0644]
src/gromacs/options/ioptionscontainer.h
src/gromacs/options/ioptionscontainerwithsections.h
src/gromacs/options/options-impl.h
src/gromacs/options/options.cpp
src/gromacs/options/options.h
src/gromacs/options/optionsection.h
src/gromacs/options/tests/optionsassigner.cpp

index 71c36c82f8c60efbe4c089f76bf6697630f4e50d..298af06134711dfa31a2b5234672f6ec02fad7d7 100644 (file)
 # To help us fund GROMACS development, we humbly ask that you cite
 # the research papers on the package. Check out http://www.gromacs.org.
 
-file(GLOB OPTIONS_SOURCES *.cpp)
-set(LIBGROMACS_SOURCES ${LIBGROMACS_SOURCES} ${OPTIONS_SOURCES} PARENT_SCOPE)
+gmx_add_libgromacs_sources(
+    abstractoption.cpp
+    abstractsection.cpp
+    basicoptions.cpp
+    behaviorcollection.cpp
+    filenameoption.cpp
+    filenameoptionmanager.cpp
+    options.cpp
+    optionsassigner.cpp
+    optionsvisitor.cpp
+    timeunitmanager.cpp
+    )
 
 gmx_install_headers(
     abstractoption.h
+    abstractsection.h
     basicoptions.h
     filenameoption.h
     filenameoptionmanager.h
diff --git a/src/gromacs/options/abstractsection.cpp b/src/gromacs/options/abstractsection.cpp
new file mode 100644 (file)
index 0000000..1e6289d
--- /dev/null
@@ -0,0 +1,67 @@
+/*
+ * 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.
+ */
+/*! \internal \file
+ * \brief
+ * Implements classes from abstractsection.h.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \ingroup module_options
+ */
+#include "gmxpre.h"
+
+#include "abstractsection.h"
+
+#include "options-impl.h"
+
+namespace gmx
+{
+
+IOptionsContainer &AbstractOptionSectionHandle::addGroup()
+{
+    return section_->addGroup();
+}
+
+internal::OptionSectionImpl *
+AbstractOptionSectionHandle::addSectionImpl(const AbstractOptionSection &section)
+{
+    return section_->addSectionImpl(section);
+}
+
+OptionInfo *AbstractOptionSectionHandle::addOptionImpl(const AbstractOption &settings)
+{
+    return section_->addOptionImpl(settings);
+}
+
+} // namespace gmx
diff --git a/src/gromacs/options/abstractsection.h b/src/gromacs/options/abstractsection.h
new file mode 100644 (file)
index 0000000..c3c5678
--- /dev/null
@@ -0,0 +1,135 @@
+/*
+ * 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.
+ */
+/*! \file
+ * \brief
+ * Declares base classes for declaring option sections.
+ *
+ * This header defines base classes for option section settings that are used
+ * with IOptionsContainerWithSections::addSection().  These classes implement
+ * the "named parameter" idiom for specifying section properties.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \ingroup module_options
+ */
+#ifndef GMX_OPTIONS_ABSTRACTSECTION_H
+#define GMX_OPTIONS_ABSTRACTSECTION_H
+
+#include "gromacs/options/ioptionscontainerwithsections.h"
+#include "gromacs/utility/classhelpers.h"
+
+namespace gmx
+{
+
+namespace internal
+{
+class OptionSectionImpl;
+}
+
+/*! \brief
+ * Base class for specifying option section properties.
+ *
+ * \ingroup module_options
+ */
+class AbstractOptionSection
+{
+    protected:
+        //! \cond libapi
+        //! Initializes option properties with the given name.
+        explicit AbstractOptionSection(const char *name) : name_(name) {}
+        //! \endcond
+
+    private:
+        const char *name_;
+
+        friend class internal::OptionSectionImpl;
+};
+
+/*! \brief
+ * Base class for handles to option sections.
+ *
+ * This class implements the common functionality for adding options and
+ * subsections to option sections.
+ *
+ * \ingroup module_options
+ */
+class AbstractOptionSectionHandle : public IOptionsContainerWithSections
+{
+    public:
+        // From IOptionsContainer
+        //! \copydoc IOptionsContainer::addGroup()
+        virtual IOptionsContainer &addGroup();
+
+    protected:
+        //! \cond libapi
+        //! Wraps a given section storage object.
+        explicit AbstractOptionSectionHandle(internal::OptionSectionImpl *section)
+            : section_(section)
+        {
+        }
+        //! \endcond
+
+    private:
+        // From IOptionsContainerWithSections
+        virtual internal::OptionSectionImpl *
+        addSectionImpl(const AbstractOptionSection &section);
+        // From IOptionsContainer
+        virtual OptionInfo *addOptionImpl(const AbstractOption &settings);
+
+        internal::OptionSectionImpl *section_;
+};
+
+class AbstractOptionSectionInfo
+{
+    public:
+        //! Wraps a given section storage object.
+        explicit AbstractOptionSectionInfo(internal::OptionSectionImpl *section)
+            : section_(*section)
+        {
+        }
+
+        //! Returns the wrapped section storage object.
+        internal::OptionSectionImpl       &section() { return section_; }
+        //! Returns the wrapped section storage object.
+        const internal::OptionSectionImpl &section() const { return section_; }
+
+    private:
+        internal::OptionSectionImpl &section_;
+
+        GMX_DISALLOW_COPY_AND_ASSIGN(AbstractOptionSectionInfo);
+};
+
+} // namespace gmx
+
+#endif
index 0cc29291c91d5ce5c99ce881e535272e2d6d8400..3cdaafc0210a4253255816dcf63301ff4401be8a 100644 (file)
@@ -92,18 +92,6 @@ class IOptionsContainer
          * output.
          */
         virtual IOptionsContainer &addGroup() = 0;
-        /*! \brief
-         * Adds a recognized option.
-         *
-         * \param[in] settings Option description.
-         * \returns   OptionInfo object for the created option (never NULL).
-         * \throws    APIError if invalid option settings are provided.
-         *
-         * This method provides the internal implementation, but in most cases
-         * the templated method is called from user code.
-         * See the templated method for more details.
-         */
-        virtual OptionInfo *addOption(const AbstractOption &settings) = 0;
         /*! \brief
          * Adds a recognized option.
          *
@@ -130,7 +118,7 @@ class IOptionsContainer
         typename OptionType::InfoType *addOption(const OptionType &settings)
         {
             OptionInfo *info
-                = addOption(static_cast<const AbstractOption &>(settings));
+                = addOptionImpl(static_cast<const AbstractOption &>(settings));
             GMX_ASSERT(info->isType<typename OptionType::InfoType>(),
                        "Mismatching option info type declaration and implementation");
             return info->toType<typename OptionType::InfoType>();
@@ -141,6 +129,19 @@ class IOptionsContainer
         // (no need for the virtual, but some compilers warn otherwise)
         virtual ~IOptionsContainer();
 
+        /*! \brief
+         * Adds a recognized option.
+         *
+         * \param[in] settings Option description.
+         * \returns   OptionInfo object for the created option (never NULL).
+         * \throws    APIError if invalid option settings are provided.
+         *
+         * This method provides the internal implementation, but the templated
+         * method is called from user code.  See the templated method for more
+         * details.
+         */
+        virtual OptionInfo *addOptionImpl(const AbstractOption &settings) = 0;
+
         GMX_DEFAULT_CONSTRUCTORS(IOptionsContainer);
 };
 
index 93f7791bddb0ede57887a2ac478f00b0ac96ff36..198e0caf8080887030e571665ef4310167992838 100644 (file)
 namespace gmx
 {
 
-class OptionSection;
-class OptionSectionInfo;
+class AbstractOptionSection;
+class AbstractOptionSectionHandle;
+
+namespace internal
+{
+class OptionSectionImpl;
+}
 
 /*! \brief
  * Interface for adding input options with sections.
@@ -68,14 +73,40 @@ class IOptionsContainerWithSections : public IOptionsContainer
         /*! \brief
          * Adds a section to this collection.
          *
-         * \param[in] section Section to add.
+         * \tparam    SectionType Type of the section description object.
+         * \param[in] section     Section description.
+         * \returns   AbstractOptionSectionHandle object for the created option.
+         * \throws    APIError if invalid option settings are provided.
+         *
+         * Options can be added to the section through the returned handle.
+         *
+         * \internal
+         * \p SectionType::HandleType must specify a type that derives from
+         * AbstractinOptionSectionHandle and has a suitable constructor.
          */
-        virtual IOptionsContainerWithSections &addSection(const OptionSection &section) = 0;
+        template <class SectionType>
+        typename SectionType::HandleType addSection(const SectionType &section)
+        {
+            internal::OptionSectionImpl *storage
+                = addSectionImpl(static_cast<const AbstractOptionSection &>(section));
+            return typename SectionType::HandleType(storage);
+        }
 
     protected:
         // Disallow deletion through the interface.
         // (no need for the virtual, but some compilers warn otherwise)
         virtual ~IOptionsContainerWithSections();
+
+        /*! \brief
+         * Adds a section to this container.
+         *
+         * \param[in] section     Section description.
+         * \returns   Pointer to the internal section representation object.
+         */
+        virtual internal::OptionSectionImpl *
+        addSectionImpl(const AbstractOptionSection &section) = 0;
+
+        GMX_DEFAULT_CONSTRUCTORS(IOptionsContainerWithSections);
 };
 
 } // namespace
index 064cafc54ec1d1f6e145dc9d264dea0e9fe045cb..6d47fefe304ea13ef14b6fef78d65838b4c68f97 100644 (file)
@@ -94,7 +94,7 @@ class OptionSectionImpl : public IOptionsContainerWithSections
 
                 // From IOptionsContainer
                 virtual IOptionsContainer &addGroup();
-                virtual OptionInfo *addOption(const AbstractOption &settings);
+                virtual OptionInfo *addOptionImpl(const AbstractOption &settings);
 
                 //! Containing options object.
                 OptionSectionImpl  *parent_;
@@ -126,11 +126,11 @@ class OptionSectionImpl : public IOptionsContainerWithSections
         }
 
         // From IOptionsContainerWithSections
-        virtual IOptionsContainerWithSections &addSection(const OptionSection &section);
+        virtual OptionSectionImpl *addSectionImpl(const AbstractOptionSection &section);
 
         // From IOptionsContainer
         virtual IOptionsContainer &addGroup();
-        virtual OptionInfo *addOption(const AbstractOption &settings);
+        virtual OptionInfo *addOptionImpl(const AbstractOption &settings);
 
         //! Returns section info object for this section.
         OptionSectionInfo       &info() { return info_; }
index 3240dce28ca2111c3301f405ff70444ac75fdbca..750dcecec030fd82a9c006fdcfa9c2f08a12a35a 100644 (file)
@@ -98,13 +98,14 @@ OptionsImpl::OptionsImpl()
  * OptionSectionImpl
  */
 
-IOptionsContainerWithSections &OptionSectionImpl::addSection(const OptionSection &section)
+OptionSectionImpl *
+OptionSectionImpl::addSectionImpl(const AbstractOptionSection &section)
 {
     const char *name = section.name_;
     // Make sure that there are no duplicate sections.
     GMX_RELEASE_ASSERT(findSection(name) == NULL, "Duplicate subsection name");
     subsections_.push_back(SectionPointer(new OptionSectionImpl(managers_, name)));
-    return *subsections_.back();
+    return subsections_.back().get();
 }
 
 IOptionsContainer &OptionSectionImpl::addGroup()
@@ -112,9 +113,9 @@ IOptionsContainer &OptionSectionImpl::addGroup()
     return rootGroup_.addGroup();
 }
 
-OptionInfo *OptionSectionImpl::addOption(const AbstractOption &settings)
+OptionInfo *OptionSectionImpl::addOptionImpl(const AbstractOption &settings)
 {
-    return rootGroup_.addOption(settings);
+    return rootGroup_.addOptionImpl(settings);
 }
 
 OptionSectionImpl *OptionSectionImpl::findSection(const char *name) const
@@ -189,7 +190,7 @@ IOptionsContainer &OptionSectionImpl::Group::addGroup()
     return subgroups_.back();
 }
 
-OptionInfo *OptionSectionImpl::Group::addOption(const AbstractOption &settings)
+OptionInfo *OptionSectionImpl::Group::addOptionImpl(const AbstractOption &settings)
 {
     OptionSectionImpl::AbstractOptionStoragePointer
          option(settings.createStorage(parent_->managers_));
@@ -237,9 +238,9 @@ void Options::addManager(IOptionManager *manager)
     impl_->managers_.add(manager);
 }
 
-IOptionsContainerWithSections &Options::addSection(const OptionSection &section)
+internal::OptionSectionImpl *Options::addSectionImpl(const AbstractOptionSection &section)
 {
-    return impl_->rootSection_.addSection(section);
+    return impl_->rootSection_.addSectionImpl(section);
 }
 
 IOptionsContainer &Options::addGroup()
@@ -247,9 +248,9 @@ IOptionsContainer &Options::addGroup()
     return impl_->rootSection_.addGroup();
 }
 
-OptionInfo *Options::addOption(const AbstractOption &settings)
+OptionInfo *Options::addOptionImpl(const AbstractOption &settings)
 {
-    return impl_->rootSection_.addOption(settings);
+    return impl_->rootSection_.addOptionImpl(settings);
 }
 
 OptionSectionInfo &Options::rootSection()
index 40944859fce593c955a12af3b59a97c28b438bb5..a4ba843633ac1889012a00e235958f30b8167844 100644 (file)
@@ -128,13 +128,8 @@ class Options : public IOptionsContainerWithSections
          */
         void addManager(IOptionManager *manager);
 
-        // From IOptionsContainerWithSections
-        virtual IOptionsContainerWithSections &addSection(const OptionSection &section);
-
         // From IOptionsContainer
         virtual IOptionsContainer &addGroup();
-        virtual OptionInfo *addOption(const AbstractOption &settings);
-        using IOptionsContainer::addOption;
 
         //! Returns a handle to the root section.
         OptionSectionInfo       &rootSection();
@@ -158,6 +153,12 @@ class Options : public IOptionsContainerWithSections
         void finish();
 
     private:
+        // From IOptionsContainerWithSections
+        virtual internal::OptionSectionImpl *
+        addSectionImpl(const AbstractOptionSection &section);
+        // From IOptionsContainer
+        virtual OptionInfo *addOptionImpl(const AbstractOption &settings);
+
         PrivateImplPointer<internal::OptionsImpl> impl_;
 
         //! Needed to be able to extend the interface of this object.
index 6f708ec01474f95f303774580a01778c2c06ac67..b811facd7deeb8b35fe21fc8ef0c73b704fcd01e 100644 (file)
 #ifndef GMX_OPTIONS_OPTIONSECTION_H
 #define GMX_OPTIONS_OPTIONSECTION_H
 
+#include "gromacs/options/abstractsection.h"
 #include "gromacs/utility/classhelpers.h"
 
 namespace gmx
 {
 
-namespace internal
-{
-class OptionSectionImpl;
-}
+class OptionSectionHandle;
 
-class OptionSection
+/*! \brief
+ * Declares a simple option section.
+ *
+ * This class declares a simple section that only provides structure for
+ * grouping the options, but does not otherwise influence the behavior of the
+ * contained options.
+ *
+ * \inpublicapi
+ * \ingroup module_options
+ */
+class OptionSection : public AbstractOptionSection
 {
     public:
-        explicit OptionSection(const char *name) : name_(name) {}
+        //! AbstractOptionSectionHandle corresponding to this option type.
+        typedef OptionSectionHandle HandleType;
 
-    private:
-        const char *name_;
+        //! Creates a section with the given name.
+        explicit OptionSection(const char *name) : AbstractOptionSection(name) {}
+};
 
-        friend class internal::OptionSectionImpl;
+/*! \brief
+ * Allows adding options to an OptionSection.
+ *
+ * An instance of this class is returned from
+ * IOptionsContainerWithSections::addSection(), and supports adding options and
+ * subsections to a section created with OptionSection.
+ *
+ * \inpublicapi
+ * \ingroup module_options
+ */
+class OptionSectionHandle : public AbstractOptionSectionHandle
+{
+    public:
+        //! Wraps a given section storage object.
+        explicit OptionSectionHandle(internal::OptionSectionImpl *section)
+            : AbstractOptionSectionHandle(section)
+        {
+        }
 };
 
-class OptionSectionInfo
+class OptionSectionInfo : public AbstractOptionSectionInfo
 {
     public:
         //! Wraps a given section storage object.
         explicit OptionSectionInfo(internal::OptionSectionImpl *section)
-            : section_(*section)
+            : AbstractOptionSectionInfo(section)
         {
         }
-
-        //! Returns the wrapped section storage object.
-        internal::OptionSectionImpl       &section() { return section_; }
-        //! Returns the wrapped section storage object.
-        const internal::OptionSectionImpl &section() const { return section_; }
-
-    private:
-        internal::OptionSectionImpl &section_;
-
-        GMX_DISALLOW_COPY_AND_ASSIGN(OptionSectionInfo);
 };
 
 } // namespace gmx
index a12201f92d3b993c3fbf375d6a718840fff023f2..2518a850fc593e1748e35b13075c1a10b3adcda6 100644 (file)
@@ -222,8 +222,8 @@ TEST(OptionsAssignerTest, HandlesSections)
 {
     using gmx::OptionSection;
     gmx::Options options;
-    auto        &sub1   = options.addSection(OptionSection("section1"));
-    auto        &sub2   = options.addSection(OptionSection("section2"));
+    auto         sub1   = options.addSection(OptionSection("section1"));
+    auto         sub2   = options.addSection(OptionSection("section2"));
     int          value  = 3;
     int          value1 = 1;
     int          value2 = 2;