Improve enum option interface.
authorTeemu Murtola <teemu.murtola@gmail.com>
Tue, 12 Feb 2013 12:44:17 +0000 (14:44 +0200)
committerTeemu Murtola <teemu.murtola@gmail.com>
Thu, 7 Mar 2013 18:57:28 +0000 (20:57 +0200)
It is no longer necessary to have explicit NULL terminating arrays
passed to StringOption::enumValue().  Instead, let the compiler deduce
the size of the array.  Reduces the potential for coding mistakes, as it
is now impossible to forget to terminate the array.

Added a separate enumValueFromNullTerminatedArray() for those rare cases
where the array is programmatically constructed or the compiler for
other reasons doesn't know its size at compile time.

Change-Id: I5e7d63db1eeea6d9d271fa299c98c781f52bd89c

src/gromacs/analysisdata/modules/plot.cpp
src/gromacs/commandline/tests/cmdlinehelpwriter.cpp
src/gromacs/options/basicoptions.cpp
src/gromacs/options/basicoptions.h
src/gromacs/options/tests/option.cpp
src/gromacs/options/tests/optionsassigner.cpp
src/gromacs/options/timeunitmanager.cpp
src/gromacs/selection/selectioncollection.cpp
src/gromacs/trajectoryanalysis/modules/angle.cpp
src/gromacs/trajectoryanalysis/modules/select.cpp

index a6e309ae8463b901e0e4a6df0af9379f30dd0817..b3db5ea4faa2964fa26dc32808e93d159030e11e 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -68,7 +68,7 @@ namespace
 
 //! Enum values for plot formats.
 const char *const g_plotFormats[] = {
-    "none", "xmgrace", "xmgr", NULL
+    "none", "xmgrace", "xmgr"
 };
 
 } // namespace
index 6b66d87873bed2d1f58e9fe2d33805316a6bb712..a211134ef7f461eeb0e57b0d2f8e25bf4ab89232 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2012, by the GROMACS development team, led by
+ * Copyright (c) 2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -117,7 +117,7 @@ TEST_F(CommandLineHelpWriterTest, HandlesOptionTypes)
                           .timeValue().defaultValue(10.0));
     options.addOption(StringOption("string").description("String option")
                           .defaultValue("test"));
-    const char * const enumValues[] = {"no", "opt1", "opt2", NULL};
+    const char * const enumValues[] = { "no", "opt1", "opt2" };
     options.addOption(StringOption("enum").description("Enum option")
                           .enumValue(enumValues).defaultEnumIndex(0));
 
index f7209545a7ce28c948a60e40f45125efb538f594..33eca9bb6101dd049fb802a6229771bf5fea38ce 100644 (file)
@@ -319,8 +319,21 @@ StringOptionStorage::StringOptionStorage(const StringOption &settings)
         enumIndexStore_ = settings.enumIndexStore_;
         const std::string *defaultValue = settings.defaultValue();
         int                match        = -1;
