]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Alternative implementation to convert iterators
authorJean-Paul Pelteret <jppelteret@gmail.com>
Fri, 15 Jul 2022 23:27:21 +0000 (01:27 +0200)
committerJean-Paul Pelteret <jppelteret@gmail.com>
Tue, 19 Jul 2022 21:33:40 +0000 (23:33 +0200)
doc/news/changes/minor/20220715Pelteret
include/deal.II/dofs/dof_handler.h
include/deal.II/grid/tria_accessor.h
source/grid/tria_accessor.cc
tests/dofs/dof_iterators_01.cc
tests/dofs/dof_iterators_02.cc

index 8cd72e453a52db34cae09f139b028b112f7ef911..fea68c0116eaadda7147b3acc3517087b379ff56 100644 (file)
@@ -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.
 <br>
 (Jean-Paul Pelteret, 2022/07/15)
index f3beda6d13881c6167ddc0d898696f46390f33ed..87926405a4c696cafe7a78bf3c01e85558c2698b 100644 (file)
@@ -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 <int dim, int spacedim>
-typename DoFHandler<dim, spacedim>::active_cell_iterator
-convert_active_cell_iterator(
-  const typename Triangulation<dim, spacedim>::active_cell_iterator &iterator,
-  const DoFHandler<dim, spacedim> &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<dim, spacedim>::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 <int dim, int spacedim>
-typename DoFHandler<dim, spacedim>::active_cell_iterator
-convert_active_cell_iterator(
-  const typename DoFHandler<dim, spacedim>::active_cell_iterator &iterator,
-  const DoFHandler<dim, spacedim> &                               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<dim, spacedim>::active_cell_iterator(
-    &dof_handler.get_triangulation(),
-    iterator->level(),
-    iterator->index(),
-    &dof_handler);
-}
-
-
 namespace internal
 {
   namespace hp
index 6d8b32f152f26e860d37467ea991132154f4e980..c3ea1f39f5050a8335d8acd1c38d6b764ecf56b5 100644 (file)
@@ -55,6 +55,12 @@ namespace parallel
   class TriangulationBase;
 }
 
+template <int dim, int spacedim>
+class DoFHandler;
+template <int dim, int spacedim, bool lda>
+class DoFCellAccessor;
+
+
 template <int dim, int spacedim>
 class Manifold;
 
@@ -3135,6 +3141,31 @@ public:
   CellAccessor<dim, spacedim> &
   operator=(CellAccessor<dim, spacedim> &&) = 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<DoFCellAccessor<dim, spacedim, false>>
+  as_dof_handler_iterator(const DoFHandler<dim, spacedim> &dof_handler) const;
+
   /**
    * @}
    */
index 2b0b436f6f606746cdbd52132f69ebcaf79b181e..dec542b95e57aa5976d8c80e814c1ef953b69a55 100644 (file)
@@ -16,6 +16,8 @@
 #include <deal.II/base/geometry_info.h>
 #include <deal.II/base/quadrature.h>
 
+#include <deal.II/dofs/dof_accessor.h>
+
 #include <deal.II/fe/fe_q.h>
 #include <deal.II/fe/mapping_q1.h>
 
@@ -2116,6 +2118,27 @@ CellAccessor<3>::point_inside(const Point<3> &p) const
 
 /*------------------- Functions: CellAccessor<dim,spacedim> -----------------*/
 
+// The return type is the same as DoFHandler<dim,spacedim>::active_cell_iterator
+template <int dim, int spacedim>
+TriaActiveIterator<DoFCellAccessor<dim, spacedim, false>>
+CellAccessor<dim, spacedim>::as_dof_handler_iterator(
+  const DoFHandler<dim, spacedim> &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<dim, spacedim>::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
index 862ab254cf8f1fcc057553d4efb4969b94a9d136..62201674886dba19e56b7f9284ab6cf69e8ba7eb 100644 (file)
@@ -26,6 +26,8 @@
 #include <deal.II/grid/grid_generator.h>
 #include <deal.II/grid/tria.h>
 
+#include <deal.II/lac/full_matrix.h>
+
 #include "../tests.h"
 
 template <int dim, int spacedim>
@@ -39,15 +41,21 @@ test()
   DoFHandler<dim, spacedim> dof_handler(triangulation);
   dof_handler.distribute_dofs(fe);
 
+  std::vector<types::global_dof_index> 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);
     }
 }
 
index 98a1609b1d89195d9284bf6b35a3e0ee6967b58e..e45e617b796687812f190d27c024e7a4d7d5fd20 100644 (file)
@@ -26,6 +26,8 @@
 #include <deal.II/grid/grid_generator.h>
 #include <deal.II/grid/tria.h>
 
+#include <deal.II/lac/full_matrix.h>
+
 #include "../tests.h"
 
 template <int dim, int spacedim>
@@ -42,15 +44,22 @@ test()
   dof_handler_1.distribute_dofs(fe_1);
   dof_handler_2.distribute_dofs(fe_2);
 
+  std::vector<types::global_dof_index> 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);
     }
 }
 

In the beginning the Universe was created. This has made a lot of people very angry and has been widely regarded as a bad move.

Douglas Adams


Typeset in Trocchi and Trocchi Bold Sans Serif.