Unit tests for parse_common_args()
authorTeemu Murtola <teemu.murtola@gmail.com>
Sat, 21 Dec 2013 18:33:44 +0000 (20:33 +0200)
committerRoland Schulz <roland@utk.edu>
Wed, 15 Jan 2014 22:53:09 +0000 (17:53 -0500)
Add unit tests for most of the basic command line parsing functionality
in parse_common_args(), including the file name handling.

To support these tests, extend arrayref.h to also include ArrayRef,
which provides a mutable wrapper for a C array.  This is complete enough
for the purpose of the tests, but may need additional work for more
general use.  Also fix some documentation typos noticed in the existing
code.

Also, fix unintuitive behavior of etTIME options that was revealed by
the tests.

Change-Id: I4b8b3ad92b0efe24fdf478cd39ba246692b036c7

src/gromacs/commandline/pargs.cpp
src/gromacs/commandline/tests/CMakeLists.txt
src/gromacs/commandline/tests/pargs.cpp [new file with mode: 0644]
src/gromacs/utility/arrayref.h

index 6b9f50ea5e0fa9236643e19710e4ca4033cd0212..c4e296a968950ef10840fb0f08c9808869d7b45f 100644 (file)
@@ -3,7 +3,7 @@
  *
  * Copyright (c) 1991-2000, University of Groningen, The Netherlands.
  * Copyright (c) 2001-2004, The GROMACS development team.