-        for (int i = 0; settings.enumValues_[i] != NULL; ++i)
+        int                count        = settings.enumValuesCount_;
+        if (count < 0)
         {
+            count = 0;
+            while (settings.enumValues_[count] != NULL)
+            {
+                ++count;
+            }
+        }
+        for (int i = 0; i < count; ++i)
+        {
+            if (settings.enumValues_[i] == NULL)
+            {
+                GMX_THROW(APIError("Enumeration value cannot be NULL"));
+            }
             if (defaultValue != NULL && settings.enumValues_[i] == *defaultValue)
             {
                 match = i;
index 8a318f14e0cad6cc16b608209a92d76e0c0af61b..356f868274ca6494f07ddb95420412fddfdabf54 100644 (file)
@@ -201,7 +201,7 @@ class DoubleOption : public OptionTemplate<double, DoubleOption>
    std::string  str;
    options.addOption(StringOption("str").store(&str));
    // Option that only accepts predefined values
-   const char * const  allowed[] = { "atom", "residue", "molecule", NULL };
+   const char * const  allowed[] = { "atom", "residue", "molecule" };
    std::string  str;
    int          type;
    options.addOption(StringOption("type").enumValue(allowed).store(&str)
@@ -220,16 +220,15 @@ class StringOption : public OptionTemplate<std::string, StringOption>
 
         //! Initializes an option with the given name.
         explicit StringOption(const char *name)
-            : MyBase(name), enumValues_(NULL), defaultEnumIndex_(-1),
-              enumIndexStore_(NULL)
+            : MyBase(name), enumValues_(NULL), enumValuesCount_(0),
+              defaultEnumIndex_(-1), enumIndexStore_(NULL)
         {
         }
 
         /*! \brief
          * Sets the option to only accept one of a fixed set of strings.
          *
-         * \param[in] values  Array of strings to accept, a NULL pointer
-         *      following the last string.
+         * \param[in] values  Array of strings to accept.
          *
          * Also accepts prefixes of the strings; if a prefix matches more than
          * one of the possible strings, the shortest one is used (in a tie, the
@@ -241,8 +240,35 @@ class StringOption : public OptionTemplate<std::string, StringOption>
          *
          * The strings are copied once the option is created.
          */
-        MyClass &enumValue(const char *const *values)
-        { enumValues_ = values; return me(); }
+        template <size_t count>
+        MyClass &enumValue(const char *const (&values)[count])
+        {
+            GMX_ASSERT(enumValues_ == NULL,
+                       "Multiple sets of enumerated values specified");
+            enumValues_      = values;
+            enumValuesCount_ = count;
+            return me();
+        }
+        /*! \brief
+         * Sets the option to only accept one of a fixed set of strings.
+         *
+         * \param[in] values  Array of strings to accept, with a NULL pointer
+         *      following the last string.
+         *
+         * Works otherwise as the array version, but accepts a pointer to
+         * an array of undetermined length.  The end of the array is indicated
+         * by a NULL pointer in the array.
+         *
+         * \see enumValue()
+         */
+        MyClass &enumValueFromNullTerminatedArray(const char *const *values)
+        {
+            GMX_ASSERT(enumValues_ == NULL,
+                       "Multiple sets of enumerated values specified");
+            enumValues_      = values;
+            enumValuesCount_ = -1;
+            return me();
+        }
         /*! \brief
          * Sets the default value using an index into the enumeration table.
          *
@@ -250,7 +276,7 @@ class StringOption : public OptionTemplate<std::string, StringOption>
          */
         MyClass &defaultEnumIndex(int index)
         {
-            GMX_RELEASE_ASSERT(index >= 0, "Invalid enumeration index");
+            GMX_ASSERT(index >= 0, "Invalid enumeration index");
             defaultEnumIndex_ = index;
             return me();
         }
@@ -264,6 +290,10 @@ class StringOption : public OptionTemplate<std::string, StringOption>
          * and there is no default value, -1 is stored.
          *
          * Cannot be specified without enumValue().
+         *
+         * \todo
+         * Implement this such that it is also possible to store the value
+         * directly into a real enum type.
          */
         MyClass &storeEnumIndex(int *store)
         { enumIndexStore_ = store; return me(); }
@@ -273,6 +303,7 @@ class StringOption : public OptionTemplate<std::string, StringOption>
         virtual AbstractOptionStoragePointer createStorage() const;
 
         const char *const      *enumValues_;
+        int                     enumValuesCount_;
         int                     defaultEnumIndex_;
         int                    *enumIndexStore_;
 
index 2cb3f3f5414db4c9031832d0bbb6b0999c90005d..bd0a3f927da795782b641cbf17707962cdb68f8f 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -68,7 +68,7 @@ TEST(OptionsTest, FailsOnIncorrectEnumDefaultValue)
 {
     gmx::Options                options(NULL, NULL);
     std::string                 value;
-    const char * const          allowed[] = { "none", "test", "value", NULL };
+    const char * const          allowed[] = { "none", "test", "value" };
     using gmx::StringOption;
     ASSERT_THROW(options.addOption(StringOption("name").store(&value)
                                        .enumValue(allowed)
index d45d636631284af70f9bf420dd55a521f2e49128..bb7cb44df4d84d64431dddd0c5a94ace74d6d151 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -699,7 +699,7 @@ TEST(OptionsAssignerStringTest, HandlesEnumValue)
 {
     gmx::Options           options(NULL, NULL);
     std::string            value;
-    const char * const     allowed[] = { "none", "test", "value", NULL };
+    const char * const     allowed[] = { "none", "test", "value" };
     int                    index     = -1;
     using gmx::StringOption;
     ASSERT_NO_THROW(options.addOption(
@@ -718,13 +718,37 @@ TEST(OptionsAssignerStringTest, HandlesEnumValue)
     EXPECT_EQ(1, index);
 }
 
-TEST(OptionsAssignerStringTest, HandlesIncorrectEnumValue)
+TEST(OptionsAssignerStringTest, HandlesEnumValueFromNullTerminatedArray)
 {
     gmx::Options           options(NULL, NULL);
     std::string            value;
     const char * const     allowed[] = { "none", "test", "value", NULL };
     int                    index     = -1;
     using gmx::StringOption;
+    ASSERT_NO_THROW(options.addOption(
+                            StringOption("p").store(&value)
+                                .enumValueFromNullTerminatedArray(allowed)
+                                .storeEnumIndex(&index)));
+
+    gmx::OptionsAssigner assigner(&options);
+    EXPECT_NO_THROW(assigner.start());
+    ASSERT_NO_THROW(assigner.startOption("p"));
+    ASSERT_NO_THROW(assigner.appendValue("value"));
+    EXPECT_NO_THROW(assigner.finishOption());
+    EXPECT_NO_THROW(assigner.finish());
+    EXPECT_NO_THROW(options.finish());
+
+    EXPECT_EQ("value", value);
+    EXPECT_EQ(2, index);
+}
+
+TEST(OptionsAssignerStringTest, HandlesIncorrectEnumValue)
+{
+    gmx::Options           options(NULL, NULL);
+    std::string            value;
+    const char * const     allowed[] = { "none", "test", "value" };
+    int                    index     = -1;
+    using gmx::StringOption;
     ASSERT_NO_THROW(options.addOption(
                             StringOption("p").store(&value)
                                 .enumValue(allowed).storeEnumIndex(&index)));
@@ -739,7 +763,7 @@ TEST(OptionsAssignerStringTest, CompletesEnumValue)
 {
     gmx::Options           options(NULL, NULL);
     std::string            value;
-    const char * const     allowed[] = { "none", "test", "value", NULL };
+    const char * const     allowed[] = { "none", "test", "value" };
     int                    index     = -1;
     using gmx::StringOption;
     ASSERT_NO_THROW(options.addOption(
@@ -762,7 +786,7 @@ TEST(OptionsAssignerStringTest, HandlesEnumWithNoValue)
 {
     gmx::Options           options(NULL, NULL);
     std::string            value;
-    const char * const     allowed[] = { "none", "test", "value", NULL };
+    const char * const     allowed[] = { "none", "test", "value" };
     int                    index     = -3;
     using gmx::StringOption;
     ASSERT_NO_THROW(options.addOption(
@@ -781,7 +805,7 @@ TEST(OptionsAssignerStringTest, HandlesEnumDefaultValue)
 {
     gmx::Options           options(NULL, NULL);
     std::string            value;
-    const char * const     allowed[] = { "none", "test", "value", NULL };
+    const char * const     allowed[] = { "none", "test", "value" };
     int                    index     = -1;
     using gmx::StringOption;
     ASSERT_NO_THROW(options.addOption(
@@ -804,7 +828,7 @@ TEST(OptionsAssignerStringTest, HandlesEnumDefaultIndex)
 {
     gmx::Options           options(NULL, NULL);
     std::string            value;
-    const char * const     allowed[] = { "none", "test", "value", NULL };
+    const char * const     allowed[] = { "none", "test", "value" };
     int                    index     = -1;
     using gmx::StringOption;
     ASSERT_NO_THROW(options.addOption(
index 8827e23cdf019334196fd85431231d8bff6a7f62..9c0baeded878c0856cf27a96170cdb8f7f0458bc 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2012, by the GROMACS development team, led by
+ * Copyright (c) 2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -55,7 +55,7 @@ namespace
  * These must correspond to the TimeUnit enum in the header!
  */
 const char *const g_timeUnits[] = {
-    "fs", "ps", "ns", "us", "ms",  "s", NULL
+    "fs", "ps", "ns", "us", "ms",  "s"
 };
 /*! \brief
  * Scaling factors from each time unit to internal units (=picoseconds).
index fad1dd62029f6009038b2d67c8f15bb915987eb2..e8e609ef31de4938018982bd518df9a17dc23cc6 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -369,8 +369,8 @@ SelectionCollection::~SelectionCollection()
 void
 SelectionCollection::initOptions(Options *options)
 {
-    static const char * const debug_levels[]
-        = {"no", "basic", "compile", "eval", "full", NULL};
+    const char * const debug_levels[]
+        = {"no", "basic", "compile", "eval", "full" };
     /*
        static const char * const desc[] = {
         "This program supports selections in addition to traditional",
@@ -382,10 +382,12 @@ SelectionCollection::initOptions(Options *options)
      */
 
     const char *const *postypes = PositionCalculationCollection::typeEnumValues;
-    options->addOption(StringOption("selrpos").enumValue(postypes)
+    options->addOption(StringOption("selrpos")
+                           .enumValueFromNullTerminatedArray(postypes)
                            .store(&impl_->rpost_).defaultValue(postypes[0])
                            .description("Selection reference positions"));
-    options->addOption(StringOption("seltype").enumValue(postypes)
+    options->addOption(StringOption("seltype")
+                           .enumValueFromNullTerminatedArray(postypes)
                            .store(&impl_->spost_).defaultValue(postypes[0])
                            .description("Default selection output positions"));
     GMX_RELEASE_ASSERT(impl_->debugLevel_ >= 0 && impl_->debugLevel_ <= 4,
index fa966dfe1f825c336b9b6133340e3131304c6761..4e74bfb7bdecbeb831330041ac7c7b704198475c 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -134,9 +134,9 @@ Angle::initOptions(Options *options, TrajectoryAnalysisSettings * /*settings*/)
          */
     };
     static const char *const cGroup1TypeEnum[] =
-    { "angle", "dihedral", "vector", "plane", NULL };
+    { "angle", "dihedral", "vector", "plane" };
     static const char *const cGroup2TypeEnum[] =
-    { "none", "vector", "plane", "t0", "z", "sphnorm", NULL };
+    { "none", "vector", "plane", "t0", "z", "sphnorm" };
 
     options->setDescription(concatenateStrings(desc));
 
index 36fac8473728ecd91ba4eef78467de3c742b5489..2f7c6d5a1298a75f09c86673cdc034baa30c03f9 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2009,2010,2011,2012, by the GROMACS development team, led by
+ * Copyright (c) 2009,2010,2011,2012,2013, by the GROMACS development team, led by
  * David van der Spoel, Berk Hess, Erik Lindahl, and including many
  * others, as listed in the AUTHORS file in the top-level source
  * directory and at http://www.gromacs.org.
@@ -377,11 +377,11 @@ Select::initOptions(Options *options, TrajectoryAnalysisSettings * /*settings*/)
                            .description("Normalize by total number of positions with -os"));
     options->addOption(BooleanOption("cfnorm").store(&bFracNorm_)
                            .description("Normalize by covered fraction with -os"));
-    const char *const cResNumberEnum[] = { "number", "index", NULL };
+    const char *const cResNumberEnum[] = { "number", "index" };
     options->addOption(StringOption("resnr").store(&resNumberType_)
                            .enumValue(cResNumberEnum).defaultEnumIndex(0)
                            .description("Residue number output type with -oi and -on"));
-    const char *const cPDBAtomsEnum[] = { "all", "maxsel", "selected", NULL };
+    const char *const cPDBAtomsEnum[] = { "all", "maxsel", "selected" };
     options->addOption(StringOption("pdbatoms").store(&pdbAtoms_)
                            .enumValue(cPDBAtomsEnum).defaultEnumIndex(0)
                            .description("Atoms to write with -ofpdb"));