Removed time unit from OptionsGlobalProperties.
authorTeemu Murtola <teemu.murtola@gmail.com>
Wed, 8 Feb 2012 17:34:57 +0000 (19:34 +0200)
committerTeemu Murtola <teemu.murtola@gmail.com>
Mon, 20 Feb 2012 18:26:43 +0000 (20:26 +0200)
Replaced the shared time unit value with a separate TimeUnitManager
class that provides basic functionality for time unit conversions, can
add an option for setting the time unit, and can scale time options
according to a time unit. Control flow should be easier to follow now
that the time options don't implicitly interact.

Part of issue #839.

Change-Id: Ib89459fae744a17e7d7c98a4d7e878b0c2dcf1ab

19 files changed:
src/gromacs/options/CMakeLists.txt
src/gromacs/options/abstractoption.cpp
src/gromacs/options/abstractoption.h
src/gromacs/options/basicoptioninfo.h
src/gromacs/options/basicoptions.cpp
src/gromacs/options/basicoptions.h
src/gromacs/options/basicoptionstorage.h
src/gromacs/options/globalproperties.cpp
src/gromacs/options/globalproperties.h
src/gromacs/options/optionflags.h
src/gromacs/options/optionstoragetemplate.h
src/gromacs/options/tests/CMakeLists.txt
src/gromacs/options/tests/timeunitmanager.cpp [new file with mode: 0644]
src/gromacs/options/timeunitmanager.cpp [new file with mode: 0644]
src/gromacs/options/timeunitmanager.h [new file with mode: 0644]
src/gromacs/trajectoryanalysis/analysissettings-impl.h
src/gromacs/trajectoryanalysis/cmdlinerunner.cpp
src/gromacs/trajectoryanalysis/runnercommon.cpp
src/gromacs/trajectoryanalysis/runnercommon.h

index 93bfef8a5ba1e0f20731386826fc14cf07f4ffc1..2242af8cbdaa48e197511165437fb30fe7ad9c78 100644 (file)
@@ -8,7 +8,8 @@ set(OPTIONS_PUBLIC_HEADERS
     optionfiletype.h
     optionflags.h
     optioninfo.h
-    options.h)
+    options.h
+    timeunitmanager.h)
 install(FILES ${OPTIONS_PUBLIC_HEADERS}
         DESTINATION ${INCL_INSTALL_DIR}/gromacs/options
         COMPONENT development)
index 627e9838fbea43440557c011e58fa78c63e14e23..2c30e741cb3aadb654923881a6ff5ca30b1ac7c5 100644 (file)
@@ -90,7 +90,7 @@ bool AbstractOptionStorage::isBoolean() const
 
 void AbstractOptionStorage::startSource()
 {
-    setFlag(efHasDefaultValue);
+    setFlag(efClearOnNextSet);
 }
 
 void AbstractOptionStorage::startSet()
@@ -99,7 +99,7 @@ void AbstractOptionStorage::startSet()
     // The last condition takes care of the situation where multiple
     // sources are used, and a later source should be able to reassign
     // the value even though the option is already set.
-    if (isSet() && !hasFlag(efMulti) && !hasFlag(efHasDefaultValue))
+    if (isSet() && !hasFlag(efMulti) && !hasFlag(efClearOnNextSet))
     {
         GMX_THROW(InvalidInputError("Option specified multiple times"));
     }
index 89c273bf5243038a299234ea09f156e36d319413..58aaaeeea54eebfe859b9a5026e917988b14e114 100644 (file)
@@ -301,7 +301,7 @@ class OptionTemplate : public AbstractOption
          * Values are added to the vector after each successful set of values
          * is parsed.  Note that for some options, the value may be changed
          * later, and is only guaranteed to be correct after Options::finish()
-         * has been called (see, e.g., DoubleOption::timeValue()).
+         * has been called.
          *
          * The pointer provided should remain valid as long as the associated
          * Options object exists.
index 5f6d6e9ed5b146f6266cb991242d215770c0f0e6..6ff9c7d4a521e894c8f8e5b86940dfb78195598f 100644 (file)
@@ -85,6 +85,23 @@ class DoubleOptionInfo : public OptionInfo
 {
     public:
         explicit DoubleOptionInfo(DoubleOptionStorage *option);
+
+        //! Whether the option provides a time value.
+        bool isTime() const;
+
+        /*! \brief
+         * Sets a scale factor for user-provided values.
+         *
+         * Any user-provided value is scaled by the provided factor.
+         * Programmatically set default values are not scaled.
+         * If called multiple times, later calls override the previously set
+         * value.  In other words, the scaling is not cumulative.
+         */
+        void setScaleFactor(double factor);
+
+    private:
+        DoubleOptionStorage &option();
+        const DoubleOptionStorage &option() const;
 };
 
 /*! \brief
index 8a8d00066796be0fc3abcecdd05976eaf9aba15e..7e6f33d8287e298c5a24bb38fe339f75a37c3110 100644 (file)
@@ -168,28 +168,18 @@ AbstractOptionStorage *IntegerOption::createDefaultStorage(Options *options) con
  */
 
 DoubleOptionStorage::DoubleOptionStorage(const DoubleOption &settings, Options *options)
