]> https://gitweb.dealii.org/ - dealii-svn.git/commitdiff
added locally_owned_dofs_with_subdomain into DoFTools and use it in the new policy...
authordenis.davydov <denis.davydov@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 21 Mar 2014 15:27:18 +0000 (15:27 +0000)
committerdenis.davydov <denis.davydov@0785d39b-7218-0410-832d-ea1e28bc413d>
Fri, 21 Mar 2014 15:27:18 +0000 (15:27 +0000)
git-svn-id: https://svn.dealii.org/branches/branch_sharedtria@32681 0785d39b-7218-0410-832d-ea1e28bc413d

deal.II/include/deal.II/dofs/dof_handler_policy.h
deal.II/include/deal.II/dofs/dof_tools.h
deal.II/include/deal.II/hp/dof_handler.h
deal.II/source/dofs/dof_handler.cc
deal.II/source/dofs/dof_handler_policy.cc
deal.II/source/dofs/dof_handler_policy.inst.in
deal.II/source/dofs/dof_tools.cc
deal.II/source/dofs/dof_tools.inst.in
deal.II/source/hp/dof_handler.cc

index 16ef98281bb8d20f4f31ae49864d4443ce554166..f588ba63f8f2c325f1ec57d492c4cf2767d7fb24 100644 (file)
@@ -22,6 +22,7 @@
 #include <deal.II/base/config.h>
 #include <deal.II/base/exceptions.h>
 #include <deal.II/base/template_constraints.h>
+#include <deal.II/dofs/dof_tools.h>
 
 #include <vector>
 #include <map>
index 22a6a74efa38bcda9043beeaf722fc440c503dbb..ab0d3ed1f03c9a3da5ad150c32432034537b2e65 100644 (file)
@@ -1372,6 +1372,11 @@ namespace DoFTools
   extract_locally_relevant_dofs (const DH &dof_handler,
                                  IndexSet &dof_set);
 
