#include <base/subscriptor.h>
+template <typename T> class vector2d;
+
+
/**
* Two-dimensional array of arbitrary data.
*
* This is an implementation of a two-dimensional array with access by
- * pairs of coordinates.
+ * pairs of coordinates. Access is either by @p{x(i,j)} to keep with
+ * matrix notation, or @p{x[i][j]} in accordance to C-style arrays.
*
* The name is chosen to be @p{vector2d} to be conformant with the
* standard C++ classes. This, although array would have been
class vector2d : public Subscriptor
{
public:
+ /**
+ * Object representing one row of
+ * a @ref{vector2d} object. It
+ * allows to access the elements
+ * through the
+ * @p{operator[]}. Since objects
+ * of this type are also
+ * generated through the
+ * @p{vector2d::operator[]}, this
+ * allows accessing the elements
+ * of @p{vector2d} objects just
+ * as those of two-dimensional
+ * C-style arrays.
+ *
+ * This class is used for
+ * constant @p{vector2d}
+ * objects. It only allows
+ * read-only access to the
+ * elements.
+ */
+ class ConstRowAccessor
+ {
+ public:
+ /**
+ * Constructor.
+ */
+ ConstRowAccessor (const vector2d<T> &parent,
+ const unsigned int row);
+
+ /**
+ * Access operator.
+ */
+ T operator [] (const unsigned int column) const;
+
+ protected:
+ /**
+ * Pointer to the parent
+ * object. Used only to check
+ * access indices for
+ * validity.
+ */
+ const vector2d<T> &parent;
+
+ /**
+ * Pointer to the start of
+ * the row pointed to by this
+ * object.
+ */
+ const T * const row_start;
+ };
+
+
+ /**
+ * Object representing one row of
+ * a @ref{vector2d} object. It
+ * allows to access the elements
+ * through the
+ * @p{operator[]}. Since objects
+ * of this type are also
+ * generated through the
+ * @p{vector2d::operator[]}, this
+ * allows accessing the elements
+ * of @p{vector2d} objects just
+ * as those of two-dimensional
+ * C-style arrays.
+ *
+ * This class is used for
+ * non-constant @p{vector2d}
+ * objects. It allows read and
+ * write access to the elements.
+ */
+ class NonConstRowAccessor
+ {
+ public:
+ /**
+ * Constructor.
+ */
+ NonConstRowAccessor (vector2d<T> &parent,
+ const unsigned int row);
+
+
+ /**
+ * Access operator.
+ */
+ T & operator [] (const unsigned int column) const;
+
+ private:
+ /**
+ * Pointer to the parent
+ * object. Used only to check
+ * access indices for
+ * validity.
+ */
+ const vector2d<T> &parent;
+
+ /**
+ * Pointer to the start of
+ * the row pointed to by this
+ * object.
+ */
+ T * const row_start;
+ };
/**
* Constructor for a quadratic
* @p{rows} by @p{rows} array. The
T & operator() (const unsigned int i,
const unsigned int j);
+ /**
+ * Return an object representing
+ * a certain row of this
+ * array. This object in turn has
+ * an @p{operator[]}, so that the
+ * elements of this array can be
+ * accessed through @p{x[i][j]}
+ * as well as through @p{x(i,j)}.
+ *
+ * This function is for constant
+ * objects. The returned row
+ * representing object only
+ * allows read access to the
+ * elements of the row pointed
+ * to.
+ */
+ ConstRowAccessor operator [] (const unsigned int row) const;
+
+ /**
+ * Return an object representing
+ * a certain row of this
+ * array. This object in turn has
+ * an @p{operator[]}, so that the
+ * elements of this array can be
+ * accessed through @p{x[i][j]}
+ * as well as through @p{x(i,j)}.
+ *
+ * This function is for
+ * non-constant objects. The
+ * returned row representing
+ * object allows read and write
+ * access to the elements of the
+ * row pointed to.
+ */
+ NonConstRowAccessor operator [] (const unsigned int row);
+
/**
* Set all entries to their
/**
* This is unfortunately needed.
*/
- template <typename T2> friend class FullMatrix;
+ template <typename T2> friend class FullMatrix;
+
+ friend class ConstRowAccessor;
+ friend class NonConstRowAccessor;
};
/* --------------------- Template and inline functions ---------------- */
+template <typename T>
+inline
+vector2d<T>::ConstRowAccessor::
+ConstRowAccessor (const vector2d<T> &parent,
+ const unsigned int row)
+ :
+ parent (parent),
+ row_start(&parent.val[row*parent.n_cols()])
+{
+ Assert (row < parent.n_rows(), ExcInternalError());
+};
+
+
+template <typename T>
+inline
+T
+vector2d<T>::ConstRowAccessor::
+operator [] (const unsigned int column) const
+{
+ Assert (column < parent.n_cols(), ExcInternalError());
+ return *(row_start+column);
+};
+
+
+
+template <typename T>
+inline
+vector2d<T>::NonConstRowAccessor::
+NonConstRowAccessor (vector2d<T> &parent,
+ const unsigned int row)
+ :
+ parent (parent),
+ row_start(&parent.val[row*parent.n_cols()])
+{
+ Assert (row < parent.n_rows(), ExcInternalError());
+};
+
+
+template <typename T>
+inline
+T &
+vector2d<T>::NonConstRowAccessor::
+operator [] (const unsigned int column) const
+{
+ Assert (column < parent.n_cols(), ExcInternalError());
+ return *(row_start+column);
+};
+
+
+
template <typename T>
inline
vector2d<T>::~vector2d ()
template <typename T>
-inline T &
+inline
+T &
vector2d<T>::operator() (const unsigned int i, const unsigned int j)
{
Assert (i<num_rows, ExcIndexRange (i, 0, num_rows));
+template <typename T>
+inline
+typename vector2d<T>::ConstRowAccessor
+vector2d<T>::operator [] (const unsigned int row) const
+{
+ // Note: check for validity of row
+ // number is done in the
+ // constructor of the created
+ // object
+ return ConstRowAccessor(*this, row);
+};
+
+
+
+template <typename T>
+inline
+typename vector2d<T>::NonConstRowAccessor
+vector2d<T>::operator [] (const unsigned int row)
+{
+ // Note: check for validity of row
+ // number is done in the
+ // constructor of the created
+ // object
+ return NonConstRowAccessor(*this, row);
+};
+
+
+
+
+
template <typename T>
inline
const T *
--- /dev/null
+//---------------------------- $RCSfile$ ---------------------------
+// $Id$
+// Version: $Name$
+//
+// Copyright (C) 1998, 1999, 2000, 2001, 2002 by the deal.II authors
+//
+// This file is subject to QPL and may not be distributed
+// without copyright and license information. Please refer
+// to the file deal.II/doc/license.html for the text and
+// further information on this license.
+//
+//---------------------------- $RCSfile$ ---------------------------
+
+
+#include <cmath>
+#include <fstream>
+#include <iostream>
+#include <iomanip>
+#include <cstdlib>
+
+#include <base/logstream.h>
+#include <base/vector2d.h>
+
+const int entries[9] = { 11,12,13,21,22,23,31,32,33 };
+
+int
+main ()
+{
+ std::ofstream logfile("vector2d.output");
+ logfile.setf(std::ios::fixed);
+ logfile.precision(3);
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+
+ vector2d<double> Td(3,3);
+ vector2d<int> Ti(3,3);
+
+ for (unsigned int i=0; i<9; ++i)
+ {
+ Td[i/3][i%3] = entries[i];
+ Ti[i/3][i%3] = entries[i];
+ };
+
+ for (unsigned int i=0; i<3; ++i)
+ for (unsigned int j=0; j<3; ++j)
+ {
+ Assert (Td[i][j] == Td(i,j), ExcInternalError());
+ Assert (Ti[i][j] == Ti(i,j), ExcInternalError());
+ Assert (Ti[i][j] == Td(i,j), ExcInternalError());
+
+ logfile << i << " " << j << " " << Td[i][j] << " ok" << std::endl;
+ };
+};
+
+
+