-    : MyBase(settings, options), _info(this), _bTime(settings._bTime)
+    : MyBase(settings, options), info_(this), bTime_(settings._bTime), factor_(1.0)
 {
-    if (_bTime)
-    {
-        options->globalProperties().request(eogpTimeScaleFactor);
-    }
 }
 
 const char *DoubleOptionStorage::typeString() const
 {
-    return hasFlag(efVector) ? "vector" : (_bTime ? "time" : "double");
+    return hasFlag(efVector) ? "vector" : (isTime() ? "time" : "double");
 }
 
 std::string DoubleOptionStorage::formatValue(int i) const
 {
-    double value = values()[i];
-    if (_bTime)
-    {
-        double factor = hostOptions().globalProperties().timeScaleFactor();
-        value /= factor;
-    }
-    return formatString("%g", value);
+    return formatString("%g", values()[i] / factor_);
 }
 
 void DoubleOptionStorage::convertValue(const std::string &value)
@@ -201,7 +191,7 @@ void DoubleOptionStorage::convertValue(const std::string &value)
     {
         GMX_THROW(InvalidInputError("Invalid value: " + value));
     }
-    addValue(dval);
+    addValue(dval * factor_);
 }
 
 void DoubleOptionStorage::processSetValues(ValueList *values)
@@ -214,16 +204,22 @@ void DoubleOptionStorage::processSetValues(ValueList *values)
 
 void DoubleOptionStorage::processAll()
 {
-    if (_bTime)
+}
+
+void DoubleOptionStorage::setScaleFactor(double factor)
+{
+    GMX_RELEASE_ASSERT(factor > 0.0, "Invalid scaling factor");
+    if (!hasFlag(efHasDefaultValue))
     {
-        double factor = hostOptions().globalProperties().timeScaleFactor();
+        double scale = factor / factor_;
         ValueList::iterator i;
         for (i = values().begin(); i != values().end(); ++i)
         {
-            (*i) *= factor;
+            (*i) *= scale;
         }
         refreshValues();
     }
+    factor_ = factor;
 }
 
 /********************************************************************
@@ -235,6 +231,26 @@ DoubleOptionInfo::DoubleOptionInfo(DoubleOptionStorage *option)
 {
 }
 
+DoubleOptionStorage &DoubleOptionInfo::option()
+{
+    return static_cast<DoubleOptionStorage &>(OptionInfo::option());
+}
+
+const DoubleOptionStorage &DoubleOptionInfo::option() const
+{
+    return static_cast<const DoubleOptionStorage &>(OptionInfo::option());
+}
+
+bool DoubleOptionInfo::isTime() const
+{
+    return option().isTime();
+}
+
+void DoubleOptionInfo::setScaleFactor(double factor)
+{
+    option().setScaleFactor(factor);
+}
+
 /********************************************************************
  * DoubleOption
  */
index 4236e301a53177c94b588505dd0f2d703f85ef8f..cd620616282f9faec598193999af9ba4845bfc28 100644 (file)
@@ -150,15 +150,16 @@ class DoubleOption : public OptionTemplate<double, DoubleOption>
         //! \copydoc IntegerOption::vector()
         MyClass &vector() { setVector(); return me(); }
         /*! \brief
-         * Sets the option to obey time conversion rules.
+         * Marks this option as providing a time value whose unit can be changed.
          *
-         * For options with this flag, the user can specify a global time unit
-         * using a global option that is added by Options::addDefaultOptions().
-         * For the program, the value is always converted to ps.
-         *
-         * When this flag is specified, the option value is available only
-         * after Options::finish() has been called instead of immediately after
-         * assignment.
+         * By itself, this option does nothing.  It marks the option as a time
+         * value such that TimeUnitManager::scaleTimeOptions() can process it.
+         * In typical cases, Gromacs scales the time options just before
+         * Options::finish() has been called, so the option value is only
+         * available after all option values have been processed.
+         * All values in the program are in ps (including any default value);
+         * user-provided values are scaled according to the time unit set in
+         * TimeUnitManager.
          */
         MyClass &timeValue() { _bTime = true; return me(); }
 
