* Pointer to the @p ith vertex bounding this object. Throw an exception if
* <code>dim=1</code>.
*/
- TriaAccessor<0,dim,spacedim> vertex_accessor (const unsigned int i) const;
+ typename dealii::internal::Triangulation::Iterators<dim,spacedim>::vertex_iterator
+ vertex_it (const unsigned int i) const;
/**
* Return the global index of i-th vertex of the current object. The
* <code>dim</code> (i.e. 1 for a triangulation of lines, 2 for a triangulation
* of quads, and 3 for a triangulation of hexes) that is embedded in a space of
* dimensionality <code>spacedim</code> (for <code>spacedim==dim</code> the
- * triangulation represents a domain in $R^{dim}$, for
+ * triangulation represents a domain in ${\mathbb R}^\text{dim}$, for
* <code>spacedim@>dim</code> the triangulation is of a manifold embedded in
* a higher dimensional space).
*
/**
* Constructor. This constructor exists in order to maintain interface
- * compatibility with the other accessor classes. However, it doesn't do
- * anything useful here and so may not actually be called.
+ * compatibility with the other accessor classes. @p index can be used to set
+ * the global index of the vertex we point to.
*/
- TriaAccessor (const Triangulation<dim,spacedim> *tria = NULL,
- const int = 0,
- const int = 0,
- const AccessorData * = 0);
+ TriaAccessor (const Triangulation<dim,spacedim> *tria = NULL,
+ const int level = 0,
+ const int index = 0,
+ const AccessorData * = 0);
/**
* Constructor. Should never be called and thus produces an error.
+template <int structdim, int dim, int spacedim>
+inline
+typename dealii::internal::Triangulation::Iterators<dim,spacedim>::vertex_iterator
+TriaAccessor<structdim,dim,spacedim>::vertex_it (const unsigned int i) const
+{
+ return typename dealii::internal::Triangulation::Iterators<dim,spacedim>::vertex_iterator
+ (this->tria, 0, vertex_index (i));
+}
+
+
+
template <int structdim, int dim, int spacedim>
unsigned int
TriaAccessor<structdim, dim, spacedim>::
const int index,
const AccessorData *)
:
- tria (NULL),
- global_vertex_index (numbers::invalid_unsigned_int)
-{
- Assert(false, ExcImpossibleInDim(0));
-}
+ tria (tria),
+ global_vertex_index (index)
+{}
* @ref Iterators
* module for more information.
*
+ * A @p vertex_iterator is typedef'd to an iterator operating on the @p
+ * vertices member variable of a <tt>Triangulation<1></tt> object.
+ *
* A @p line_iterator is typedef'd to an iterator operating on the @p
* lines member variable of a <tt>Triangulation<1></tt> object. An @p
* active_line_iterator only operates on the active lines. @p
template <int spacedim>
struct Iterators<1,spacedim>
{
+ typedef TriaIterator <dealii::TriaAccessor<0, 1, spacedim> > vertex_iterator;
+
typedef TriaRawIterator <dealii::CellAccessor<1,spacedim> > raw_line_iterator;
typedef TriaIterator <dealii::CellAccessor<1,spacedim> > line_iterator;
typedef TriaActiveIterator<dealii::CellAccessor<1,spacedim> > active_line_iterator;
* @ref Iterators
* module for more information.
*
+ * A @p vertex_iterator is typedef'd to an iterator operating on the @p
+ * vertices member variable of a <tt>Triangulation<2></tt> object.
+ *
* A @p line_iterator is typedef'd to an iterator operating on the @p
* lines member variable of a <tt>Triangulation<2></tt> object. An @p
* active_line_iterator only operates on the active lines. @p
template <int spacedim>
struct Iterators<2,spacedim>
{
+ typedef TriaIterator <dealii::TriaAccessor<0, 2, spacedim> > vertex_iterator;
+
typedef TriaRawIterator <dealii::TriaAccessor<1, 2, spacedim> > raw_line_iterator;
typedef TriaIterator <dealii::TriaAccessor<1, 2, spacedim> > line_iterator;
typedef TriaActiveIterator<dealii::TriaAccessor<1, 2, spacedim> > active_line_iterator;
template <int spacedim>
struct Iterators<3,spacedim>
{
+ typedef TriaIterator <dealii::TriaAccessor<0, 3, spacedim> > vertex_iterator;
+
typedef TriaRawIterator <dealii::TriaAccessor<1, 3, spacedim> > raw_line_iterator;
typedef TriaIterator <dealii::TriaAccessor<1, 3, spacedim> > line_iterator;
typedef TriaActiveIterator<dealii::TriaAccessor<1, 3, spacedim> > active_line_iterator;
-template <>
-TriaAccessor<0,1,1>
-TriaAccessor<1,1,1>::vertex_accessor (const unsigned int i) const
-{
- // Silence a warning
- (void) i;
- return TriaAccessor<0,1,1>();
-}
-
-
-
-template <>
-TriaAccessor<0,1,2>
-TriaAccessor<1,1,2>::vertex_accessor (const unsigned int i) const
-{
- // Silence a warning
- (void) i;
- return TriaAccessor<0,1,2>();
-}
-
-
-
-template <>
-TriaAccessor<0,1,3>
-TriaAccessor<1,1,3>::vertex_accessor (const unsigned int i) const
-{
- // Silence a warning
- (void) i;
- return TriaAccessor<0,1,3>();
-}
-
-
-
-template <>
-TriaAccessor<0,2,2>
-TriaAccessor<1,2,2>::vertex_accessor (const unsigned int i) const
-{
- return TriaAccessor<0,2,2>(this->tria, vertex_index (i));
-}
-
-
-
-template <>
-TriaAccessor<0,2,3>
-TriaAccessor<1,2,3>::vertex_accessor (const unsigned int i) const
-{
- return TriaAccessor<0,2,3>(this->tria, vertex_index (i));
-}
-
-
-
-template <>
-TriaAccessor<0,3,3>
-TriaAccessor<1,3,3>::vertex_accessor (const unsigned int i) const
-{
- return TriaAccessor<0,3,3>(this->tria, vertex_index (i));
-}
-
-
-
-template <>
-TriaAccessor<0,2,2>
-TriaAccessor<2,2,2>::vertex_accessor (const unsigned int i) const
-{
- return TriaAccessor<0,2,2>(this->tria, vertex_index (i));
-}
-
-
-
-template <>
-TriaAccessor<0,2,3>
-TriaAccessor<2,2,3>::vertex_accessor (const unsigned int i) const
-{
- return TriaAccessor<0,2,3>(this->tria, vertex_index (i));
-}
-
-
-
-template <>
-TriaAccessor<0,3,3>
-TriaAccessor<2,3,3>::vertex_accessor (const unsigned int i) const
-{
- return TriaAccessor<0,3,3>(this->tria, vertex_index (i));
-}
-
-
-
-template <>
-TriaAccessor<0,3,3>
-TriaAccessor<3,3,3>::vertex_accessor (const unsigned int i) const
-{
- return TriaAccessor<0,3,3>(this->tria, vertex_index (i));
-}
-
-
-
template <int structdim, int dim, int spacedim>
Point<spacedim>
TriaAccessor<structdim, dim, spacedim>::
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2015 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, 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.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+
+// test TriaAccessor<0,dim,spacedim>
+
+
+#include "../tests.h"
+#include <deal.II/grid/tria.h>
+#include <deal.II/grid/tria_iterator.h>
+#include <deal.II/grid/tria_accessor.h>
+#include <deal.II/grid/grid_generator.h>
+
+#include <fstream>
+
+
+void test ()
+{
+ Triangulation<2> tria;
+ GridGenerator::hyper_cube (tria);
+ tria.refine_global (2);
+
+ for (typename Triangulation<2>::active_cell_iterator
+ cell = tria.begin_active(); cell != tria.end(); ++cell)
+ {
+ for (unsigned int i = 0; i<4; ++i)
+ deallog << cell->vertex_it(i)->center() <<std::endl;
+ }
+}
+
+
+
+int main ()
+{
+ std::ofstream logfile("output");
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+
+ test ();
+
+ return 0;
+}
--- /dev/null
+
+DEAL::0.00000 0.00000
+DEAL::0.250000 0.00000
+DEAL::0.00000 0.250000
+DEAL::0.250000 0.250000
+DEAL::0.250000 0.00000
+DEAL::0.500000 0.00000
+DEAL::0.250000 0.250000
+DEAL::0.500000 0.250000
+DEAL::0.00000 0.250000
+DEAL::0.250000 0.250000
+DEAL::0.00000 0.500000
+DEAL::0.250000 0.500000
+DEAL::0.250000 0.250000
+DEAL::0.500000 0.250000
+DEAL::0.250000 0.500000
+DEAL::0.500000 0.500000
+DEAL::0.500000 0.00000
+DEAL::0.750000 0.00000
+DEAL::0.500000 0.250000
+DEAL::0.750000 0.250000
+DEAL::0.750000 0.00000
+DEAL::1.00000 0.00000
+DEAL::0.750000 0.250000
+DEAL::1.00000 0.250000
+DEAL::0.500000 0.250000
+DEAL::0.750000 0.250000
+DEAL::0.500000 0.500000
+DEAL::0.750000 0.500000
+DEAL::0.750000 0.250000
+DEAL::1.00000 0.250000
+DEAL::0.750000 0.500000
+DEAL::1.00000 0.500000
+DEAL::0.00000 0.500000
+DEAL::0.250000 0.500000
+DEAL::0.00000 0.750000
+DEAL::0.250000 0.750000
+DEAL::0.250000 0.500000
+DEAL::0.500000 0.500000
+DEAL::0.250000 0.750000
+DEAL::0.500000 0.750000
+DEAL::0.00000 0.750000
+DEAL::0.250000 0.750000
+DEAL::0.00000 1.00000
+DEAL::0.250000 1.00000
+DEAL::0.250000 0.750000
+DEAL::0.500000 0.750000
+DEAL::0.250000 1.00000
+DEAL::0.500000 1.00000
+DEAL::0.500000 0.500000
+DEAL::0.750000 0.500000
+DEAL::0.500000 0.750000
+DEAL::0.750000 0.750000
+DEAL::0.750000 0.500000
+DEAL::1.00000 0.500000
+DEAL::0.750000 0.750000
+DEAL::1.00000 0.750000
+DEAL::0.500000 0.750000
+DEAL::0.750000 0.750000
+DEAL::0.500000 1.00000
+DEAL::0.750000 1.00000
+DEAL::0.750000 0.750000
+DEAL::1.00000 0.750000
+DEAL::0.750000 1.00000
+DEAL::1.00000 1.00000