From 8a88c590fbe1ea1a049ae016c597999ff549ab36 Mon Sep 17 00:00:00 2001 From: Wolfgang Bangerth Date: Thu, 7 Aug 2014 10:56:53 -0500 Subject: [PATCH] Provide functions for C++11 range-based for loops. Document C++11 interaction. This commit consists of the following individual commits: Specifically, add Triangulation::cells(), Triangulation::active_cells() and similar for the ones on a single level. Fix up code to make it compile. Fix one compile error. Also provide similar functions for the two DoFHandler classes. Augment documentation. Update the specification of how we implement C++11 range-based for loops. Link to the new documentation module (still coming up). Document C++11 features and interaction. --- doc/doxygen/headers/c++11.h | 103 ++++++ doc/doxygen/headers/std_cxx1x.h | 11 +- doc/news/changes.h | 21 +- include/deal.II/base/iterator_range.h | 334 ++++++++++++++++++ include/deal.II/dofs/dof_handler.h | 123 ++++++- include/deal.II/grid/tria.h | 89 +++++ include/deal.II/hp/dof_handler.h | 93 +++++ source/dofs/dof_handler.cc | 68 +++- source/grid/tria.cc | 43 +++ source/hp/dof_handler.cc | 43 +++ tests/deal.II/range_based_for_tria.cc | 86 +++++ ...ange_based_for_tria.with_cxx11=true.output | 3 + 12 files changed, 1010 insertions(+), 7 deletions(-) create mode 100644 doc/doxygen/headers/c++11.h create mode 100644 include/deal.II/base/iterator_range.h create mode 100644 tests/deal.II/range_based_for_tria.cc create mode 100644 tests/deal.II/range_based_for_tria.with_cxx11=true.output diff --git a/doc/doxygen/headers/c++11.h b/doc/doxygen/headers/c++11.h new file mode 100644 index 0000000000..c543616d6f --- /dev/null +++ b/doc/doxygen/headers/c++11.h @@ -0,0 +1,103 @@ +// --------------------------------------------------------------------- +// $Id$ +// +// Copyright (C) 2014 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. +// +// --------------------------------------------------------------------- + + +/** + * @defgroup CPP11 deal.II and the C++11 standard + * + * At present, deal.II only requires a compiler that conforms to the + * C++98 + * standard and does not rely on compilers to either + * provide the features introduced in + * C++03 or + * C++11 + * + * That said, deal.II interfaces with C++11 in several ways as + * outlined below. + * + * + *

Use of C++11 classes and substitution by BOOST

+ * + * deal.II makes use of many of the classes that were only + * added as part of C++11. This includes std::shared_ptr, + * std::function, std::bind, std::tuple and a number of others. + * Because we do not assume that the compiler actually supports + * C++11, there needs to be a way to ensure that these classes + * are available also for pre-C++11 compilers. This is done using + * the following approach: + * + * - We create a namespace std_cxx1x. + * - If the compiler supports C++11, we import the relevant classes + * and functions into this namespace using statements such as + * @code + * namespace std_cxx1x { using std::shared_ptr; } + * @endcode + * - If the compiler does not support C++11, if its support for + * C++11 is incomplete, or if it is buggy, then we use as a fallback + * the corresponding classes and functions provided by the + * BOOST library through + * statements such as + * @code + * namespace std_cxx1x { using boost::shared_ptr; } + * @endcode + * + * Consequently, namespace std_cxx1x contains all of the symbols + * we require. The classes that can be used this way are obviously + * a subset of the intersection between C++11 and what BOOST provides. + * + * + *

Support for C++11 range-based for loops

