Introduce StringToEnumValue
authorMark Abraham <mark.j.abraham@gmail.com>
Fri, 5 Mar 2021 08:02:59 +0000 (08:02 +0000)
committerMark Abraham <mark.j.abraham@gmail.com>
Fri, 5 Mar 2021 08:02:59 +0000 (08:02 +0000)
This helper functor seems useful in a few places. Adds tests and
improves documentation.

Fixes #3963

docs/dev-manual/naming.rst
src/gromacs/fileio/pdbio.cpp
src/gromacs/gmxpreprocess/toppush.cpp
src/gromacs/utility/stringtoenumvalueconverter.h [new file with mode: 0644]
src/gromacs/utility/tests/CMakeLists.txt
src/gromacs/utility/tests/stringtoenumvalueconverter.cpp [new file with mode: 0644]

index 20df0d8d2b6bfbceadf9ae42c9195ef4c33c149a..89985e7db2aeb87d736d383c2bb87e3d9c3f21d4 100644 (file)
@@ -106,6 +106,7 @@ C++ code
   typically including a verb for boolean variables, like ``foundAtom``.
 * Prefer class enums over regular ones, so that unexpected conversions to
   int do not happen.
+* Name functions to convert class enum values to strings as ``enumValueToString``.
 * When using a non-class enum, prefer to include the name of the enumeration type
   as a base in the name of enum values, e.g., ``HelpOutputFormat_Console``,
   in particular for settings exposed to other modules.
index c46a2ccb2ef8aafe898cc0f7b6c60f88e8688b0a..2e5d1ec39f89b2eae176300556a9ff94ae49cf17 100644 (file)
@@ -63,6 +63,7 @@
 #include "gromacs/utility/futil.h"
 #include "gromacs/utility/smalloc.h"
 #include "gromacs/utility/snprintf.h"
