]> https://gitweb.dealii.org/ - dealii.git/commitdiff
Implementation of DoFTools::locally_relevant_dofs_per_subdomain.
authorJean-Paul Pelteret <jppelteret@gmail.com>
Thu, 24 Sep 2015 13:44:44 +0000 (15:44 +0200)
committerJean-Paul Pelteret <jppelteret@gmail.com>
Fri, 25 Sep 2015 10:48:07 +0000 (12:48 +0200)
Also added Denis' authorship for DofTools::locally_owned_dofs_per_subdomain

doc/news/changes.h
include/deal.II/dofs/dof_tools.h
source/dofs/dof_tools.cc
source/dofs/dof_tools.inst.in

index 673990896e0d9d7642369dae442bf8abbf030bda..5984fc3f5932d042a7ce3493a5d101357a8d9d2e 100644 (file)
@@ -356,6 +356,13 @@ inconvenience this causes.
 
 <ol>
 
+  <li> Added: DoFTools::locally_relevant_dofs_per_subdomain can be used
+  to extract an IndexSet of locally relevant DoFs for a Triangulation
+  partitioned using METIS or with a parallel::shared::Triangulation .
+  <br>
+  (Jean-Paul Pelteret, 2015/09/24)
+  </li>
+
   <li> Fixed: hp::SolutionTransfer could get confused when dealing with
   FE_Nothing elements. This is now fixed.
   <br>
index f4922d3b7914fcaf6bb0545ce01913ce889fbcf3..729dcb0420e041982c18589b54031fb5096e45d6 100644 (file)
@@ -1427,11 +1427,31 @@ namespace DoFTools
    * we do not have information about all cells of the triangulation available locally,
    * and consequently can not say anything definitive about the degrees of freedom active on other
    * processors' locally owned cells.
+   *
+   * @author Denis Davydov, 2015
    */
   template <class DH>
   std::vector<IndexSet>
   locally_owned_dofs_per_subdomain (const DH   &dof_handler);
 