+  template <class DH>
+  IndexSet
+  locally_owned_dofs_with_subdomain (const DH                  &dof_handler,
+                                            const types::subdomain_id  subdomain   );
+
   /**
    * For each DoF, return in the output array to which subdomain (as
    * given by the <tt>cell->subdomain_id()</tt> function) it
index 8b0c2b1ace78a288f534d71d405a0194316e2a9c..ef0e929767d96b6ad5ed0bd0d67130a125a58fe6 100644 (file)
@@ -26,6 +26,7 @@
 #include <deal.II/dofs/function_map.h>
 #include <deal.II/dofs/dof_iterator_selector.h>
 #include <deal.II/dofs/number_cache.h>
+//#include <deal.II/dofs/dof_tools.h>
 #include <deal.II/hp/fe_collection.h>
 
 #include <vector>
index 503e2357d5b630f2495b2aa661235b1973b053f6..a87e4fcb04866844335bc234822c50e5df2bbe60 100644 (file)
@@ -73,6 +73,8 @@ namespace internal
       policy_name = "Policy::Sequential<";
     else if (dynamic_cast<const typename dealii::internal::DoFHandler::Policy::ParallelDistributed<dim,spacedim>*>(&policy))
       policy_name = "Policy::ParallelDistributed<";
+    else if (dynamic_cast<const typename dealii::internal::DoFHandler::Policy::ParallelShared<dim,spacedim>*>(&policy))
+      policy_name = "Policy::ParallelShared<";
     else
       AssertThrow(false, ExcNotImplemented());
     policy_name += Utilities::int_to_string(dim)+
@@ -763,6 +765,11 @@ DoFHandler<dim,spacedim>::DoFHandler (const Triangulation<dim,spacedim> &tria)
   // decide whether we need a
   // sequential or a parallel
   // distributed policy
+  if (dynamic_cast<const parallel::shared::Triangulation< dim, spacedim>*>
+         (&tria)
+         != 0)
+       policy.reset (new internal::DoFHandler::Policy::ParallelShared<dim,spacedim>());
+  else
   if (dynamic_cast<const parallel::distributed::Triangulation< dim, spacedim >*>
       (&tria)
       == 0)
index 30926ebeba0c12eae3a76091bcba815ef879c72f..da667b0400bc6f364c1e908334c3d6996e1ecfc1 100644 (file)
@@ -971,8 +971,7 @@ namespace internal
       {
          NumberCache number_cache = Sequential<dim,spacedim>::distribute_dofs (dof_handler);
          //correct number_cache:
-
-
+         number_cache.locally_owned_dofs = dealii::DoFTools::locally_owned_dofs_with_subdomain(dof_handler,dof_handler.get_tria().locally_owned_subdomain() );
          return number_cache;
       }
 
@@ -982,9 +981,8 @@ namespace internal
       distribute_mg_dofs (DoFHandler<dim,spacedim> &dof_handler,
                           std::vector<NumberCache> &number_caches) const
       {
-         NumberCache number_cache = Sequential<dim,spacedim>:: distribute_mg_dofs (dof_handler, number_caches);
+         Sequential<dim,spacedim>:: distribute_mg_dofs (dof_handler, number_caches);
          Assert(false,ExcMessage("Not implemented"));
-         return number_cache;
       }
 
       template <int dim, int spacedim>
@@ -995,7 +993,7 @@ namespace internal
       {
          NumberCache number_cache = Sequential<dim,spacedim>::renumber_dofs (new_numbers,dof_handler);
          //correct number_cache:
-
+         number_cache.locally_owned_dofs = dealii::DoFTools::locally_owned_dofs_with_subdomain(dof_handler,dof_handler.get_tria().locally_owned_subdomain() );
          return number_cache;
       }
 
index 9b86beb03f213fecafa74be4337e3a66dc0f293f..b1fbf341f6ebbdc8aa679133096e86a5a4029eab 100644 (file)
@@ -25,17 +25,20 @@ namespace internal
     \{
       template class PolicyBase<deal_II_dimension,deal_II_dimension>;
       template class Sequential<deal_II_dimension,deal_II_dimension>;
+      template class ParallelShared<deal_II_dimension,deal_II_dimension>;
       template class ParallelDistributed<deal_II_dimension,deal_II_dimension>;
 
 #if deal_II_dimension==1 || deal_II_dimension==2
       template class PolicyBase<deal_II_dimension,deal_II_dimension+1>;
       template class Sequential<deal_II_dimension,deal_II_dimension+1>;
+      template class ParallelShared<deal_II_dimension,deal_II_dimension+1>;
       template class ParallelDistributed<deal_II_dimension,deal_II_dimension+1>;
 #endif
 
 #if deal_II_dimension==3
       template class PolicyBase<1,3>;
       template class Sequential<1,3>;
+      template class ParallelShared<1,3>;
       template class ParallelDistributed<1,3>;
 #endif
     \}
index 3ea549933447ed05ab9bf7388230f450e8851844..a2ac1ad4b62bafec386469d6db47e0a8dcc0b13c 100644 (file)
@@ -1058,7 +1058,44 @@ namespace DoFTools
       active_fe_indices[index] = cell->active_fe_index();
   }
 
-
+  template <class DH>
+  IndexSet
+  locally_owned_dofs_with_subdomain (const DH                  &dof_handler,
+                                            const types::subdomain_id  subdomain   )
+  {
+         IndexSet index_set (dof_handler.n_dofs());
+         std::vector< types::subdomain_id > subdomain_association(dof_handler.n_dofs());
+         DoFTools::get_subdomain_association(dof_handler,subdomain_association);
+         types::global_dof_index i_min=0, i_max_plus_1=0;
+
+         //loop over subdomain_association
+         //and populate IndexSet when
+         //a given subdomain is found
+         bool found_start = false;
+         for (types::global_dof_index ind = 0; ind < subdomain_association.size(); ind++)
+         {
+                 if (!found_start && (subdomain_association[ind]== subdomain) )
+                 {
+                         found_start = true;
+                         i_min = ind;
+                 }
+                 if (found_start && (subdomain_association[ind] != subdomain) )
+                 {
+                         found_start = false;
+                         i_max_plus_1 = ind;
+                         index_set.add_range(i_min,i_max_plus_1);
+                 }
+         }
+
+         // we have found the starting element, but not the end.
+         // the last element must have the same subdomain id,
+         // thus use it to add the range
+         if (found_start) {
+                 index_set.add_range(i_min,subdomain_association.size());
+         }
+         index_set.compress();
+         return index_set;
+  }
 
   template <class DH>
   void
index 9309ae6bc2e8a07d28da25dbf422782257cb83c4..4e16d7b6d878e8268d36785d8c0512c7fe81c6c0 100644 (file)
@@ -269,7 +269,17 @@ void
 DoFTools::get_subdomain_association<hp::DoFHandler<deal_II_dimension> >
 (const hp::DoFHandler<deal_II_dimension> &dof_handler,
  std::vector<types::subdomain_id>           &subdomain_association);
-
+template
+IndexSet
+DoFTools::locally_owned_dofs_with_subdomain<DoFHandler<deal_II_dimension> >
+(const DoFHandler<deal_II_dimension> &dof_handler,
+ const types::subdomain_id            subdomain);
+template
+IndexSet
+DoFTools::locally_owned_dofs_with_subdomain<hp::DoFHandler<deal_II_dimension> >
+(const hp::DoFHandler<deal_II_dimension> &dof_handler,
+ const types::subdomain_id                subdomain);
 
 template
 unsigned int
@@ -352,6 +362,17 @@ void
 DoFTools::get_subdomain_association<hp::DoFHandler<deal_II_dimension,deal_II_dimension+1> >
 (const hp::DoFHandler<deal_II_dimension,deal_II_dimension+1> &dof_handler,
  std::vector<types::subdomain_id>           &subdomain_association);
+template
+IndexSet
+DoFTools::locally_owned_dofs_with_subdomain<DoFHandler<deal_II_dimension,deal_II_dimension+1> >
+(const DoFHandler<deal_II_dimension,deal_II_dimension+1> &dof_handler,
+ const types::subdomain_id                                subdomain);
+template
+IndexSet
+DoFTools::locally_owned_dofs_with_subdomain<hp::DoFHandler<deal_II_dimension,deal_II_dimension+1> >
+(const hp::DoFHandler<deal_II_dimension,deal_II_dimension+1> &dof_handler,
+ const types::subdomain_id                                    subdomain);
 
 template
 void
@@ -380,6 +401,17 @@ void
 DoFTools::get_subdomain_association<hp::DoFHandler<1,3> >
 (const hp::DoFHandler<1,3> &dof_handler,
  std::vector<types::subdomain_id>           &subdomain_association);
+template
+IndexSet
+DoFTools::locally_owned_dofs_with_subdomain<DoFHandler<1,3> >
+(const DoFHandler<1,3>     &dof_handler,
+ const types::subdomain_id  subdomain);
+template
+IndexSet
+DoFTools::locally_owned_dofs_with_subdomain<hp::DoFHandler<1,3> >
+(const hp::DoFHandler<1,3> &dof_handler,
+ const types::subdomain_id  subdomain);
 
 template
 unsigned int
index 2f271224d7dd89320fa48683269137af903eccd7..d73508e06a4ab4f76f33ad543a0193f7755d79cc 100644 (file)
@@ -2690,10 +2690,17 @@ namespace hp
     number_cache.n_global_dofs        = next_free_dof;
     number_cache.n_locally_owned_dofs = number_cache.n_global_dofs;
 
-    number_cache.locally_owned_dofs
-      = IndexSet (number_cache.n_global_dofs);
-    number_cache.locally_owned_dofs.add_range (0,
-                                               number_cache.n_global_dofs);
+    if (dynamic_cast<const parallel::shared::Triangulation< dim, spacedim >*>
+        (&this->get_tria())
+        == 0) {
+      number_cache.locally_owned_dofs
+        = IndexSet (number_cache.n_global_dofs);
+      number_cache.locally_owned_dofs.add_range (0,
+                                                 number_cache.n_global_dofs);
+    } else {
+       Assert(false, ExcMessage("not implemented"));
+       //number_cache.locally_owned_dofs = dealii::DoFTools::locally_owned_dofs_with_subdomain(this,tria->locally_owned_subdomain() );
+    }
     Assert (number_cache.n_global_dofs < std::numeric_limits<unsigned int>::max (),
             ExcMessage ("Global number of degrees of freedom is too large."));
     number_cache.n_locally_owned_dofs_per_processor

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.