+#include "gromacs/utility/stringtoenumvalueconverter.h"
 
 typedef struct
 {
@@ -463,25 +464,6 @@ void write_pdbfile(FILE*          out,
     sfree(index);
 }
 
-static PdbRecordType line2type(const char* line)
-{
-    int  k;
-    char type[8];
-
-    for (k = 0; (k < 6); k++)
-    {
-        type[k] = line[k];
-    }
-    type[k] = '\0';
-
-    using PdbRecordArray = gmx::EnumerationArray<PdbRecordType, bool>;
-    auto entry           = std::find_if(
-            PdbRecordArray::keys().begin(), PdbRecordArray::keys().end(), [type](const auto& key) {
-                return std::strncmp(type, enumValueToString(key), strlen(enumValueToString(key))) == 0;
-            });
-    return (entry != PdbRecordArray::keys().end()) ? *entry : PdbRecordType::Count;
-}
-
 static void read_anisou(char line[], int natom, t_atoms* atoms)
 {
     int  i, j, k, atomnr;
@@ -891,15 +873,25 @@ int read_pdbfile(FILE*      in,
     title[0] = '\0';
     natom    = 0;
     chainnum = 0;
+    static const gmx::StringToEnumValueConverter<PdbRecordType, enumValueToString, gmx::StringCompareType::Exact, gmx::StripStrings::Yes> pdbRecordTypeIdentifier;
     while (!bStop && (fgets2(line, STRLEN, in) != nullptr))
     {
-        PdbRecordType line_type = line2type(line);
+        // Convert line to a null-terminated string, then take a substring of at most 6 chars
+        std::string                  lineAsString(line);
+        std::optional<PdbRecordType> line_type =
+                pdbRecordTypeIdentifier.valueFrom(lineAsString.substr(0, 6));
+        if (!line_type)
+        {
+            // Skip this line because it does not start with a
+            // recognized leading string describing a PDB record type.
+            continue;
+        }
 
-        switch (line_type)
+        switch (line_type.value())
         {
             case PdbRecordType::Atom:
             case PdbRecordType::Hetatm:
-                natom = read_atom(symtab, line, line_type, natom, atoms, x, chainnum);
+                natom = read_atom(symtab, line, line_type.value(), natom, atoms, x, chainnum);
                 break;
 
             case PdbRecordType::Anisou:
index ef613221efb7821bdd7a6d2b199d09e88f5ece07..42fd101db5d585326613a8d3574a51c45123308b 100644 (file)
 #include "gromacs/topology/symtab.h"
 #include "gromacs/utility/arrayref.h"
 #include "gromacs/utility/cstringutil.h"
+#include "gromacs/utility/enumerationhelpers.h"
 #include "gromacs/utility/fatalerror.h"
 #include "gromacs/utility/gmxassert.h"
 #include "gromacs/utility/smalloc.h"
+#include "gromacs/utility/stringtoenumvalueconverter.h"
 #include "gromacs/utility/stringutil.h"
 
 void generate_nbparams(CombinationRule         comb,
@@ -272,6 +274,15 @@ static void copy_B_from_A(int ftype, double* c)
     }
 }
 
+//! Local definition that supersedes the central one, as we only want the leading letter
+static const char* enumValueToLetterAsString(ParticleType enumValue)
+{
+    static constexpr gmx::EnumerationArray<ParticleType, const char*> particleTypeLetters = {
+        "A", "N", "S", "B", "V"
+    };
+    return particleTypeLetters[enumValue];
+}
+
 void push_at(t_symtab*                  symtab,
              PreprocessingAtomTypes*    at,
              PreprocessingBondAtomType* bondAtomType,
@@ -281,18 +292,6 @@ void push_at(t_symtab*                  symtab,
              t_nbparam***               pair,
              warninp*                   wi)
 {
-    struct t_xlate
-    {
-        const char*  entry;
-        ParticleType ptype;
-    };
-    gmx::EnumerationArray<ParticleType, t_xlate> xl;
-    xl[ParticleType::Atom]    = { "A", ParticleType::Atom };
-    xl[ParticleType::Nucleus] = { "N", ParticleType::Nucleus };
-    xl[ParticleType::Shell]   = { "S", ParticleType::Shell };
-    xl[ParticleType::Bond]    = { "B", ParticleType::Bond };
-    xl[ParticleType::VSite]   = { "V", ParticleType::VSite };
-
     int     nfields, nfp0 = -1;
     int     nread;
     char    type[STRLEN], btype[STRLEN], ptype[STRLEN];
@@ -531,19 +530,18 @@ void push_at(t_symtab*                  symtab,
     {
         sprintf(ptype, "V");
     }
-    auto entry = std::find_if(xl.begin(), xl.end(), [ptype](const t_xlate& type) {
-        return gmx_strcasecmp(ptype, type.entry) == 0;
-    });
-    if (entry == xl.end())
+    static const gmx::StringToEnumValueConverter<ParticleType, enumValueToLetterAsString, gmx::StringCompareType::CaseInsensitive, gmx::StripStrings::No>
+                                s_stringToParticleType;
+    std::optional<ParticleType> pt = s_stringToParticleType.valueFrom(ptype);
+    if (!pt)
     {
         auto message = gmx::formatString("Invalid particle type %s", ptype);
         warning_error_and_exit(wi, message, FARGS);
     }
-    ParticleType pt = entry->ptype;
 
     atom->q     = q;
     atom->m     = m;
-    atom->ptype = pt;
+    atom->ptype = pt.value();
     for (int i = 0; i < MAXFORCEPARAM; i++)
     {
         forceParam[i] = c[i];
diff --git a/src/gromacs/utility/stringtoenumvalueconverter.h b/src/gromacs/utility/stringtoenumvalueconverter.h
new file mode 100644 (file)
index 0000000..bc45039
--- /dev/null
@@ -0,0 +1,150 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2021, 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 Defines helper function object for class enumerations
+ *
+ * This helper type facilitates efficient lookup of values from
+ * an enumeration identified by a string, given a pre-declared mapping of
+ * values to such strings.
+ *
+ * It is separated from enumerationhelpers.h because it is not as
+ * widely used and brings in several extra dependencies not needed by
+ * that header.
+ *
+ * \author Mark Abraham <mark.j.abraham@gmail.com>
+ * \inlibraryapi
+ * \ingroup module_utility
+ */
+#ifndef GMX_UTILITY_STRINGTOENUMVALUECONVERTER_H
+#define GMX_UTILITY_STRINGTOENUMVALUECONVERTER_H
+
+#include <map>
+#include <optional>
+#include <string>
+
+#include "enumerationhelpers.h"
+#include "stringcompare.h"
+#include "stringutil.h"
+
+namespace gmx
+{
+
+/*! \brief Enum class for whether StringToEnumValueConverter will strip strings
+ * of leading and trailing whitespace before comparison. */
+enum class StripStrings : int
+{
+    //! Do not strip strings
+    No,
+    //! Strip strings
+    Yes
+};
+
+/*! \brief A class to convert a string to an enum value of type \c EnumType.
+ *
+ * It can be configured:
+ *  - to match different enumerations,
+ *  - to convert enumerations to strings in a custom way,
+ *  - either strip strings of leading and trailing whitespace before
+ *    comparison or not,
+ *  - with different handling of how the string characters are compared.
+ *
+ * Usage example:
+ *
+ *    enum class Foo : int {
+ *       Fizz, Buzz, Count, Default = Fizz
+ *    };
+ *    StringToEnumValueConverter<Foo, enumValueToString> converter;
+ *    Foo type = converter.valueFrom(theString);
+ *
+ * \tparam EnumType                    A class enum for which enumValueToString
+ *                                     is defined and maps all values
+ *                                     (except EnumType::Count) to a string.
+ * \tparam enumValueToStringFunction   Function to convert EnumValue to string, which
+ *                                     is typically enumValueToString, per convention
+ * \tparam stringCompareType           Indicates how the string should be compared
+ *                                     with respect to case, hyphens, underscores, etc.
+ * \tparam stripStrings                Indicates whether strings should have leading
+ *                                     and trailing whitespace removed before comparison
+ */
+template<typename EnumType,
+         const char*(enumValueToStringFunction)(EnumType),
+         StringCompareType stringCompareType = StringCompareType::Exact,
+         StripStrings      stripStrings      = StripStrings::No>
+class StringToEnumValueConverter
+{
+public:
+    StringToEnumValueConverter() : stringToEnumValue_(stringCompareType)
+    {
+        for (const auto type : EnumerationWrapper<EnumType>{})
+        {
+            GMX_RELEASE_ASSERT(type != EnumType::Count,
+                               "EnumerationWrapper<EnumType> should never return EnumType::Count");
+            std::string stringFromType = enumValueToStringFunction(type);
+            if (stripStrings == StripStrings::Yes)
+            {
+                // Ensure leading and trailing spaces are removed
+                stringFromType = stripString(stringFromType);
+            }
+            stringToEnumValue_[stringFromType] = type;
+        }
+    }
+
+    //! Return an optional enum value identified from the \c s (which is never EnumType::Count)
+    std::optional<EnumType> valueFrom(const std::string& s) const
+    {
+        typename decltype(stringToEnumValue_)::const_iterator typeIt;
+        if (stripStrings == StripStrings::Yes)
+        {
+            // Ensure leading and trailing spaces are removed
+            typeIt = stringToEnumValue_.find(stripString(s));
+        }
+        else
+        {
+            typeIt = stringToEnumValue_.find(s);
+        }
+        return (typeIt != stringToEnumValue_.end()) ? std::make_optional(typeIt->second) : std::nullopt;
+    }
+
+private:
+    /*! \brief Map of strings to enumeration values.
+     *
+     * By construction, those values never include
+     * EnumType::Count. */
+    std::map<std::string, EnumType, StringCompare> stringToEnumValue_;
+};
+
+} // namespace gmx
+
+#endif
index 88506a393a3e55b7af0549aed54ed5472ae61480..856a949540ac8f4da493f465103102dcc2784ff0 100644 (file)
@@ -52,6 +52,7 @@ gmx_add_unit_test(UtilityUnitTests utility-test
         physicalnodecommunicator.cpp
         range.cpp
         strconvert.cpp
+        stringtoenumvalueconverter.cpp
         stringutil.cpp
         template_mp.cpp
         textreader.cpp
diff --git a/src/gromacs/utility/tests/stringtoenumvalueconverter.cpp b/src/gromacs/utility/tests/stringtoenumvalueconverter.cpp
new file mode 100644 (file)
index 0000000..7a984d9
--- /dev/null
@@ -0,0 +1,230 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2021, 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 Tests for string-to-enum-value helper functor
+ *
+ * \author Mark Abraham <mark.j.abraham@gmail.com>
+ * \ingroup module_utility
+ */
+#include "gmxpre.h"
+
+#include "gromacs/utility/stringtoenumvalueconverter.h"
+
+#include <gtest/gtest.h>
+
+namespace gmx
+{
+
+namespace
+{
+
+//! Type to use in testing
+enum class Foo : int
+{
+    Bar,
+    Ugh,
+    Fooz,
+    Weird,
+    Empty,
+    Count,
+    Default = Ugh
+};
+
+//! Declare the conventional conversion function
+const char* enumValueToString(Foo f)
+{
+    static constexpr gmx::EnumerationArray<Foo, const char*> toString = {
+        "Bar", "Ugh ", "Foo z", "We-i_rd", ""
+    };
+    return toString[f];
+}
+
+//! Declare an atypical conversion function
+const char* enumValueToLetterAsString(Foo f)
+{
+    static constexpr gmx::EnumerationArray<Foo, const char*> toString = { "B", "U", "F", "W", "" };
+    return toString[f];
+}
+
+//! Pretty-printer for GoogleTest to use
+::std::ostream& operator<<(::std::ostream& os, const std::optional<Foo>& f)
+{
+    if (f)
+    {
+        return os << enumValueToString(f.value());
+    }
+    else
+    {
+        return os << "Out-of-range Foo value";
+    }
+}
+
+TEST(StringToEnumValueConverterTest, ExactStringComparisonWorksWithoutStripping)
+{
+    StringToEnumValueConverter<Foo, enumValueToString, StringCompareType::Exact, StripStrings::No> converter;
+    EXPECT_EQ(converter.valueFrom("Bar"), Foo::Bar);
+    EXPECT_FALSE(converter.valueFrom("bar"));
+    EXPECT_EQ(converter.valueFrom("Ugh "), Foo::Ugh);
+    EXPECT_FALSE(converter.valueFrom("ugh "));
+    EXPECT_EQ(converter.valueFrom("Foo z"), Foo::Fooz);
+    EXPECT_FALSE(converter.valueFrom("foo z"));
+    EXPECT_EQ(converter.valueFrom("We-i_rd"), Foo::Weird);
+    EXPECT_FALSE(converter.valueFrom("we-i_rd"));
+    EXPECT_EQ(converter.valueFrom(""), Foo::Empty);
+    EXPECT_FALSE(converter.valueFrom("count"));
+    EXPECT_FALSE(converter.valueFrom("Unknown"));
+    EXPECT_FALSE(converter.valueFrom("BarFoo z"));
+    EXPECT_FALSE(converter.valueFrom("Ugh"));
+    EXPECT_FALSE(converter.valueFrom("Default"));
+}
+
+TEST(StringToEnumValueConverterTest, CaseInsensitiveStringComparisonWorksWithoutStripping)
+{
+    StringToEnumValueConverter<Foo, enumValueToString, StringCompareType::CaseInsensitive, StripStrings::No> converter;
+    EXPECT_EQ(converter.valueFrom("Bar"), Foo::Bar);
+    EXPECT_EQ(converter.valueFrom("bar"), Foo::Bar);
+    EXPECT_EQ(converter.valueFrom("Ugh "), Foo::Ugh);
+    EXPECT_EQ(converter.valueFrom("ugh "), Foo::Ugh);
+    EXPECT_EQ(converter.valueFrom("Foo z"), Foo::Fooz);
+    EXPECT_EQ(converter.valueFrom("foo z"), Foo::Fooz);
+    EXPECT_EQ(converter.valueFrom("We-i_rd"), Foo::Weird);
+    EXPECT_EQ(converter.valueFrom("we-i_rd"), Foo::Weird);
+    EXPECT_EQ(converter.valueFrom(""), Foo::Empty);
+    EXPECT_FALSE(converter.valueFrom("Count"));
+    EXPECT_FALSE(converter.valueFrom("count"));
+    EXPECT_FALSE(converter.valueFrom("Unknown"));
+    EXPECT_FALSE(converter.valueFrom("barfoo z"));
+    EXPECT_FALSE(converter.valueFrom("Ugh"));
+    EXPECT_FALSE(converter.valueFrom("Default"));
+}
+
+TEST(StringToEnumValueConverterTest, CaseAndDashInsensitiveStringComparisonWorksWithoutStripping)
+{
+    StringToEnumValueConverter<Foo, enumValueToString, StringCompareType::CaseAndDashInsensitive, StripStrings::No> converter;
+    EXPECT_EQ(converter.valueFrom("Bar"), Foo::Bar);
+    EXPECT_EQ(converter.valueFrom("b-ar"), Foo::Bar);
+    EXPECT_EQ(converter.valueFrom("Ugh "), Foo::Ugh);
+    EXPECT_EQ(converter.valueFrom("_ugh "), Foo::Ugh);
+    EXPECT_EQ(converter.valueFrom("Foo z"), Foo::Fooz);
+    EXPECT_EQ(converter.valueFrom("fo_o z"), Foo::Fooz);
+    EXPECT_EQ(converter.valueFrom("We-i_rd"), Foo::Weird);
+    EXPECT_EQ(converter.valueFrom("we-i_rd"), Foo::Weird);
+    EXPECT_EQ(converter.valueFrom(""), Foo::Empty);
+    EXPECT_FALSE(converter.valueFrom("Count"));
+    EXPECT_FALSE(converter.valueFrom("count"));
+    EXPECT_FALSE(converter.valueFrom("Unknown"));
+    EXPECT_FALSE(converter.valueFrom("Bar-Foo z"));
+    EXPECT_FALSE(converter.valueFrom("Ugh"));
+    EXPECT_FALSE(converter.valueFrom("Default"));
+}
+
+TEST(StringToEnumValueConverterTest, ExactStringComparisonWorksWithStripping)
+{
+    StringToEnumValueConverter<Foo, enumValueToString, StringCompareType::Exact, StripStrings::Yes> converter;
+    EXPECT_EQ(converter.valueFrom("Bar "), Foo::Bar);
+    EXPECT_FALSE(converter.valueFrom("Ba r"));
+    EXPECT_EQ(converter.valueFrom("Ugh"), Foo::Ugh);
+    EXPECT_FALSE(converter.valueFrom("ugh"));
+    EXPECT_EQ(converter.valueFrom("  Foo z "), Foo::Fooz);
+    EXPECT_FALSE(converter.valueFrom(" foo z"));
+    EXPECT_EQ(converter.valueFrom("We-i_rd  "), Foo::Weird);
+    EXPECT_FALSE(converter.valueFrom("  we-i_rd "));
+    EXPECT_EQ(converter.valueFrom(""), Foo::Empty);
+    EXPECT_FALSE(converter.valueFrom(" Count"));
+    EXPECT_FALSE(converter.valueFrom("count  "));
+    EXPECT_FALSE(converter.valueFrom("Unknown"));
+    EXPECT_FALSE(converter.valueFrom("Bar Foo z"));
+    EXPECT_EQ(converter.valueFrom(" Ugh"), Foo::Ugh);
+    EXPECT_FALSE(converter.valueFrom("Default"));
+}
+
+TEST(StringToEnumValueConverterTest, CaseInsensitiveStringComparisonWorksWithStripping)
+{
+    StringToEnumValueConverter<Foo, enumValueToString, StringCompareType::CaseInsensitive, StripStrings::Yes> converter;
+    EXPECT_EQ(converter.valueFrom("Bar "), Foo::Bar);
+    EXPECT_FALSE(converter.valueFrom("ba r"));
+    EXPECT_EQ(converter.valueFrom("Ugh "), Foo::Ugh);
+    EXPECT_EQ(converter.valueFrom("  ugh "), Foo::Ugh);
+    EXPECT_EQ(converter.valueFrom("  Foo z "), Foo::Fooz);
+    EXPECT_EQ(converter.valueFrom("foo z   "), Foo::Fooz);
+    EXPECT_EQ(converter.valueFrom("We-i_rd  "), Foo::Weird);
+    EXPECT_EQ(converter.valueFrom("  we-i_rd "), Foo::Weird);
+    EXPECT_EQ(converter.valueFrom(""), Foo::Empty);
+    EXPECT_FALSE(converter.valueFrom("  Count"));
+    EXPECT_FALSE(converter.valueFrom("count "));
+    EXPECT_FALSE(converter.valueFrom("Unknown"));
+    EXPECT_FALSE(converter.valueFrom("  barfoo z"));
+    EXPECT_EQ(converter.valueFrom(" Ugh"), Foo::Ugh);
+    EXPECT_FALSE(converter.valueFrom("Default"));
+}
+
+TEST(StringToEnumValueConverterTest, CaseAndDashInsensitiveStringComparisonWorksWithStripping)
+{
+    StringToEnumValueConverter<Foo, enumValueToString, StringCompareType::CaseAndDashInsensitive, StripStrings::Yes> converter;
+    EXPECT_EQ(converter.valueFrom("Bar "), Foo::Bar);
+    EXPECT_FALSE(converter.valueFrom("b-a r"));
+    EXPECT_EQ(converter.valueFrom("Ugh "), Foo::Ugh);
+    EXPECT_EQ(converter.valueFrom(" _ugh "), Foo::Ugh);
+    EXPECT_EQ(converter.valueFrom(" Foo z "), Foo::Fooz);
+    EXPECT_EQ(converter.valueFrom("fo_o z  "), Foo::Fooz);
+    EXPECT_EQ(converter.valueFrom("We-i_rd  "), Foo::Weird);
+    EXPECT_EQ(converter.valueFrom("  we-i_rd "), Foo::Weird);
+    EXPECT_EQ(converter.valueFrom(""), Foo::Empty);
+    EXPECT_FALSE(converter.valueFrom("  Count"));
+    EXPECT_FALSE(converter.valueFrom("coun-t "));
+    EXPECT_FALSE(converter.valueFrom("Unknown"));
+    EXPECT_FALSE(converter.valueFrom("  Bar-Foo z   "));
+    EXPECT_EQ(converter.valueFrom("Ugh  "), Foo::Ugh);
+    EXPECT_FALSE(converter.valueFrom("Default"));
+}
+
+TEST(StringToEnumValueConverterTest, CustomConverterWorks)
+{
+    StringToEnumValueConverter<Foo, enumValueToLetterAsString> converter;
+    EXPECT_EQ(converter.valueFrom("B"), Foo::Bar);
+    EXPECT_FALSE(converter.valueFrom("b"));
+    EXPECT_EQ(converter.valueFrom("U"), Foo::Ugh);
+    EXPECT_FALSE(converter.valueFrom("u"));
+    EXPECT_EQ(converter.valueFrom("F"), Foo::Fooz);
+    EXPECT_FALSE(converter.valueFrom("f"));
+    EXPECT_EQ(converter.valueFrom(""), Foo::Empty);
+    EXPECT_FALSE(converter.valueFrom("C"));
+    EXPECT_FALSE(converter.valueFrom("c"));
+    EXPECT_FALSE(converter.valueFrom("X"));
+    EXPECT_FALSE(converter.valueFrom("Default"));
+}
+
+} // namespace
+} // namespace gmx