index d95e449c6fe6043be336c8e4973c1d4cc1d27163..c605dd0e01da7017cf5c2854bcf9a31ad53cbaed 100644 (file)
@@ -119,17 +119,21 @@ class DoubleOptionStorage : public OptionStorageTemplate<double>
         //! \copydoc IntegerOptionStorage::IntegerOptionStorage()
         DoubleOptionStorage(const DoubleOption &settings, Options *options);
 
-        virtual OptionInfo &optionInfo() { return _info; }
+        virtual OptionInfo &optionInfo() { return info_; }
         virtual const char *typeString() const;
         virtual std::string formatValue(int i) const;
 
+        bool isTime() const { return bTime_; }
+        void setScaleFactor(double factor);
+
     private:
         virtual void convertValue(const std::string &value);
         virtual void processSetValues(ValueList *values);
         virtual void processAll();
 
-        DoubleOptionInfo        _info;
-        bool                    _bTime;
+        DoubleOptionInfo        info_;
+        bool                    bTime_;
+        double                  factor_;
 };
 
 /*! \internal \brief
index 54f28d6d02bad430ee886b95b7ab83bb4c6b547e..6338a550036f86ed080b71d825a5386b0f12503d 100644 (file)
 #include <smalloc.h>
 #include <statutil.h>
 
-#include "gromacs/fatalerror/gmxassert.h"
 #include "gromacs/options/basicoptions.h"
 #include "gromacs/options/options.h"
 
 namespace gmx
 {
 
-static const char *const timeUnits[] = {
-    "fs", "ps", "ns", "us", "ms",  "s", NULL
-};
 static const char *const plotFormats[] = {
     "none", "xmgrace", "xmgr", NULL
 };
-static const double timeScaleFactors[] = {
-    1e-3,    1,  1e3,  1e6,  1e9, 1e12
-};
 
 
 OptionsGlobalProperties::OptionsGlobalProperties()
-    : _usedProperties(0), _timeUnit(1), _plotFormat(1),
+    : _usedProperties(0), _plotFormat(1),
       _oenv(NULL)
 {
     // TODO: If/when this is refactored, exception safety should be considered
@@ -79,24 +72,8 @@ OptionsGlobalProperties::~OptionsGlobalProperties()
 }
 
 
-double OptionsGlobalProperties::timeScaleFactor() const
-{
-    GMX_RELEASE_ASSERT(_timeUnit >= 0
-        && (size_t)_timeUnit < sizeof(timeScaleFactors)/sizeof(timeScaleFactors[0]),
-        "Time unit index has become out of range");
-    return timeScaleFactors[_timeUnit];
-}
-
-
 void OptionsGlobalProperties::addDefaultOptions(Options *options)
 {
-    if (isPropertyUsed(eogpTimeScaleFactor))
-    {
-        options->addOption(StringOption("tu").enumValue(timeUnits)
-                               .defaultValue("ps")
-                               .storeEnumIndex(&_timeUnit)
-                               .description("Unit for time values"));
-    }
     if (isPropertyUsed(eogpPlotFormat))
     {
         options->addOption(StringOption("xvg").enumValue(plotFormats)
@@ -109,10 +86,6 @@ void OptionsGlobalProperties::addDefaultOptions(Options *options)
 
 void OptionsGlobalProperties::finish()
 {
-    if (isPropertyUsed(eogpTimeScaleFactor))
-    {
-        _oenv->time_unit = static_cast<time_unit_t>(_timeUnit + 1);
-    }
     if (isPropertyUsed(eogpPlotFormat))
     {
         if (_plotFormat == 0)
index 28ca7f4466c873f4737ce6f32ba9f15a735e674e..ad297be16f7f5db5ef34889874eaa5b8391474cd 100644 (file)
@@ -51,7 +51,6 @@ class Options;
  */
 enum OptionGlobalPropertyId
 {
-    eogpTimeScaleFactor,
     eogpPlotFormat,
 };
 
@@ -84,8 +83,6 @@ class OptionsGlobalProperties
             _usedProperties |= (1<<id);
         }
 
-        //! Returns the scaling factor to get times in ps.
-        double timeScaleFactor() const;
         /*! \brief
          * Returns an output environment structure for interfacing with old
          * code.
@@ -129,7 +126,6 @@ class OptionsGlobalProperties
         void finish();
 
         unsigned long           _usedProperties;
-        int                     _timeUnit;
         int                     _plotFormat;
         output_env_t            _oenv;
 
index abc88a376c387b34ea14e6cfe0800b5e725ec269..83fb51c2bf2bbdb4da1e3d10600be0934ebb5ae8 100644 (file)
@@ -57,13 +57,15 @@ enum OptionFlag
 {
     //! %Option has been set.
     efSet                 = 1<<0,
+    //! The current value of the option is a programmatic default value.
+    efHasDefaultValue     = 1<<1,
     /*! \brief
-     * The current value of the option is a default value.
+     * Next assignment to the option clears old values.
      *
-     * This flag is also set when a new option source starts, such that values
+     * This flag is set when a new option source starts, such that values
      * from the new source will overwrite old ones.
      */