- * Copyright (c) 2013, by the GROMACS development team, led by
+ * Copyright (c) 2013,2014, 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.
@@ -866,7 +866,7 @@ gmx_bool parse_common_args(int *argc, char *argv[], unsigned long Flags,
 
     for (i = 0; i < npall; i++)
     {
-        if ((all_pa[i].type == etTIME) && (*all_pa[i].u.r >= 0))
+        if (all_pa[i].type == etTIME && all_pa[i].bSet)
         {
             *all_pa[i].u.r *= output_env_get_time_invfactor(*oenv);
         }
index 58b0021fe7eb0d83fd675d0f7dc31f4aa563c7a9..5d6404d2f2d645d1d0b127d90c86141f15b74943 100644 (file)
@@ -1,7 +1,7 @@
 #
 # This file is part of the GROMACS molecular simulation package.
 #
-# Copyright (c) 2012,2014, by the GROMACS development team, led by
+# Copyright (c) 2012,2013,2014, 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.
@@ -36,4 +36,5 @@ gmx_add_unit_test(CommandLineUnitTests commandline-test
                   cmdlinehelpwriter.cpp
                   cmdlinemodulemanager.cpp
                   cmdlineparser.cpp
+                  pargs.cpp
                   $<TARGET_OBJECTS:onlinehelp-test-shared>)
diff --git a/src/gromacs/commandline/tests/pargs.cpp b/src/gromacs/commandline/tests/pargs.cpp
new file mode 100644 (file)
index 0000000..24e253f
--- /dev/null
@@ -0,0 +1,369 @@
+/*
+ * This file is part of the GROMACS molecular simulation package.
+ *
+ * Copyright (c) 2013,2014, 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 parse_common_args().
+ *
+ * Currently, negative tests are not possible, because those trigger
+ * gmx_fatal() and terminate the whole test binary.
+ *
+ * \todo
+ * Add tests that exercise the machinery triggered by ffREAD to detect the
+ * extension for certain types of files.  Also some other paths in the file
+ * name logic may not get tested by the current set.
+ *
+ * \author Teemu Murtola <teemu.murtola@gmail.com>
+ * \ingroup module_commandline
+ */
+#include <gtest/gtest.h>
+
+#include "gromacs/commandline/pargs.h"
+#include "gromacs/utility/arrayref.h"
+
+#include "testutils/cmdlinetest.h"
+
+namespace
+{
+
+using gmx::test::CommandLine;
+
+class ParseCommonArgsTest : public ::testing::Test
+{
+    public:
+        ParseCommonArgsTest()
+            : oenv_(NULL), fileCount_(0)
+        {
+        }
+        virtual ~ParseCommonArgsTest()
+        {
+            output_env_done(oenv_);
+        }
+
+        int nfile() const { return fileCount_; }
+
+        void parse(gmx::ConstArrayRef<const char *> cmdline,
+                   unsigned long                    flags,
+                   gmx::ArrayRef<t_filenm>          fnm,
+                   gmx::ArrayRef<t_pargs>           pa)
+        {
+            args_.initFromArray(cmdline);
+            fileCount_ = fnm.size();
+            bool bOk = parse_common_args(&args_.argc(), args_.argv(), flags,
+                                         fnm.size(), fnm.data(),
+                                         pa.size(), pa.data(),
+                                         0, NULL, 0, NULL, &oenv_);
+            EXPECT_TRUE(bOk);
+        }
+
+        CommandLine              args_;
+        output_env_t             oenv_;
+
+    private:
+        size_t                   fileCount_;
+};
+
+/********************************************************************
+ * Tests for different types of options
+ */
+
+TEST_F(ParseCommonArgsTest, ParsesIntegerArgs)
+{
+    int               value1 = 0, value2 = 0, value3 = 3;
+    t_pargs           pa[]   = {
+        { "-i1", FALSE, etINT, {&value1}, "Description" },
+        { "-i2", FALSE, etINT, {&value2}, "Description" },
+        { "-i3", FALSE, etINT, {&value3}, "Description" }
+    };
+    const char *const cmdline[] = {
+        "test", "-i1", "2", "-i2", "-3"
+    };
+    parse(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    EXPECT_EQ( 2, value1);
+    EXPECT_EQ(-3, value2);
+    EXPECT_EQ( 3, value3);
+}
+
+TEST_F(ParseCommonArgsTest, ParsesInt64Args)
+{
+    gmx_int64_t       value1 = 0, value2 = 0, value3 = 3;
+    t_pargs           pa[]   = {
+        { "-i1", FALSE, etINT64, {&value1}, "Description" },
+        { "-i2", FALSE, etINT64, {&value2}, "Description" },
+        { "-i3", FALSE, etINT64, {&value3}, "Description" }
+    };
+    const char *const cmdline[] = {
+        "test", "-i1", "2", "-i2", "-3"
+    };
+    parse(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    EXPECT_EQ( 2, value1);
+    EXPECT_EQ(-3, value2);
+    EXPECT_EQ( 3, value3);
+}
+
+TEST_F(ParseCommonArgsTest, ParsesRealArgs)
+{
+    real              value1 = 0.0, value2 = 0.0, value3 = 2.5;
+    t_pargs           pa[]   = {
+        { "-r1", FALSE, etREAL, {&value1}, "Description" },
+        { "-r2", FALSE, etREAL, {&value2}, "Description" },
+        { "-r3", FALSE, etREAL, {&value3}, "Description" }
+    };
+    const char *const cmdline[] = {
+        "test", "-r1", "2", "-r2", "-.5"
+    };
+    parse(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    EXPECT_EQ( 2.0, value1);
+    EXPECT_EQ(-0.5, value2);
+    EXPECT_EQ( 2.5, value3);
+}
+
+TEST_F(ParseCommonArgsTest, ParsesStringArgs)
+{
+    const char       *value1 = "def", *value2 = "", *value3 = "default";
+    t_pargs           pa[]   = {
+        { "-s1", FALSE, etSTR, {&value1}, "Description" },
+        { "-s2", FALSE, etSTR, {&value2}, "Description" },
+        { "-s3", FALSE, etSTR, {&value3}, "Description" }
+    };
+    const char *const cmdline[] = {
+        "test", "-s1", "", "-s2", "test"
+    };
+    parse(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    EXPECT_STREQ("", value1);
+    EXPECT_STREQ("test", value2);
+    EXPECT_STREQ("default", value3);
+}
+
+TEST_F(ParseCommonArgsTest, ParsesBooleanArgs)
+{
+    gmx_bool          value1 = TRUE, value2 = FALSE, value3 = TRUE;
+    t_pargs           pa[]   = {
+        { "-b1", FALSE, etBOOL, {&value1}, "Description" },
+        { "-b2", FALSE, etBOOL, {&value2}, "Description" },
+        { "-b3", FALSE, etBOOL, {&value3}, "Description" }
+    };
+    const char *const cmdline[] = {
+        "test", "-nob1", "-b2"
+    };
+    parse(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    EXPECT_FALSE(value1);
+    EXPECT_TRUE(value2);
+    EXPECT_TRUE(value3);
+}
+
+TEST_F(ParseCommonArgsTest, ParsesVectorArgs)
+{
+    rvec              value1 = {0, 0, 0}, value2 = {0, 0, 0}, value3 = {1, 2, 3};
+    t_pargs           pa[]   = {
+        { "-v1", FALSE, etRVEC, {&value1}, "Description" },
+        { "-v2", FALSE, etRVEC, {&value2}, "Description" },
+        { "-v3", FALSE, etRVEC, {&value3}, "Description" }
+    };
+    const char *const cmdline[] = {
+        "test", "-v1", "2", "1", "3", "-v2", "1"
+    };
+    parse(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    EXPECT_EQ(2.0, value1[XX]);
+    EXPECT_EQ(1.0, value1[YY]);
+    EXPECT_EQ(3.0, value1[ZZ]);
+    EXPECT_EQ(1.0, value2[XX]);
+    EXPECT_EQ(1.0, value2[YY]);
+    EXPECT_EQ(1.0, value2[ZZ]);
+    EXPECT_EQ(1.0, value3[XX]);
+    EXPECT_EQ(2.0, value3[YY]);
+    EXPECT_EQ(3.0, value3[ZZ]);
+}
+
+TEST_F(ParseCommonArgsTest, ParsesTimeArgs)
+{
+    real              value1 = 1.0, value2 = 2.0, value3 = 2.5;
+    t_pargs           pa[]   = {
+        { "-t1", FALSE, etTIME, {&value1}, "Description" },
+        { "-t2", FALSE, etTIME, {&value2}, "Description" },
+        { "-t3", FALSE, etTIME, {&value3}, "Description" }
+    };
+    const char *const cmdline[] = {
+        "test", "-t1", "2", "-t2", "-.5"
+    };
+    parse(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    EXPECT_EQ( 2.0, value1);
+    EXPECT_EQ(-0.5, value2);
+    EXPECT_EQ( 2.5, value3);
+}
+
+TEST_F(ParseCommonArgsTest, ParsesTimeArgsWithTimeUnit)
+{
+    real              value1 = 1.0, value2 = 2.0, value3 = 2.5;
+    t_pargs           pa[]   = {
+        { "-t1", FALSE, etTIME, {&value1}, "Description" },
+        { "-t2", FALSE, etTIME, {&value2}, "Description" },
+        { "-t3", FALSE, etTIME, {&value3}, "Description" }
+    };
+    const char *const cmdline[] = {
+        "test", "-t1", "2", "-t2", "-.5", "-tu", "ns"
+    };
+    parse(cmdline, PCA_TIME_UNIT, gmx::EmptyArrayRef(), pa);
+    EXPECT_EQ( 2000.0, value1);
+    EXPECT_EQ(-500.0, value2);
+    EXPECT_EQ( 2.5, value3);
+}
+
+TEST_F(ParseCommonArgsTest, ParsesEnumArgs)
+{
+    const char       *value1[] = {NULL, "none", "on", "off", NULL };
+    const char       *value2[] = {NULL, "def", "value", "value_other", NULL };
+    const char       *value3[] = {NULL, "default", "value", NULL };
+    t_pargs           pa[]     = {
+        { "-s1", FALSE, etENUM, {value1}, "Description" },
+        { "-s2", FALSE, etENUM, {value2}, "Description" },
+        { "-s3", FALSE, etENUM, {value3}, "Description" }
+    };
+    const char *const cmdline[] = {
+        "test", "-s1", "off", "-s2", "val"
+    };
+    parse(cmdline, 0, gmx::EmptyArrayRef(), pa);
+    EXPECT_STREQ("off", value1[0]);
+    EXPECT_STREQ("value", value2[0]);
+    EXPECT_STREQ("default", value3[0]);
+    EXPECT_EQ(value1[nenum(value1)], value1[0]);
+    EXPECT_EQ(value2[nenum(value2)], value2[0]);
+    EXPECT_EQ(value3[nenum(value3)], value3[0]);
+}
+
+/********************************************************************
+ * Tests for file name options
+ */
+
+TEST_F(ParseCommonArgsTest, ParsesFileArgs)
+{
+    t_filenm          fnm[] = {
+        { efXVG, "-o1", "out1", ffOPTWR },
+        { efXVG, "-o2", "out2", ffOPTWR },
+        { efXVG, "-om", "outm", ffWRMULT },
+        { efXVG, "-om2", "outm2", ffWRMULT },
+    };
+    const char *const cmdline[] = {
+        "test", "-o1", "-o2", "test", "-om", "test1", "test2.xvg", "-om2"
+    };
+    parse(cmdline, 0, fnm, gmx::EmptyArrayRef());
+    EXPECT_STREQ("out1.xvg", opt2fn_null("-o1", nfile(), fnm));
+    EXPECT_STREQ("test.xvg", opt2fn_null("-o2", nfile(), fnm));
+    char **files;
+    EXPECT_EQ(2, opt2fns(&files, "-om", nfile(), fnm));
+    EXPECT_STREQ("test1.xvg", files[0]);
+    EXPECT_STREQ("test2.xvg", files[1]);
+    EXPECT_STREQ("outm2.xvg", opt2fn("-om2", nfile(), fnm));
+    done_filenms(nfile(), fnm);
+}
+
+TEST_F(ParseCommonArgsTest, ParsesFileArgsWithDefaults)
+{
+    t_filenm          fnm[] = {
+        { efTPS, NULL,  NULL,   ffWRITE },
+        { efTRX, "-f2", NULL,   ffOPTWR },
+        { efTRX, "-f3", "trj3", ffWRITE },
+        { efXVG, "-o",  "out",  ffWRITE },
+        { efXVG, "-om", "outm", ffWRMULT },
+    };
+    const char *const cmdline[] = {
+        "test"
+    };
+    parse(cmdline, 0, fnm, gmx::EmptyArrayRef());
+    EXPECT_STREQ("topol.tpr", ftp2fn(efTPS, nfile(), fnm));
+    EXPECT_STREQ("traj.xtc", opt2fn("-f2", nfile(), fnm));
+    EXPECT_EQ(NULL, opt2fn_null("-f2", nfile(), fnm));
+    EXPECT_STREQ("trj3.xtc", opt2fn("-f3", nfile(), fnm));
+    EXPECT_STREQ("out.xvg", opt2fn("-o", nfile(), fnm));
+    EXPECT_STREQ("outm.xvg", opt2fn("-om", nfile(), fnm));
+    done_filenms(nfile(), fnm);
+}
+
+TEST_F(ParseCommonArgsTest, ParsesFileArgsWithDefaultFileName)
+{
+    t_filenm          fnm[] = {
+        { efTPS, "-s",  NULL,   ffWRITE },
+        { efTRX, "-f2", NULL,   ffWRITE },
+        { efTRX, "-f3", "trj3", ffWRITE },
+        { efXVG, "-o",  "out",  ffWRITE },
+        { efXVG, "-om", "outm", ffWRMULT },
+    };
+    const char *const cmdline[] = {
+        "test", "-deffnm", "def", "-f2", "other", "-o"
+    };
+    parse(cmdline, PCA_CAN_SET_DEFFNM, fnm, gmx::EmptyArrayRef());
+    EXPECT_STREQ("def.tpr", ftp2fn(efTPS, nfile(), fnm));
+    EXPECT_STREQ("other.xtc", opt2fn("-f2", nfile(), fnm));
+    EXPECT_STREQ("def.xtc", opt2fn("-f3", nfile(), fnm));
+    EXPECT_STREQ("def.xvg", opt2fn("-o", nfile(), fnm));
+    EXPECT_STREQ("def.xvg", opt2fn("-om", nfile(), fnm));
+    done_filenms(nfile(), fnm);
+}
+
+/********************************************************************
+ * Tests for general behavior
+ */
+
+TEST_F(ParseCommonArgsTest, CanKeepUnknownArgs)
+{
+    int               ivalue = 0;
+    gmx_bool          bvalue = FALSE;
+    t_pargs           pa[]   = {
+        { "-i", FALSE, etINT, {&ivalue}, "Description" },
+        { "-b", FALSE, etBOOL, {&bvalue}, "Description" },
+    };
+    t_filenm          fnm[] = {
+        { efXVG, "-o1", "out1", ffOPTWR },
+        { efXVG, "-o2", "out2", ffOPTWR }
+    };
+    const char *const cmdline[] = {
+        "test", "foo", "-unk", "-o1", "-unk2", "-o2", "test",
+        "-i", "2", "-unk3", "-b", "-unk4"
+    };
+    parse(cmdline, PCA_NOEXIT_ON_ARGS, fnm, pa);
+    EXPECT_EQ(2, ivalue);
+    EXPECT_TRUE(bvalue);
+    EXPECT_STREQ("out1.xvg", opt2fn_null("-o1", nfile(), fnm));
+    EXPECT_STREQ("test.xvg", opt2fn_null("-o2", nfile(), fnm));
+    EXPECT_EQ(6, args_.argc());
+    EXPECT_STREQ(cmdline[0],  args_.arg(0));
+    EXPECT_STREQ(cmdline[1],  args_.arg(1));
+    EXPECT_STREQ(cmdline[2],  args_.arg(2));
+    EXPECT_STREQ(cmdline[4],  args_.arg(3));
+    EXPECT_STREQ(cmdline[9],  args_.arg(4));
+    EXPECT_STREQ(cmdline[11], args_.arg(5));
+    done_filenms(nfile(), fnm);
+}
+
+} // namespace
index e87f31569408093172672575b0a2c0c09af7eba6..c02f120ab0570bb1c86deb033fd4dd3f72cae886 100644 (file)
@@ -34,7 +34,7 @@
  */
 /*! \file
  * \brief
- * Declares gmx::ConstArrayRef.
+ * Declares gmx::ArrayRef and gmx::ConstArrayRef.
  *
  * \author Teemu Murtola <teemu.murtola@gmail.com>
  * \inpublicapi
@@ -56,7 +56,234 @@ namespace gmx
 {
 
 /*! \brief
- * STL container non-mutable interface for a C array (or part of a std::vector).
+ * Tag type to initialize empty array references.
+ *
+ * This type (together with appropriate constructors in ArrayRef and
+ * ConstArrayRef) allows initializing any array reference to an empty value
+ * without explicitly specifying its type.  This is convenient when calling
+ * a function that takes an array reference, where constructing an empty
+ * reference explicitly would otherwise require specifying the full array
+ * reference type, including the template parameter.
+ */
+struct EmptyArrayRef {};
+
+/*! \brief
+ * STL-like container for an interface to a C array (or part of a std::vector).
+ *
+ * \tparam T  Value type of elements.
+ *
+ * This class provides an interface similar to \c std::vector<T>, with the
+ * following main differences:
+ *  - This class does not have its own storage.  Instead, it references an
+ *    existing array of values (either a C-style array or part of an existing
+ *    std::vector<T>).
+ *  - It is only possible to modify the values themselves through ArrayRef;
+ *    it is not possible to add or remove values.
+ *  - Copying objects of this type is cheap, and the copies behave identically
+ *    to the original object: the copy references the same set of values.
+ *
+ * This class is useful for writing wrappers that expose a different view of
+ * the internal data stored as a single vector/array.
+ *
+ * Methods in this class do not throw, except where indicated.
+ *
+ * Note that due to a Doxygen limitation, the constructor that takes a C array
+ * whose size is known at compile time does not appear in the documentation.
+ *
+ * \todo
+ * This class is not complete.  At least, it should be possible to convert an
+ * ArrayRef to a ConstArrayRef.  There are likely also methods missing (not
+ * required for current usage).
+ *
+ * \inpublicapi
+ * \ingroup module_utility
+ */
+template <typename T>
+class ArrayRef
+{
+    public:
+        //! Type of values stored in the container.
+        typedef T         value_type;
+        //! Type for representing size of the container.
+        typedef size_t    size_type;
+        //! Type for representing difference between two container indices.
+        typedef ptrdiff_t difference_type;
+        //! Const reference to a container element.
+        typedef const T  &const_reference;
+        //! Const pointer to a container element.
+        typedef const T  *const_pointer;
+        //! Const iterator type for the container.
+        typedef const T  *const_iterator;
+        //! Reference to a container element.
+        typedef T        &reference;
+        //! Pointer to a container element.
+        typedef T        *pointer;
+        //! Iterator type for the container.
+        typedef T        *iterator;
+        //! Standard reverse iterator.
+        typedef std::reverse_iterator<iterator>       reverse_iterator;
+        //! Standard reverse iterator.
+        typedef std::reverse_iterator<const_iterator> const_reverse_iterator;
+
+        /*! \brief
+         * Constructs an empty reference.
+         */
+        ArrayRef() : begin_(NULL), end_(NULL) {}
+        /*! \brief
+         * Constructs an empty reference.
+         *
+         * This is provided for convenience, such that EmptyArrayRef can be
+         * used to initialize any ArrayRef, without specifying the template
+         * type.  It is not explicit to enable that usage.
+         */
+        ArrayRef(const EmptyArrayRef &) : begin_(NULL), end_(NULL) {}
+        /*! \brief
+         * Constructs a reference to a particular range.
+         *
+         * \param[in] begin  Pointer to the beginning of a range.
+         * \param[in] end    Pointer to the end of a range.
+         *
+         * Passed pointers must remain valid for the lifetime of this object.
+         */
+        ArrayRef(pointer begin, pointer end)
+            : begin_(begin), end_(end)
+        {
+            GMX_ASSERT(end >= begin, "Invalid range");
+        }
+        /*! \brief
+         * Constructs a reference to a particular range in a std::vector.
+         *
+         * \param[in] begin  Iterator to the beginning of a range.
+         * \param[in] end    Iterator to the end of a range.
+         *
+         * The referenced vector must remain valid and not be reallocated for
+         * the lifetime of this object.
+         */
+        ArrayRef(typename std::vector<value_type>::iterator begin,
+                 typename std::vector<value_type>::iterator end)
+            : begin_((begin != end) ? &*begin : NULL),
+              end_(begin_+(end-begin))
+        {
+            GMX_ASSERT(end >= begin, "Invalid range");
+        }
+        /*! \brief
+         * Constructs a reference to an array.
+         *
+         * \param[in] begin  Pointer to the beginning of the array.
+         *      May be NULL if \p size is zero.
+         * \param[in] size   Number of elements in the array.
+         *
+         * Passed pointer must remain valid for the lifetime of this object.
+         */
+        ArrayRef(pointer begin, size_type size)
+            : begin_(begin), end_(begin + size)
+        {
+        }
+        //! \cond
+        // Doxygen 1.8.5 doesn't parse the declaration correctly...
+        /*! \brief
+         * Constructs a reference to a C array.
+         *
+         * \param[in] array  C array to reference.
+         * \tparam    count  Deduced number of elements in \p array.
+         *
+         * This constructor can only be used with a real array (not with a
+         * pointer).  It constructs a reference to the whole array, without
+         * a need to pass the number of elements explicitly.  The compiler
+         * must be able to deduce the array size.
+         *
+         * Passed array must remain valid for the lifetime of this object.
+         *
+         * This constructor is not explicit to allow directly passing
+         * a C array to a function that takes an ArrayRef parameter.
+         */
+        template <size_t count>
+        ArrayRef(value_type (&array)[count])
+            : begin_(array), end_(array + count)
+        {
+        }
+        //! \endcond
+
+        //! Returns an iterator to the beginning of the container.
+        iterator begin() { return begin_; }
+        //! Returns an iterator to the beginning of the container.
+        const_iterator begin() const { return begin_; }
+        //! Returns an iterator to the end of the container.
+        iterator end() { return end_; }
+        //! Returns an iterator to the end of the container.
+        const_iterator end() const { return end_; }
+        //! Returns an iterator to the reverse beginning of the container.
+        iterator rbegin() { return reverse_iterator(end()); }
+        //! Returns an iterator to the reverse beginning of the container.
+        const_iterator rbegin() const { return reverse_iterator(end()); }
+        //! Returns an iterator to the reverse end of the container.
+        iterator rend() { return reverse_iterator(begin()); }
+        //! Returns an iterator to the reverse end of the container.
+        const_iterator rend() const { return reverse_iterator(begin()); }
+
+        //! Returns the size of the container.
+        size_type size() const { return end_ - begin_; }
+        //! Identical to size().
+        size_type capacity() const { return end_ - begin_; }
+        //! Whether the container is empty.
+        bool empty() const { return begin_ == end_; }
+
+        //! Access container element.
+        reference operator[](size_type n) { return begin_[n]; }
+        //! Access container element.
+        const_reference operator[](size_type n) const { return begin_[n]; }
+        //! Access container element (throws on out-of-range error).
+        reference at(size_type n)
+        {
+            if (n >= size())
+            {
+                throw std::out_of_range("Vector index out of range");
+            }
+            return begin_[n];
+        }
+        //! Access container element (throws on out-of-range error).
+        const_reference at(size_type n) const
+        {
+            if (n >= size())
+            {
+                throw std::out_of_range("Vector index out of range");
+            }
+            return begin_[n];
+        }
+        //! Returns the first element in the container.
+        reference front() { return *begin_; }
+        //! Returns the first element in the container.
+        const_reference front() const { return *begin_; }
+        //! Returns the last element in the container.
+        reference back() { return *(end_ - 1); }
+        //! Returns the last element in the container.
+        const_reference back() const { return *(end_ - 1); }
+
+        //! Returns a raw pointer to the contents of the array.
+        pointer data() { return begin_; }
+        //! Returns a raw pointer to the contents of the array.
+        const_pointer data() const { return begin_; }
+
+        /*! \brief
+         * Swaps referenced memory with the other object.
+         *
+         * The actual memory areas are not modified, only the references are
+         * swapped.
+         */
+        void swap(ArrayRef<T> &other)
+        {
+            std::swap(begin_, other.begin_);
+            std::swap(end_, other.end_);
+        }
+
+    private:
+        pointer           begin_;
+        pointer           end_;
+};
+
+/*! \brief
+ * STL-like container for non-mutable interface to a C array (or part of a
+ * std::vector).
  *
  * \tparam T  Value type of elements.
  *
@@ -112,6 +339,14 @@ class ConstArrayRef
          * Constructs an empty reference.
          */
         ConstArrayRef() : begin_(NULL), end_(NULL) {}
+        /*! \brief
+         * Constructs an empty reference.
+         *
+         * This is provided for convenience, such that EmptyArrayRef can be
+         * used to initialize any Const ArrayRef, without specifying the
+         * template type.  It is not explicit to enable that usage.
+         */
+        ConstArrayRef(const EmptyArrayRef &) : begin_(NULL), end_(NULL) {}
         /*! \brief
          * Constructs a reference to a particular range.
          *
@@ -178,13 +413,13 @@ class ConstArrayRef
         }
         //! \endcond
 
-        //! Returns an interator to the beginning of the container.
+        //! Returns an iterator to the beginning of the container.
         const_iterator begin() const { return begin_; }
-        //! Returns an interator to the end of the container.
+        //! Returns an iterator to the end of the container.
         const_iterator end() const { return end_; }
-        //! Returns an interator to the reverse beginning of the container.
+        //! Returns an iterator to the reverse beginning of the container.
         const_iterator rbegin() const { return reverse_iterator(end()); }
-        //! Returns an interator to the reverse end of the container.
+        //! Returns an iterator to the reverse end of the container.
         const_iterator rend() const { return reverse_iterator(begin()); }
 
         //! Returns the size of the container.
@@ -230,6 +465,19 @@ class ConstArrayRef
         const_pointer           end_;
 };
 
+/*! \brief
+ * Simple swap method for ArrayRef objects.
+ *
+ * \see ArrayRef::swap()
+ *
+ * \ingroup module_utility
+ */
+template <typename T>
+void swap(ArrayRef<T> &a, ArrayRef<T> &b)
+{
+    a.swap(b);
+}
+
 /*! \brief
  * Simple swap method for ConstArrayRef objects.
  *