+  /**
+   *
+   * For each processor, determine the set of locally relevant degrees of freedom as an IndexSet.
+   * This function then returns a vector of index sets, where the vector has size equal to the
+   * number of MPI processes that participate in the DoF handler object.
+   *
+   * The function can be used for objects of type dealii::Triangulation or parallel::shared::Triangulation.
+   * It will not work for objects of type parallel::distributed::Triangulation since for such triangulations
+   * we do not have information about all cells of the triangulation available locally,
+   * and consequently can not say anything definitive about the degrees of freedom active on other
+   * processors' locally owned cells.
+   *
+   * @author Jean-Paul Pelteret, 2015
+   */
+  template <class DH>
+  std::vector<IndexSet>
+  locally_relevant_dofs_per_subdomain (const DH   &dof_handler);
+
   /**
    * For each DoF, return in the output array to which subdomain (as given by
    * the <tt>cell->subdomain_id()</tt> function) it belongs. The output array
index 12665d0565c68900ca8b52810b856af040511d2e..16c0e42cd981d50c1fbd58270275efceb69f99ce 100644 (file)
@@ -27,6 +27,7 @@
 #include <deal.II/grid/tria.h>
 #include <deal.II/grid/tria_iterator.h>
 #include <deal.II/grid/intergrid_map.h>
+#include <deal.II/grid/filtered_iterator.h>
 #include <deal.II/grid/grid_tools.h>
 #include <deal.II/dofs/dof_handler.h>
 #include <deal.II/dofs/dof_accessor.h>
@@ -1123,6 +1124,60 @@ namespace DoFTools
     return index_sets;
   }
 
+  template <class DH>
+  std::vector<IndexSet>
+  locally_relevant_dofs_per_subdomain (const DH  &dof_handler)
+  {
+    // Collect all the locally owned DoFs
+    // Note: Even though the distribution of DoFs by the locally_owned_dofs_per_subdomain
+    // function is pseudo-random, we will collect all the DoFs on the subdomain
+    // and its layer cell. Therefore, the random nature of this function does
+    // not play a role in the extraction of the locally relevant DoFs
+    std::vector<IndexSet> dof_set = locally_owned_dofs_per_subdomain(dof_handler);
+    const dealii::types::subdomain_id n_subdomains = dof_set.size();
+
+    // Add the DoFs on the adjacent (equivalent ghost) cells to the IndexSet,
+    // cache them in a set. Need to check each DoF manually because we can't
+    // be sure that the DoF range of locally_owned_dofs is really contiguous.
+    for (dealii::types::subdomain_id subdomain_id = 0;
+         subdomain_id < n_subdomains; ++subdomain_id)
+      {
+        // Extract the layer of cells around this subdomain
+        std_cxx11::function<bool (const typename DH::active_cell_iterator &)> predicate
+          = IteratorFilters::SubdomainEqualTo(subdomain_id);
+        const std::vector<typename DH::active_cell_iterator>
+        active_halo_layer = GridTools::compute_active_cell_halo_layer (dof_handler,
+                            predicate);
+
+        // Extract DoFs associated with halo layer
+        std::vector<types::global_dof_index> local_dof_indices;
+        std::set<types::global_dof_index> subdomain_halo_global_dof_indices;
+        for (typename std::vector<typename DH::active_cell_iterator>::const_iterator
+             it_cell = active_halo_layer.begin(); it_cell!=active_halo_layer.end(); ++it_cell)
+          {
+            const typename DH::active_cell_iterator &cell = *it_cell;
+            Assert(cell->subdomain_id() != subdomain_id,
+                   ExcMessage("The subdomain ID of the halo cell should not match that of the vector entry."));
+
+            local_dof_indices.resize(cell->get_fe().dofs_per_cell);
+            cell->get_dof_indices(local_dof_indices);
+
+            for (std::vector<types::global_dof_index>::iterator it=local_dof_indices.begin();
+                 it!=local_dof_indices.end();
+                 ++it)
+              if (!dof_set[subdomain_id].is_element(*it))
+                subdomain_halo_global_dof_indices.insert(*it);
+          }
+
+        dof_set[subdomain_id].add_indices(subdomain_halo_global_dof_indices.begin(),
+                                          subdomain_halo_global_dof_indices.end());
+
+        dof_set[subdomain_id].compress();
+      }
+
+    return dof_set;
+  }
+
   template <class DH>
   void
   get_subdomain_association (const DH                  &dof_handler,
index 0ac9e83efbf49f3ba0e115539f0de1a3934444b8..607b040036f7ee2c2c39ad73846f77845cdc945f 100644 (file)
@@ -279,6 +279,15 @@ std::vector<IndexSet>
 DoFTools::locally_owned_dofs_per_subdomain<hp::DoFHandler<deal_II_dimension> >
 (const hp::DoFHandler<deal_II_dimension> &dof_handler);
 
+template
+std::vector<IndexSet>
+DoFTools::locally_relevant_dofs_per_subdomain<DoFHandler<deal_II_dimension> >
+(const DoFHandler<deal_II_dimension> &dof_handler);
+template
+std::vector<IndexSet>
+DoFTools::locally_relevant_dofs_per_subdomain<hp::DoFHandler<deal_II_dimension> >
+(const hp::DoFHandler<deal_II_dimension> &dof_handler);
+
 template
 unsigned int
 DoFTools::count_dofs_with_subdomain_association<DoFHandler<deal_II_dimension> >
@@ -380,6 +389,15 @@ std::vector<IndexSet>
 DoFTools::locally_owned_dofs_per_subdomain<hp::DoFHandler<deal_II_dimension,deal_II_dimension+1> >
 (const hp::DoFHandler<deal_II_dimension,deal_II_dimension+1> &dof_handler);
 
+template
+std::vector<IndexSet>
+DoFTools::locally_relevant_dofs_per_subdomain<DoFHandler<deal_II_dimension,deal_II_dimension+1> >
+(const DoFHandler<deal_II_dimension,deal_II_dimension+1> &dof_handler);
+template
+std::vector<IndexSet>
+DoFTools::locally_relevant_dofs_per_subdomain<hp::DoFHandler<deal_II_dimension,deal_II_dimension+1> >
+(const hp::DoFHandler<deal_II_dimension,deal_II_dimension+1> &dof_handler);
+
 template
 void
 DoFTools::count_dofs_with_subdomain_association<hp::DoFHandler<deal_II_dimension,deal_II_dimension+1> >
@@ -425,6 +443,15 @@ std::vector<IndexSet>
 DoFTools::locally_owned_dofs_per_subdomain<hp::DoFHandler<1,3> >
 (const hp::DoFHandler<1,3> &dof_handler);
 
+template
+std::vector<IndexSet>
+DoFTools::locally_relevant_dofs_per_subdomain<DoFHandler<1,3> >
+(const DoFHandler<1,3>     &dof_handler);
+template
+std::vector<IndexSet>
+DoFTools::locally_relevant_dofs_per_subdomain<hp::DoFHandler<1,3> >
+(const hp::DoFHandler<1,3> &dof_handler);
+
 template
 unsigned int
 DoFTools::count_dofs_with_subdomain_association<DoFHandler<1,3> >

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.