-    efHasDefaultValue     = 1<<1,
+    efClearOnNextSet      = 1<<5,
     //! %Option is required to be set.
     efRequired            = 1<<2,
     //! %Option can be specified multiple times.
index 3ad175709a53e075cec72c5b81a9f30113cecf41..c4815b5a5901cbf678cd69e364fc3cbe326cc9d7 100644 (file)
@@ -277,11 +277,11 @@ OptionStorageTemplate<T>::OptionStorageTemplate(const OptionTemplate<T, U> &sett
     }
     if (!hasFlag(efNoDefaultValue))
     {
+        setFlag(efHasDefaultValue);
         if (settings._defaultValue != NULL)
         {
             _values->clear();
             addValue(*settings._defaultValue);
-            setFlag(efHasDefaultValue);
             // TODO: This is a bit hairy, as it indirectly calls a virtual function.
             commitValues();
         }
@@ -294,13 +294,17 @@ OptionStorageTemplate<T>::OptionStorageTemplate(const OptionTemplate<T, U> &sett
             {
                 _values->push_back(_store[i]);
             }
-            setFlag(efHasDefaultValue);
         }
         if (settings._defaultValueIfSet != NULL)
         {
+            if (hasFlag(efMulti))
+            {
+                GMX_THROW(APIError("defaultValueIfSet() is not supported with allowMultiple()"));
+            }
             _defaultValueIfSet.reset(new T(*settings._defaultValueIfSet));
         }
     }
+    setFlag(efClearOnNextSet);
     valueGuard.release();
 }
 
@@ -329,6 +333,11 @@ void OptionStorageTemplate<T>::processSet()
     if (_setValues.empty() && _defaultValueIfSet.get() != NULL)
     {
         addValue(*_defaultValueIfSet);
+        setFlag(efHasDefaultValue);
+    }
+    else
+    {
+        clearFlag(efHasDefaultValue);
     }
     if (!hasFlag(efDontCheckMinimumCount)
         && _setValues.size() < static_cast<size_t>(minValueCount()))
