]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Allow hp::DoFHandler with parallel::DistributedTriangulationBase 10036/head
authorSebastian Stark <stark.sebastian@gmx.de>
Mon, 16 Sep 2019 07:56:55 +0000 (09:56 +0200)
committerSebastian Stark <stark.sebastian@gmx.de>
Wed, 6 May 2020 20:33:08 +0000 (22:33 +0200)
This allows the distribution of dofs for an hp::DoFHandler combined with a triangulation class derived from parallel::DistributedTriangulationBase other than parallel::distributed::Triangulation.
This makes it possible to handle the triangulation without p4est (e.g. manual partitioning).

doc/news/changes/minor/20200506SebastianStark [new file with mode: 0644]
source/hp/dof_handler.cc

diff --git a/doc/news/changes/minor/20200506SebastianStark b/doc/news/changes/minor/20200506SebastianStark
new file mode 100644 (file)
index 0000000..e1358ad
--- /dev/null
@@ -0,0 +1,3 @@
+Changed: Allow hp::DoFHandler to distribute dofs based on a triangulation class derived from parallel::DistributedTriangulationBase other than parallel::distributed::Triangulation.
+<br>
+(Sebastian Stark, 2020/05/06)
index 2068ac9c1e48bd9d6ce2d0b3888b5cc22c586c28..e6187b0e9b0abc8f70a07fb8407ce871c589065e 100644 (file)
@@ -968,7 +968,7 @@ namespace internal
          * For parallel::shared::Triangulation objects,
          * this information is distributed on both ghost and artificial cells.
          *
-         * In case a parallel::distributed::Triangulation is used,
+         * In case a parallel::DistributedTriangulationBase is used,
          * indices are communicated only to ghost cells.
          */
         template <int dim, int spacedim>
@@ -1014,11 +1014,12 @@ namespace internal
                     cell->index(),
                     active_fe_indices[cell->active_cell_index()]);
             }
-          else if (const dealii::parallel::distributed::Triangulation<dim,
-                                                                      spacedim>
-                     *tr = dynamic_cast<const dealii::parallel::distributed::
-                                          Triangulation<dim, spacedim> *>(
-                       &dof_handler.get_triangulation()))
+          else if (const dealii::parallel::
+                     DistributedTriangulationBase<dim, spacedim> *tr =
+                       dynamic_cast<
+                         const dealii::parallel::
+                           DistributedTriangulationBase<dim, spacedim> *>(
+                         &dof_handler.get_triangulation()))
             {
               // For completely distributed meshes, use the function that is
               // able to move data from locally owned cells on one processor to
@@ -1217,7 +1218,7 @@ namespace internal
 
         /**
          * Coarsening strategy for the CellDataTransfer object responsible for
-         * tranferring the active_fe_index of each cell on
+         * transferring the active_fe_index of each cell on
          * parallel::distributed::Triangulation objects that have been refined.
          *
          * A finite element index needs to be determined for the (not yet
@@ -1716,7 +1717,7 @@ namespace hp
     // decide whether we need a sequential or a parallel shared/distributed
     // policy and attach corresponding callback functions dealing with the
     // transfer of active_fe_indices
-    if (dynamic_cast<const parallel::distributed::Triangulation<dim, spacedim>
+    if (dynamic_cast<const parallel::DistributedTriangulationBase<dim, spacedim>
                        *>(&this->get_triangulation()))
       {
         policy = std_cxx14::make_unique<
@@ -2024,8 +2025,18 @@ namespace hp
   DoFHandler<dim, spacedim>::pre_distributed_active_fe_index_transfer()
   {
 #ifndef DEAL_II_WITH_P4EST
-    Assert(false, ExcInternalError());
+    Assert(false,
+           ExcMessage(
+             "You are attempting to use a functionality that is only available "
+             "if deal.II was configured to use p4est, but cmake did not find a "
+             "valid p4est library."));
 #else
+    // the implementation below requires a p:d:T currently
+    Assert((dynamic_cast<
+              const parallel::distributed::Triangulation<dim, spacedim> *>(
+              &this->get_triangulation()) != nullptr),
+           ExcNotImplemented());
+
     // Finite elements need to be assigned to each cell by calling
     // distribute_dofs() first to make this functionality available.
     if (fe_collection.size() > 0)
@@ -2107,8 +2118,18 @@ namespace hp
   DoFHandler<dim, spacedim>::post_distributed_active_fe_index_transfer()
   {
 #ifndef DEAL_II_WITH_P4EST
-    Assert(false, ExcInternalError());
+    Assert(false,
+           ExcMessage(
+             "You are attempting to use a functionality that is only available "
+             "if deal.II was configured to use p4est, but cmake did not find a "
+             "valid p4est library."));
 #else
+    // the implementation below requires a p:d:T currently
+    Assert((dynamic_cast<
+              const parallel::distributed::Triangulation<dim, spacedim> *>(
+              &this->get_triangulation()) != nullptr),
+           ExcNotImplemented());
+
     // Finite elements need to be assigned to each cell by calling
     // distribute_dofs() first to make this functionality available.
     if (fe_collection.size() > 0)
@@ -2147,11 +2168,11 @@ namespace hp
              "if deal.II was configured to use p4est, but cmake did not find a "
              "valid p4est library."));
 #else
-    Assert(
-      (dynamic_cast<const parallel::distributed::Triangulation<dim, spacedim>
-                      *>(&this->get_triangulation()) != nullptr),
-      ExcMessage(
-        "This functionality requires a parallel::distributed::Triangulation object."));
+    // the implementation below requires a p:d:T currently
+    Assert((dynamic_cast<
+              const parallel::distributed::Triangulation<dim, spacedim> *>(
+              &this->get_triangulation()) != nullptr),
+           ExcNotImplemented());
 
     // Finite elements need to be assigned to each cell by calling
     // distribute_dofs() first to make this functionality available.
@@ -2204,6 +2225,12 @@ namespace hp
              "if deal.II was configured to use p4est, but cmake did not find a "
              "valid p4est library."));
 #else
+    // the implementation below requires a p:d:T currently
+    Assert((dynamic_cast<
+              const parallel::distributed::Triangulation<dim, spacedim> *>(
+              &this->get_triangulation()) != nullptr),
+           ExcNotImplemented());
+
     if (fe_collection.size() > 0)
       {
         Assert(active_fe_index_transfer != nullptr, ExcInternalError());
@@ -2227,11 +2254,11 @@ namespace hp
              "if deal.II was configured to use p4est, but cmake did not find a "
              "valid p4est library."));
 #else
-    Assert(
-      (dynamic_cast<const parallel::distributed::Triangulation<dim, spacedim>
-                      *>(&this->get_triangulation()) != nullptr),
-      ExcMessage(
-        "This functionality requires a parallel::distributed::Triangulation object."));
+    // the implementation below requires a p:d:T currently
+    Assert((dynamic_cast<
+              const parallel::distributed::Triangulation<dim, spacedim> *>(
+              &this->get_triangulation()) != nullptr),
+           ExcNotImplemented());
 
     // Finite elements need to be assigned to each cell by calling
     // distribute_dofs() first to make this functionality available.

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.