Add a simple C++ file wrapper class.
authorTeemu Murtola <teemu.murtola@gmail.com>
Sun, 22 Apr 2012 07:52:43 +0000 (10:52 +0300)
committerTeemu Murtola <teemu.murtola@gmail.com>
Mon, 23 Apr 2012 13:43:20 +0000 (16:43 +0300)
This helps in writing exception-safe code, and also takes care of error
reporting (using exceptions).  Currently only a minimal implementation,
but more methods can be added when needed.

Also made a minor adjustmend to the internal boost to avoid compiler
warnings.

Change-Id: I751f1ad242e7f94d06a9165dcca57525b6a99c96

src/external/boost/README
src/external/boost/boost/exception/detail/error_info_impl.hpp
src/gromacs/selection/selectioncollection.cpp
src/gromacs/utility/CMakeLists.txt
src/gromacs/utility/exceptions.cpp
src/gromacs/utility/exceptions.h
src/gromacs/utility/file.cpp [new file with mode: 0644]
src/gromacs/utility/file.h [new file with mode: 0644]

index b69a87543aee03d7dc92d79e3443bad6dee226cc..39bf8d155980dcdce8172491f6f25d749cc5e908 100644 (file)
@@ -10,4 +10,5 @@ Steps to produce minimal version of BOOST:
    Make sure that they are really not needed and that activating the feature they implement gives a compiler
    not linker error (by adding an appropriate #error directive). If any source files are added make sure 
    to add them to cmake.
-
+6) Check that boost/exception/detail/error_info_impl.hpp has virtual
+   destructors to avoid compiler warnings.
index f4fe8baf63b1cb020e733b5df15cacfeaaa44ed3..b730790578fd4db423353af927422ecd1f74c4c9 100644 (file)
@@ -30,7 +30,7 @@ boost
 \r
             protected:\r
 \r
-            ~error_info_base() throw()\r
+            virtual ~error_info_base() throw()\r
                 {\r
                 }\r
             };\r
@@ -46,7 +46,7 @@ boost
         typedef T value_type;\r
 \r
         error_info( value_type const & value );\r
-        ~error_info() throw();\r
+        virtual ~error_info() throw();\r
 \r
         value_type const &\r
         value() const\r
index 32dd0144afd4e550e27c32a7c58aa7e805e85b28..2f5d77a734ac845aa06eebcc8705f0cb0caf93bb 100644 (file)
@@ -51,6 +51,7 @@
 #include "gromacs/selection/selection.h"
 #include "gromacs/selection/selectioncollection.h"
 #include "gromacs/utility/exceptions.h"
+#include "gromacs/utility/file.h"
 #include "gromacs/utility/gmxassert.h"
 #include "gromacs/utility/messagestringcollector.h"
 
@@ -491,24 +492,14 @@ SelectionCollection::parseFromFile(const std::string &filename,
                                    SelectionList *output)
 {
     yyscan_t scanner;
-    FILE *fp;
 
     _gmx_sel_init_lexer(&scanner, &_impl->_sc, false, -1,
                         _impl->_bExternalGroupsSet,
                         _impl->_grps);
-    fp = ffopen(filename.c_str(), "r");
-    _gmx_sel_set_lex_input_file(scanner, fp);
-    // TODO: Use RAII
-    try
-    {
-        _impl->runParser(scanner, -1, output);
-    }
-    catch (...)
-    {
-        ffclose(fp);
-        throw;
-    }
-    ffclose(fp);
+    File file(filename, "r");
+    _gmx_sel_set_lex_input_file(scanner, file.handle());
+    _impl->runParser(scanner, -1, output);
+    file.close();
 }
 
 
index 0e25ccfe64061448c5194c2517d5da9a405fcdfd..0a56432656d9d8b398fd4297438fea43dc36b583 100644 (file)
@@ -6,6 +6,7 @@ set(UTILITY_PUBLIC_HEADERS
     common.h
     errorcodes.h
     exceptions.h
+    file.h
     flags.h
     format.h
     gmxassert.h
index 510297162febcdd24b84a4a5e972be64da150da2..9074acf66bb1fb73ccc97710d173ce0681ae6428 100644 (file)
@@ -120,6 +120,7 @@ std::string formatErrorMessage(const std::exception &ex)
         file = *boost::get_error_info<boost::throw_file>(*gmxEx);
         line = *boost::get_error_info<boost::throw_line>(*gmxEx);
     }
