]> https://gitweb.dealii.org/ - dealii.git/commitdiff
MGTransferGlobalCoarseningTools::create_geometric_coarsening_sequence load forest... 11715/head
authorPeter Munch <peterrmuench@gmail.com>
Tue, 9 Feb 2021 09:30:48 +0000 (10:30 +0100)
committerPeter Munch <peterrmuench@gmail.com>
Tue, 9 Feb 2021 17:05:52 +0000 (18:05 +0100)
doc/news/changes/minor/20210209Munch [new file with mode: 0644]
include/deal.II/distributed/p4est_wrappers.h
include/deal.II/distributed/tria.h
include/deal.II/multigrid/mg_transfer_global_coarsening.templates.h
source/distributed/p4est_wrappers.cc
source/distributed/tria.cc

diff --git a/doc/news/changes/minor/20210209Munch b/doc/news/changes/minor/20210209Munch
new file mode 100644 (file)
index 0000000..e8e7e83
--- /dev/null
@@ -0,0 +1,6 @@
+New: The method parallel::distributed::Triagnulation::load() 
+can now also accept a p4est/p8est forest, which can be queried
+from an existing triangulation via parallel::distributed::Triagnulation::get_p4est().
+<br>
+(Marc Fehling, Peter Munch, 2021/02/09)
+
index 0d48339d3ad021a187576de5f14c7131a3c13348..27894a5a7407877390a67c3a96aa8af5c92d1866 100644 (file)
@@ -180,6 +180,9 @@ namespace internal
         p4est_init_t            init_fn,
         void *                  user_pointer);
 
+      static types<2>::forest *(&copy_forest)(types<2>::forest *input,
+                                              int               copy_data);
+
       static void (&destroy)(types<2>::forest *p4est);
 
       static void (&refine)(types<2>::forest *p4est,
@@ -350,6 +353,9 @@ namespace internal
         p8est_init_t            init_fn,
         void *                  user_pointer);
 