+ * + * C++11 provides many new core language features, such as + * rvalue references and move semantics, initialized lists, tuples, + * variadic templates and + * others. For a complete list, see http://en.wikipedia.org/wiki/C++11 . + * We can not use most of these in deal.II itself because we cannot rely + * on compilers supporting them. + * + * However, this does not preclude users from using such features in their + * own applications if they can be reasonably sure that the compilers on + * all of the systems they will work on do support C++11. An example are + * automatically + * typed variables. + * + * deal.II does provide some features that make programming simpler when using + * C++11. This is true, in particular, for + * range-based + * for loops. deal.II-based codes often have many loops of the kind + * @code + * Triangulation triangulation; + * ... + * typename Triangulation::active_cell_iterator + * cell = triangulation.begin_active(), + * endc = triangulation.end(); + * for (; cell!=endc; ++cell) + * cell->set_refine_flag(); + * @endcode + * Using C++11's range-based for loops, you can now write this as follows: + * @code + * Triangulation triangulation; + * ... + * for (auto cell : triangulation.active_cell_iterators()) + * cell->set_refine_flag(); + * @endcode + * This relies on functions such as Triangulation::active_cell_iterators(), + * and equivalents in the DoF handler classes, + * DoFHandler::active_cell_iterators(), hp::DoFHandler::active_cell_iterators(). + * There are variants of these functions that provide iterator ranges + * for all cells (not just the active ones) and for cells on individual + * levels. + */ diff --git a/doc/doxygen/headers/std_cxx1x.h b/doc/doxygen/headers/std_cxx1x.h index 5cc44b5fc5..02f038e2ea 100644 --- a/doc/doxygen/headers/std_cxx1x.h +++ b/doc/doxygen/headers/std_cxx1x.h @@ -16,12 +16,17 @@ /** - * A namespace that contains a selection of classes and functions that will be part of the - * next C++ standard (tentatively called C++1x) and that are currently - * available via the BOOST library. The + * A namespace that contains a selection of classes and functions that are part + * of the C++11 standard and that are also provided by the + * BOOST library. The * elements that are available through the current namespace are either * imported from namespace std (if a compiler's library provides * them) or from namespace boost. + * + * For more information on the topic, + * see also @ref CPP11 "C++11 standard" + * + * @ingroup CPP11 */ namespace std_cxx1x { diff --git a/doc/news/changes.h b/doc/news/changes.h index 469c185d1f..7631e8d60e 100644 --- a/doc/news/changes.h +++ b/doc/news/changes.h @@ -97,7 +97,12 @@ inconvenience this causes.
    - +
  1. New: There is now a documentation module that describes + deal.II's support for and interaction with the + @ref CPP11 "C++11 standard". +
    + (Wolfgang Bangerth, 2014/08/14) +
  2. New: Added FunctionManifold descritpion.
    @@ -110,7 +115,6 @@ inconvenience this causes. (Luca Heltai, 2014/08/07)
  3. -
  4. New: Added CylindricalManifold descritpion.
    This class allows refinement of cylindrical manifolds. It is a good @@ -246,6 +250,19 @@ inconvenience this causes.

    Specific improvements

      +
    1. New: To better support applications that want to use C++11's + range-based + for loops, there are now functions Triangulation::cell_iterators(), + Triangulation::all_cell_iterators() and similarly in classes DoFHandler + and hp::DoFHandler + that return a range object that can then be used in range-based for loops. + The underlying implementation uses the new IteratorRange class. +
      + See the new @ref CPP11 "C++11" page for more information. +
      + (Wolfgang Bangerth, 2014/08/07) +
    2. +
    3. New: TrilinosWrappers::PreconditionAMG can now be initialized from an object of type Epetra_RowMatrix, which allows using it with more arbitrary matrix objects, including matrix-free methods. diff --git a/include/deal.II/base/iterator_range.h b/include/deal.II/base/iterator_range.h new file mode 100644 index 0000000000..c752cc2954 --- /dev/null +++ b/include/deal.II/base/iterator_range.h @@ -0,0 +1,334 @@ +// --------------------------------------------------------------------- +// $Id$ +// +// Copyright (C) 2014 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. +// +// --------------------------------------------------------------------- + +#ifndef __deal2__iterator_range_h +#define __deal2__iterator_range_h + + +#include + +#include + + +DEAL_II_NAMESPACE_OPEN + + +/** + * A class that is used to denote a collection of iterators that can + * be expressed in terms of a range of iterators characterized by a begin + * and an end iterator. As is common in C++, these ranges are specified as + * half open intervals defined by a begin iterator and a one-past-the-end + * iterator. + * + * The purpose of this class is so that classes such as Triangulation and + * DoFHandler can return ranges of cell iterators using an object of the + * current type from functions such as Triangulation::cells() and that such an + * object can then be used in a range-based for loop as supported by C++11, + * see also @ref CPP11 "C++11 standard". + * + * For example, such a loop could look like this if the goal is to set + * the user flag on every active cell: + * @code + * Triangulation triangulation; + * ... + * for (auto cell : triangulation.active_cell_iterators()) + * cell->set_user_flag(); + * @endcode + * In other words, the cell objects are iterators, and the + * range object returned by Triangulation::active_cell_iterators() and + * similar functions are conceptually thought of as collections of + * iterators. + * + * Of course, the class may also be used to denote other iterator + * ranges using different kinds of iterators into other containers. + * + * + *

      Class design: Motivation

      + * + * Informally, the way the C++11 standard describes + * range-based + * for loops works as follows: A range-based for loop of the form + * @code + * Container c; + * for (auto v : c) + * statement; + * @endcode + * where c is a container or collection, is equivalent to the following + * loop: + * @code + * Container c; + * for (auto tmp=c.begin(); tmp!=c.end(); ++tmp) + * { + * auto v = *tmp; + * statement; + * } + * @endcode + * In other words, the compiler introduces a temporary variable that iterates + * over the elements of the container or collection, and the original variable + * v that appeared in the range-based for loop represents the + * dereferenced state of these iterators -- in other words, the + * elements of the collection. + * + * In the context of loops over cells, we typically want to retain the fact that + * the loop variable is an iterator, not a value. This is because in deal.II, + * we never actually use the dereferenced state of a cell iterator: + * conceptually, it would represent a cell, and technically it is implemented + * by classes such as CellAccessor and DoFCellAccessor, but these classes are + * never used explicitly. Consequently, what we would like is that a call + * such as Triangulation::active_cell_iterators() returns an object that + * represents a collection of iterators of the kind + * {begin, begin+1, ..., end-1}. This is conveniently expressed + * as the half open interval [begin,end). The loop variable in the + * range-based for loop would then take on each of these iterators in turn. + * + * + *

      Class design: Implementation

      + * + * To represent the desired semantics as outlined above, this class + * stores a half-open range of iterators [b,e) of + * the given template type. Secondly, the class needs to provide begin() + * and end() functions in such a way that if you dereference the + * result of IteratorRange::begin(), you get the b iterator. + * Furthermore, you must be able to increment the object returned by + * IteratorRange::begin() so that *(++begin()) == b+1. + * In other words, IteratorRange::begin() must return an iterator that + * when dereferenced returns an iterator of the template type + * Iterator: It is an iterator over iterators in the same + * sense as if you had a pointer into an array of pointers. + * + * This is implemented in the form of the IteratorRange::IteratorOverIterators + * class. + * + * @ingroup CPP11 + * @author Wolfgang Bangerth, 2014 + */ +template +class IteratorRange +{ +public: + /** + * A class that implements the semantics of iterators over iterators + * as discussed in the design sections of the IteratorRange class. + */ + class IteratorOverIterators : public std::iterator + { + public: + /** + * Typedef the elements of the collection to give them a name that is + * more distinct. + */ + typedef Iterator BaseIterator; + + /** + * Constructor. Initialize this iterator-over-iterator in + * such a way that it points to the given argument. + * + * @param iterator An iterator to which this object + * is supposed to point. + */ + IteratorOverIterators (const BaseIterator &iterator); + + /** + * Dereferencing operator. + * @return The iterator within the collection currently pointed to. + */ + BaseIterator operator* () const; + + /** + * Dereferencing operator. + * @return The iterator within the collection currently pointed to. + */ + const BaseIterator * operator-> () const; + + /** + * Prefix increment operator. Move the current iterator to the next + * element of the collection and return the new value. + */ + IteratorOverIterators & operator ++ (); + + /** + * Postfix increment operator. Move the current iterator to the next + * element of the collection, but return the previous value of the + * iterator. + */ + IteratorOverIterators operator ++ (int); + + /** + * Comparison operator + * @param i_o_i Another iterator over iterators. + * @return Returns whether the current iterator points to a + * different object than the iterator represented by the + * argument. + */ + bool operator != (const IteratorOverIterators &i_o_i); + + private: + /** + * The object this iterator currently points to. + */ + BaseIterator element_of_iterator_collection; + }; + + + /** + * Typedef for the iterator type represent by this class. + */ + typedef Iterator iterator; + + /** + * Default constructor. Create a range represented by two + * default constructed iterators. This range is likely (depending + * on the type of the iterators) empty. + */ + IteratorRange(); + + /** + * Constructor. Constructs a range given the begin and end iterators. + * + * @param[in] begin An iterator pointing to the first element of the range + * @param[in] end An iterator pointing past the last element represented + * by this range. + */ + IteratorRange (const iterator begin, + const iterator end); + + /** + * Return the iterator pointing to the first element of this range. + */ + IteratorOverIterators begin(); + + /** + * Return the iterator pointing to the element past the last + * element of this range. + */ + IteratorOverIterators end(); + +private: + /** + * Iterators characterizing the begin and end of the range. + */ + const iterator it_begin; + const iterator it_end; +}; + + +// ------------------- template member functions + + +template +inline +IteratorRange::IteratorOverIterators:: +IteratorOverIterators (const BaseIterator &iterator) + : + element_of_iterator_collection (iterator) +{} + + + +template +inline +typename IteratorRange::IteratorOverIterators::BaseIterator +IteratorRange::IteratorOverIterators::operator* () const +{ + return element_of_iterator_collection; +} + + + +template +inline +const typename IteratorRange::IteratorOverIterators::BaseIterator * +IteratorRange::IteratorOverIterators::operator-> () const +{ + return &element_of_iterator_collection; +} + + + +template +inline +typename IteratorRange::IteratorOverIterators & +IteratorRange::IteratorOverIterators::operator ++ () +{ + ++element_of_iterator_collection; + return *this; +} + + + +template +inline +typename IteratorRange::IteratorOverIterators +IteratorRange::IteratorOverIterators::operator ++ (int) +{ + const IteratorOverIterators old_value = *this; + ++element_of_iterator_collection; + return *old_value; +} + + + +template +inline +bool +IteratorRange::IteratorOverIterators::operator != (const IteratorOverIterators &i_o_i) +{ + return element_of_iterator_collection != i_o_i.element_of_iterator_collection; +} + + +template +inline +IteratorRange::IteratorRange () + : + it_begin(), + it_end() +{} + + + +template +inline +IteratorRange::IteratorRange (const iterator b, + const iterator e) + : + it_begin(b), + it_end(e) +{} + + +template +inline +typename IteratorRange::IteratorOverIterators +IteratorRange::begin() +{ + return IteratorOverIterators(it_begin); +} + + +template +inline +typename IteratorRange::IteratorOverIterators +IteratorRange::end() +{ + return IteratorOverIterators(it_end); +} + + +DEAL_II_NAMESPACE_CLOSE + +#endif diff --git a/include/deal.II/dofs/dof_handler.h b/include/deal.II/dofs/dof_handler.h index e5b155ead0..df011c759c 100644 --- a/include/deal.II/dofs/dof_handler.h +++ b/include/deal.II/dofs/dof_handler.h @@ -1,7 +1,7 @@ // --------------------------------------------------------------------- // $Id$ // -// Copyright (C) 1998 - 2013 by the deal.II authors +// Copyright (C) 1998 - 2014 by the deal.II authors // // This file is part of the deal.II library. // @@ -23,6 +23,7 @@ #include #include #include +#include #include #include #include @@ -632,6 +633,126 @@ public: */ level_cell_iterator end_mg () const; + /*@}*/ + + /** + * @name Cell iterator functions returning ranges of iterators + */ + + /** + * Return an iterator range that contains all cells (active or not) + * that make up this DoFHandler. Such a range is useful to + * initialize range-based for loops as supported by C++11. See the + * example in the documentation of active_cell_iterators(). + * + * @return The half open range [this->begin(), this->end()) + * + * @ingroup CPP11 + */ + IteratorRange cell_iterators () const; + + /** + * Return an iterator range that contains all active cells + * that make up this DoFHandler. Such a range is useful to + * initialize range-based for loops as supported by C++11, + * see also @ref CPP11 "C++11 standard". + * + * Range-based for loops are useful in that they require much less + * code than traditional loops (see + * here + * for a discussion of how they work). An example is that without + * range-based for loops, one often writes code such as the following: + * @code + * DoFHandler dof_handler; + * ... + * typename DoFHandler::active_cell_iterator + * cell = dof_handler.begin_active(), + * endc = dof_handler.end(); + * for (; cell!=endc; ++cell) + * { + * fe_values.reinit (cell); + * ...do the local integration on 'cell'...; + * } + * @endcode + * Using C++11's range-based for loops, this is now entirely + * equivalent to the following: + * @code + * DoFHandler dof_handler; + * ... + * for (auto cell : dof_handler.active_cell_iterators()) + * { + * fe_values.reinit (cell); + * ...do the local integration on 'cell'...; + * } + * @endcode + * To use this feature, you need a compiler that supports C++11. + * + * @return The half open range [this->begin_active(), this->end()) + * + * @ingroup CPP11 + */ + IteratorRange active_cell_iterators () const; + + /** + * Return an iterator range that contains all cells (active or not) + * that make up this DoFHandler in their level-cell form. Such a range is useful to + * initialize range-based for loops as supported by C++11. See the + * example in the documentation of active_cell_iterators(). + * + * @return The half open range [this->begin_mg(), this->end_mg()) + * + * @ingroup CPP11 + */ + IteratorRange mg_cell_iterators () const; + + /** + * Return an iterator range that contains all cells (active or not) + * that make up the given level of this DoFHandler. Such a range is useful to + * initialize range-based for loops as supported by C++11. See the + * example in the documentation of active_cell_iterators(). + * + * @param[in] level A given level in the refinement hierarchy of this + * triangulation. + * @return The half open range [this->begin(level), this->end(level)) + * + * @pre level must be less than this->n_levels(). + * + * @ingroup CPP11 + */ + IteratorRange cell_iterators_on_level (const unsigned int level) const; + + /** + * Return an iterator range that contains all active cells + * that make up the given level of this DoFHandler. Such a range is useful to + * initialize range-based for loops as supported by C++11. See the + * example in the documentation of active_cell_iterators(). + * + * @param[in] level A given level in the refinement hierarchy of this + * triangulation. + * @return The half open range [this->begin_active(level), this->end(level)) + * + * @pre level must be less than this->n_levels(). + * + * @ingroup CPP11 + */ + IteratorRange active_cell_iterators_on_level (const unsigned int level) const; + + /** + * Return an iterator range that contains all cells (active or not) + * that make up the given level of this DoFHandler in their level-cell form. + * Such a range is useful to + * initialize range-based for loops as supported by C++11. See the + * example in the documentation of active_cell_iterators(). + * + * @param[in] level A given level in the refinement hierarchy of this + * triangulation. + * @return The half open range [this->begin_mg(level), this->end_mg(level)) + * + * @pre level must be less than this->n_levels(). + * + * @ingroup CPP11 + */ + IteratorRange mg_cell_iterators_on_level (const unsigned int level) const; //@} diff --git a/include/deal.II/grid/tria.h b/include/deal.II/grid/tria.h index ed43e3f220..1eacc53dc0 100644 --- a/include/deal.II/grid/tria.h +++ b/include/deal.II/grid/tria.h @@ -23,6 +23,7 @@ #include #include #include +#include #include #include #include @@ -2544,6 +2545,94 @@ public: active_cell_iterator last_active () const; /*@}*/ + /** + * @name Cell iterator functions returning ranges of iterators + */ + + /** + * Return an iterator range that contains all cells (active or not) + * that make up this triangulation. Such a range is useful to + * initialize range-based for loops as supported by C++11. See the + * example in the documentation of active_cell_iterators(). + * + * @return The half open range [this->begin(), this->end()) + * + * @ingroup CPP11 + */ + IteratorRange cell_iterators () const; + + /** + * Return an iterator range that contains all active cells + * that make up this triangulation. Such a range is useful to + * initialize range-based for loops as supported by C++11, + * see also @ref CPP11 "C++11 standard". + * + * Range-based for loops are useful in that they require much less + * code than traditional loops (see + * here + * for a discussion of how they work). An example is that without + * range-based for loops, one often writes code such as the following + * (assuming for a moment that our goal is setting the user flag + * on every active cell): + * @code + * Triangulation triangulation; + * ... + * typename Triangulation::active_cell_iterator + * cell = triangulation.begin_active(), + * endc = triangulation.end(); + * for (; cell!=endc; ++cell) + * cell->set_user_flag(); + * @endcode + * Using C++11's range-based for loops, this is now entirely + * equivalent to the following: + * @code + * Triangulation triangulation; + * ... + * for (auto cell : triangulation.active_cell_iterators()) + * cell->set_user_flag(); + * @endcode + * To use this feature, you need a compiler that supports C++11. + * + * @return The half open range [this->begin_active(), this->end()) + * + * @ingroup CPP11 + */ + IteratorRange active_cell_iterators () const; + + /** + * Return an iterator range that contains all cells (active or not) + * that make up the given level of this triangulation. Such a range is useful to + * initialize range-based for loops as supported by C++11. See the + * example in the documentation of active_cell_iterators(). + * + * @param[in] level A given level in the refinement hierarchy of this + * triangulation. + * @return The half open range [this->begin(level), this->end(level)) + * + * @pre level must be less than this->n_levels(). + * + * @ingroup CPP11 + */ + IteratorRange cell_iterators_on_level (const unsigned int level) const; + + /** + * Return an iterator range that contains all active cells + * that make up the given level of this triangulation. Such a range is useful to + * initialize range-based for loops as supported by C++11. See the + * example in the documentation of active_cell_iterators(). + * + * @param[in] level A given level in the refinement hierarchy of this + * triangulation. + * @return The half open range [this->begin_active(level), this->end(level)) + * + * @pre level must be less than this->n_levels(). + * + * @ingroup CPP11 + */ + IteratorRange active_cell_iterators_on_level (const unsigned int level) const; + + /*@}*/ + /*---------------------------------------*/ /*---------------------------------------*/ diff --git a/include/deal.II/hp/dof_handler.h b/include/deal.II/hp/dof_handler.h index 10f0522892..7103b9273a 100644 --- a/include/deal.II/hp/dof_handler.h +++ b/include/deal.II/hp/dof_handler.h @@ -23,6 +23,7 @@ #include #include #include +#include #include #include #include @@ -404,6 +405,98 @@ namespace hp /*@}*/ + /** + * @name Cell iterator functions returning ranges of iterators + */ + + /** + * Return an iterator range that contains all cells (active or not) + * that make up this DoFHandler. Such a range is useful to + * initialize range-based for loops as supported by C++11. See the + * example in the documentation of active_cell_iterators(). + * + * @return The half open range [this->begin(), this->end()) + * + * @ingroup CPP11 + */ + IteratorRange cell_iterators () const; + + /** + * Return an iterator range that contains all active cells + * that make up this DoFHandler. Such a range is useful to + * initialize range-based for loops as supported by C++11, + * see also @ref CPP11 "C++11 standard". + * + * Range-based for loops are useful in that they require much less + * code than traditional loops (see + * here + * for a discussion of how they work). An example is that without + * range-based for loops, one often writes code such as the following: + * @code + * DoFHandler dof_handler; + * ... + * typename DoFHandler::active_cell_iterator + * cell = dof_handler.begin_active(), + * endc = dof_handler.end(); + * for (; cell!=endc; ++cell) + * { + * fe_values.reinit (cell); + * ...do the local integration on 'cell'...; + * } + * @endcode + * Using C++11's range-based for loops, this is now entirely + * equivalent to the following: + * @code + * DoFHandler dof_handler; + * ... + * for (auto cell : dof_handler.active_cell_iterators()) + * { + * fe_values.reinit (cell); + * ...do the local integration on 'cell'...; + * } + * @endcode + * To use this feature, you need a compiler that supports C++11. + * + * @return The half open range [this->begin_active(), this->end()) + * + * @ingroup CPP11 + */ + IteratorRange active_cell_iterators () const; + + /** + * Return an iterator range that contains all cells (active or not) + * that make up the given level of this DoFHandler. Such a range is useful to + * initialize range-based for loops as supported by C++11. See the + * example in the documentation of active_cell_iterators(). + * + * @param[in] level A given level in the refinement hierarchy of this + * triangulation. + * @return The half open range [this->begin(level), this->end(level)) + * + * @pre level must be less than this->n_levels(). + * + * @ingroup CPP11 + */ + IteratorRange cell_iterators_on_level (const unsigned int level) const; + + /** + * Return an iterator range that contains all active cells + * that make up the given level of this DoFHandler. Such a range is useful to + * initialize range-based for loops as supported by C++11. See the + * example in the documentation of active_cell_iterators(). + * + * @param[in] level A given level in the refinement hierarchy of this + * triangulation. + * @return The half open range [this->begin_active(level), this->end(level)) + * + * @pre level must be less than this->n_levels(). + * + * @ingroup CPP11 + */ + IteratorRange active_cell_iterators_on_level (const unsigned int level) const; + + /*@}*/ + /*---------------------------------------*/ diff --git a/source/dofs/dof_handler.cc b/source/dofs/dof_handler.cc index 503e2357d5..0756879b4d 100644 --- a/source/dofs/dof_handler.cc +++ b/source/dofs/dof_handler.cc @@ -1,7 +1,7 @@ // --------------------------------------------------------------------- // $Id$ // -// Copyright (C) 1998 - 2013 by the deal.II authors +// Copyright (C) 1998 - 2014 by the deal.II authors // // This file is part of the deal.II library. // @@ -914,6 +914,72 @@ DoFHandler::end_mg () const +template +IteratorRange::cell_iterator> +DoFHandler::cell_iterators () const +{ + return + IteratorRange::cell_iterator> + (begin(), end()); +} + + +template +IteratorRange::active_cell_iterator> +DoFHandler::active_cell_iterators () const +{ + return + IteratorRange::active_cell_iterator> + (begin_active(), end()); +} + + + +template +IteratorRange::level_cell_iterator> +DoFHandler::mg_cell_iterators () const +{ + return + IteratorRange::level_cell_iterator> + (begin_mg(), end_mg()); +} + + + + +template +IteratorRange::cell_iterator> +DoFHandler::cell_iterators_on_level (const unsigned int level) const +{ + return + IteratorRange::cell_iterator> + (begin(level), end(level)); +} + + + +template +IteratorRange::active_cell_iterator> +DoFHandler::active_cell_iterators_on_level (const unsigned int level) const +{ + return + IteratorRange::active_cell_iterator> + (begin_active(level), end_active(level)); +} + + + +template +IteratorRange::level_cell_iterator> +DoFHandler::mg_cell_iterators_on_level (const unsigned int level) const +{ + return + IteratorRange::level_cell_iterator> + (begin_mg(level), end_mg(level)); +} + + + //--------------------------------------------------------------------------- diff --git a/source/grid/tria.cc b/source/grid/tria.cc index 14da3cc2b9..ca640d4ed7 100644 --- a/source/grid/tria.cc +++ b/source/grid/tria.cc @@ -10567,6 +10567,49 @@ Triangulation::end_active (const unsigned int level) const } + +template +IteratorRange::cell_iterator> +Triangulation::cell_iterators () const +{ + return + IteratorRange::cell_iterator> + (begin(), end()); +} + + +template +IteratorRange::active_cell_iterator> +Triangulation::active_cell_iterators () const +{ + return + IteratorRange::active_cell_iterator> + (begin_active(), end()); +} + + + +template +IteratorRange::cell_iterator> +Triangulation::cell_iterators_on_level (const unsigned int level) const +{ + return + IteratorRange::cell_iterator> + (begin(level), end(level)); +} + + + +template +IteratorRange::active_cell_iterator> +Triangulation::active_cell_iterators_on_level (const unsigned int level) const +{ + return + IteratorRange::active_cell_iterator> + (begin_active(level), end_active(level)); +} + + /*------------------------ Face iterator functions ------------------------*/ diff --git a/source/hp/dof_handler.cc b/source/hp/dof_handler.cc index 2f271224d7..8c821f0143 100644 --- a/source/hp/dof_handler.cc +++ b/source/hp/dof_handler.cc @@ -1740,6 +1740,49 @@ namespace hp + template + IteratorRange::cell_iterator> + DoFHandler::cell_iterators () const + { + return + IteratorRange::cell_iterator> + (begin(), end()); + } + + + template + IteratorRange::active_cell_iterator> + DoFHandler::active_cell_iterators () const + { + return + IteratorRange::active_cell_iterator> + (begin_active(), end()); + } + + + + template + IteratorRange::cell_iterator> + DoFHandler::cell_iterators_on_level (const unsigned int level) const + { + return + IteratorRange::cell_iterator> + (begin(level), end(level)); + } + + + + template + IteratorRange::active_cell_iterator> + DoFHandler::active_cell_iterators_on_level (const unsigned int level) const + { + return + IteratorRange::active_cell_iterator> + (begin_active(level), end_active(level)); + } + + + //------------------------------------------------------------------ diff --git a/tests/deal.II/range_based_for_tria.cc b/tests/deal.II/range_based_for_tria.cc new file mode 100644 index 0000000000..8c924ffd08 --- /dev/null +++ b/tests/deal.II/range_based_for_tria.cc @@ -0,0 +1,86 @@ +// --------------------------------------------------------------------- +// $Id$ +// +// Copyright (C) 2014 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. +// +// --------------------------------------------------------------------- + + + +// Check range-based for loops for triangulations + +#include "../tests.h" +#include +#include +#include + +#include +#include + +#include +#include + + +template +void check() +{ + Triangulation tr; + GridGenerator::hyper_cube(tr); + tr.refine_global(2); + + + { + // set flags on active cells + tr.clear_user_flags (); + for (auto cell : tr.active_cell_iterators()) + cell->set_user_flag(); + + // now verify that it is really only the active cells + for (auto cell : tr.cell_iterators()) + Assert (cell->user_flag_set() == !cell->has_children(), + ExcInternalError()); + } + + // now do the same again for all levels of the triangulation + for (unsigned int l=0; lset_user_flag(); + + for (auto cell : tr.cell_iterators_on_level(l)) + Assert (cell->user_flag_set() == !cell->has_children(), + ExcInternalError()); + + for (auto cell : tr.cell_iterators()) + Assert ((cell->user_flag_set() == !cell->has_children()) + || + (l != cell->level()), + ExcInternalError()); + } + + deallog << "OK" << std::endl; +} + + +int main() +{ + deal_II_exceptions::disable_abort_on_exception(); + + std::ofstream logfile("output"); + deallog.attach(logfile); + deallog.depth_console(0); + deallog.threshold_double(1.e-10); + + check<2>(); + check<3>(); +} diff --git a/tests/deal.II/range_based_for_tria.with_cxx11=true.output b/tests/deal.II/range_based_for_tria.with_cxx11=true.output new file mode 100644 index 0000000000..8b3b075900 --- /dev/null +++ b/tests/deal.II/range_based_for_tria.with_cxx11=true.output @@ -0,0 +1,3 @@ + +DEAL::OK +DEAL::OK -- 2.39.5