@@ -355,10 +364,10 @@ void OptionStorageTemplate<T>::addValue(const T &value)
 template <typename T>
 void OptionStorageTemplate<T>::commitValues()
 {
-    if (hasFlag(efHasDefaultValue))
+    if (hasFlag(efClearOnNextSet))
     {
         _values->swap(_setValues);
-        clearFlag(efHasDefaultValue);
+        clearFlag(efClearOnNextSet);
     }
     else
     {
index 7c3f3df74dd7a478c897f09e4780a32d73f21ffc..1e0aaaed3b94c422c3cc2d127f35559491989298 100644 (file)
@@ -1,3 +1,4 @@
 add_gtest_or_gmock_test(OptionsUnitTests options-test
                         cmdlineparser.cpp option.cpp optionsassigner.cpp
+                        timeunitmanager.cpp
                         GMOCK_ONLY abstractoptionstorage.cpp)
diff --git a/src/gromacs/options/tests/timeunitmanager.cpp b/src/gromacs/options/tests/timeunitmanager.cpp
new file mode 100644 (file)
index 0000000..73c0495
--- /dev/null
@@ -0,0 +1,176 @@
+/*
+ *
+ *                This source code is part of
+ *
+ *                 G   R   O   M   A   C   S
+ *
+ *          GROningen MAchine for Chemical Simulations
+ *
+ * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
+ * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
+ * Copyright (c) 2001-2009, The GROMACS development team,
+ * check out http://www.gromacs.org for more information.
+
+ * This program is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU General Public License
+ * as published by the Free Software Foundation; either version 2
+ * of the License, or (at your option) any later version.
+ *
+ * If you want to redistribute modifications, 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 www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the papers on the package - you can find them in the top README file.
+ *
+ * For more info, check our website at http://www.gromacs.org
+ */
+/*! \internal \file
+ * \brief
+ * Tests handling of time units with gmx::TimeUnitManager.
+ *
+ * Also related functionality in gmx::DoubleOptionStorage is tested.
+ *
+ * \author Teemu Murtola <teemu.murtola@cbr.su.se>
+ * \ingroup module_options
+ */
+#include <gtest/gtest.h>
+
+#include "gromacs/options/basicoptions.h"
+#include "gromacs/options/options.h"
+#include "gromacs/options/optionsassigner.h"
+#include "gromacs/options/timeunitmanager.h"
+
+namespace
+{
+
+TEST(TimeUnitManagerTest, BasicOperations)
+{
+    gmx::TimeUnitManager manager;
+    EXPECT_EQ(gmx::eTimeUnit_ps, manager.timeUnit());
+    EXPECT_DOUBLE_EQ(1.0, manager.timeScaleFactor());
+    manager.setTimeUnit(gmx::eTimeUnit_ns);
+    EXPECT_EQ(gmx::eTimeUnit_ns, manager.timeUnit());
+    EXPECT_DOUBLE_EQ(1e3, manager.timeScaleFactor());
+    EXPECT_DOUBLE_EQ(1e-3, manager.inverseTimeScaleFactor());
+}
+
+TEST(TimeUnitManagerTest, ScalesAssignedOptionValue)
+{
+    gmx::TimeUnitManager manager;
+
+    gmx::Options options(NULL, NULL);
+    double value = 0.0;
+    using gmx::DoubleOption;
+    ASSERT_NO_THROW(options.addOption(DoubleOption("p").store(&value).timeValue()));
+
+    gmx::OptionsAssigner assigner(&options);
+    EXPECT_NO_THROW(assigner.start());
+    ASSERT_NO_THROW(assigner.startOption("p"));
+    ASSERT_NO_THROW(assigner.appendValue("1.5"));
+    EXPECT_NO_THROW(assigner.finishOption());
+    EXPECT_NO_THROW(assigner.finish());
+
+    EXPECT_DOUBLE_EQ(1.5, value);
+    manager.setTimeUnit(gmx::eTimeUnit_ns);
+    manager.scaleTimeOptions(&options);
+    EXPECT_DOUBLE_EQ(1500, value);
+
+    EXPECT_NO_THROW(options.finish());
+
+    manager.setTimeUnit(gmx::eTimeUnit_us);
+    manager.scaleTimeOptions(&options);
+    EXPECT_DOUBLE_EQ(1500000, value);
+
+    manager.setTimeUnit(gmx::eTimeUnit_fs);
+    manager.scaleTimeOptions(&options);
+    EXPECT_DOUBLE_EQ(0.0015, value);
+
+    manager.setTimeUnit(gmx::eTimeUnit_ps);
+    manager.scaleTimeOptions(&options);
+    EXPECT_DOUBLE_EQ(1.5, value);
+}
+
+TEST(TimeUnitManagerTest, DoesNotScaleDefaultValues)
+{
+    gmx::TimeUnitManager manager;
+
+    gmx::Options options(NULL, NULL);
+    double value = 1.5, value2 = 0.0;
+    using gmx::DoubleOption;
+    ASSERT_NO_THROW(options.addOption(DoubleOption("p").store(&value).timeValue()));
+    ASSERT_NO_THROW(options.addOption(DoubleOption("q").store(&value2).timeValue()
+                        .defaultValueIfSet(2.5)));
+
+    gmx::OptionsAssigner assigner(&options);
+    EXPECT_NO_THROW(assigner.start());
+    ASSERT_NO_THROW(assigner.startOption("q"));
+    EXPECT_NO_THROW(assigner.finishOption());
+    EXPECT_NO_THROW(assigner.finish());
+    EXPECT_NO_THROW(options.finish());
+
+    EXPECT_DOUBLE_EQ(2.5, value2);
+    manager.setTimeUnit(gmx::eTimeUnit_ns);
+    manager.scaleTimeOptions(&options);
+    EXPECT_DOUBLE_EQ(1.5, value);
+    EXPECT_DOUBLE_EQ(2.5, value2);
+}
+
+TEST(TimeUnitManagerTest, ScalesUserInputWithMultipleSources)
+{
+    gmx::TimeUnitManager manager;
+
+    gmx::Options options(NULL, NULL);
+    double value = 0.0;
+    using gmx::DoubleOption;
+    ASSERT_NO_THROW(options.addOption(DoubleOption("p").store(&value).timeValue()));
+
+    gmx::OptionsAssigner assigner(&options);
+    EXPECT_NO_THROW(assigner.start());
+    ASSERT_NO_THROW(assigner.startOption("p"));
+    ASSERT_NO_THROW(assigner.appendValue("1.5"));
+    EXPECT_NO_THROW(assigner.finishOption());
+    EXPECT_NO_THROW(assigner.finish());
+    gmx::OptionsAssigner assigner2(&options);
+    EXPECT_NO_THROW(assigner2.start());
+    EXPECT_NO_THROW(assigner2.finish());
+    EXPECT_NO_THROW(options.finish());
+
+    EXPECT_DOUBLE_EQ(1.5, value);
+    manager.setTimeUnit(gmx::eTimeUnit_ns);
+    manager.scaleTimeOptions(&options);
+    EXPECT_DOUBLE_EQ(1500, value);
+}
+
+TEST(TimeUnitManagerTest, TimeUnitOptionWorks)
+{
+    gmx::TimeUnitManager manager;
+
+    gmx::Options options(NULL, NULL);
+    double value = 0.0;
+    using gmx::DoubleOption;
+    ASSERT_NO_THROW(options.addOption(DoubleOption("p").store(&value).timeValue()));
+    ASSERT_NO_THROW(manager.addTimeUnitOption(&options, "tu"));
+
+    gmx::OptionsAssigner assigner(&options);
+    EXPECT_NO_THROW(assigner.start());
+    ASSERT_NO_THROW(assigner.startOption("p"));
+    ASSERT_NO_THROW(assigner.appendValue("1.5"));
+    EXPECT_NO_THROW(assigner.finishOption());
+    ASSERT_NO_THROW(assigner.startOption("tu"));
+    ASSERT_NO_THROW(assigner.appendValue("ns"));
+    EXPECT_NO_THROW(assigner.finishOption());
+    EXPECT_NO_THROW(assigner.finish());
+
+    EXPECT_DOUBLE_EQ(1.5, value);
+    EXPECT_EQ(gmx::eTimeUnit_ns, manager.timeUnit());
+    manager.scaleTimeOptions(&options);
+    EXPECT_DOUBLE_EQ(1500, value);
+
+    EXPECT_NO_THROW(options.finish());
+}
+
+} // namespace
diff --git a/src/gromacs/options/timeunitmanager.cpp b/src/gromacs/options/timeunitmanager.cpp
new file mode 100644 (file)
index 0000000..a86f72f
--- /dev/null
@@ -0,0 +1,143 @@
+/*
+ *
+ *                This source code is part of
+ *
+ *                 G   R   O   M   A   C   S
+ *
+ *          GROningen MAchine for Chemical Simulations
+ *
+ * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
+ * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
+ * Copyright (c) 2001-2009, The GROMACS development team,
+ * check out http://www.gromacs.org for more information.
+
+ * This program is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU General Public License
+ * as published by the Free Software Foundation; either version 2
+ * of the License, or (at your option) any later version.
+ *
+ * If you want to redistribute modifications, 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 www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the papers on the package - you can find them in the top README file.
+ *
+ * For more info, check our website at http://www.gromacs.org
+ */
+/*! \internal \file
+ * \brief
+ * Implements gmx::TimeUnitManager.
+ *
+ * \author Teemu Murtola <teemu.murtola@cbr.su.se>
+ * \ingroup module_options
+ */
+#include "gromacs/options/timeunitmanager.h"
+
+#include "gromacs/fatalerror/gmxassert.h"
+#include "gromacs/options/basicoptioninfo.h"
+#include "gromacs/options/basicoptions.h"
+#include "gromacs/options/options.h"
+#include "gromacs/options/optionsvisitor.h"
+
+namespace gmx
+{
+
+// These must correspond to the TimeUnit enum in the header!
+static const char *const g_timeUnits[] = {
+    "fs", "ps", "ns", "us", "ms",  "s", NULL
+};
+static const double g_timeScaleFactors[] = {
+    1e-3,    1,  1e3,  1e6,  1e9, 1e12
+};
+
+TimeUnitManager::TimeUnitManager()
+    : timeUnit_(eTimeUnit_ps)
+{
+}
+
+TimeUnitManager::TimeUnitManager(TimeUnit unit)
+{
+    setTimeUnit(unit);
+}
+
+void TimeUnitManager::setTimeUnit(TimeUnit unit)
+{
+    GMX_RELEASE_ASSERT(unit >= 0 && unit <= eTimeUnit_s,
+                       "Invalid time unit");
+    timeUnit_ = unit;
+}
+
+const char *TimeUnitManager::timeUnitAsString() const
+{
+    GMX_RELEASE_ASSERT(timeUnit_ >= 0 && timeUnit_ <= eTimeUnit_s,
+                       "Invalid time unit");
+    return g_timeUnits[timeUnit_];
+}
+
+double TimeUnitManager::timeScaleFactor() const
+{
+    GMX_RELEASE_ASSERT(timeUnit_ >= 0
+        && (size_t)timeUnit_ < sizeof(g_timeScaleFactors)/sizeof(g_timeScaleFactors[0]),
+        "Time unit index has become out-of-range");
+    return g_timeScaleFactors[timeUnit_];
+}
+
+double TimeUnitManager::inverseTimeScaleFactor() const
+{
+    return 1.0 / timeScaleFactor();
+}
+
+void TimeUnitManager::addTimeUnitOption(Options *options, const char *name)
+{
+    options->addOption(StringOption(name).enumValue(g_timeUnits)
+                           .defaultValue(g_timeUnits[timeUnit()])
+                           .storeEnumIndex(&timeUnit_)
+                           .description("Unit for time values"));
+}
+
+namespace
+{
+
+/*! \internal \brief
+ * Option visitor that scales time options.
+ *
+ * \ingroup module_options
+ */
+class TimeOptionScaler : public OptionsModifyingTypeVisitor<DoubleOptionInfo>
+{
+    public:
+        //! Initializes a scaler with the given factor.
+        explicit TimeOptionScaler(double factor) : factor_(factor) {}
+
+        void visitSubSection(Options *section)
+        {
+            OptionsModifyingIterator iterator(section);
+            iterator.acceptSubSections(this);
+            iterator.acceptOptions(this);
+        }
+
+        void visitOptionType(DoubleOptionInfo *option)
+        {
+            if (option->isTime())
+            {
+                option->setScaleFactor(factor_);
+            }
+        }
+
+    private:
+        double                  factor_;
+};
+
+} // namespace
+
+void TimeUnitManager::scaleTimeOptions(Options *options) const
+{
+    double factor = timeScaleFactor();
+    TimeOptionScaler(factor).visitSubSection(options);
+}
+
+} // namespace gmx
diff --git a/src/gromacs/options/timeunitmanager.h b/src/gromacs/options/timeunitmanager.h
new file mode 100644 (file)
index 0000000..53b0869
--- /dev/null
@@ -0,0 +1,152 @@
+/*
+ *
+ *                This source code is part of
+ *
+ *                 G   R   O   M   A   C   S
+ *
+ *          GROningen MAchine for Chemical Simulations
+ *
+ * Written by David van der Spoel, Erik Lindahl, Berk Hess, and others.
+ * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
+ * Copyright (c) 2001-2009, The GROMACS development team,
+ * check out http://www.gromacs.org for more information.
+
+ * This program is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU General Public License
+ * as published by the Free Software Foundation; either version 2
+ * of the License, or (at your option) any later version.
+ *
+ * If you want to redistribute modifications, 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 www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the papers on the package - you can find them in the top README file.
+ *
+ * For more info, check our website at http://www.gromacs.org
+ */
+/*! \file
+ * \brief
+ * Declares gmx::TimeUnitManager.
+ *
+ * \author Teemu Murtola <teemu.murtola@cbr.su.se>
+ * \inpublicapi
+ * \ingroup module_options
+ */
+#ifndef GMX_OPTIONS_TIMEUNITMANAGER_H
+#define GMX_OPTIONS_TIMEUNITMANAGER_H
+
+#include "../fatalerror/gmxassert.h"
+
+namespace gmx
+{
+
+class Options;
+
+/*! \brief
+ * Time values for TimeUnitManager.
+ *
+ * \if internal
+ * Currently, this should match with the ::time_unit_t enum defined in oenv.h
+ * except that there is no NULL first item in this enum.
+ * \endif
+ *
+ * \inpublicapi
+ */
+enum TimeUnit
+{
+    eTimeUnit_fs, //!< Femtoseconds.
+    eTimeUnit_ps, //!< Picoseconds.
+    eTimeUnit_ns, //!< Nanoseconds.
+    eTimeUnit_us, //!< Microseconds.
+    eTimeUnit_ms, //!< Milliseconds.
+    eTimeUnit_s   //!< Seconds.
+};
+
+/*! \brief
+ * Provides common functionality for time unit conversions.
+ *
+ * Methods/objects that need to deal with time units can either take a
+ * TimeUnitManager object, or they can take a TimeUnit value and construct a
+ * TimeUnitManager object internally.
+ *
+ * Default copy constructor and assignment are used: the copy is an independent
+ * object that is initialized with the same time unit as the original.
+ *
+ * \if internal
+ * \todo
+ * Most of this class is independent of the options implementation.
+ * To ease reuse, it could be split such that the generic part is moved to the
+ * utility module, and only the options-specific parts left in the options
+ * module.
+ * \endif
+ *
+ * \inpublicapi
+ * \ingroup module_options
+ */
+class TimeUnitManager
+{
+    public:
+        //! Creates a time unit manager with the default (ps) time unit.
+        TimeUnitManager();
+        //! Creates a time unit manager with the given time unit.
+        explicit TimeUnitManager(TimeUnit unit);
+
+        //! Returns the currently selected time unit.
+        TimeUnit timeUnit() const
+        {
+            GMX_ASSERT(timeUnit_ >= 0 && timeUnit_ <= eTimeUnit_s,
+                       "Time unit index has become out-of-range");
+            return static_cast<TimeUnit>(timeUnit_);
+        }
+        //! Set a new time unit for the manager.
+        void setTimeUnit(TimeUnit unit);
+
+        //! Returns a string constant corresponding to the current time unit.
+        const char *timeUnitAsString() const;
+
+        //! Returns the scaling factor to convert times to ps.
+        double timeScaleFactor() const;
+        //! Returns the scaling factor to convert times from ps.
+        double inverseTimeScaleFactor() const;
+
+        /*! \brief
+         * Adds a common option for selecting the time unit.
+         *
+         * \param[in,out] options Options to which the common option is added.
+         * \param[in]     name    Name of the option to add.
+         *
+         * Adds an enum option to \p options to select the time unit for this
+         * manager.
+         */
+        void addTimeUnitOption(Options *options, const char *name);
+        /*! \brief
+         * Scales user input values given to time options.
+         *
+         * \param[in,out] options Options in which to scale times.
+         *
+         * Scales each time option (see DoubleOption::timeValue()) in
+         * \p options such that any user-given values are interpreted as given
+         * in the time unit specified by this manager, and scaled to
+         * picoseconds.  Programmatically given values (e.g., as default values
+         * for the options) are not scaled.
+         */
+        void scaleTimeOptions(Options *options) const;
+
+    private:
+        /*! \brief
+         * Currently set time unit for this manager.
+         *
+         * Type is int to make it possible to use it with
+         * StringOption::storeEnumIndex(), but it should always one of the
+         * allowed values for TimeUnit.
+         */
+        int                     timeUnit_;
+};
+
+} // namespace gmx
+
+#endif
index 3c59c30b4f592fabb708a0ef2ba77d48f8e283b7..c9eaece10d8d610b34e39e5a6900a54dc491496a 100644 (file)
@@ -40,6 +40,8 @@
 
 #include "analysissettings.h"
 
+#include "../options/timeunitmanager.h"
+
 namespace gmx
 {
 
@@ -54,6 +56,7 @@ class TrajectoryAnalysisSettings::Impl
         //! Initializes the default values for the settings object.
         Impl() : flags(0), frflags(0), bRmPBC(true), bPBC(true) {}
 
+        TimeUnitManager      timeUnitManager;
         unsigned long        flags;
         int                  frflags;
 
index 37a32d97909b558c8250b5988b6416f11ab36263..a87c335ce646733f0cdf7d6ec14055924fab636d 100644 (file)
@@ -121,8 +121,6 @@ TrajectoryAnalysisCommandLineRunner::Impl::parseOptions(
         Options *options,
         int *argc, char *argv[])
 {
-    int rc;
-
     Options *moduleOptions = _module->initOptions(settings);
     GMX_RELEASE_ASSERT(moduleOptions != NULL, "Module returned NULL options");
     Options *commonOptions = common->initOptions();
@@ -147,6 +145,7 @@ TrajectoryAnalysisCommandLineRunner::Impl::parseOptions(
             throw;
         }
         printHelp(*options, *common);
+        common->scaleTimeOptions(options);
         options->finish();
     }
 
index cb0030d8880b58a924ec3d77f39acfd45fc95252..e30135bb3b9b5dc9567d3bd3b78927c921f573b8 100644 (file)
@@ -197,6 +197,9 @@ TrajectoryAnalysisRunnerCommon::initOptions()
     options.addOption(DoubleOption("dt").store(&_impl->_deltaTime).timeValue()
                           .description("Only use frame if t MOD dt == first time (%t)"));
 
+    // Add time unit option.
+    settings._impl->timeUnitManager.addTimeUnitOption(&options, "tu");
+
     // Add common options for trajectory processing.
     if (!settings.hasFlag(TrajectoryAnalysisSettings::efNoUserRmPBC))
     {
@@ -213,6 +216,13 @@ TrajectoryAnalysisRunnerCommon::initOptions()
 }
 
 
+void
+TrajectoryAnalysisRunnerCommon::scaleTimeOptions(Options *options)
+{
+    _impl->_settings._impl->timeUnitManager.scaleTimeOptions(options);
+}
+
+
 bool
 TrajectoryAnalysisRunnerCommon::initOptionsDone()
 {
index bf5317865fcdc532c6efcb8c4815750a72d08496..3a4575f02a2305a881d996ea05242fbf742d68ce 100644 (file)
@@ -75,6 +75,7 @@ class TrajectoryAnalysisRunnerCommon
         ~TrajectoryAnalysisRunnerCommon();
 
         Options *initOptions();
+        void scaleTimeOptions(Options *options);
         bool initOptionsDone();
         void initIndexGroups(SelectionCollection *selections);
         void doneIndexGroups(SelectionCollection *selections);