+      static types<3>::forest *(&copy_forest)(types<3>::forest *input,
+                                              int               copy_data);
+
       static void (&destroy)(types<3>::forest *p8est);
 
       static void (&refine)(types<3>::forest *p8est,
index 272ed414a6d83068a197bed3bdaeebdb7ad0002b..238db08b55c627375bd8552e59d94ba86a7ebf38 100644 (file)
@@ -604,6 +604,14 @@ namespace parallel
       load(const std::string &filename,
            const bool         autopartition = true) override;
 
+      /**
+       * Load the refinement information from a given parallel forest. This
+       * forest might be obtained from the function call to
+       * parallel::distributed::Triangulation::get_p4est();
+       */
+      void
+      load(const typename dealii::internal::p4est::types<dim>::forest *forest);
+
       /**
        * Return a permutation vector for the order the coarse cells are handed
        * off to p4est. For example the value of the $i$th element in this
index 4ba56a2f64c8fe396e49c9ea37e0e87009d42c51..9c6bf048e8992c87db46a2a69c6f12086439bcc8 100644 (file)
@@ -1694,6 +1694,30 @@ namespace MGTransferGlobalCoarseningTools
   }
 
 
+  namespace internal
+  {
+    template <int dim, int spacedim>
+    void
+    load_forest(const parallel::distributed::Triangulation<dim, spacedim> &in,
+                parallel::distributed::Triangulation<dim, spacedim> &      out)
+    {
+#ifndef DEAL_II_WITH_P4EST
+      (void)in;
+      (void)out;
+#else
+      out.load(in.get_p4est());
+#endif
+    }
+
+    template <int spacedim>
+    void
+    load_forest(const parallel::distributed::Triangulation<1, spacedim> &in,
+                parallel::distributed::Triangulation<1, spacedim> &      out)
+    {
+      (void)in;
+      (void)out;
+    }
+  } // namespace internal
 
   template <int dim, int spacedim>
   std::vector<std::shared_ptr<Triangulation<dim, spacedim>>>
@@ -1730,13 +1754,6 @@ namespace MGTransferGlobalCoarseningTools
     // create coarse meshes
     for (unsigned int l = max_level; l > 0; --l)
       {
-        // store finer mesh to file
-        dynamic_cast<
-          const parallel::distributed::Triangulation<dim, spacedim> *>(
-          coarse_grid_triangulations[l].get())
-          ->save("mesh");
-        MPI_Barrier(fine_triangulation->get_communicator());
-
         // create empty triangulation
         auto new_tria =
           std::make_shared<parallel::distributed::Triangulation<dim, spacedim>>(
@@ -1753,8 +1770,11 @@ namespace MGTransferGlobalCoarseningTools
             new_tria->set_manifold(i, fine_triangulation->get_manifold(i));
 
         // create refinement hierarchy (by loading stored mesh)
-        new_tria->load("mesh", false);
-        MPI_Barrier(fine_triangulation->get_communicator());
+        internal::load_forest(
+          *dynamic_cast<
+            const parallel::distributed::Triangulation<dim, spacedim> *>(
+            coarse_grid_triangulations[l].get()),
+          *new_tria);
 
         // coarsen mesh
         new_tria->coarsen_global();
index c1de26ef33500e4e41a2b5885437d3036a834379..07b15cdb744a86e6571a2ff006930228ad10740f 100644 (file)
@@ -401,6 +401,9 @@ namespace internal
       p4est_init_t            init_fn,
       void *                  user_pointer) = p4est_new_ext;
 
+    types<2>::forest *(&functions<2>::copy_forest)(types<2>::forest *input,
+                                                   int copy_data) = p4est_copy;
+
     void (&functions<2>::destroy)(types<2>::forest *p4est) = p4est_destroy;
 
     void (&functions<2>::refine)(types<2>::forest *p4est,
@@ -588,6 +591,9 @@ namespace internal
       p8est_init_t            init_fn,
       void *                  user_pointer) = p8est_new_ext;
 
+    types<3>::forest *(&functions<3>::copy_forest)(types<3>::forest *input,
+                                                   int copy_data) = p8est_copy;
+
     void (&functions<3>::destroy)(types<3>::forest *p8est) = p8est_destroy;
 
     void (&functions<3>::refine)(types<3>::forest *p8est,
index 60898ada3540e1f62a2794989f9fbc9df6b9178c..cf2b51eae344df8074dafb893665599bf7ba9eb6 100644 (file)
@@ -1608,6 +1608,49 @@ namespace parallel
     }
 
 
+    template <int dim, int spacedim>
+    void
+    Triangulation<dim, spacedim>::load(
+      const typename dealii::internal::p4est::types<dim>::forest *forest)
+    {
+      if (parallel_ghost != nullptr)
+        {
+          dealii::internal::p4est::functions<dim>::ghost_destroy(
+            parallel_ghost);
+          parallel_ghost = nullptr;
+        }
+      dealii::internal::p4est::functions<dim>::destroy(parallel_forest);
+      parallel_forest = nullptr;
+
+      // note: we can keep the connectivity, since the coarse grid does not
+      // change
+
+      // create deep copy
+      typename dealii::internal::p4est::types<dim>::forest *temp =
+        const_cast<typename dealii::internal::p4est::types<dim>::forest *>(
+          forest);
+      parallel_forest =
+        dealii::internal::p4est::functions<dim>::copy_forest(temp, false);
+      parallel_forest->user_pointer = this;
+
+      try
+        {
+          copy_local_forest_to_triangulation();
+        }
+      catch (const typename Triangulation<dim>::DistortedCellList &)
+        {
+          // the underlying
+          // triangulation should not
+          // be checking for
+          // distorted cells
+          Assert(false, ExcInternalError());
+        }
+
+      this->update_periodic_face_map();
+      this->update_number_cache();
+    }
+
+
 
     template <int dim, int spacedim>
     unsigned int

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.