Extract string parsing routines to utility
authorTeemu Murtola <teemu.murtola@gmail.com>
Thu, 30 Jun 2016 18:55:09 +0000 (21:55 +0300)
committerTeemu Murtola <teemu.murtola@gmail.com>
Thu, 30 Jun 2016 18:57:23 +0000 (21:57 +0300)
Extract some string parsing utility functions from the options module to
utility/ for reuse elsewhere.

Change-Id: Ie5d779c99644f25faffc877a75c55d4c67da9fb2

src/gromacs/options/basicoptions.cpp
src/gromacs/utility.h
src/gromacs/utility/strconvert.cpp [new file with mode: 0644]
src/gromacs/utility/strconvert.h [new file with mode: 0644]

index 1a09d40569fbd53474e2d370ef55982a5956613f..9679510d601fe0b3e3534b4d804500fa03f5dbeb 100644 (file)
@@ -1,7 +1,7 @@
 /*
  * This file is part of the GROMACS molecular simulation package.
  *
- * Copyright (c) 2010,2011,2012,2013,2014,2015, by the GROMACS development team, led by
+ * Copyright (c) 2010,2011,2012,2013,2014,2015,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.
@@ -54,6 +54,7 @@
 
 #include "gromacs/utility/cstringutil.h"
 #include "gromacs/utility/exceptions.h"
+#include "gromacs/utility/strconvert.h"
 #include "gromacs/utility/stringutil.h"
 
 #include "basicoptionstorage.h"
@@ -197,23 +198,7 @@ std::string IntegerOptionStorage::formatSingleValue(const int &value) const
 
 void IntegerOptionStorage::convertValue(const std::string &value)
 {
-    const char *ptr = value.c_str();
-    char       *endptr;
-    errno = 0;
-    long int    ival = std::strtol(ptr, &endptr, 10);
-    if (errno == ERANGE
-        || ival < std::numeric_limits<int>::min()
-        || ival > std::numeric_limits<int>::max())
-    {
-        GMX_THROW(InvalidInputError("Invalid value: '" + value
-                                    + "'; it causes an integer overflow"));
-    }
-    if (*ptr == '\0' || *endptr != '\0')
-    {
-        GMX_THROW(InvalidInputError("Invalid value: '" + value
-                                    + "'; expected an integer"));
-    }
-    addValue(ival);
+    addValue(fromString<int>(value));
 }
 
 void IntegerOptionStorage::processSetValues(ValueList *values)
@@ -255,21 +240,7 @@ std::string Int64OptionStorage::formatSingleValue(const gmx_int64_t &value) cons
 
 void Int64OptionStorage::convertValue(const std::string &value)
 {
-    const char       *ptr = value.c_str();
-    char             *endptr;
-    errno = 0;
-    const gmx_int64_t ival = str_to_int64_t(ptr, &endptr);
-    if (errno == ERANGE)
-    {
-        GMX_THROW(InvalidInputError("Invalid value: '" + value
-                                    + "'; it causes an integer overflow"));
-    }
-    if (*ptr == '\0' || *endptr != '\0')
-    {
-        GMX_THROW(InvalidInputError("Invalid value: '" + value
-                                    + "'; expected an integer"));
-    }
-    addValue(ival);
+    addValue(fromString<gmx_int64_t>(value));
 }
 
 /********************************************************************
@@ -313,21 +284,7 @@ std::string DoubleOptionStorage::formatSingleValue(const double &value) const
 
 void DoubleOptionStorage::convertValue(const std::string &value)
 {
-    const char *ptr = value.c_str();
-    char       *endptr;
-    errno = 0;
-    double      dval = std::strtod(ptr, &endptr);
-    if (errno == ERANGE)
-    {
-        GMX_THROW(InvalidInputError("Invalid value: '" + value
-                                    + "'; it causes an overflow/underflow"));
-    }
-    if (*ptr == '\0' || *endptr != '\0')
-    {
-        GMX_THROW(InvalidInputError("Invalid value: '" + value
-                                    + "'; expected a number"));
-    }
-    addValue(dval * factor_);
+    addValue(fromString<double>(value) * factor_);
 }
 
 void DoubleOptionStorage::processSetValues(ValueList *values)
@@ -415,23 +372,9 @@ std::string FloatOptionStorage::formatSingleValue(const float &value) const
 
 void FloatOptionStorage::convertValue(const std::string &value)
 {
-    const char *ptr = value.c_str();
-    char       *endptr;
-    errno = 0;
-    double      dval = std::strtod(ptr, &endptr);
-    if (errno == ERANGE
-        || dval * factor_ < -std::numeric_limits<float>::max()
-        || dval * factor_ >  std::numeric_limits<float>::max())
-    {
-        GMX_THROW(InvalidInputError("Invalid value: '" + value
-                                    + "'; it causes an overflow/underflow"));
-    }
-    if (*ptr == '\0' || *endptr != '\0')
-    {
-        GMX_THROW(InvalidInputError("Invalid value: '" + value
-                                    + "'; expected a number"));
-    }
-    addValue(dval * factor_);
+    // TODO: Consider testing for overflow when scaling with factor_ (also for
+    // the double precision case).
+    addValue(fromString<float>(value) * factor_);
 }
 
 void FloatOptionStorage::processSetValues(ValueList *values)
index 6af9b3752b181f7c11a22b4c027e0c622a70d34f..609b906cbeaa7b23473ff8a0881e9f39e659338c 100644 (file)
  *
  * \if libapi
  *
+ * The header strconvert.h declares string parsing routines.
+ *
  * The header gmxmpi.h abstracts away the mechanism of including either MPI or
  * thread-MPI headers depending on compilation options.
  * Similarly, gmxomp.h removes the need to use conditional compilation for code
diff --git a/src/gromacs/utility/strconvert.cpp b/src/gromacs/utility/strconvert.cpp
new file mode 100644 (file)
index 0000000..40fe2aa
--- /dev/null
@@ -0,0 +1,138 @@
+/*
+ * 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 functions in strconvert.h.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \ingroup module_utility
+ */
+#include "gmxpre.h"
+
+#include "strconvert.h"
+
+#include <cerrno>
+#include <cstdlib>
+
+#include <limits>
+#include <string>
+
+#include "gromacs/utility/cstringutil.h"
+#include "gromacs/utility/exceptions.h"
+
+namespace gmx
+{
+
+//! \cond libapi
+
+int intFromString(const char *str)
+{
+    errno = 0;
+    char           *endptr;
+    const long int  value = std::strtol(str, &endptr, 10);
+    if (errno == ERANGE
+        || value < std::numeric_limits<int>::min()
+        || value > std::numeric_limits<int>::max())
+    {
+        GMX_THROW(InvalidInputError("Invalid value: '" + std::string(str)
+                                    + "'; it causes an integer overflow"));
+    }
+    if (str[0] == '\0' || *endptr != '\0')
+    {
+        GMX_THROW(InvalidInputError("Invalid value: '" + std::string(str)
+                                    + "'; expected an integer"));
+    }
+    return value;
+}
+
+gmx_int64_t int64FromString(const char *str)
+{
+    errno = 0;
+    char              *endptr;
+    const gmx_int64_t  value = str_to_int64_t(str, &endptr);
+    if (errno == ERANGE)
+    {
+        GMX_THROW(InvalidInputError("Invalid value: '" + std::string(str)
+                                    + "'; it causes an integer overflow"));
+    }
+    if (str[0] == '\0' || *endptr != '\0')
+    {
+        GMX_THROW(InvalidInputError("Invalid value: '" + std::string(str)
+                                    + "'; expected an integer"));
+    }
+    return value;
+}
+
+float floatFromString(const char *str)
+{
+    errno = 0;
+    char         *endptr;
+    const double  value = std::strtod(str, &endptr);
+    if (errno == ERANGE
+        || value < -std::numeric_limits<float>::max()
+        || value >  std::numeric_limits<float>::max())
+    {
+        GMX_THROW(InvalidInputError("Invalid value: '" + std::string(str)
+                                    + "'; it causes an overflow/underflow"));
+    }
+    if (str[0] == '\0' || *endptr != '\0')
+    {
+        GMX_THROW(InvalidInputError("Invalid value: '" + std::string(str)
+                                    + "'; expected a number"));
+    }
+    return value;
+}
+
+double doubleFromString(const char *str)
+{
+    errno = 0;
+    char         *endptr;
+    const double  value = std::strtod(str, &endptr);
+    if (errno == ERANGE)
+    {
+        GMX_THROW(InvalidInputError("Invalid value: '" + std::string(str)
+                                    + "'; it causes an overflow/underflow"));
+    }
+    if (str[0] == '\0' || *endptr != '\0')
+    {
+        GMX_THROW(InvalidInputError("Invalid value: '" + std::string(str)
+                                    + "'; expected a number"));
+    }
+    return value;
+}
+
+//! \endcond
+
+} // namespace gmx
diff --git a/src/gromacs/utility/strconvert.h b/src/gromacs/utility/strconvert.h
new file mode 100644 (file)
index 0000000..2a720b2
--- /dev/null
@@ -0,0 +1,124 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2016, by the GROMACS development team, led by
+ * Mark Abraham, David van der Spoel, Berk Hess, and Erik Lindahl,
+ * and including many others, as listed in the AUTHORS file in the
+ * top-level source directory and at http://www.gromacs.org.
+ *
+ * GROMACS is free software; you can redistribute it and/or
+ * modify it under the terms of the GNU Lesser General Public License
+ * as published by the Free Software Foundation; either version 2.1
+ * of the License, or (at your option) any later version.
+ *
+ * GROMACS is distributed in the hope that it will be useful,
+ * but WITHOUT ANY WARRANTY; without even the implied warranty of
+ * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
+ * Lesser General Public License for more details.
+ *
+ * You should have received a copy of the GNU Lesser General Public
+ * License along with GROMACS; if not, see
+ * http://www.gnu.org/licenses, or write to the Free Software Foundation,
+ * Inc., 51 Franklin Street, Fifth Floor, Boston, MA  02110-1301  USA.
+ *
+ * If you want to redistribute modifications to GROMACS, please
+ * consider that scientific software is very special. Version
+ * control is crucial - bugs must be traceable. We will be happy to
+ * consider code for inclusion in the official distribution, but
+ * derived work must not be called official GROMACS. Details are found
+ * in the README & COPYING files - if they are missing, get the
+ * official version at http://www.gromacs.org.
+ *
+ * To help us fund GROMACS development, we humbly ask that you cite
+ * the research papers on the package. Check out http://www.gromacs.org.
+ */
+/*! \libinternal \file
+ * \brief
+ * Declares common utility functions for conversions from strings.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \inlibraryapi
+ * \ingroup module_utility
+ */
+#ifndef GMX_UTILITY_STRCONVERT_H
+#define GMX_UTILITY_STRCONVERT_H
+
+#include <string>
+
+#include "gromacs/utility/basedefinitions.h"
+
+namespace gmx
+{
+
+//! \cond libapi
+//! \addtogroup module_utility
+//! \{
+
+/*! \brief
+ * Parses an integer from a string.
+ *
+ * \throws  InvalidInputError if `str` is not a valid integer.
+ *
+ * Also checks for overflow.
+ */
+int intFromString(const char *str);
+/*! \brief
+ * Parses a 64-bit integer from a string.
+ *
+ * \throws  InvalidInputError if `str` is not a valid integer.
+ *
+ * Also checks for overflow.
+ */
+gmx_int64_t int64FromString(const char *str);
+/*! \brief
+ * Parses a float value from a string.
+ *
+ * \throws  InvalidInputError if `str` is not a valid number.
+ *
+ * Also checks for overflow.
+ */
+float floatFromString(const char *str);
+/*! \brief
+ * Parses a double value from a string.
+ *
+ * \throws  InvalidInputError if `str` is not a valid number.
+ *
+ * Also checks for overflow.
+ */
+double doubleFromString(const char *str);
+
+/*! \brief
+ * Parses a value from a string to a given type.
+ *
+ * \tparam T Type of value to parse.
+ *
+ * `T` can only be one of the types that is explicity supported.
+ * The main use for this function is to write `fromString<real>(value)`,
+ * but it can also be used for other types for consistency.
+ */
+template <typename T> static inline T fromString(const char *str);
+//! \copydoc fromString(const char *)
+template <typename T> static inline T fromString(const std::string &str)
+{
+    return fromString<T>(str.c_str());
+}
+
+//! Implementation for integer values.
+template <> inline
+int fromString<int>(const char *str) { return intFromString(str); }
+//! Implementation for 64-bit integer values.
+template <> inline
+gmx_int64_t fromString<gmx_int64_t>(const char *str) { return int64FromString(str); }
+//! Implementation for float values.
+template <> inline
+float fromString<float>(const char *str) { return floatFromString(str); }
+//! Implementation for double values.
+template <> inline
+double fromString<double>(const char *str) { return doubleFromString(str); }
+
+//! \}
+//! \endcond
+
+} // namespace gmx
+
+#endif