+    // TODO: Treat errno information in boost exceptions
     return internal::formatFatalError(title, ex.what(), func, file, line);
 }
 
index 2cf437663c46767c20e4baee294282ae5d7e12b0..4415f1ab1dedf9620ea3ac73effdc88abad873e2 100644 (file)
@@ -42,6 +42,8 @@
 #include <exception>
 #include <string>
 
+#include <boost/exception/errinfo_api_function.hpp>
+#include <boost/exception/errinfo_errno.hpp>
 #include <boost/exception/exception.hpp>
 #include <boost/exception/info.hpp>
 #include <boost/throw_exception.hpp>
@@ -249,6 +251,37 @@ if (value < 0)
 #define GMX_THROW(e) \
     BOOST_THROW_EXCEPTION((e))
 
+/*! \brief
+ * Macro for throwing an exception based on errno.
+ *
+ * \param[in] e       Exception object to throw.
+ * \param[in] syscall Name of the syscall that returned the error.
+ * \param[in] err     errno value returned by the syscall.
+ *
+ * This macro provides a convenience interface for throwing an exception to
+ * report an error based on a errno value.  In addition to adding the necessary
+ * information to the exception object, the macro also ensures that \p errno is
+ * evaluated before, e.g., the constructor of \p e may call other functions
+ * that could overwrite the errno value.
+ * \p e should evaluate to an instance of an object derived from
+ * GromacsException.
+ *
+ * Typical usage (note that gmx::File wraps this particular case):
+ * \code
+FILE *fp = fopen("filename.txt", "r");
+if (fp == NULL)
+{
+    GMX_THROW(FileIOError("Could not open file"), "fopen", errno);
+}
+ * \endcode
+ */
+#define GMX_THROW_WITH_ERRNO(e, syscall, err) \
+    do { \
+        int stored_errno_ = (err); \
+        GMX_THROW((e) << boost::errinfo_errno(stored_errno_) \
+                      << boost::errinfo_api_function(syscall)); \
+    } while(0)
+
 /*! \brief
  * Formats a standard error message for reporting an error.
  *
diff --git a/src/gromacs/utility/file.cpp b/src/gromacs/utility/file.cpp
new file mode 100644 (file)
index 0000000..39209f0
--- /dev/null
@@ -0,0 +1,161 @@
+/*
+ *
+ *                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::File.
+ *
+ * \author Teemu Murtola <teemu.murtola@cbr.su.se>
+ * \ingroup module_utility
+ */
+#include "file.h"
+
+#include <cerrno>
+#include <cstdio>
+
+#include <string>
+#include <vector>
+
+#include "gromacs/utility/exceptions.h"
+#include "gromacs/utility/gmxassert.h"
+#include "gromacs/utility/format.h"
+
+namespace gmx
+{
+
+File::File(const char *filename, const char *mode)
+    : fp_(NULL)
+{
+    open(filename, mode);
+}
+
+File::File(const std::string &filename, const char *mode)
+    : fp_(NULL)
+{
+    open(filename, mode);
+}
+
+File::~File()
+{
+    if (fp_ != NULL)
+    {
+        if (fclose(fp_) != 0)
+        {
+            // TODO: Log the error somewhere
+        }
+    }
+}
+
+void File::open(const char *filename, const char *mode)
+{
+    GMX_RELEASE_ASSERT(fp_ == NULL,
+                       "Attempted to open the same file object twice");
+    // TODO: Port all necessary functionality from ffopen() here.
+    fp_ = fopen(filename, mode);
+    if (fp_ == NULL)
+    {
+        GMX_THROW_WITH_ERRNO(
+                FileIOError(formatString("Could not open file '%s'", filename)),
+                "fopen", errno);
+    }
+}
+
+void File::open(const std::string &filename, const char *mode)
+{
+    open(filename.c_str(), mode);
+}
+
+void File::close()
+{
+    GMX_RELEASE_ASSERT(fp_ != NULL,
+                       "Attempted to close a file object that is not open");
+    bool bOk = (fclose(fp_) == 0);
+    fp_ = NULL;
+    if (!bOk)
+    {
+        GMX_THROW_WITH_ERRNO(
+                FileIOError("Error while closing file"), "fclose", errno);
+    }
+}
+
+FILE *File::handle()
+{
+    GMX_RELEASE_ASSERT(fp_ != NULL,
+                       "Attempted to access a file object that is not open");
+    return fp_;
+}
+
+void File::readBytes(void *buffer, size_t bytes)
+{
+    GMX_RELEASE_ASSERT(fp_ != NULL,
+                       "Attempted to access a file object that is not open");
+    errno = 0;
+    // TODO: Retry based on errno or something else?
+    size_t bytesRead = std::fread(buffer, 1, bytes, fp_);
+    if (bytesRead != bytes)
+    {
+        if (feof(fp_))
+        {
+            GMX_THROW(FileIOError("Premature end of file"));
+        }
+        else
+        {
+            GMX_THROW_WITH_ERRNO(FileIOError("Error while reading file"),
+                                 "fread", errno);
+        }
+    }
+}
+
+// static
+std::string File::readToString(const char *filename)
+{
+    File file(filename, "r");
+    FILE *fp = file.handle();
+
+    // TODO: Full error checking.
+    std::fseek(fp, 0L, SEEK_END);
+    long len = std::ftell(fp);
+    std::fseek(fp, 0L, SEEK_SET);
+
+    std::vector<char> data(len);
+    file.readBytes(&data[0], len);
+    std::string result(&data[0], len);
+
+    file.close();
+    return result;
+}
+
+// static
+std::string File::readToString(const std::string &filename)
+{
+    return readToString(filename.c_str());
+}
+
+} // namespace gmx
diff --git a/src/gromacs/utility/file.h b/src/gromacs/utility/file.h
new file mode 100644 (file)
index 0000000..f158e44
--- /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
+ */
+/*! \file
+ * \brief
+ * Declares functions for file handling.
+ *
+ * \author Teemu Murtola <teemu.murtola@cbr.su.se>
+ * \inpublicapi
+ * \ingroup module_utility
+ */
+#ifndef GMX_UTILITY_FILE_H
+#define GMX_UTILITY_FILE_H
+
+#include <cstdio>
+
+#include <string>
+
+#include "common.h"
+
+namespace gmx
+{
+
+/*! \brief
+ * Basic file object.
+ *
+ * This class provides basic file I/O functionality and uses exceptions
+ * (FileIOError) for error reporting.
+ *
+ * \inpublicapi
+ * \ingroup module_utility
+ */
+class File
+{
+    public:
+        /*! \brief
+         * Creates a file object and opens a file.
+         *
+         * \param[in] filename  Path of the file to open.
+         * \param[in] mode      Mode to open the file in (for fopen()).
+         * \throws    FileIOError on any I/O error.
+         *
+         * \see open(const char *, const char *)
+         */
+        File(const char *filename, const char *mode);
+        //! \copydoc File(const char *, const char *)
+        File(const std::string &filename, const char *mode);
+        /*! \brief
+         * Destroys the file object.
+         *
+         * If the file is still open, it is closed.
+         * Any error conditions will be ignored.
+         */
+        ~File();
+
+        /*! \brief
+         * Opens a file.
+         *
+         * \param[in] filename  Path of the file to open.
+         * \param[in] mode      Mode to open the file in (for fopen()).
+         * \throws    FileIOError on any I/O error.
+         *
+         * The file object must not be open.
+         */
+        void open(const char *filename, const char *mode);
+        //! \copydoc open(const char *, const char *)
+        void open(const std::string &filename, const char *mode);
+        /*! \brief
+         * Closes the file object.
+         *
+         * \throws  FileIOError on any I/O error.
+         *
+         * The file must be open.
+         */
+        void close();
+
+        /*! \brief
+         * Returns a file handle for interfacing with C functions.
+         *
+         * The file must be open.
+         * Does not throw.
+         */
+        FILE *handle();
+
+        /*! \brief
+         * Reads given number of bytes from the file.
+         *
+         * \param[out] buffer  Pointer to buffer that receives the bytes.
+         * \param[in]  bytes   Number of bytes to read.
+         * \throws     FileIOError on any I/O error.
+         *
+         * The file must be open.
+         */
+        void readBytes(void *buffer, size_t bytes);
+
+        /*! \brief
+         * Reads contents of a file to a std::string.
+         *
+         * \param[in] filename  File to read.
+         * \returns   The contents of \p filename.
+         * \throws    std::bad_alloc if out of memory.
+         * \throws    FileIOError on any I/O error.
+         */
+        static std::string readToString(const char *filename);
+        //! \copydoc readToString(const char *)
+        static std::string readToString(const std::string &filename);
+
+    private:
+        FILE                   *fp_;
+
+        GMX_DISALLOW_COPY_AND_ASSIGN(File);
+};
+
+} // namespace gmx
+
+#endif