From: Jean-Paul Pelteret Date: Fri, 15 Jul 2022 23:27:21 +0000 (+0200) Subject: Alternative implementation to convert iterators X-Git-Tag: v9.5.0-rc1~1071^2~1 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=482dd5bf8faf95eb5d32870ea8572a5a38a52e2a;p=dealii.git Alternative implementation to convert iterators --- diff --git a/doc/news/changes/minor/20220715Pelteret b/doc/news/changes/minor/20220715Pelteret index 8cd72e453a..fea68c0116 100644 --- a/doc/news/changes/minor/20220715Pelteret +++ b/doc/news/changes/minor/20220715Pelteret @@ -1,5 +1,6 @@ -New: The convert_active_cell_iterator() function can be used to convert from a -Triangulation active cell iterator to a DoFHandler active cell iterator, or -from an active cell iterator of one DoFHandler to that of another DoFHandler. +New: The CellAccessor::as_dof_handler_iterator() function can be used +to convert from a Triangulation active cell iterator to a DoFHandler active cell +iterator, or from an active cell iterator of one DoFHandler to that of another +DoFHandler.
(Jean-Paul Pelteret, 2022/07/15) diff --git a/include/deal.II/dofs/dof_handler.h b/include/deal.II/dofs/dof_handler.h index f3beda6d13..87926405a4 100644 --- a/include/deal.II/dofs/dof_handler.h +++ b/include/deal.II/dofs/dof_handler.h @@ -1674,66 +1674,6 @@ private: #endif }; - -/** - * A function that converts a Triangulation active cell iterator to a - * DoFHandler active cell iterator. The triangulation @p iterator must be - * associated with the triangulation of the @p dof_handler. - * - * @param iterator The input iterator to be converted. - * @param dof_handler The DoFHandler for the output active cell iterator. - * @return An active cell iterator for the @p dof_handler, matching the cell - * referenced by the input triangulation @p iterator. - */ -template -typename DoFHandler::active_cell_iterator -convert_active_cell_iterator( - const typename Triangulation::active_cell_iterator &iterator, - const DoFHandler &dof_handler) -{ - Assert( - &iterator->get_triangulation() == &dof_handler.get_triangulation(), - ExcMessage( - "The triangulation associated with the iterator does not match that of the dof_handler.")); - - return typename DoFHandler::active_cell_iterator( - &dof_handler.get_triangulation(), - iterator->level(), - iterator->index(), - &dof_handler); -} - - -/** - * A function that converts a DoFHandler active cell iterator to an active - * cell iterator of another DoFHandler. The degree-of-freedom handler to which - * @p iterator is associated must have the same triangulation as the - * second @p dof_handler. - * - * @param iterator The input iterator to be converted. - * @param dof_handler The DoFHandler for the output active cell iterator. - * @return An active cell iterator for the @p dof_handler, matching the cell - * referenced by the input @p iterator for another DoFHandler. - */ -template -typename DoFHandler::active_cell_iterator -convert_active_cell_iterator( - const typename DoFHandler::active_cell_iterator &iterator, - const DoFHandler & dof_handler) -{ - Assert( - &iterator->get_triangulation() == &dof_handler.get_triangulation(), - ExcMessage( - "The triangulation associated with the iterator does not match that of the dof_handler.")); - - return typename DoFHandler::active_cell_iterator( - &dof_handler.get_triangulation(), - iterator->level(), - iterator->index(), - &dof_handler); -} - - namespace internal { namespace hp diff --git a/include/deal.II/grid/tria_accessor.h b/include/deal.II/grid/tria_accessor.h index 6d8b32f152..c3ea1f39f5 100644 --- a/include/deal.II/grid/tria_accessor.h +++ b/include/deal.II/grid/tria_accessor.h @@ -55,6 +55,12 @@ namespace parallel class TriangulationBase; } +template +class DoFHandler; +template +class DoFCellAccessor; + + template class Manifold; @@ -3135,6 +3141,31 @@ public: CellAccessor & operator=(CellAccessor &&) = default; // NOLINT + /** + * @} + */ + + /** + * @name Converting iterators + */ + /** + * @{ + */ + + /** + * A function that converts a Triangulation active cell iterator to a + * DoFHandler active cell iterator, or a DoFHandler active cell iterator + * to an active cell iterator of another DoFHandler. The @p iterator must be + * associated with the triangulation of the @p dof_handler. + * + * @param dof_handler The DoFHandler for the output active cell iterator. + * @return An active cell iterator for the @p dof_handler, matching the cell + * referenced by the input @p iterator. The type of the + * returned object is a DoFHandler::active_cell_iterator. + */ + TriaActiveIterator> + as_dof_handler_iterator(const DoFHandler &dof_handler) const; + /** * @} */ diff --git a/source/grid/tria_accessor.cc b/source/grid/tria_accessor.cc index 2b0b436f6f..dec542b95e 100644 --- a/source/grid/tria_accessor.cc +++ b/source/grid/tria_accessor.cc @@ -16,6 +16,8 @@ #include #include +#include + #include #include @@ -2116,6 +2118,27 @@ CellAccessor<3>::point_inside(const Point<3> &p) const /*------------------- Functions: CellAccessor -----------------*/ +// The return type is the same as DoFHandler::active_cell_iterator +template +TriaActiveIterator> +CellAccessor::as_dof_handler_iterator( + const DoFHandler &dof_handler) const +{ + Assert(is_active(), + ExcMessage("The current iterator points to an inactive cell. " + "You cannot convert it to an iterator to an active cell.")); + Assert(&this->get_triangulation() == &dof_handler.get_triangulation(), + ExcMessage("The triangulation associated with the iterator does not " + "match that of the DoFHandler.")); + + return typename DoFHandler::active_cell_iterator( + &dof_handler.get_triangulation(), + this->level(), + this->index(), + &dof_handler); +} + + // For codim>0 we proceed as follows: // 1) project point onto manifold and // 2) transform to the unit cell with a Q1 mapping diff --git a/tests/dofs/dof_iterators_01.cc b/tests/dofs/dof_iterators_01.cc index 862ab254cf..6220167488 100644 --- a/tests/dofs/dof_iterators_01.cc +++ b/tests/dofs/dof_iterators_01.cc @@ -26,6 +26,8 @@ #include #include +#include + #include "../tests.h" template @@ -39,15 +41,21 @@ test() DoFHandler dof_handler(triangulation); dof_handler.distribute_dofs(fe); + std::vector local_dof_indices(fe.n_dofs_per_cell()); + for (const auto &it_tria : triangulation.active_cell_iterators()) { - const auto it_dh = convert_active_cell_iterator(it_tria, dof_handler); + const auto it_dh = it_tria->as_dof_handler_iterator(dof_handler); Assert(it_tria->level() == it_dh->level(), ExcMessage("Iterator conversion failed: Level.")); Assert(it_tria->index() == it_dh->index(), ExcMessage("Iterator conversion failed: Index.")); Assert(it_tria->id() == it_dh->id(), ExcMessage("Iterator conversion failed: Id.")); + + // Check that some basic features work (i.e. that we have the right + // accessor type) + it_dh->get_dof_indices(local_dof_indices); } } diff --git a/tests/dofs/dof_iterators_02.cc b/tests/dofs/dof_iterators_02.cc index 98a1609b1d..e45e617b79 100644 --- a/tests/dofs/dof_iterators_02.cc +++ b/tests/dofs/dof_iterators_02.cc @@ -26,6 +26,8 @@ #include #include +#include + #include "../tests.h" template @@ -42,15 +44,22 @@ test() dof_handler_1.distribute_dofs(fe_1); dof_handler_2.distribute_dofs(fe_2); + std::vector local_dof_indices( + fe_2.n_dofs_per_cell()); + for (const auto &it_dh_1 : dof_handler_1.active_cell_iterators()) { - const auto it_dh_2 = convert_active_cell_iterator(it_dh_1, dof_handler_2); + const auto it_dh_2 = it_dh_1->as_dof_handler_iterator(dof_handler_2); Assert(it_dh_1->level() == it_dh_2->level(), ExcMessage("Iterator conversion failed: Level.")); Assert(it_dh_1->index() == it_dh_2->index(), ExcMessage("Iterator conversion failed: Index.")); Assert(it_dh_1->id() == it_dh_2->id(), ExcMessage("Iterator conversion failed: Id.")); + + // Check that some basic features work (i.e. that we have the right + // accessor type) + it_dh_2->get_dof_indices(local_dof_indices); } }